MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mfem::NURBSExtension Class Reference

NURBSExtension generally contains multiple NURBSPatch objects spanning an entire Mesh. It also defines and manages DOFs in NURBS finite element spaces. More...

#include <nurbs.hpp>

Inheritance diagram for mfem::NURBSExtension:
[legend]
Collaboration diagram for mfem::NURBSExtension:
[legend]

Public Member Functions

 NURBSExtension (const NURBSExtension &orig)
 Copy constructor: deep copy.
 
 NURBSExtension (std::istream &input, bool spacing=false)
 Read-in a NURBSExtension from a stream input..
 
 NURBSExtension (NURBSExtension *parent, int newOrder)
 Create a NURBSExtension with elevated order by repeating the endpoints of the KnotVectors and using uniform weights of 1.
 
 NURBSExtension (NURBSExtension *parent, const Array< int > &newOrders, Mode mode=Mode::H_1)
 Create a NURBSExtension with elevated KnotVector orders (by repeating the endpoints of the KnotVectors and using uniform weights of 1) as given by the array newOrders.
 
 NURBSExtension (Mesh *mesh_array[], int num_pieces)
 Construct a NURBSExtension by merging a partitioned NURBS mesh.
 
 NURBSExtension (const Mesh *patch_topology, const Array< const NURBSPatch * > p)
 
NURBSExtensionoperator= (const NURBSExtension &)=delete
 Copy assignment not supported.
 
void ConnectBoundaries (Array< int > &master, Array< int > &slave)
 Generate connections between boundaries, such as periodic BCs.
 
const Array< int > & GetMaster () const
 
Array< int > & GetMaster ()
 
const Array< int > & GetSlave () const
 
Array< int > & GetSlave ()
 
void MergeGridFunctions (GridFunction *gf_array[], int num_pieces, GridFunction &merged)
 Set the DOFs of merged to values from active elements in num_pieces of Gridfunctions gf_array.
 
virtual ~NURBSExtension ()
 Destroy a NURBSExtension.
 
void Print (std::ostream &os, const std::string &comments="") const
 Writes all patch data to the stream os.
 
void PrintCharacteristics (std::ostream &os) const
 Print various mesh characteristics to the stream os.
 
void PrintFunctions (const char *basename, int samples=11) const
 Call KnotVector::PrintFunctions for all KnotVectors, using a separate, newly created ofstream with filename "basename_i.dat" for KnotVector i.
 
int Dimension () const
 Return the dimension of the reference space (not physical space).
 
int GetNP () const
 Return the number of patches.
 
int GetNBP () const
 Return the number of boundary patches.
 
const Array< int > & GetOrders () const
 Read-only access to the orders of all KnotVectors.
 
int GetOrder () const
 If all KnotVector orders are identical, return that number. Otherwise, return NURBSFECollection::VariableOrder.
 
int GetNKV () const
 Return the number of KnotVectors.
 
int GetGNV () const
 Return the global number of vertices.
 
int GetNV () const
 Return the local number of active vertices.
 
int GetGNE () const
 Return the global number of elements.
 
int GetNE () const
 Return the number of active elements.
 
int GetGNBE () const
 Return the global number of boundary elements.
 
int GetNBE () const
 Return the number of active boundary elements.
 
int GetNTotalDof () const
 Return the total number of DOFs.
 
int GetNDof () const
 Return the number of active DOFs.
 
int GetActiveDof (int glob) const
 Return the local DOF number for a given global DOF number glob.
 
int DofMap (int dof) const
 Return the dof index whilst accounting for periodic boundaries.
 
void GetPatchKnotVectors (int p, Array< const KnotVector * > &kv) const
 Return KnotVectors in kv in each dimension for patch p.
 
void GetBdrPatchKnotVectors (int bp, Array< const KnotVector * > &kv) const
 Return KnotVectors in kv in each dimension for boundary patch bp.
 
const KnotVectorGetKnotVector (int i) const
 KnotVector read-only access function.
 
void GetElementTopo (Array< Element * > &elements) const
 Generate the active mesh elements and return them in elements.
 
void GetBdrElementTopo (Array< Element * > &boundary) const
 Generate the active mesh boundary elements and return them in boundary.
 
bool HavePatches () const
 Return true if at least 1 patch is defined, false otherwise.
 
TableGetElementDofTable ()
 
TableGetBdrElementDofTable ()
 
void GetVertexLocalToGlobal (Array< int > &lvert_vert)
 Get the local to global vertex index map lvert_vert.
 
void GetElementLocalToGlobal (Array< int > &lelem_elem)
 Get the local to global element index map lelem_elem.
 
void SetPatchAttribute (int i, int attr)
 Set the attribute for patch i, which is set to all elements in the patch.
 
int GetPatchAttribute (int i) const
 Get the attribute for patch i, which is set to all elements in the patch.
 
void SetPatchBdrAttribute (int i, int attr)
 Set the attribute for patch boundary element i to attr, which is set to all boundary elements in the patch.
 
int GetPatchBdrAttribute (int i) const
 Get the attribute for boundary patch element i, which is set to all boundary elements in the patch.
 
void LoadFE (int i, const FiniteElement *FE) const
 Load element i into FE.
 
void LoadBE (int i, const FiniteElement *BE) const
 Load boundary element i into BE.
 
const VectorGetWeights () const
 Access function to the vector of weights weights.
 
VectorGetWeights ()
 
void ConvertToPatches (const Vector &Nodes)
 Define patches in IKJ (B-net) format, using FE coordinates in Nodes.
 
void SetKnotsFromPatches ()
 Set KnotVectors from patches and construct mesh and space data.
 
void SetCoordsFromPatches (Vector &Nodes)
 Set FE coordinates in Nodes, using data from patches, and erase patches.
 
void LoadSolution (std::istream &input, GridFunction &sol) const
 Read a GridFunction sol from stream input, written patch-by-patch, e.g. with PrintSolution().
 
void PrintSolution (const GridFunction &sol, std::ostream &os) const
 Write a GridFunction sol patch-by-patch to stream os.
 
void DegreeElevate (int rel_degree, int degree=16)
 Call DegreeElevate for all KnotVectors of all patches. For each KnotVector, the new degree is max(old_degree, min(old_degree + rel_degree, degree)).
 
void UniformRefinement (int rf=2)
 Refine with optional refinement factor rf. Uniform means refinement is done everywhere by the same factor, although nonuniform spacing functions may be used.
 
void UniformRefinement (Array< int > const &rf)
 
void Coarsen (int cf=2, real_t tol=1.0e-12)
 
void Coarsen (Array< int > const &cf, real_t tol=1.0e-12)
 
void KnotInsert (Array< KnotVector * > &kv)
 Insert knots from kv into all KnotVectors in all patches. The size of kv should be the same as knotVectors.
 
void KnotInsert (Array< Vector * > &kv)
 
NURBSExtensionGetDivExtension (int component)
 
NURBSExtensionGetCurlExtension (int component)
 
void KnotRemove (Array< Vector * > &kv, real_t tol=1.0e-12)
 
void GetCoarseningFactors (Array< int > &f) const
 
int GetElementPatch (int elem) const
 Returns the index of the patch containing element elem.
 
void GetElementIJK (int elem, Array< int > &ijk)
 Return Cartesian indices (i,j) in 2D or (i,j,k) in 3D of element elem, in the knot-span tensor product ordering for its patch.
 
void GetPatchDofs (const int patch, Array< int > &dofs)
 Return the degrees of freedom in dofs on patch patch, in Cartesian order.
 
const Array< int > & GetPatchElements (int patch)
 Return the array of indices of all elements in patch patch.
 
const Array< int > & GetPatchBdrElements (int patch)
 Return the array of indices of all boundary elements in patch patch.
 

Protected Types

enum class  Mode { H_1 , H_DIV , H_CURL }
 Flag for indicating what type of NURBS fespace this extension is used for. More...
 

Protected Member Functions

int KnotInd (int edge) const
 Return the unsigned index of the KnotVector for edge edge.
 
KnotVectorKnotVec (int edge)
 
const KnotVectorKnotVec (int edge) const
 
const KnotVectorKnotVec (int edge, int oedge, int *okv) const
 
void CheckPatches ()
 Throw an error if any patch has an inconsistent edge-to-knot mapping.
 
void CheckBdrPatches ()
 Throw an error if any boundary patch has invalid KnotVector orientation.
 
void CheckKVDirection (int p, Array< int > &kvdir)
 Return the directions in kvdir of the KnotVectors in patch p based on the patch edge orientations. Each entry of kvdir is -1 if the KnotVector direction is flipped, +1 otherwise.
 
void CreateComprehensiveKV ()
 Create the comprehensive set of KnotVectors. In 1D, this set is identical to the unique set of KnotVectors.
 
void UpdateUniqueKV ()
 
bool ConsistentKVSets ()
 Check if the comprehensive array of KnotVectors agrees with the unique set of KnotVectors, on each patch. Return false if there is a difference, true otherwise. This function throws an error in 1D.
 
void GetPatchKnotVectors (int p, Array< KnotVector * > &kv)
 Return KnotVectors in kv in each dimension for patch p.
 
void GetBdrPatchKnotVectors (int bp, Array< KnotVector * > &kv)
 Return KnotVectors in kv in each dimension for boundary patch bp.
 
void SetOrderFromOrders ()
 Set overall order mOrder based on KnotVector orders.
 
void SetOrdersFromKnotVectors ()
 Set orders from KnotVector orders.
 
void InitDofMap ()
 Set DOF map sizes to 0.
 
void ConnectBoundaries ()
 Set DOF maps for periodic BC.
 
void ConnectBoundaries1D (int bnd0, int bnd1)
 
void ConnectBoundaries2D (int bnd0, int bnd1)
 
void ConnectBoundaries3D (int bnd0, int bnd1)
 
void GenerateOffsets ()
 Set the mesh and space offsets, and also count the global NumOfVertices and the global NumOfDofs.
 
void CountElements ()
 Count the global NumOfElements.
 
void CountBdrElements ()
 Count the global NumOfBdrElements.
 
void Get1DElementTopo (Array< Element * > &elements) const
 Generate the active mesh elements and return them in elements.
 
void Get2DElementTopo (Array< Element * > &elements) const
 
void Get3DElementTopo (Array< Element * > &elements) const
 
void Get1DBdrElementTopo (Array< Element * > &boundary) const
 Generate the active mesh boundary elements and return them in boundary.
 
void Get2DBdrElementTopo (Array< Element * > &boundary) const
 
void Get3DBdrElementTopo (Array< Element * > &boundary) const
 
void GenerateElementDofTable ()
 Based on activeElem, count NumOfActiveDofs and generate el_dof, el_to_patch, el_to_IJK, activeDof map (global-to-local).
 
void Generate1DElementDofTable ()
 Generate elem_to_global-dof table for the active elements, and define el_to_patch, el_to_IJK, activeDof (as bool).
 
void Generate2DElementDofTable ()
 
void Generate3DElementDofTable ()
 
void GenerateBdrElementDofTable ()
 Call after GenerateElementDofTable to set boundary element DOF table.
 
void Generate1DBdrElementDofTable ()
 Generate the table of global DOFs for active boundary elements, and define bel_to_patch, bel_to_IJK.
 
void Generate2DBdrElementDofTable ()
 
void Generate3DBdrElementDofTable ()
 
void GetPatchNets (const Vector &coords, int vdim)
 Set the B-NET on each patch using values from coords.
 
void Get1DPatchNets (const Vector &coords, int vdim)
 
void Get2DPatchNets (const Vector &coords, int vdim)
 
void Get3DPatchNets (const Vector &coords, int vdim)
 
void SetSolutionVector (Vector &coords, int vdim)
 Return in coords the coordinates from each patch. Side effects: delete the patches and update the weights from the patches.
 
void Set1DSolutionVector (Vector &coords, int vdim)
 
void Set2DSolutionVector (Vector &coords, int vdim)
 
void Set3DSolutionVector (Vector &coords, int vdim)
 
void GenerateActiveVertices ()
 Determine activeVert, NumOfActiveVertices from the activeElem array.
 
void GenerateActiveBdrElems ()
 Determine activeBdrElem, NumOfActiveBdrElems.
 
void MergeWeights (Mesh *mesh_array[], int num_pieces)
 Set the weights in this object to values from active elements in num_pieces meshes in mesh_array.
 
void SetPatchToElements ()
 Set patch_to_el.
 
void SetPatchToBdrElements ()
 Set patch_to_bel.
 
 NURBSExtension ()
 To be used by ParNURBSExtension constructor(s)
 

Protected Attributes

Mode mode = Mode::H_1
 
int mOrder
 Order of KnotVectors, see GetOrder() for description.
 
Array< int > mOrders
 Orders of all KnotVectors.
 
int NumOfKnotVectors
 Number of KnotVectors.
 
int NumOfVertices
 Global entity counts.
 
int NumOfElements
 
int NumOfBdrElements
 
int NumOfDofs
 
int NumOfActiveVertices
 Local entity counts.
 
int NumOfActiveElems
 
int NumOfActiveBdrElems
 
int NumOfActiveDofs
 
Array< int > activeVert
 
Array< bool > activeElem
 
Array< bool > activeBdrElem
 
Array< int > activeDof
 
MeshpatchTopo
 Patch topology mesh.
 
bool own_topo
 Whether this object owns patchTopo.
 
Array< int > edge_to_knot
 Map from edge indices to KnotVector indices.
 
Array< KnotVector * > knotVectors
 Set of unique KnotVectors.
 
Array< KnotVector * > knotVectorsCompr
 Comprehensive set of all KnotVectors, one for every edge.
 
Vector weights
 Weights for each control point or DOF.
 
Array< int > d_to_d
 Periodic BC info:
 
Array< int > master
 
Array< int > slave
 
Array< int > v_meshOffsets
 Global mesh offsets, meshOffsets == meshVertexOffsets.
 
Array< int > e_meshOffsets
 
Array< int > f_meshOffsets
 
Array< int > p_meshOffsets
 
Array< int > v_spaceOffsets
 Global space offsets, spaceOffsets == dofOffsets.
 
Array< int > e_spaceOffsets
 
Array< int > f_spaceOffsets
 
Array< int > p_spaceOffsets
 
Tableel_dof
 Table of DOFs for each element (el_dof) or boundary element (bel_dof).
 
Tablebel_dof
 
Array< int > el_to_patch
 Map from element indices to patch indices.
 
Array< int > bel_to_patch
 Map from boundary element indices to patch indices.
 
Array2D< int > el_to_IJK
 Map from element indices to IJK knot span indices.
 
Array2D< int > bel_to_IJK
 
std::vector< Array< int > > patch_to_el
 For each patch p, patch_to_el[p] lists all elements in the patch.
 
std::vector< Array< int > > patch_to_bel
 For each patch p, patch_to_bel[p] lists all boundary elements in the patch.
 
Array< NURBSPatch * > patches
 Array of all patches in the mesh.
 

Friends

class ParNURBSExtension
 
class NURBSPatchMap
 

Detailed Description

NURBSExtension generally contains multiple NURBSPatch objects spanning an entire Mesh. It also defines and manages DOFs in NURBS finite element spaces.

Definition at line 448 of file nurbs.hpp.

Member Enumeration Documentation

◆ Mode

enum class mfem::NURBSExtension::Mode
strongprotected

Flag for indicating what type of NURBS fespace this extension is used for.

Enumerator
H_1 
H_DIV 

‍Extension for a standard scalar-valued space

H_CURL 

‍Extension for a divergence conforming vector-valued space

Definition at line 458 of file nurbs.hpp.

Constructor & Destructor Documentation

◆ NURBSExtension() [1/7]

mfem::NURBSExtension::NURBSExtension ( )
inlineprotected

To be used by ParNURBSExtension constructor(s)

Definition at line 677 of file nurbs.hpp.

◆ NURBSExtension() [2/7]

mfem::NURBSExtension::NURBSExtension ( const NURBSExtension & orig)

Copy constructor: deep copy.

Definition at line 1970 of file nurbs.cpp.

◆ NURBSExtension() [3/7]

mfem::NURBSExtension::NURBSExtension ( std::istream & input,
bool spacing = false )

Read-in a NURBSExtension from a stream input..

Definition at line 2024 of file nurbs.cpp.

◆ NURBSExtension() [4/7]

mfem::NURBSExtension::NURBSExtension ( NURBSExtension * parent,
int newOrder )

Create a NURBSExtension with elevated order by repeating the endpoints of the KnotVectors and using uniform weights of 1.

Note
If a KnotVector in parent already has order greater than or equal to newOrder, it will be used unmodified.

Definition at line 2225 of file nurbs.cpp.

◆ NURBSExtension() [5/7]

mfem::NURBSExtension::NURBSExtension ( NURBSExtension * parent,
const Array< int > & newOrders,
Mode mode = Mode::H_1 )

Create a NURBSExtension with elevated KnotVector orders (by repeating the endpoints of the KnotVectors and using uniform weights of 1) as given by the array newOrders.

note If a KnotVector in parent already has order greater than or equal to the corresponding entry in newOrder, it will be used unmodified.

Definition at line 2278 of file nurbs.cpp.

◆ NURBSExtension() [6/7]

mfem::NURBSExtension::NURBSExtension ( Mesh * mesh_array[],
int num_pieces )

Construct a NURBSExtension by merging a partitioned NURBS mesh.

Definition at line 2334 of file nurbs.cpp.

◆ NURBSExtension() [7/7]

mfem::NURBSExtension::NURBSExtension ( const Mesh * patch_topology,
const Array< const NURBSPatch * > p )

Definition at line 2379 of file nurbs.cpp.

◆ ~NURBSExtension()

mfem::NURBSExtension::~NURBSExtension ( )
virtual

Destroy a NURBSExtension.

Definition at line 2446 of file nurbs.cpp.

Member Function Documentation

◆ CheckBdrPatches()

void mfem::NURBSExtension::CheckBdrPatches ( )
protected

Throw an error if any boundary patch has invalid KnotVector orientation.

Definition at line 2963 of file nurbs.cpp.

◆ CheckKVDirection()

void mfem::NURBSExtension::CheckKVDirection ( int p,
Array< int > & kvdir )
protected

Return the directions in kvdir of the KnotVectors in patch p based on the patch edge orientations. Each entry of kvdir is -1 if the KnotVector direction is flipped, +1 otherwise.

Definition at line 2991 of file nurbs.cpp.

◆ CheckPatches()

void mfem::NURBSExtension::CheckPatches ( )
protected

Throw an error if any patch has an inconsistent edge-to-knot mapping.

Definition at line 2926 of file nurbs.cpp.

◆ Coarsen() [1/2]

void mfem::NURBSExtension::Coarsen ( Array< int > const & cf,
real_t tol = 1.0e-12 )

Definition at line 4454 of file nurbs.cpp.

◆ Coarsen() [2/2]

void mfem::NURBSExtension::Coarsen ( int cf = 2,
real_t tol = 1.0e-12 )

Definition at line 4469 of file nurbs.cpp.

◆ ConnectBoundaries() [1/2]

void mfem::NURBSExtension::ConnectBoundaries ( )
protected

Set DOF maps for periodic BC.

Definition at line 2584 of file nurbs.cpp.

◆ ConnectBoundaries() [2/2]

void mfem::NURBSExtension::ConnectBoundaries ( Array< int > & master,
Array< int > & slave )

Generate connections between boundaries, such as periodic BCs.

Definition at line 2577 of file nurbs.cpp.

◆ ConnectBoundaries1D()

void mfem::NURBSExtension::ConnectBoundaries1D ( int bnd0,
int bnd1 )
protected

Definition at line 2650 of file nurbs.cpp.

◆ ConnectBoundaries2D()

void mfem::NURBSExtension::ConnectBoundaries2D ( int bnd0,
int bnd1 )
protected

Definition at line 2664 of file nurbs.cpp.

◆ ConnectBoundaries3D()

void mfem::NURBSExtension::ConnectBoundaries3D ( int bnd0,
int bnd1 )
protected

Definition at line 2710 of file nurbs.cpp.

◆ ConsistentKVSets()

bool mfem::NURBSExtension::ConsistentKVSets ( )
protected

Check if the comprehensive array of KnotVectors agrees with the unique set of KnotVectors, on each patch. Return false if there is a difference, true otherwise. This function throws an error in 1D.

Definition at line 3186 of file nurbs.cpp.

◆ ConvertToPatches()

void mfem::NURBSExtension::ConvertToPatches ( const Vector & Nodes)

Define patches in IKJ (B-net) format, using FE coordinates in Nodes.

Definition at line 4256 of file nurbs.cpp.

◆ CountBdrElements()

void mfem::NURBSExtension::CountBdrElements ( )
protected

Count the global NumOfBdrElements.

Definition at line 3477 of file nurbs.cpp.

◆ CountElements()

void mfem::NURBSExtension::CountElements ( )
protected

Count the global NumOfElements.

Definition at line 3457 of file nurbs.cpp.

◆ CreateComprehensiveKV()

void mfem::NURBSExtension::CreateComprehensiveKV ( )
protected

Create the comprehensive set of KnotVectors. In 1D, this set is identical to the unique set of KnotVectors.

Definition at line 3058 of file nurbs.cpp.

◆ DegreeElevate()

void mfem::NURBSExtension::DegreeElevate ( int rel_degree,
int degree = 16 )

Call DegreeElevate for all KnotVectors of all patches. For each KnotVector, the new degree is max(old_degree, min(old_degree + rel_degree, degree)).

Definition at line 4391 of file nurbs.cpp.

◆ Dimension()

int mfem::NURBSExtension::Dimension ( ) const
inline

Return the dimension of the reference space (not physical space).

Definition at line 742 of file nurbs.hpp.

◆ DofMap()

int mfem::NURBSExtension::DofMap ( int dof) const
inline

Return the dof index whilst accounting for periodic boundaries.

Definition at line 782 of file nurbs.hpp.

◆ Generate1DBdrElementDofTable()

void mfem::NURBSExtension::Generate1DBdrElementDofTable ( )
protected

Generate the table of global DOFs for active boundary elements, and define bel_to_patch, bel_to_IJK.

Definition at line 4018 of file nurbs.cpp.

◆ Generate1DElementDofTable()

void mfem::NURBSExtension::Generate1DElementDofTable ( )
protected

Generate elem_to_global-dof table for the active elements, and define el_to_patch, el_to_IJK, activeDof (as bool).

Definition at line 3772 of file nurbs.cpp.

◆ Generate2DBdrElementDofTable()

void mfem::NURBSExtension::Generate2DBdrElementDofTable ( )
protected

Definition at line 4048 of file nurbs.cpp.

◆ Generate2DElementDofTable()

void mfem::NURBSExtension::Generate2DElementDofTable ( )
protected

Definition at line 3815 of file nurbs.cpp.

◆ Generate3DBdrElementDofTable()

void mfem::NURBSExtension::Generate3DBdrElementDofTable ( )
protected

Definition at line 4111 of file nurbs.cpp.

◆ Generate3DElementDofTable()

void mfem::NURBSExtension::Generate3DElementDofTable ( )
protected

Definition at line 3869 of file nurbs.cpp.

◆ GenerateActiveBdrElems()

void mfem::NURBSExtension::GenerateActiveBdrElems ( )
protected

Determine activeBdrElem, NumOfActiveBdrElems.

Definition at line 2854 of file nurbs.cpp.

◆ GenerateActiveVertices()

void mfem::NURBSExtension::GenerateActiveVertices ( )
protected

Determine activeVert, NumOfActiveVertices from the activeElem array.

Definition at line 2781 of file nurbs.cpp.

◆ GenerateBdrElementDofTable()

void mfem::NURBSExtension::GenerateBdrElementDofTable ( )
protected

Call after GenerateElementDofTable to set boundary element DOF table.

Definition at line 3984 of file nurbs.cpp.

◆ GenerateElementDofTable()

void mfem::NURBSExtension::GenerateElementDofTable ( )
protected

Based on activeElem, count NumOfActiveDofs and generate el_dof, el_to_patch, el_to_IJK, activeDof map (global-to-local).

Definition at line 3736 of file nurbs.cpp.

◆ GenerateOffsets()

void mfem::NURBSExtension::GenerateOffsets ( )
protected

Set the mesh and space offsets, and also count the global NumOfVertices and the global NumOfDofs.

Definition at line 3365 of file nurbs.cpp.

◆ Get1DBdrElementTopo()

void mfem::NURBSExtension::Get1DBdrElementTopo ( Array< Element * > & boundary) const
protected

Generate the active mesh boundary elements and return them in boundary.

Definition at line 3644 of file nurbs.cpp.

◆ Get1DElementTopo()

void mfem::NURBSExtension::Get1DElementTopo ( Array< Element * > & elements) const
protected

Generate the active mesh elements and return them in elements.

Definition at line 3515 of file nurbs.cpp.

◆ Get1DPatchNets()

void mfem::NURBSExtension::Get1DPatchNets ( const Vector & coords,
int vdim )
protected

Definition at line 4691 of file nurbs.cpp.

◆ Get2DBdrElementTopo()

void mfem::NURBSExtension::Get2DBdrElementTopo ( Array< Element * > & boundary) const
protected

Definition at line 3667 of file nurbs.cpp.

◆ Get2DElementTopo()

void mfem::NURBSExtension::Get2DElementTopo ( Array< Element * > & elements) const
protected

Definition at line 3545 of file nurbs.cpp.

◆ Get2DPatchNets()

void mfem::NURBSExtension::Get2DPatchNets ( const Vector & coords,
int vdim )
protected

Definition at line 4716 of file nurbs.cpp.

◆ Get3DBdrElementTopo()

void mfem::NURBSExtension::Get3DBdrElementTopo ( Array< Element * > & boundary) const
protected

Definition at line 3698 of file nurbs.cpp.

◆ Get3DElementTopo()

void mfem::NURBSExtension::Get3DElementTopo ( Array< Element * > & elements) const
protected

Definition at line 3581 of file nurbs.cpp.

◆ Get3DPatchNets()

void mfem::NURBSExtension::Get3DPatchNets ( const Vector & coords,
int vdim )
protected

Definition at line 4743 of file nurbs.cpp.

◆ GetActiveDof()

int mfem::NURBSExtension::GetActiveDof ( int glob) const
inline

Return the local DOF number for a given global DOF number glob.

Definition at line 779 of file nurbs.hpp.

◆ GetBdrElementDofTable()

Table * mfem::NURBSExtension::GetBdrElementDofTable ( )
inline

Access function for the boundary element DOF table bel_dof.

Note
The returned object should NOT be deleted by the caller.

Definition at line 812 of file nurbs.hpp.

◆ GetBdrElementTopo()

void mfem::NURBSExtension::GetBdrElementTopo ( Array< Element * > & boundary) const

Generate the active mesh boundary elements and return them in boundary.

Definition at line 3626 of file nurbs.cpp.

◆ GetBdrPatchKnotVectors() [1/2]

void mfem::NURBSExtension::GetBdrPatchKnotVectors ( int bp,
Array< const KnotVector * > & kv ) const

Return KnotVectors in kv in each dimension for boundary patch bp.

Definition at line 3320 of file nurbs.cpp.

◆ GetBdrPatchKnotVectors() [2/2]

void mfem::NURBSExtension::GetBdrPatchKnotVectors ( int bp,
Array< KnotVector * > & kv )
protected

Return KnotVectors in kv in each dimension for boundary patch bp.

Definition at line 3300 of file nurbs.cpp.

◆ GetCoarseningFactors()

void mfem::NURBSExtension::GetCoarseningFactors ( Array< int > & f) const

Calls GetCoarseningFactors for each patch and finds the minimum factor for each direction that ensures refinement will work in the case of non-nested spacing functions.

Definition at line 4476 of file nurbs.cpp.

◆ GetCurlExtension()

NURBSExtension * mfem::NURBSExtension::GetCurlExtension ( int component)

Returns the NURBSExtension to be used for component of an H(curl) conforming NURBS space. Caller gets ownership of the returned object, and is responsible for deletion.

Definition at line 4422 of file nurbs.cpp.

◆ GetDivExtension()

NURBSExtension * mfem::NURBSExtension::GetDivExtension ( int component)

Returns the NURBSExtension to be used for component of an H(div) conforming NURBS space. Caller gets ownership of the returned object, and is responsible for deletion.

Definition at line 4407 of file nurbs.cpp.

◆ GetElementDofTable()

Table * mfem::NURBSExtension::GetElementDofTable ( )
inline

Access function for the element DOF table el_dof.

Note
The returned object should NOT be deleted by the caller.

Definition at line 808 of file nurbs.hpp.

◆ GetElementIJK()

void mfem::NURBSExtension::GetElementIJK ( int elem,
Array< int > & ijk )

Return Cartesian indices (i,j) in 2D or (i,j,k) in 3D of element elem, in the knot-span tensor product ordering for its patch.

Definition at line 4875 of file nurbs.cpp.

◆ GetElementLocalToGlobal()

void mfem::NURBSExtension::GetElementLocalToGlobal ( Array< int > & lelem_elem)

Get the local to global element index map lelem_elem.

Definition at line 4202 of file nurbs.cpp.

◆ GetElementPatch()

int mfem::NURBSExtension::GetElementPatch ( int elem) const
inline

Returns the index of the patch containing element elem.

Definition at line 904 of file nurbs.hpp.

◆ GetElementTopo()

void mfem::NURBSExtension::GetElementTopo ( Array< Element * > & elements) const

Generate the active mesh elements and return them in elements.

Definition at line 3497 of file nurbs.cpp.

◆ GetGNBE()

int mfem::NURBSExtension::GetGNBE ( ) const
inline

Return the global number of boundary elements.

Definition at line 769 of file nurbs.hpp.

◆ GetGNE()

int mfem::NURBSExtension::GetGNE ( ) const
inline

Return the global number of elements.

Definition at line 765 of file nurbs.hpp.

◆ GetGNV()

int mfem::NURBSExtension::GetGNV ( ) const
inline

Return the global number of vertices.

Definition at line 761 of file nurbs.hpp.

◆ GetKnotVector()

const KnotVector * mfem::NURBSExtension::GetKnotVector ( int i) const
inline

KnotVector read-only access function.

Definition at line 794 of file nurbs.hpp.

◆ GetMaster() [1/2]

Array< int > & mfem::NURBSExtension::GetMaster ( )
inline

Definition at line 709 of file nurbs.hpp.

◆ GetMaster() [2/2]

const Array< int > & mfem::NURBSExtension::GetMaster ( ) const
inline

Definition at line 708 of file nurbs.hpp.

◆ GetNBE()

int mfem::NURBSExtension::GetNBE ( ) const
inline

Return the number of active boundary elements.

Definition at line 771 of file nurbs.hpp.

◆ GetNBP()

int mfem::NURBSExtension::GetNBP ( ) const
inline

Return the number of boundary patches.

Definition at line 748 of file nurbs.hpp.

◆ GetNDof()

int mfem::NURBSExtension::GetNDof ( ) const
inline

Return the number of active DOFs.

Definition at line 776 of file nurbs.hpp.

◆ GetNE()

int mfem::NURBSExtension::GetNE ( ) const
inline

Return the number of active elements.

Definition at line 767 of file nurbs.hpp.

◆ GetNKV()

int mfem::NURBSExtension::GetNKV ( ) const
inline

Return the number of KnotVectors.

Definition at line 758 of file nurbs.hpp.

◆ GetNP()

int mfem::NURBSExtension::GetNP ( ) const
inline

Return the number of patches.

Definition at line 745 of file nurbs.hpp.

◆ GetNTotalDof()

int mfem::NURBSExtension::GetNTotalDof ( ) const
inline

Return the total number of DOFs.

Definition at line 774 of file nurbs.hpp.

◆ GetNV()

int mfem::NURBSExtension::GetNV ( ) const
inline

Return the local number of active vertices.

Definition at line 763 of file nurbs.hpp.

◆ GetOrder()

int mfem::NURBSExtension::GetOrder ( ) const
inline

If all KnotVector orders are identical, return that number. Otherwise, return NURBSFECollection::VariableOrder.

Definition at line 755 of file nurbs.hpp.

◆ GetOrders()

const Array< int > & mfem::NURBSExtension::GetOrders ( ) const
inline

Read-only access to the orders of all KnotVectors.

Definition at line 751 of file nurbs.hpp.

◆ GetPatchAttribute()

int mfem::NURBSExtension::GetPatchAttribute ( int i) const
inline

Get the attribute for patch i, which is set to all elements in the patch.

Definition at line 825 of file nurbs.hpp.

◆ GetPatchBdrAttribute()

int mfem::NURBSExtension::GetPatchBdrAttribute ( int i) const
inline

Get the attribute for boundary patch element i, which is set to all boundary elements in the patch.

Definition at line 834 of file nurbs.hpp.

◆ GetPatchBdrElements()

const Array< int > & mfem::NURBSExtension::GetPatchBdrElements ( int patch)

Return the array of indices of all boundary elements in patch patch.

Definition at line 4910 of file nurbs.cpp.

◆ GetPatchDofs()

void mfem::NURBSExtension::GetPatchDofs ( const int patch,
Array< int > & dofs )

Return the degrees of freedom in dofs on patch patch, in Cartesian order.

Definition at line 3935 of file nurbs.cpp.

◆ GetPatchElements()

const Array< int > & mfem::NURBSExtension::GetPatchElements ( int patch)

Return the array of indices of all elements in patch patch.

Definition at line 4903 of file nurbs.cpp.

◆ GetPatchKnotVectors() [1/2]

void mfem::NURBSExtension::GetPatchKnotVectors ( int p,
Array< const KnotVector * > & kv ) const

Return KnotVectors in kv in each dimension for patch p.

Definition at line 3276 of file nurbs.cpp.

◆ GetPatchKnotVectors() [2/2]

void mfem::NURBSExtension::GetPatchKnotVectors ( int p,
Array< KnotVector * > & kv )
protected

Return KnotVectors in kv in each dimension for patch p.

Definition at line 3253 of file nurbs.cpp.

◆ GetPatchNets()

void mfem::NURBSExtension::GetPatchNets ( const Vector & coords,
int vdim )
protected

Set the B-NET on each patch using values from coords.

Definition at line 4675 of file nurbs.cpp.

◆ GetSlave() [1/2]

Array< int > & mfem::NURBSExtension::GetSlave ( )
inline

Definition at line 711 of file nurbs.hpp.

◆ GetSlave() [2/2]

const Array< int > & mfem::NURBSExtension::GetSlave ( ) const
inline

Definition at line 710 of file nurbs.hpp.

◆ GetVertexLocalToGlobal()

void mfem::NURBSExtension::GetVertexLocalToGlobal ( Array< int > & lvert_vert)

Get the local to global vertex index map lvert_vert.

Definition at line 4192 of file nurbs.cpp.

◆ GetWeights() [1/2]

Vector & mfem::NURBSExtension::GetWeights ( )
inline

Definition at line 846 of file nurbs.hpp.

◆ GetWeights() [2/2]

const Vector & mfem::NURBSExtension::GetWeights ( ) const
inline

Access function to the vector of weights weights.

Definition at line 845 of file nurbs.hpp.

◆ HavePatches()

bool mfem::NURBSExtension::HavePatches ( ) const
inline

Return true if at least 1 patch is defined, false otherwise.

Definition at line 804 of file nurbs.hpp.

◆ InitDofMap()

void mfem::NURBSExtension::InitDofMap ( )
protected

Set DOF map sizes to 0.

Definition at line 2570 of file nurbs.cpp.

◆ KnotInd()

int mfem::NURBSExtension::KnotInd ( int edge) const
inlineprotected

Return the unsigned index of the KnotVector for edge edge.

Definition at line 1139 of file nurbs.hpp.

◆ KnotInsert() [1/2]

void mfem::NURBSExtension::KnotInsert ( Array< KnotVector * > & kv)

Insert knots from kv into all KnotVectors in all patches. The size of kv should be the same as knotVectors.

Definition at line 4503 of file nurbs.cpp.

◆ KnotInsert() [2/2]

void mfem::NURBSExtension::KnotInsert ( Array< Vector * > & kv)

Definition at line 4553 of file nurbs.cpp.

◆ KnotRemove()

void mfem::NURBSExtension::KnotRemove ( Array< Vector * > & kv,
real_t tol = 1.0e-12 )

Definition at line 4615 of file nurbs.cpp.

◆ KnotVec() [1/3]

KnotVector * mfem::NURBSExtension::KnotVec ( int edge)
inlineprotected

Access function for the KnotVector associated with edge edge.

Note
The returned object should NOT be deleted by the caller.

Definition at line 1145 of file nurbs.hpp.

◆ KnotVec() [2/3]

const KnotVector * mfem::NURBSExtension::KnotVec ( int edge) const
inlineprotected

Const access function for the KnotVector associated with edge edge.

Note
The returned object should NOT be deleted by the caller.

Definition at line 1150 of file nurbs.hpp.

◆ KnotVec() [3/3]

const KnotVector * mfem::NURBSExtension::KnotVec ( int edge,
int oedge,
int * okv ) const
inlineprotected

Definition at line 1155 of file nurbs.hpp.

◆ LoadBE()

void mfem::NURBSExtension::LoadBE ( int i,
const FiniteElement * BE ) const

Load boundary element i into BE.

Definition at line 4233 of file nurbs.cpp.

◆ LoadFE()

void mfem::NURBSExtension::LoadFE ( int i,
const FiniteElement * FE ) const

Load element i into FE.

Definition at line 4212 of file nurbs.cpp.

◆ LoadSolution()

void mfem::NURBSExtension::LoadSolution ( std::istream & input,
GridFunction & sol ) const

Read a GridFunction sol from stream input, written patch-by-patch, e.g. with PrintSolution().

Definition at line 4316 of file nurbs.cpp.

◆ MergeGridFunctions()

void mfem::NURBSExtension::MergeGridFunctions ( GridFunction * gf_array[],
int num_pieces,
GridFunction & merged )

Set the DOFs of merged to values from active elements in num_pieces of Gridfunctions gf_array.

Definition at line 2901 of file nurbs.cpp.

◆ MergeWeights()

void mfem::NURBSExtension::MergeWeights ( Mesh * mesh_array[],
int num_pieces )
protected

Set the weights in this object to values from active elements in num_pieces meshes in mesh_array.

Definition at line 2876 of file nurbs.cpp.

◆ operator=()

NURBSExtension & mfem::NURBSExtension::operator= ( const NURBSExtension & )
delete

Copy assignment not supported.

◆ Print()

void mfem::NURBSExtension::Print ( std::ostream & os,
const std::string & comments = "" ) const

Writes all patch data to the stream os.

The optional input argument comments is a string of comments to be printed after the first line (containing version number) of a mesh file. The output is formatted for writing a mesh to file. This function is called by Mesh::Printer.

Definition at line 2472 of file nurbs.cpp.

◆ PrintCharacteristics()

void mfem::NURBSExtension::PrintCharacteristics ( std::ostream & os) const

Print various mesh characteristics to the stream os.

Definition at line 2527 of file nurbs.cpp.

◆ PrintFunctions()

void mfem::NURBSExtension::PrintFunctions ( const char * basename,
int samples = 11 ) const

Call KnotVector::PrintFunctions for all KnotVectors, using a separate, newly created ofstream with filename "basename_i.dat" for KnotVector i.

Definition at line 2557 of file nurbs.cpp.

◆ PrintSolution()

void mfem::NURBSExtension::PrintSolution ( const GridFunction & sol,
std::ostream & os ) const

Write a GridFunction sol patch-by-patch to stream os.

Definition at line 4353 of file nurbs.cpp.

◆ Set1DSolutionVector()

void mfem::NURBSExtension::Set1DSolutionVector ( Vector & coords,
int vdim )
protected

Definition at line 4789 of file nurbs.cpp.

◆ Set2DSolutionVector()

void mfem::NURBSExtension::Set2DSolutionVector ( Vector & coords,
int vdim )
protected

Definition at line 4816 of file nurbs.cpp.

◆ Set3DSolutionVector()

void mfem::NURBSExtension::Set3DSolutionVector ( Vector & coords,
int vdim )
protected

Definition at line 4844 of file nurbs.cpp.

◆ SetCoordsFromPatches()

void mfem::NURBSExtension::SetCoordsFromPatches ( Vector & Nodes)

Set FE coordinates in Nodes, using data from patches, and erase patches.

Definition at line 4267 of file nurbs.cpp.

◆ SetKnotsFromPatches()

void mfem::NURBSExtension::SetKnotsFromPatches ( )

Set KnotVectors from patches and construct mesh and space data.

Definition at line 4275 of file nurbs.cpp.

◆ SetOrderFromOrders()

void mfem::NURBSExtension::SetOrderFromOrders ( )
protected

Set overall order mOrder based on KnotVector orders.

Definition at line 3341 of file nurbs.cpp.

◆ SetOrdersFromKnotVectors()

void mfem::NURBSExtension::SetOrdersFromKnotVectors ( )
protected

Set orders from KnotVector orders.

Definition at line 3355 of file nurbs.cpp.

◆ SetPatchAttribute()

void mfem::NURBSExtension::SetPatchAttribute ( int i,
int attr )
inline

Set the attribute for patch i, which is set to all elements in the patch.

Definition at line 821 of file nurbs.hpp.

◆ SetPatchBdrAttribute()

void mfem::NURBSExtension::SetPatchBdrAttribute ( int i,
int attr )
inline

Set the attribute for patch boundary element i to attr, which is set to all boundary elements in the patch.

Definition at line 829 of file nurbs.hpp.

◆ SetPatchToBdrElements()

void mfem::NURBSExtension::SetPatchToBdrElements ( )
protected

Set patch_to_bel.

Definition at line 4892 of file nurbs.cpp.

◆ SetPatchToElements()

void mfem::NURBSExtension::SetPatchToElements ( )
protected

Set patch_to_el.

Definition at line 4881 of file nurbs.cpp.

◆ SetSolutionVector()

void mfem::NURBSExtension::SetSolutionVector ( Vector & coords,
int vdim )
protected

Return in coords the coordinates from each patch. Side effects: delete the patches and update the weights from the patches.

Definition at line 4773 of file nurbs.cpp.

◆ UniformRefinement() [1/2]

void mfem::NURBSExtension::UniformRefinement ( Array< int > const & rf)

Definition at line 4439 of file nurbs.cpp.

◆ UniformRefinement() [2/2]

void mfem::NURBSExtension::UniformRefinement ( int rf = 2)

Refine with optional refinement factor rf. Uniform means refinement is done everywhere by the same factor, although nonuniform spacing functions may be used.

Definition at line 4447 of file nurbs.cpp.

◆ UpdateUniqueKV()

void mfem::NURBSExtension::UpdateUniqueKV ( )
protected

Update the unique set of KnotVectors. In 1D, this set is identical to the comprehensive set of KnotVectors.

Definition at line 3108 of file nurbs.cpp.

Friends And Related Symbol Documentation

◆ NURBSPatchMap

friend class NURBSPatchMap
friend

Definition at line 453 of file nurbs.hpp.

◆ ParNURBSExtension

friend class ParNURBSExtension
friend

Definition at line 451 of file nurbs.hpp.

Member Data Documentation

◆ activeBdrElem

Array<bool> mfem::NURBSExtension::activeBdrElem
protected

Definition at line 484 of file nurbs.hpp.

◆ activeDof

Array<int> mfem::NURBSExtension::activeDof
protected

Definition at line 485 of file nurbs.hpp.

◆ activeElem

Array<bool> mfem::NURBSExtension::activeElem
protected

Definition at line 483 of file nurbs.hpp.

◆ activeVert

Array<int> mfem::NURBSExtension::activeVert
protected

Definition at line 482 of file nurbs.hpp.

◆ bel_dof

Table * mfem::NURBSExtension::bel_dof
protected

Definition at line 525 of file nurbs.hpp.

◆ bel_to_IJK

Array2D<int> mfem::NURBSExtension::bel_to_IJK
protected

Definition at line 534 of file nurbs.hpp.

◆ bel_to_patch

Array<int> mfem::NURBSExtension::bel_to_patch
protected

Map from boundary element indices to patch indices.

Definition at line 530 of file nurbs.hpp.

◆ d_to_d

Array<int> mfem::NURBSExtension::d_to_d
protected

Periodic BC info:

  • dof 2 dof map
  • master and slave boundary indices

Definition at line 508 of file nurbs.hpp.

◆ e_meshOffsets

Array<int> mfem::NURBSExtension::e_meshOffsets
protected

Definition at line 514 of file nurbs.hpp.

◆ e_spaceOffsets

Array<int> mfem::NURBSExtension::e_spaceOffsets
protected

Definition at line 520 of file nurbs.hpp.

◆ edge_to_knot

Array<int> mfem::NURBSExtension::edge_to_knot
protected

Map from edge indices to KnotVector indices.

Definition at line 494 of file nurbs.hpp.

◆ el_dof

Table* mfem::NURBSExtension::el_dof
protected

Table of DOFs for each element (el_dof) or boundary element (bel_dof).

Definition at line 525 of file nurbs.hpp.

◆ el_to_IJK

Array2D<int> mfem::NURBSExtension::el_to_IJK
protected

Map from element indices to IJK knot span indices.

Definition at line 533 of file nurbs.hpp.

◆ el_to_patch

Array<int> mfem::NURBSExtension::el_to_patch
protected

Map from element indices to patch indices.

Definition at line 528 of file nurbs.hpp.

◆ f_meshOffsets

Array<int> mfem::NURBSExtension::f_meshOffsets
protected

Definition at line 515 of file nurbs.hpp.

◆ f_spaceOffsets

Array<int> mfem::NURBSExtension::f_spaceOffsets
protected

Definition at line 521 of file nurbs.hpp.

◆ knotVectors

Array<KnotVector *> mfem::NURBSExtension::knotVectors
protected

Set of unique KnotVectors.

Definition at line 497 of file nurbs.hpp.

◆ knotVectorsCompr

Array<KnotVector *> mfem::NURBSExtension::knotVectorsCompr
protected

Comprehensive set of all KnotVectors, one for every edge.

Definition at line 500 of file nurbs.hpp.

◆ master

Array<int> mfem::NURBSExtension::master
protected

Definition at line 509 of file nurbs.hpp.

◆ mode

Mode mfem::NURBSExtension::mode = Mode::H_1
protected

Definition at line 464 of file nurbs.hpp.

◆ mOrder

int mfem::NURBSExtension::mOrder
protected

Order of KnotVectors, see GetOrder() for description.

Definition at line 467 of file nurbs.hpp.

◆ mOrders

Array<int> mfem::NURBSExtension::mOrders
protected

Orders of all KnotVectors.

Definition at line 470 of file nurbs.hpp.

◆ NumOfActiveBdrElems

int mfem::NURBSExtension::NumOfActiveBdrElems
protected

Definition at line 479 of file nurbs.hpp.

◆ NumOfActiveDofs

int mfem::NURBSExtension::NumOfActiveDofs
protected

Definition at line 480 of file nurbs.hpp.

◆ NumOfActiveElems

int mfem::NURBSExtension::NumOfActiveElems
protected

Definition at line 479 of file nurbs.hpp.

◆ NumOfActiveVertices

int mfem::NURBSExtension::NumOfActiveVertices
protected

Local entity counts.

Definition at line 479 of file nurbs.hpp.

◆ NumOfBdrElements

int mfem::NURBSExtension::NumOfBdrElements
protected

Definition at line 476 of file nurbs.hpp.

◆ NumOfDofs

int mfem::NURBSExtension::NumOfDofs
protected

Definition at line 476 of file nurbs.hpp.

◆ NumOfElements

int mfem::NURBSExtension::NumOfElements
protected

Definition at line 476 of file nurbs.hpp.

◆ NumOfKnotVectors

int mfem::NURBSExtension::NumOfKnotVectors
protected

Number of KnotVectors.

Definition at line 473 of file nurbs.hpp.

◆ NumOfVertices

int mfem::NURBSExtension::NumOfVertices
protected

Global entity counts.

Definition at line 476 of file nurbs.hpp.

◆ own_topo

bool mfem::NURBSExtension::own_topo
protected

Whether this object owns patchTopo.

Definition at line 491 of file nurbs.hpp.

◆ p_meshOffsets

Array<int> mfem::NURBSExtension::p_meshOffsets
protected

Definition at line 516 of file nurbs.hpp.

◆ p_spaceOffsets

Array<int> mfem::NURBSExtension::p_spaceOffsets
protected

Definition at line 522 of file nurbs.hpp.

◆ patch_to_bel

std::vector<Array<int> > mfem::NURBSExtension::patch_to_bel
protected

For each patch p, patch_to_bel[p] lists all boundary elements in the patch.

Definition at line 539 of file nurbs.hpp.

◆ patch_to_el

std::vector<Array<int> > mfem::NURBSExtension::patch_to_el
protected

For each patch p, patch_to_el[p] lists all elements in the patch.

Definition at line 537 of file nurbs.hpp.

◆ patches

Array<NURBSPatch*> mfem::NURBSExtension::patches
protected

Array of all patches in the mesh.

Definition at line 542 of file nurbs.hpp.

◆ patchTopo

Mesh* mfem::NURBSExtension::patchTopo
protected

Patch topology mesh.

Definition at line 488 of file nurbs.hpp.

◆ slave

Array<int> mfem::NURBSExtension::slave
protected

Definition at line 510 of file nurbs.hpp.

◆ v_meshOffsets

Array<int> mfem::NURBSExtension::v_meshOffsets
protected

Global mesh offsets, meshOffsets == meshVertexOffsets.

Definition at line 513 of file nurbs.hpp.

◆ v_spaceOffsets

Array<int> mfem::NURBSExtension::v_spaceOffsets
protected

Global space offsets, spaceOffsets == dofOffsets.

Definition at line 519 of file nurbs.hpp.

◆ weights

Vector mfem::NURBSExtension::weights
protected

Weights for each control point or DOF.

Definition at line 503 of file nurbs.hpp.


The documentation for this class was generated from the following files: