MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Classes | Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Protected Attributes | Friends | List of all members
mfem::FiniteElementSpace Class Reference

Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of degrees of freedom. More...

#include <fespace.hpp>

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

Classes

class  DerefinementOperator
 Derefinement operator, used by the friend class InterpolationGridTransfer. More...
 
struct  key_hash
 
class  RefinementOperator
 GridFunction interpolation operator applicable after mesh refinement. More...
 

Public Member Functions

 FiniteElementSpace ()
 Default constructor: the object is invalid until initialized using the method Load(). More...
 
 FiniteElementSpace (const FiniteElementSpace &orig, Mesh *mesh=NULL, const FiniteElementCollection *fec=NULL)
 Copy constructor: deep copy all data from orig except the Mesh, the FiniteElementCollection, ans some derived data. More...
 
 FiniteElementSpace (Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
 
 FiniteElementSpace (Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
 Construct a NURBS FE space based on the given NURBSExtension, ext. More...
 
FiniteElementSpaceoperator= (const FiniteElementSpace &)=delete
 Copy assignment not supported. More...
 
MeshGetMesh () const
 Returns the mesh. More...
 
const NURBSExtensionGetNURBSext () const
 
NURBSExtensionGetNURBSext ()
 
NURBSExtensionStealNURBSext ()
 
bool Conforming () const
 
bool Nonconforming () const
 
void SetElementOrder (int i, int p)
 Sets the order of the i'th finite element. More...
 
int GetElementOrder (int i) const
 Returns the order of the i'th finite element. More...
 
int GetMaxElementOrder () const
 Return the maximum polynomial order. More...
 
bool IsVariableOrder () const
 Returns true if the space contains elements of varying polynomial orders. More...
 
const SparseMatrixGetConformingProlongation () const
 The returned SparseMatrix is owned by the FiniteElementSpace. More...
 
const SparseMatrixGetConformingRestriction () const
 The returned SparseMatrix is owned by the FiniteElementSpace. More...
 
const SparseMatrixGetHpConformingRestriction () const
 The returned SparseMatrix is owned by the FiniteElementSpace. More...
 
virtual const OperatorGetProlongationMatrix () const
 The returned Operator is owned by the FiniteElementSpace. More...
 
virtual const OperatorGetRestrictionTransposeOperator () const
 Return an operator that performs the transpose of GetRestrictionOperator. More...
 
virtual const OperatorGetRestrictionOperator () const
 An abstract operator that performs the same action as GetRestrictionMatrix. More...
 
virtual const SparseMatrixGetRestrictionMatrix () const
 The returned SparseMatrix is owned by the FiniteElementSpace. More...
 
virtual const SparseMatrixGetHpRestrictionMatrix () const
 The returned SparseMatrix is owned by the FiniteElementSpace. More...
 
const OperatorGetElementRestriction (ElementDofOrdering e_ordering) const
 Return an Operator that converts L-vectors to E-vectors. More...
 
virtual const FaceRestrictionGetFaceRestriction (ElementDofOrdering e_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
 Return an Operator that converts L-vectors to E-vectors on each face. More...
 
const QuadratureInterpolatorGetQuadratureInterpolator (const IntegrationRule &ir) const
 Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). More...
 
const QuadratureInterpolatorGetQuadratureInterpolator (const QuadratureSpace &qs) const
 Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). More...
 
const FaceQuadratureInterpolatorGetFaceQuadratureInterpolator (const IntegrationRule &ir, FaceType type) const
 Return a FaceQuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). More...
 
int GetOrder (int i) const
 Returns the polynomial degree of the i'th finite element. More...
 
int GetEdgeOrder (int edge, int variant=0) const
 
int GetFaceOrder (int face, int variant=0) const
 Returns the polynomial degree of the i'th face finite element. More...
 
int GetVDim () const
 Returns vector dimension. More...
 
int GetNDofs () const
 Returns number of degrees of freedom. More...
 
int GetVSize () const
 Return the number of vector dofs, i.e. GetNDofs() x GetVDim(). More...
 
virtual int GetTrueVSize () const
 Return the number of vector true (conforming) dofs. More...
 
int GetNConformingDofs () const
 
int GetConformingVSize () const
 
Ordering::Type GetOrdering () const
 Return the ordering method. More...
 
const FiniteElementCollectionFEColl () const
 
int GetNVDofs () const
 Number of all scalar vertex dofs. More...
 
int GetNEDofs () const
 Number of all scalar edge-interior dofs. More...
 
int GetNFDofs () const
 Number of all scalar face-interior dofs. More...
 
int GetNV () const
 Returns number of vertices in the mesh. More...
 
int GetNE () const
 Returns number of elements in the mesh. More...
 
int GetNF () const
 Returns number of faces (i.e. co-dimension 1 entities) in the mesh. More...
 
int GetNBE () const
 Returns number of boundary elements in the mesh. More...
 
int GetNFbyType (FaceType type) const
 Returns the number of faces according to the requested type. More...
 
int GetElementType (int i) const
 Returns the type of element i. More...
 
void GetElementVertices (int i, Array< int > &vertices) const
 Returns the vertices of element i. More...
 
int GetBdrElementType (int i) const
 Returns the type of boundary element i. More...
 
ElementTransformationGetElementTransformation (int i) const
 Returns ElementTransformation for the i-th element. More...
 
void GetElementTransformation (int i, IsoparametricTransformation *ElTr)
 Returns the transformation defining the i-th element in the user-defined variable ElTr. More...
 
ElementTransformationGetBdrElementTransformation (int i) const
 Returns ElementTransformation for the i-th boundary element. More...
 
int GetAttribute (int i) const
 
int GetBdrAttribute (int i) const
 
virtual DofTransformationGetElementDofs (int elem, Array< int > &dofs) const
 Returns indices of degrees of freedom of element 'elem'. More...
 
virtual DofTransformationGetBdrElementDofs (int bel, Array< int > &dofs) const
 Returns indices of degrees of freedom for boundary element 'bel'. More...
 
virtual int GetFaceDofs (int face, Array< int > &dofs, int variant=0) const
 Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edges and the vertices of the face. More...
 
int GetEdgeDofs (int edge, Array< int > &dofs, int variant=0) const
 Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vertices of the edge. More...
 
void GetVertexDofs (int i, Array< int > &dofs) const
 
void GetElementInteriorDofs (int i, Array< int > &dofs) const
 
void GetFaceInteriorDofs (int i, Array< int > &dofs) const
 
int GetNumElementInteriorDofs (int i) const
 
void GetEdgeInteriorDofs (int i, Array< int > &dofs) const
 
void GetVDofs (int vd, Array< int > &dofs, int ndofs=-1) const
 Returns the indices of all of the VDofs for the specified dimension 'vd'. More...
 
void DofsToVDofs (Array< int > &dofs, int ndofs=-1) const
 
void DofsToVDofs (int vd, Array< int > &dofs, int ndofs=-1) const
 
int DofToVDof (int dof, int vd, int ndofs=-1) const
 
int VDofToDof (int vdof) const
 
DofTransformationGetElementVDofs (int i, Array< int > &vdofs) const
 Returns indexes of degrees of freedom in array dofs for i'th element. More...
 
DofTransformationGetBdrElementVDofs (int i, Array< int > &vdofs) const
 Returns indexes of degrees of freedom for i'th boundary element. More...
 
void GetFaceVDofs (int i, Array< int > &vdofs) const
 Returns indexes of degrees of freedom for i'th face element (2D and 3D). More...
 
void GetEdgeVDofs (int i, Array< int > &vdofs) const
 Returns indexes of degrees of freedom for i'th edge. More...
 
void GetVertexVDofs (int i, Array< int > &vdofs) const
 
void GetElementInteriorVDofs (int i, Array< int > &vdofs) const
 
void GetEdgeInteriorVDofs (int i, Array< int > &vdofs) const
 
MFEM_DEPRECATED void RebuildElementToDofTable ()
 ( More...
 
void ReorderElementToDofTable ()
 Reorder the scalar DOFs based on the element ordering. More...
 
const TableGetElementToFaceOrientationTable () const
 
const TableGetElementToDofTable () const
 Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element, as returned by GetElementDofs(). More...
 
const TableGetBdrElementToDofTable () const
 Return a reference to the internal Table that stores the lists of scalar dofs, for each boundary mesh element, as returned by GetBdrElementDofs(). More...
 
const TableGetFaceToDofTable () const
 Return a reference to the internal Table that stores the lists of scalar dofs, for each face in the mesh, as returned by GetFaceDofs(). In this context, "face" refers to a (dim-1)-dimensional mesh entity. More...
 
void BuildDofToArrays ()
 Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof(). More...
 
int GetElementForDof (int i) const
 Return the index of the first element that contains dof i. More...
 
int GetLocalDofForDof (int i) const
 Return the local dof index in the first element that contains dof i. More...
 
virtual const FiniteElementGetFE (int i) const
 Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in the mesh object. More...
 
const FiniteElementGetBE (int i) const
 Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary face in the mesh object. More...
 
const FiniteElementGetFaceElement (int i) const
 Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the mesh object. Faces in this case refer to the MESHDIM-1 primitive so in 2D they are segments and in 1D they are points. More...
 
const FiniteElementGetEdgeElement (int i, int variant=0) const
 Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th edge in the mesh object. More...
 
const FiniteElementGetTraceElement (int i, Geometry::Type geom_type) const
 Return the trace element from element 'i' to the given 'geom_type'. More...
 
virtual void GetEssentialVDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
 Mark degrees of freedom associated with boundary elements with the specified boundary attributes (marked in 'bdr_attr_is_ess'). For spaces with 'vdim' > 1, the 'component' parameter can be used to restricts the marked vDOFs to the specified component. More...
 
virtual void GetEssentialTrueDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
 Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in the array bdr_attr_is_ess. For spaces with 'vdim' > 1, the 'component' parameter can be used to restricts the marked tDOFs to the specified component. More...
 
void GetBoundaryTrueDofs (Array< int > &boundary_dofs, int component=-1)
 Get a list of all boundary true dofs, boundary_dofs. For spaces with 'vdim' > 1, the 'component' parameter can be used to restricts the marked tDOFs to the specified component. Equivalent to FiniteElementSpace::GetEssentialTrueDofs with all boundary attributes marked as essential. More...
 
void ConvertToConformingVDofs (const Array< int > &dofs, Array< int > &cdofs)
 For a partially conforming FE space, convert a marker array (nonzero entries are true) on the partially conforming dofs to a marker array on the conforming dofs. A conforming dofs is marked iff at least one of its dependent dofs is marked. More...
 
void ConvertFromConformingVDofs (const Array< int > &cdofs, Array< int > &dofs)
 For a partially conforming FE space, convert a marker array (nonzero entries are true) on the conforming dofs to a marker array on the (partially conforming) dofs. A dof is marked iff it depends on a marked conforming dofs, where dependency is defined by the ConformingRestriction matrix; in other words, a dof is marked iff it corresponds to a marked conforming dof. More...
 
SparseMatrixD2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes)
 Generate the global restriction matrix from a discontinuous FE space to the continuous FE space of the same polynomial degree. More...
 
SparseMatrixD2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes)
 Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE space. More...
 
SparseMatrixH2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes)
 Construct the restriction matrix from the FE space given by (*this) to the lower degree FE space given by (*lfes) which is defined on the same mesh. More...
 
void GetTransferOperator (const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
 Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh. More...
 
virtual void GetTrueTransferOperator (const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
 Construct and return an Operator that can be used to transfer true-dof data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh. More...
 
virtual void Update (bool want_transform=true)
 Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation operator (unless want_transform is false). Safe to call multiple times, does nothing if space already up to date. More...
 
const OperatorGetUpdateOperator ()
 Get the GridFunction update operator. More...
 
void GetUpdateOperator (OperatorHandle &T)
 Return the update operator in the given OperatorHandle, T. More...
 
void SetUpdateOperatorOwner (bool own)
 Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator() must be deleted outside the FiniteElementSpace. More...
 
void SetUpdateOperatorType (Operator::Type tid)
 Specify the Operator::Type to be used by the update operators. More...
 
virtual void UpdatesFinished ()
 Free the GridFunction update operator (if any), to save memory. More...
 
long GetSequence () const
 
bool IsDGSpace () const
 Return whether or not the space is discontinuous (L2) More...
 
void SetRelaxedHpConformity (bool relaxed=true)
 
void Save (std::ostream &out) const
 Save finite element space to output stream out. More...
 
FiniteElementCollectionLoad (Mesh *m, std::istream &input)
 Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller. More...
 
virtual ~FiniteElementSpace ()
 

Static Public Member Functions

static void AdjustVDofs (Array< int > &vdofs)
 
static void MarkerToList (const Array< int > &marker, Array< int > &list)
 Convert a Boolean marker array to a list containing all marked indices. More...
 
static void ListToMarker (const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
 Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked with the given value and the rest are set to zero. More...
 

Protected Types

using key_face = std::tuple< bool, ElementDofOrdering, FaceType, L2FaceValues >
 The face restriction operators, see GetFaceRestriction(). More...
 
using map_L2F = std::unordered_map< const key_face, FaceRestriction *, key_hash >
 
typedef std::uint64_t VarOrderBits
 Bit-mask representing a set of orders needed by an edge/face. More...
 

Protected Member Functions

void UpdateNURBS ()
 
void Construct ()
 
void Destroy ()
 
void ConstructDoFTrans ()
 
void DestroyDoFTrans ()
 
void BuildElementToDofTable () const
 
void BuildBdrElementToDofTable () const
 
void BuildFaceToDofTable () const
 
void BuildNURBSFaceToDofTable () const
 Generates partial face_dof table for a NURBS space. More...
 
int GetElementOrderImpl (int i) const
 Return element order: internal version of GetElementOrder without checks. More...
 
void CalcEdgeFaceVarOrders (Array< VarOrderBits > &edge_orders, Array< VarOrderBits > &face_orders) const
 
int MakeDofTable (int ent_dim, const Array< int > &entity_orders, Table &entity_dofs, Array< char > *var_ent_order)
 
int FindDofs (const Table &var_dof_table, int row, int ndof) const
 Search row of a DOF table for a DOF set of size 'ndof', return first DOF. More...
 
int FindEdgeDof (int edge, int ndof) const
 
int FindFaceDof (int face, int ndof) const
 Similar to FindEdgeDof, but used for mixed meshes too. More...
 
int FirstFaceDof (int face, int variant=0) const
 
int GetNVariants (int entity, int index) const
 Return number of possible DOF variants for edge/face (var. order spaces). More...
 
int GetEntityDofs (int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID, int variant=0) const
 Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.). More...
 
int GetDegenerateFaceDofs (int index, Array< int > &dofs, Geometry::Type master_geom, int variant) const
 
int GetNumBorderDofs (Geometry::Type geom, int order) const
 
void BuildConformingInterpolation () const
 Calculate the cP and cR matrices for a nonconforming mesh. More...
 
void AddEdgeFaceDependencies (SparseMatrix &deps, Array< int > &master_dofs, const FiniteElement *master_fe, Array< int > &slave_dofs, int slave_face, const DenseMatrix *pm) const
 
void MakeVDimMatrix (SparseMatrix &mat) const
 Replicate 'mat' in the vector dimension, according to vdim ordering mode. More...
 
SparseMatrixRefinementMatrix_main (const int coarse_ndofs, const Table &coarse_elem_dof, const Table *coarse_elem_fos, const DenseTensor localP[]) const
 
void GetLocalRefinementMatrices (Geometry::Type geom, DenseTensor &localP) const
 
void GetLocalDerefinementMatrices (Geometry::Type geom, DenseTensor &localR) const
 
SparseMatrixRefinementMatrix (int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos)
 
SparseMatrixDerefinementMatrix (int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos)
 Calculate GridFunction restriction matrix after mesh derefinement. More...
 
void GetLocalRefinementMatrices (const FiniteElementSpace &coarse_fes, Geometry::Type geom, DenseTensor &localP) const
 Return in localP the local refinement matrices that map between fespaces after mesh refinement. More...
 
void Constructor (Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
 Help function for constructors + Load(). More...
 
virtual void UpdateMeshPointer (Mesh *new_mesh)
 
void UpdateElementOrders ()
 Resize the elem_order array on mesh change. More...
 
virtual void CopyProlongationAndRestriction (const FiniteElementSpace &fes, const Array< int > *perm)
 Copies the prolongation and restriction matrices from fes. More...
 

Static Protected Member Functions

static int MinOrder (VarOrderBits bits)
 Return the minimum order (least significant bit set) in the bit mask. More...
 
static int EncodeDof (int entity_base, int idx)
 Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes). More...
 
static int DecodeDof (int dof)
 Helpers to remove encoded sign from a DOF. More...
 
static int DecodeDof (int dof, double &sign)
 
static void AddDependencies (SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I, int skipfirst=0)
 
static bool DofFinalizable (int dof, const Array< bool > &finalized, const SparseMatrix &deps)
 

Protected Attributes

Meshmesh
 The mesh that FE space lives on (not owned). More...
 
const FiniteElementCollectionfec
 Associated FE collection (not owned). More...
 
int vdim
 Vector dimension (number of unknowns per degree of freedom). More...
 
Ordering::Type ordering
 
int ndofs
 Number of degrees of freedom. Number of unknowns is ndofs * vdim. More...
 
Array< char > elem_order
 
int nvdofs
 
int nedofs
 
int nfdofs
 
int nbdofs
 
int uni_fdof
 

of single face DOFs if all faces uniform; -1 otherwise

More...
 
int * bdofs
 internal DOFs of elements if mixed/var-order; NULL otherwise More...
 
Table var_edge_dofs
 
Table var_face_dofs
 NOTE: also used for spaces with mixed faces. More...
 
Array< char > var_edge_orders
 
Array< char > var_face_orders
 
Tableelem_dof
 
Tableelem_fos
 
Tablebdr_elem_dof
 
Tablebdr_elem_fos
 
Tableface_dof
 
Array< int > dof_elem_array
 
Array< int > dof_ldof_array
 
NURBSExtensionNURBSext
 
int own_ext
 
Array< int > face_to_be
 
Array< DofTransformation * > DoFTrans
 
VDofTransformation VDoFTrans
 
SparseMatrixcP
 
SparseMatrixcR
 Conforming restriction matrix such that cR.cP=I. More...
 
SparseMatrixcR_hp
 A version of the conforming restriction matrix for variable-order spaces. More...
 
bool cP_is_set
 
OperatorHandle Th
 Transformation to apply to GridFunctions after space Update(). More...
 
OperatorHandle L2E_nat
 The element restriction operators, see GetElementRestriction(). More...
 
OperatorHandle L2E_lex
 
map_L2F L2F
 
Array< QuadratureInterpolator * > E2Q_array
 
Array
< FaceQuadratureInterpolator * > 
E2IFQ_array
 
Array
< FaceQuadratureInterpolator * > 
E2BFQ_array
 
long sequence
 
long mesh_sequence
 
bool orders_changed
 True if at least one element order changed (variable-order space only). More...
 
bool relaxed_hp
 

Static Protected Attributes

static constexpr int MaxVarOrder = 8*sizeof(VarOrderBits) - 1
 

Friends

class InterpolationGridTransfer
 
class PRefinementTransferOperator
 
class LORBase
 
void Mesh::Swap (Mesh &, bool)
 

Detailed Description

Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of degrees of freedom.

Definition at line 88 of file fespace.hpp.

Member Typedef Documentation

The face restriction operators, see GetFaceRestriction().

Definition at line 162 of file fespace.hpp.

using mfem::FiniteElementSpace::map_L2F = std::unordered_map<const key_face,FaceRestriction*,key_hash>
protected

Definition at line 173 of file fespace.hpp.

typedef std::uint64_t mfem::FiniteElementSpace::VarOrderBits
protected

Bit-mask representing a set of orders needed by an edge/face.

Definition at line 212 of file fespace.hpp.

Constructor & Destructor Documentation

mfem::FiniteElementSpace::FiniteElementSpace ( )

Default constructor: the object is invalid until initialized using the method Load().

Definition at line 59 of file fespace.cpp.

mfem::FiniteElementSpace::FiniteElementSpace ( const FiniteElementSpace orig,
Mesh mesh = NULL,
const FiniteElementCollection fec = NULL 
)

Copy constructor: deep copy all data from orig except the Mesh, the FiniteElementCollection, ans some derived data.

If the mesh or fec pointers are NULL (default), then the new FiniteElementSpace will reuse the respective pointers from orig. If any of these pointers is not NULL, the given pointer will be used instead of the one used by orig.

Note
The objects pointed to by the mesh and fec parameters must be either the same objects as the ones used by orig, or copies of them. Otherwise, the behavior is undefined.
Derived data objects, such as the conforming prolongation and restriction matrices, and the update operator, will not be copied, even if they are created in the orig object.

Definition at line 72 of file fespace.cpp.

mfem::FiniteElementSpace::FiniteElementSpace ( Mesh mesh,
const FiniteElementCollection fec,
int  vdim = 1,
int  ordering = Ordering::byNODES 
)
inline

Definition at line 415 of file fespace.hpp.

mfem::FiniteElementSpace::FiniteElementSpace ( Mesh mesh,
NURBSExtension ext,
const FiniteElementCollection fec,
int  vdim = 1,
int  ordering = Ordering::byNODES 
)
inline

Construct a NURBS FE space based on the given NURBSExtension, ext.

Note
If the pointer ext is NULL, this constructor is equivalent to the standard constructor with the same arguments minus the NURBSExtension, ext.

Definition at line 424 of file fespace.hpp.

mfem::FiniteElementSpace::~FiniteElementSpace ( )
virtual

Definition at line 3180 of file fespace.cpp.

Member Function Documentation

void mfem::FiniteElementSpace::AddDependencies ( SparseMatrix deps,
Array< int > &  master_dofs,
Array< int > &  slave_dofs,
DenseMatrix I,
int  skipfirst = 0 
)
staticprotected

Definition at line 758 of file fespace.cpp.

void mfem::FiniteElementSpace::AddEdgeFaceDependencies ( SparseMatrix deps,
Array< int > &  master_dofs,
const FiniteElement master_fe,
Array< int > &  slave_dofs,
int  slave_face,
const DenseMatrix pm 
) const
protected

Definition at line 783 of file fespace.cpp.

void mfem::FiniteElementSpace::AdjustVDofs ( Array< int > &  vdofs)
static

Definition at line 267 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildBdrElementToDofTable ( ) const
protected

Definition at line 383 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildConformingInterpolation ( ) const
protected

Calculate the cP and cR matrices for a nonconforming mesh.

Definition at line 928 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildDofToArrays ( )

Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof().

Definition at line 460 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildElementToDofTable ( ) const
protected

Definition at line 342 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildFaceToDofTable ( ) const
protected

Definition at line 405 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildNURBSFaceToDofTable ( ) const
protected

Generates partial face_dof table for a NURBS space.

The table is only defined for exterior faces that coincide with a boundary.

Definition at line 2266 of file fespace.cpp.

void mfem::FiniteElementSpace::CalcEdgeFaceVarOrders ( Array< VarOrderBits > &  edge_orders,
Array< VarOrderBits > &  face_orders 
) const
protected

In a variable order space, calculate a bitmask of polynomial orders that need to be represented on each edge and face.

Definition at line 2452 of file fespace.cpp.

bool mfem::FiniteElementSpace::Conforming ( ) const
inline

Definition at line 439 of file fespace.hpp.

void mfem::FiniteElementSpace::Construct ( )
protected

Definition at line 2315 of file fespace.cpp.

void mfem::FiniteElementSpace::ConstructDoFTrans ( )
protected

Definition at line 2192 of file fespace.cpp.

void mfem::FiniteElementSpace::Constructor ( Mesh mesh,
NURBSExtension ext,
const FiniteElementCollection fec,
int  vdim = 1,
int  ordering = Ordering::byNODES 
)
protected

Help function for constructors + Load().

Definition at line 2141 of file fespace.cpp.

void mfem::FiniteElementSpace::ConvertFromConformingVDofs ( const Array< int > &  cdofs,
Array< int > &  dofs 
)

For a partially conforming FE space, convert a marker array (nonzero entries are true) on the conforming dofs to a marker array on the (partially conforming) dofs. A dof is marked iff it depends on a marked conforming dofs, where dependency is defined by the ConformingRestriction matrix; in other words, a dof is marked iff it corresponds to a marked conforming dof.

Definition at line 638 of file fespace.cpp.

void mfem::FiniteElementSpace::ConvertToConformingVDofs ( const Array< int > &  dofs,
Array< int > &  cdofs 
)

For a partially conforming FE space, convert a marker array (nonzero entries are true) on the partially conforming dofs to a marker array on the conforming dofs. A conforming dofs is marked iff at least one of its dependent dofs is marked.

Definition at line 630 of file fespace.cpp.

void mfem::FiniteElementSpace::CopyProlongationAndRestriction ( const FiniteElementSpace fes,
const Array< int > *  perm 
)
protectedvirtual

Copies the prolongation and restriction matrices from fes.

Used for low order preconditioning on non-conforming meshes. If the DOFs require a permutation, it will be supplied by non-NULL perm. NULL perm indicates that no permutation is required.

Definition at line 100 of file fespace.cpp.

SparseMatrix * mfem::FiniteElementSpace::D2C_GlobalRestrictionMatrix ( FiniteElementSpace cfes)

Generate the global restriction matrix from a discontinuous FE space to the continuous FE space of the same polynomial degree.

Definition at line 647 of file fespace.cpp.

SparseMatrix * mfem::FiniteElementSpace::D2Const_GlobalRestrictionMatrix ( FiniteElementSpace cfes)

Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE space.

Definition at line 679 of file fespace.cpp.

static int mfem::FiniteElementSpace::DecodeDof ( int  dof)
inlinestaticprotected

Helpers to remove encoded sign from a DOF.

Definition at line 254 of file fespace.hpp.

static int mfem::FiniteElementSpace::DecodeDof ( int  dof,
double &  sign 
)
inlinestaticprotected

Definition at line 257 of file fespace.hpp.

SparseMatrix * mfem::FiniteElementSpace::DerefinementMatrix ( int  old_ndofs,
const Table old_elem_dof,
const Table old_elem_fos 
)
protected

Calculate GridFunction restriction matrix after mesh derefinement.

TODO: Implement DofTransformation support

Definition at line 2040 of file fespace.cpp.

void mfem::FiniteElementSpace::Destroy ( )
protected

Definition at line 3185 of file fespace.cpp.

void mfem::FiniteElementSpace::DestroyDoFTrans ( )
protected

Definition at line 3237 of file fespace.cpp.

bool mfem::FiniteElementSpace::DofFinalizable ( int  dof,
const Array< bool > &  finalized,
const SparseMatrix deps 
)
staticprotected

Definition at line 838 of file fespace.cpp.

void mfem::FiniteElementSpace::DofsToVDofs ( Array< int > &  dofs,
int  ndofs = -1 
) const

Definition at line 215 of file fespace.cpp.

void mfem::FiniteElementSpace::DofsToVDofs ( int  vd,
Array< int > &  dofs,
int  ndofs = -1 
) const

Definition at line 230 of file fespace.cpp.

int mfem::FiniteElementSpace::DofToVDof ( int  dof,
int  vd,
int  ndofs = -1 
) const

Definition at line 251 of file fespace.cpp.

static int mfem::FiniteElementSpace::EncodeDof ( int  entity_base,
int  idx 
)
inlinestaticprotected

Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).

Definition at line 250 of file fespace.hpp.

const FiniteElementCollection* mfem::FiniteElementSpace::FEColl ( ) const
inline

Definition at line 577 of file fespace.hpp.

int mfem::FiniteElementSpace::FindDofs ( const Table var_dof_table,
int  row,
int  ndof 
) const
protected

Search row of a DOF table for a DOF set of size 'ndof', return first DOF.

Definition at line 2623 of file fespace.cpp.

int mfem::FiniteElementSpace::FindEdgeDof ( int  edge,
int  ndof 
) const
inlineprotected

In a variable order space, return edge DOFs associated with a polynomial order that has 'ndof' degrees of freedom.

Definition at line 236 of file fespace.hpp.

int mfem::FiniteElementSpace::FindFaceDof ( int  face,
int  ndof 
) const
inlineprotected

Similar to FindEdgeDof, but used for mixed meshes too.

Definition at line 240 of file fespace.hpp.

int mfem::FiniteElementSpace::FirstFaceDof ( int  face,
int  variant = 0 
) const
inlineprotected

Definition at line 243 of file fespace.hpp.

int mfem::FiniteElementSpace::GetAttribute ( int  i) const
inline

Definition at line 635 of file fespace.hpp.

int mfem::FiniteElementSpace::GetBdrAttribute ( int  i) const
inline

Definition at line 637 of file fespace.hpp.

DofTransformation * mfem::FiniteElementSpace::GetBdrElementDofs ( int  bel,
Array< int > &  dofs 
) const
virtual

Returns indices of degrees of freedom for boundary element 'bel'.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 2814 of file fespace.cpp.

const Table& mfem::FiniteElementSpace::GetBdrElementToDofTable ( ) const
inline

Return a reference to the internal Table that stores the lists of scalar dofs, for each boundary mesh element, as returned by GetBdrElementDofs().

Definition at line 731 of file fespace.hpp.

ElementTransformation* mfem::FiniteElementSpace::GetBdrElementTransformation ( int  i) const
inline

Returns ElementTransformation for the i-th boundary element.

Definition at line 632 of file fespace.hpp.

int mfem::FiniteElementSpace::GetBdrElementType ( int  i) const
inline

Returns the type of boundary element i.

Definition at line 619 of file fespace.hpp.

DofTransformation * mfem::FiniteElementSpace::GetBdrElementVDofs ( int  i,
Array< int > &  vdofs 
) const

Returns indexes of degrees of freedom for i'th boundary element.

Definition at line 297 of file fespace.cpp.

const FiniteElement * mfem::FiniteElementSpace::GetBE ( int  i) const

Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary face in the mesh object.

Definition at line 3102 of file fespace.cpp.

void mfem::FiniteElementSpace::GetBoundaryTrueDofs ( Array< int > &  boundary_dofs,
int  component = -1 
)

Get a list of all boundary true dofs, boundary_dofs. For spaces with 'vdim' > 1, the 'component' parameter can be used to restricts the marked tDOFs to the specified component. Equivalent to FiniteElementSpace::GetEssentialTrueDofs with all boundary attributes marked as essential.

Definition at line 582 of file fespace.cpp.

const SparseMatrix * mfem::FiniteElementSpace::GetConformingProlongation ( ) const

The returned SparseMatrix is owned by the FiniteElementSpace.

Definition at line 1234 of file fespace.cpp.

const SparseMatrix * mfem::FiniteElementSpace::GetConformingRestriction ( ) const

The returned SparseMatrix is owned by the FiniteElementSpace.

Definition at line 1241 of file fespace.cpp.

int mfem::FiniteElementSpace::GetConformingVSize ( ) const
inline

Definition at line 572 of file fespace.hpp.

int mfem::FiniteElementSpace::GetDegenerateFaceDofs ( int  index,
Array< int > &  dofs,
Geometry::Type  master_geom,
int  variant 
) const
protected

Definition at line 852 of file fespace.cpp.

int mfem::FiniteElementSpace::GetEdgeDofs ( int  edge,
Array< int > &  dofs,
int  variant = 0 
) const

Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vertices of the edge.

In variable order spaces, multiple sets of DOFs may exist on an edge, corresponding to the different polynomial orders of incident elements. The 'variant' parameter is the zero-based index of the desired DOF set. The variants are ordered from lowest polynomial degree to the highest.

Returns
Order of the selected variant, or -1 if there are no more variants.

Definition at line 2987 of file fespace.cpp.

const FiniteElement * mfem::FiniteElementSpace::GetEdgeElement ( int  i,
int  variant = 0 
) const

Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th edge in the mesh object.

Definition at line 3165 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEdgeInteriorDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 3066 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEdgeInteriorVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 336 of file fespace.cpp.

int mfem::FiniteElementSpace::GetEdgeOrder ( int  edge,
int  variant = 0 
) const

Return the order of an edge. In a variable order space, return the order of a specific variant, or -1 if there are no more variants.

Definition at line 2640 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEdgeVDofs ( int  i,
Array< int > &  vdofs 
) const

Returns indexes of degrees of freedom for i'th edge.

Definition at line 318 of file fespace.cpp.

DofTransformation * mfem::FiniteElementSpace::GetElementDofs ( int  elem,
Array< int > &  dofs 
) const
virtual

Returns indices of degrees of freedom of element 'elem'.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 2680 of file fespace.cpp.

int mfem::FiniteElementSpace::GetElementForDof ( int  i) const
inline

Return the index of the first element that contains dof i.

This method can be called only after setup is performed using the method BuildDofToArrays().

Definition at line 749 of file fespace.hpp.

void mfem::FiniteElementSpace::GetElementInteriorDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 3045 of file fespace.cpp.

void mfem::FiniteElementSpace::GetElementInteriorVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 330 of file fespace.cpp.

int mfem::FiniteElementSpace::GetElementOrder ( int  i) const

Returns the order of the i'th finite element.

Definition at line 178 of file fespace.cpp.

int mfem::FiniteElementSpace::GetElementOrderImpl ( int  i) const
protected

Return element order: internal version of GetElementOrder without checks.

Definition at line 189 of file fespace.cpp.

const Operator * mfem::FiniteElementSpace::GetElementRestriction ( ElementDofOrdering  e_ordering) const

Return an Operator that converts L-vectors to E-vectors.

An L-vector is a vector of size GetVSize() which is the same size as a GridFunction. An E-vector represents the element-wise discontinuous version of the FE space.

The layout of the E-vector is: ND x VDIM x NE, where ND is the number of degrees of freedom, VDIM is the vector dimension of the FE space, and NE is the number of the mesh elements.

The parameter e_ordering describes how the local DOFs in each element should be ordered, see ElementDofOrdering.

For discontinuous spaces, the element restriction corresponds to a permutation of the degrees of freedom, implemented by the L2ElementRestriction class.

The returned Operator is owned by the FiniteElementSpace.

Definition at line 1261 of file fespace.cpp.

const Table& mfem::FiniteElementSpace::GetElementToDofTable ( ) const
inline

Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element, as returned by GetElementDofs().

Definition at line 726 of file fespace.hpp.

const Table* mfem::FiniteElementSpace::GetElementToFaceOrientationTable ( ) const
inline

Definition at line 722 of file fespace.hpp.

ElementTransformation* mfem::FiniteElementSpace::GetElementTransformation ( int  i) const
inline

Returns ElementTransformation for the i-th element.

Definition at line 623 of file fespace.hpp.

void mfem::FiniteElementSpace::GetElementTransformation ( int  i,
IsoparametricTransformation ElTr 
)
inline

Returns the transformation defining the i-th element in the user-defined variable ElTr.

Definition at line 628 of file fespace.hpp.

int mfem::FiniteElementSpace::GetElementType ( int  i) const
inline

Returns the type of element i.

Definition at line 611 of file fespace.hpp.

DofTransformation * mfem::FiniteElementSpace::GetElementVDofs ( int  i,
Array< int > &  vdofs 
) const

Returns indexes of degrees of freedom in array dofs for i'th element.

Definition at line 281 of file fespace.cpp.

void mfem::FiniteElementSpace::GetElementVertices ( int  i,
Array< int > &  vertices 
) const
inline

Returns the vertices of element i.

Definition at line 615 of file fespace.hpp.

int mfem::FiniteElementSpace::GetEntityDofs ( int  entity,
int  index,
Array< int > &  dofs,
Geometry::Type  master_geom = Geometry::INVALID,
int  variant = 0 
) const
protected

Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).

Definition at line 903 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEssentialTrueDofs ( const Array< int > &  bdr_attr_is_ess,
Array< int > &  ess_tdof_list,
int  component = -1 
)
virtual

Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in the array bdr_attr_is_ess. For spaces with 'vdim' > 1, the 'component' parameter can be used to restricts the marked tDOFs to the specified component.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 564 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEssentialVDofs ( const Array< int > &  bdr_attr_is_ess,
Array< int > &  ess_vdofs,
int  component = -1 
) const
virtual

Mark degrees of freedom associated with boundary elements with the specified boundary attributes (marked in 'bdr_attr_is_ess'). For spaces with 'vdim' > 1, the 'component' parameter can be used to restricts the marked vDOFs to the specified component.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 495 of file fespace.cpp.

int mfem::FiniteElementSpace::GetFaceDofs ( int  face,
Array< int > &  dofs,
int  variant = 0 
) const
virtual

Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edges and the vertices of the face.

In variable order spaces, multiple variants of DOFs can be returned. See GetEdgeDofs for more details.

Returns
Order of the selected variant, or -1 if there are no more variants.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 2906 of file fespace.cpp.

const FiniteElement * mfem::FiniteElementSpace::GetFaceElement ( int  i) const

Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the mesh object. Faces in this case refer to the MESHDIM-1 primitive so in 2D they are segments and in 1D they are points.

Definition at line 3135 of file fespace.cpp.

void mfem::FiniteElementSpace::GetFaceInteriorDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 3078 of file fespace.cpp.

int mfem::FiniteElementSpace::GetFaceOrder ( int  face,
int  variant = 0 
) const

Returns the polynomial degree of the i'th face finite element.

Definition at line 2651 of file fespace.cpp.

const FaceQuadratureInterpolator * mfem::FiniteElementSpace::GetFaceQuadratureInterpolator ( const IntegrationRule ir,
FaceType  type 
) const

Return a FaceQuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors).

Definition at line 1358 of file fespace.cpp.

const FaceRestriction * mfem::FiniteElementSpace::GetFaceRestriction ( ElementDofOrdering  e_ordering,
FaceType  type,
L2FaceValues  mul = L2FaceValues::DoubleValued 
) const
virtual

Return an Operator that converts L-vectors to E-vectors on each face.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 1294 of file fespace.cpp.

const Table& mfem::FiniteElementSpace::GetFaceToDofTable ( ) const
inline

Return a reference to the internal Table that stores the lists of scalar dofs, for each face in the mesh, as returned by GetFaceDofs(). In this context, "face" refers to a (dim-1)-dimensional mesh entity.

Note
In the case of a NURBS space, the rows corresponding to interior faces will be empty.

Definition at line 739 of file fespace.hpp.

void mfem::FiniteElementSpace::GetFaceVDofs ( int  i,
Array< int > &  vdofs 
) const

Returns indexes of degrees of freedom for i'th face element (2D and 3D).

Definition at line 312 of file fespace.cpp.

const FiniteElement * mfem::FiniteElementSpace::GetFE ( int  i) const
virtual

Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in the mesh object.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 2783 of file fespace.cpp.

const SparseMatrix * mfem::FiniteElementSpace::GetHpConformingRestriction ( ) const

The returned SparseMatrix is owned by the FiniteElementSpace.

Return a version of the conforming restriction matrix for variable-order spaces with complex hp interfaces, where some true DOFs are not owned by any elements and need to be interpolated from higher order edge/face variants (see also SetRelaxedHpConformity()).

Definition at line 1248 of file fespace.cpp.

virtual const SparseMatrix* mfem::FiniteElementSpace::GetHpRestrictionMatrix ( ) const
inlinevirtual

The returned SparseMatrix is owned by the FiniteElementSpace.

Definition at line 491 of file fespace.hpp.

void mfem::FiniteElementSpace::GetLocalDerefinementMatrices ( Geometry::Type  geom,
DenseTensor localR 
) const
protected

Definition at line 2016 of file fespace.cpp.

int mfem::FiniteElementSpace::GetLocalDofForDof ( int  i) const
inline

Return the local dof index in the first element that contains dof i.

This method can be called only after setup is performed using the method BuildDofToArrays().

Definition at line 753 of file fespace.hpp.

void mfem::FiniteElementSpace::GetLocalRefinementMatrices ( Geometry::Type  geom,
DenseTensor localP 
) const
protected

Definition at line 1453 of file fespace.cpp.

void mfem::FiniteElementSpace::GetLocalRefinementMatrices ( const FiniteElementSpace coarse_fes,
Geometry::Type  geom,
DenseTensor localP 
) const
protected

Return in localP the local refinement matrices that map between fespaces after mesh refinement.

This method assumes that this->mesh is a refinement of coarse_fes->mesh and that the CoarseFineTransformations of this->mesh are set accordingly. Another assumption is that the FEs of this use the same MapType as the FEs of coarse_fes. Finally, it assumes that the spaces this and coarse_fes are NOT variable-order spaces.

Definition at line 2114 of file fespace.cpp.

int mfem::FiniteElementSpace::GetMaxElementOrder ( ) const
inline

Return the maximum polynomial order.

Definition at line 451 of file fespace.hpp.

Mesh* mfem::FiniteElementSpace::GetMesh ( ) const
inline

Returns the mesh.

Definition at line 433 of file fespace.hpp.

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

Returns number of boundary elements in the mesh.

Definition at line 599 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNConformingDofs ( ) const

Returns the number of conforming ("true") degrees of freedom (if the space is on a nonconforming mesh with hanging nodes).

Definition at line 1255 of file fespace.cpp.

int mfem::FiniteElementSpace::GetNDofs ( ) const
inline

Returns number of degrees of freedom.

Definition at line 560 of file fespace.hpp.

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

Returns number of elements in the mesh.

Definition at line 590 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNEDofs ( ) const
inline

Number of all scalar edge-interior dofs.

Definition at line 582 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNF ( ) const
inline

Returns number of faces (i.e. co-dimension 1 entities) in the mesh.

The co-dimension 1 entities are those that have dimension 1 less than the mesh dimension, e.g. for a 2D mesh, the faces are the 1D entities, i.e. the edges.

Definition at line 596 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNFbyType ( FaceType  type) const
inline

Returns the number of faces according to the requested type.

If type==Boundary returns only the "true" number of boundary faces contrary to GetNBE() that returns "fake" boundary faces associated to visualization for GLVis. Similarly, if type==Interior, the "fake" boundary faces associated to visualization are counted as interior faces.

Definition at line 607 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNFDofs ( ) const
inline

Number of all scalar face-interior dofs.

Definition at line 584 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNumBorderDofs ( Geometry::Type  geom,
int  order 
) const
protected

Definition at line 894 of file fespace.cpp.

int mfem::FiniteElementSpace::GetNumElementInteriorDofs ( int  i) const

Definition at line 3060 of file fespace.cpp.

const NURBSExtension* mfem::FiniteElementSpace::GetNURBSext ( ) const
inline

Definition at line 435 of file fespace.hpp.

NURBSExtension* mfem::FiniteElementSpace::GetNURBSext ( )
inline

Definition at line 436 of file fespace.hpp.

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

Returns number of vertices in the mesh.

Definition at line 587 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNVariants ( int  entity,
int  index 
) const
protected

Return number of possible DOF variants for edge/face (var. order spaces).

Definition at line 2667 of file fespace.cpp.

int mfem::FiniteElementSpace::GetNVDofs ( ) const
inline

Number of all scalar vertex dofs.

Definition at line 580 of file fespace.hpp.

int mfem::FiniteElementSpace::GetOrder ( int  i) const
inline

Returns the polynomial degree of the i'th finite element.

NOTE: it is recommended to use GetElementOrder in new code.

Definition at line 547 of file fespace.hpp.

Ordering::Type mfem::FiniteElementSpace::GetOrdering ( ) const
inline

Return the ordering method.

Definition at line 575 of file fespace.hpp.

virtual const Operator* mfem::FiniteElementSpace::GetProlongationMatrix ( ) const
inlinevirtual

The returned Operator is owned by the FiniteElementSpace.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 471 of file fespace.hpp.

const QuadratureInterpolator * mfem::FiniteElementSpace::GetQuadratureInterpolator ( const IntegrationRule ir) const

Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors).

An E-vector represents the element-wise discontinuous version of the FE space and can be obtained, for example, from a GridFunction using the Operator returned by GetElementRestriction().

All elements will use the same IntegrationRule, ir as the target quadrature points.

Definition at line 1329 of file fespace.cpp.

const QuadratureInterpolator * mfem::FiniteElementSpace::GetQuadratureInterpolator ( const QuadratureSpace qs) const

Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors).

An E-vector represents the element-wise discontinuous version of the FE space and can be obtained, for example, from a GridFunction using the Operator returned by GetElementRestriction().

The target quadrature points in the elements are described by the given QuadratureSpace, qs.

Definition at line 1343 of file fespace.cpp.

virtual const SparseMatrix* mfem::FiniteElementSpace::GetRestrictionMatrix ( ) const
inlinevirtual

The returned SparseMatrix is owned by the FiniteElementSpace.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 487 of file fespace.hpp.

virtual const Operator* mfem::FiniteElementSpace::GetRestrictionOperator ( ) const
inlinevirtual

An abstract operator that performs the same action as GetRestrictionMatrix.

In some cases this is an optimized matrix-free implementation. The returned operator is owned by the FiniteElementSpace.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 483 of file fespace.hpp.

virtual const Operator* mfem::FiniteElementSpace::GetRestrictionTransposeOperator ( ) const
inlinevirtual

Return an operator that performs the transpose of GetRestrictionOperator.

The returned operator is owned by the FiniteElementSpace. In serial this is the same as GetProlongationMatrix()

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 477 of file fespace.hpp.

long mfem::FiniteElementSpace::GetSequence ( ) const
inline

Return update counter, similar to Mesh::GetSequence(). Used by GridFunction to check if it is up to date with the space.

Definition at line 898 of file fespace.hpp.

const FiniteElement * mfem::FiniteElementSpace::GetTraceElement ( int  i,
Geometry::Type  geom_type 
) const

Return the trace element from element 'i' to the given 'geom_type'.

Definition at line 3175 of file fespace.cpp.

void mfem::FiniteElementSpace::GetTransferOperator ( const FiniteElementSpace coarse_fes,
OperatorHandle T 
) const

Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh.

It is assumed that the mesh of this FE space is a refinement of the mesh of coarse_fes and the CoarseFineTransformations returned by the method Mesh::GetRefinementTransforms() of the refined mesh are set accordingly. The Operator::Type of T can be set to request an Operator of the set type. Currently, only Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE (matrix-free) are supported. When Operator::ANY_TYPE is requested, the choice of the particular Operator sub-class is left to the method. This method also works in parallel because the transfer operator is local to the MPI task when the input is a synchronized ParGridFunction.

Definition at line 3246 of file fespace.cpp.

void mfem::FiniteElementSpace::GetTrueTransferOperator ( const FiniteElementSpace coarse_fes,
OperatorHandle T 
) const
virtual

Construct and return an Operator that can be used to transfer true-dof data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh.

This method calls GetTransferOperator() and multiplies the result by the prolongation operator of coarse_fes on the right, and by the restriction operator of this FE space on the left.

The Operator::Type of T can be set to request an Operator of the set type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE (matrix-free). In parallel, the supported types are: Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated as Operator::ANY_TYPE: the operator representation choice is made by this method.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 3273 of file fespace.cpp.

virtual int mfem::FiniteElementSpace::GetTrueVSize ( ) const
inlinevirtual

Return the number of vector true (conforming) dofs.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 566 of file fespace.hpp.

const Operator* mfem::FiniteElementSpace::GetUpdateOperator ( )
inline

Get the GridFunction update operator.

Definition at line 873 of file fespace.hpp.

void mfem::FiniteElementSpace::GetUpdateOperator ( OperatorHandle T)
inline

Return the update operator in the given OperatorHandle, T.

Definition at line 876 of file fespace.hpp.

int mfem::FiniteElementSpace::GetVDim ( ) const
inline

Returns vector dimension.

Definition at line 557 of file fespace.hpp.

void mfem::FiniteElementSpace::GetVDofs ( int  vd,
Array< int > &  dofs,
int  ndofs = -1 
) const

Returns the indices of all of the VDofs for the specified dimension 'vd'.

The 'ndofs' parameter defines the number of Dofs in the FiniteElementSpace. If 'ndofs' is -1 (the default value), then the number of Dofs is determined by the FiniteElementSpace.

Definition at line 195 of file fespace.cpp.

void mfem::FiniteElementSpace::GetVertexDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 3035 of file fespace.cpp.

void mfem::FiniteElementSpace::GetVertexVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 324 of file fespace.cpp.

int mfem::FiniteElementSpace::GetVSize ( ) const
inline

Return the number of vector dofs, i.e. GetNDofs() x GetVDim().

Definition at line 563 of file fespace.hpp.

SparseMatrix * mfem::FiniteElementSpace::H2L_GlobalRestrictionMatrix ( FiniteElementSpace lfes)

Construct the restriction matrix from the FE space given by (*this) to the lower degree FE space given by (*lfes) which is defined on the same mesh.

Definition at line 710 of file fespace.cpp.

bool mfem::FiniteElementSpace::IsDGSpace ( ) const
inline

Return whether or not the space is discontinuous (L2)

Definition at line 901 of file fespace.hpp.

bool mfem::FiniteElementSpace::IsVariableOrder ( ) const
inline

Returns true if the space contains elements of varying polynomial orders.

Definition at line 455 of file fespace.hpp.

void mfem::FiniteElementSpace::ListToMarker ( const Array< int > &  list,
int  marker_size,
Array< int > &  marker,
int  mark_val = -1 
)
static

Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked with the given value and the rest are set to zero.

Definition at line 617 of file fespace.cpp.

FiniteElementCollection * mfem::FiniteElementSpace::Load ( Mesh m,
std::istream &  input 
)

Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller.

Definition at line 3515 of file fespace.cpp.

int mfem::FiniteElementSpace::MakeDofTable ( int  ent_dim,
const Array< int > &  entity_orders,
Table entity_dofs,
Array< char > *  var_ent_order 
)
protected

Build the table var_edge_dofs (or var_face_dofs) in a variable order space; return total edge/face DOFs.

Definition at line 2568 of file fespace.cpp.

void mfem::FiniteElementSpace::MakeVDimMatrix ( SparseMatrix mat) const
protected

Replicate 'mat' in the vector dimension, according to vdim ordering mode.

Definition at line 1206 of file fespace.cpp.

void mfem::FiniteElementSpace::MarkerToList ( const Array< int > &  marker,
Array< int > &  list 
)
static

Convert a Boolean marker array to a list containing all marked indices.

Definition at line 598 of file fespace.cpp.

int mfem::FiniteElementSpace::MinOrder ( VarOrderBits  bits)
staticprotected

Return the minimum order (least significant bit set) in the bit mask.

Definition at line 2441 of file fespace.cpp.

bool mfem::FiniteElementSpace::Nonconforming ( ) const
inline

Definition at line 440 of file fespace.hpp.

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

Copy assignment not supported.

void mfem::FiniteElementSpace::RebuildElementToDofTable ( )

(

Deprecated:
) Use the Update() method if the space or mesh changed.

Definition at line 431 of file fespace.cpp.

SparseMatrix * mfem::FiniteElementSpace::RefinementMatrix ( int  old_ndofs,
const Table old_elem_dof,
const Table old_elem_fos 
)
protected

Calculate explicit GridFunction interpolation matrix (after mesh refinement). NOTE: consider using the RefinementOperator class instead of the fully assembled matrix, which can take a lot of memory.

Definition at line 1476 of file fespace.cpp.

SparseMatrix * mfem::FiniteElementSpace::RefinementMatrix_main ( const int  coarse_ndofs,
const Table coarse_elem_dof,
const Table coarse_elem_fos,
const DenseTensor  localP[] 
) const
protected

This method makes the same assumptions as the method: void GetLocalRefinementMatrices( const FiniteElementSpace &coarse_fes, Geometry::Type geom, DenseTensor &localP) const which is defined below. It also assumes that the coarse fes and this have the same vector dimension, vdim.

TODO: Implement DofTransformation support

Definition at line 1389 of file fespace.cpp.

void mfem::FiniteElementSpace::ReorderElementToDofTable ( )

Reorder the scalar DOFs based on the element ordering.

The new ordering is constructed as follows: 1) loop over all elements as ordered in the Mesh; 2) for each element, assign new indices to all of its current DOFs that are still unassigned; the new indices we assign are simply the sequence 0,1,2,...; if there are any signed DOFs their sign is preserved.

Definition at line 440 of file fespace.cpp.

void mfem::FiniteElementSpace::Save ( std::ostream &  out) const

Save finite element space to output stream out.

Definition at line 3446 of file fespace.cpp.

void mfem::FiniteElementSpace::SetElementOrder ( int  i,
int  p 
)

Sets the order of the i'th finite element.

By default, all elements are assumed to be of fec->GetOrder(). Once SetElementOrder is called, the space becomes a variable order space.

Definition at line 151 of file fespace.cpp.

void mfem::FiniteElementSpace::SetRelaxedHpConformity ( bool  relaxed = true)
inline

In variable order spaces on nonconforming (NC) meshes, this function controls whether strict conformity is enforced in cases where coarse edges/faces have higher polynomial order than their fine NC neighbors. In the default (strict) case, the coarse side polynomial order is reduced to that of the lowest order fine edge/face, so all fine neighbors can interpolate the coarse side exactly. If relaxed == true, some discontinuities in the solution in such cases are allowed and the coarse side is not restricted. For an example, see https://github.com/mfem/mfem/pull/1423#issuecomment-621340392

Definition at line 915 of file fespace.hpp.

void mfem::FiniteElementSpace::SetUpdateOperatorOwner ( bool  own)
inline

Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator() must be deleted outside the FiniteElementSpace.

The update operator ownership is automatically reset to true when a new update operator is created by the Update() method.

Definition at line 883 of file fespace.hpp.

void mfem::FiniteElementSpace::SetUpdateOperatorType ( Operator::Type  tid)
inline

Specify the Operator::Type to be used by the update operators.

The default type is Operator::ANY_TYPE which leaves the choice to this class. The other currently supported option is Operator::MFEM_SPARSEMAT which is only guaranteed to be honored for a refinement update operator. Any other type will be treated as Operator::ANY_TYPE.

Note
This operation destroys the current update operator (if owned).

Definition at line 891 of file fespace.hpp.

NURBSExtension * mfem::FiniteElementSpace::StealNURBSext ( )

Definition at line 2231 of file fespace.cpp.

void mfem::FiniteElementSpace::Update ( bool  want_transform = true)
virtual

Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation operator (unless want_transform is false). Safe to call multiple times, does nothing if space already up to date.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 3336 of file fespace.cpp.

void mfem::FiniteElementSpace::UpdateElementOrders ( )
protected

Resize the elem_order array on mesh change.

Definition at line 3314 of file fespace.cpp.

void mfem::FiniteElementSpace::UpdateMeshPointer ( Mesh new_mesh)
protectedvirtual

Updates the internal mesh pointer.

Warning
new_mesh must be topologically identical to the existing mesh. Used if the address of the Mesh object has changed, e.g. in Mesh::Swap.

Definition at line 3441 of file fespace.cpp.

void mfem::FiniteElementSpace::UpdateNURBS ( )
protected

Definition at line 2242 of file fespace.cpp.

virtual void mfem::FiniteElementSpace::UpdatesFinished ( )
inlinevirtual

Free the GridFunction update operator (if any), to save memory.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 894 of file fespace.hpp.

int mfem::FiniteElementSpace::VDofToDof ( int  vdof) const
inline

Definition at line 687 of file fespace.hpp.

Friends And Related Function Documentation

friend class InterpolationGridTransfer
friend

Definition at line 90 of file fespace.hpp.

friend class LORBase
friend

Definition at line 93 of file fespace.hpp.

void Mesh::Swap ( Mesh ,
bool   
)
friend
friend class PRefinementTransferOperator
friend

Definition at line 91 of file fespace.hpp.

Member Data Documentation

int* mfem::FiniteElementSpace::bdofs
protected

internal DOFs of elements if mixed/var-order; NULL otherwise

Definition at line 119 of file fespace.hpp.

Table* mfem::FiniteElementSpace::bdr_elem_dof
mutableprotected

Definition at line 133 of file fespace.hpp.

Table* mfem::FiniteElementSpace::bdr_elem_fos
mutableprotected

Definition at line 134 of file fespace.hpp.

SparseMatrix* mfem::FiniteElementSpace::cP
mutableprotected

Matrix representing the prolongation from the global conforming dofs to a set of intermediate partially conforming dofs, e.g. the dofs associated with a "cut" space on a non-conforming mesh.

Definition at line 149 of file fespace.hpp.

bool mfem::FiniteElementSpace::cP_is_set
mutableprotected

Definition at line 154 of file fespace.hpp.

SparseMatrix* mfem::FiniteElementSpace::cR
mutableprotected

Conforming restriction matrix such that cR.cP=I.

Definition at line 151 of file fespace.hpp.

SparseMatrix* mfem::FiniteElementSpace::cR_hp
mutableprotected

A version of the conforming restriction matrix for variable-order spaces.

Definition at line 153 of file fespace.hpp.

Array<int> mfem::FiniteElementSpace::dof_elem_array
protected

Definition at line 137 of file fespace.hpp.

Array<int> mfem::FiniteElementSpace::dof_ldof_array
protected

Definition at line 137 of file fespace.hpp.

Array<DofTransformation*> mfem::FiniteElementSpace::DoFTrans
protected

Definition at line 143 of file fespace.hpp.

Array<FaceQuadratureInterpolator*> mfem::FiniteElementSpace::E2BFQ_array
mutableprotected

Definition at line 178 of file fespace.hpp.

Array<FaceQuadratureInterpolator*> mfem::FiniteElementSpace::E2IFQ_array
mutableprotected

Definition at line 177 of file fespace.hpp.

Array<QuadratureInterpolator*> mfem::FiniteElementSpace::E2Q_array
mutableprotected

Definition at line 176 of file fespace.hpp.

Table* mfem::FiniteElementSpace::elem_dof
mutableprotected

Definition at line 131 of file fespace.hpp.

Table* mfem::FiniteElementSpace::elem_fos
mutableprotected

Definition at line 132 of file fespace.hpp.

Array<char> mfem::FiniteElementSpace::elem_order
protected

Polynomial order for each element. If empty, all elements are assumed to be of the default order (fec->GetOrder()).

Definition at line 115 of file fespace.hpp.

Table* mfem::FiniteElementSpace::face_dof
mutableprotected

Definition at line 135 of file fespace.hpp.

Array<int> mfem::FiniteElementSpace::face_to_be
mutableprotected

Definition at line 141 of file fespace.hpp.

const FiniteElementCollection* mfem::FiniteElementSpace::fec
protected

Associated FE collection (not owned).

Definition at line 100 of file fespace.hpp.

OperatorHandle mfem::FiniteElementSpace::L2E_lex
mutableprotected

Definition at line 160 of file fespace.hpp.

OperatorHandle mfem::FiniteElementSpace::L2E_nat
mutableprotected

The element restriction operators, see GetElementRestriction().

Definition at line 160 of file fespace.hpp.

map_L2F mfem::FiniteElementSpace::L2F
mutableprotected

Definition at line 174 of file fespace.hpp.

constexpr int mfem::FiniteElementSpace::MaxVarOrder = 8*sizeof(VarOrderBits) - 1
staticprotected

Definition at line 213 of file fespace.hpp.

Mesh* mfem::FiniteElementSpace::mesh
protected

The mesh that FE space lives on (not owned).

Definition at line 97 of file fespace.hpp.

long mfem::FiniteElementSpace::mesh_sequence
protected

Mesh sequence number last seen when constructing the space. The space needs updating if Mesh::GetSequence() is larger than this.

Definition at line 186 of file fespace.hpp.

int mfem::FiniteElementSpace::nbdofs
protected

Definition at line 117 of file fespace.hpp.

int mfem::FiniteElementSpace::ndofs
protected

Number of degrees of freedom. Number of unknowns is ndofs * vdim.

Definition at line 111 of file fespace.hpp.

int mfem::FiniteElementSpace::nedofs
protected

Definition at line 117 of file fespace.hpp.

int mfem::FiniteElementSpace::nfdofs
protected

Definition at line 117 of file fespace.hpp.

NURBSExtension* mfem::FiniteElementSpace::NURBSext
protected

Definition at line 139 of file fespace.hpp.

int mfem::FiniteElementSpace::nvdofs
protected

Definition at line 117 of file fespace.hpp.

Ordering::Type mfem::FiniteElementSpace::ordering
protected

Type of ordering of the vector dofs when vdim > 1.

Definition at line 108 of file fespace.hpp.

bool mfem::FiniteElementSpace::orders_changed
protected

True if at least one element order changed (variable-order space only).

Definition at line 189 of file fespace.hpp.

int mfem::FiniteElementSpace::own_ext
protected

Definition at line 140 of file fespace.hpp.

bool mfem::FiniteElementSpace::relaxed_hp
protected

Definition at line 191 of file fespace.hpp.

long mfem::FiniteElementSpace::sequence
protected

Update counter, incremented every time the space is constructed/updated. Used by GridFunctions to check if they are up to date with the space.

Definition at line 182 of file fespace.hpp.

OperatorHandle mfem::FiniteElementSpace::Th
protected

Transformation to apply to GridFunctions after space Update().

Definition at line 157 of file fespace.hpp.

int mfem::FiniteElementSpace::uni_fdof
protected

of single face DOFs if all faces uniform; -1 otherwise

Definition at line 118 of file fespace.hpp.

Table mfem::FiniteElementSpace::var_edge_dofs
protected

Variable order spaces only: DOF assignments for edges and faces, see docs in MakeDofTable. For constant order spaces the tables are empty.

Definition at line 123 of file fespace.hpp.

Array<char> mfem::FiniteElementSpace::var_edge_orders
protected

Additional data for the var_*_dofs tables: individual variant orders (these are basically alternate J arrays for var_edge/face_dofs).

Definition at line 128 of file fespace.hpp.

Table mfem::FiniteElementSpace::var_face_dofs
protected

NOTE: also used for spaces with mixed faces.

Definition at line 124 of file fespace.hpp.

Array<char> mfem::FiniteElementSpace::var_face_orders
protected

Definition at line 128 of file fespace.hpp.

int mfem::FiniteElementSpace::vdim
protected

Vector dimension (number of unknowns per degree of freedom).

Definition at line 103 of file fespace.hpp.

VDofTransformation mfem::FiniteElementSpace::VDoFTrans
mutableprotected

Definition at line 144 of file fespace.hpp.


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