MFEM
v4.3.0
Finite element discretization library
|
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of degrees of freedom. More...
#include <fespace.hpp>
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... | |
Mesh * | GetMesh () const |
Returns the mesh. More... | |
const NURBSExtension * | GetNURBSext () const |
NURBSExtension * | GetNURBSext () |
NURBSExtension * | StealNURBSext () |
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 SparseMatrix * | GetConformingProlongation () const |
The returned SparseMatrix is owned by the FiniteElementSpace. More... | |
const SparseMatrix * | GetConformingRestriction () const |
The returned SparseMatrix is owned by the FiniteElementSpace. More... | |
const SparseMatrix * | GetHpConformingRestriction () const |
The returned SparseMatrix is owned by the FiniteElementSpace. More... | |
virtual const Operator * | GetProlongationMatrix () const |
The returned Operator is owned by the FiniteElementSpace. More... | |
virtual const Operator * | GetRestrictionTransposeOperator () const |
Return an operator that performs the transpose of GetRestrictionOperator. More... | |
virtual const Operator * | GetRestrictionOperator () const |
An abstract operator that performs the same action as GetRestrictionMatrix. More... | |
virtual const SparseMatrix * | GetRestrictionMatrix () const |
The returned SparseMatrix is owned by the FiniteElementSpace. More... | |
virtual const SparseMatrix * | GetHpRestrictionMatrix () const |
The returned SparseMatrix is owned by the FiniteElementSpace. More... | |
const Operator * | GetElementRestriction (ElementDofOrdering e_ordering) const |
Return an Operator that converts L-vectors to E-vectors. More... | |
virtual const FaceRestriction * | GetFaceRestriction (ElementDofOrdering e_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const |
Return an Operator that converts L-vectors to E-vectors on each face. More... | |
const QuadratureInterpolator * | GetQuadratureInterpolator (const IntegrationRule &ir) const |
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). More... | |
const QuadratureInterpolator * | GetQuadratureInterpolator (const QuadratureSpace &qs) const |
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). More... | |
const FaceQuadratureInterpolator * | GetFaceQuadratureInterpolator (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 FiniteElementCollection * | FEColl () 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... | |
ElementTransformation * | GetElementTransformation (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... | |
ElementTransformation * | GetBdrElementTransformation (int i) const |
Returns ElementTransformation for the i-th boundary element. More... | |
int | GetAttribute (int i) const |
int | GetBdrAttribute (int i) const |
virtual void | GetElementDofs (int elem, Array< int > &dofs) const |
Returns indices of degrees of freedom of element 'elem'. More... | |
virtual void | GetBdrElementDofs (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 |
void | GetElementVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom in array dofs for i'th element. More... | |
void | GetBdrElementVDofs (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 Table & | GetElementToDofTable () 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 Table & | GetBdrElementToDofTable () 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 Table & | GetFaceToDofTable () 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 FiniteElement * | GetFE (int i) const |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in the mesh object. More... | |
const FiniteElement * | GetBE (int i) const |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary face in the mesh object. More... | |
const FiniteElement * | 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. More... | |
const FiniteElement * | GetEdgeElement (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 FiniteElement * | GetTraceElement (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... | |
SparseMatrix * | D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
Generate the global restriction matrix from a discontinuous FE space to the continuous FE space of the same polynomial degree. More... | |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE space. More... | |
SparseMatrix * | 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. 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 Operator * | GetUpdateOperator () |
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... | |
FiniteElementCollection * | Load (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 | 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... | |
SparseMatrix * | RefinementMatrix_main (const int coarse_ndofs, const Table &coarse_elem_dof, const DenseTensor localP[]) const |
void | GetLocalRefinementMatrices (Geometry::Type geom, DenseTensor &localP) const |
void | GetLocalDerefinementMatrices (Geometry::Type geom, DenseTensor &localR) const |
SparseMatrix * | RefinementMatrix (int old_ndofs, const Table *old_elem_dof) |
SparseMatrix * | DerefinementMatrix (int old_ndofs, const Table *old_elem_dof) |
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 | |
Mesh * | mesh |
The mesh that FE space lives on (not owned). More... | |
const FiniteElementCollection * | fec |
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 otherwiseMore... | |
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 |
Table * | elem_dof |
Table * | bdr_elem_dof |
Table * | face_dof |
Array< int > | dof_elem_array |
Array< int > | dof_ldof_array |
NURBSExtension * | NURBSext |
int | own_ext |
Array< int > | face_to_be |
SparseMatrix * | cP |
SparseMatrix * | cR |
Conforming restriction matrix such that cR.cP=I. More... | |
SparseMatrix * | cR_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) |
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of degrees of freedom.
Definition at line 87 of file fespace.hpp.
|
protected |
The face restriction operators, see GetFaceRestriction().
Definition at line 156 of file fespace.hpp.
|
protected |
Definition at line 167 of file fespace.hpp.
|
protected |
Bit-mask representing a set of orders needed by an edge/face.
Definition at line 203 of file fespace.hpp.
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.
Definition at line 69 of file fespace.cpp.
|
inline |
Definition at line 395 of file fespace.hpp.
|
inline |
Construct a NURBS FE space based on the given NURBSExtension, ext.
Definition at line 404 of file fespace.hpp.
|
virtual |
Definition at line 2764 of file fespace.cpp.
|
staticprotected |
Definition at line 700 of file fespace.cpp.
|
protected |
Definition at line 725 of file fespace.cpp.
|
static |
Definition at line 249 of file fespace.cpp.
|
protected |
Definition at line 327 of file fespace.cpp.
|
protected |
Calculate the cP and cR matrices for a nonconforming mesh.
Definition at line 869 of file fespace.cpp.
void mfem::FiniteElementSpace::BuildDofToArrays | ( | ) |
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof().
Definition at line 402 of file fespace.cpp.
|
protected |
Definition at line 304 of file fespace.cpp.
|
protected |
Definition at line 349 of file fespace.cpp.
|
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 1889 of file fespace.cpp.
|
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 2071 of file fespace.cpp.
|
inline |
Definition at line 416 of file fespace.hpp.
|
protected |
Definition at line 1938 of file fespace.cpp.
|
protected |
Help function for constructors + Load().
Definition at line 1807 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 580 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 572 of file fespace.cpp.
|
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 96 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 589 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 621 of file fespace.cpp.
|
inlinestaticprotected |
Helpers to remove encoded sign from a DOF.
Definition at line 245 of file fespace.hpp.
|
inlinestaticprotected |
Definition at line 248 of file fespace.hpp.
|
protected |
Calculate GridFunction restriction matrix after mesh derefinement.
Definition at line 1709 of file fespace.cpp.
|
protected |
Definition at line 2769 of file fespace.cpp.
|
staticprotected |
Definition at line 780 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | Array< int > & | dofs, |
int | ndofs = -1 |
||
) | const |
Definition at line 197 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | int | vd, |
Array< int > & | dofs, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 212 of file fespace.cpp.
int mfem::FiniteElementSpace::DofToVDof | ( | int | dof, |
int | vd, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 233 of file fespace.cpp.
|
inlinestaticprotected |
Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
Definition at line 241 of file fespace.hpp.
|
inline |
Definition at line 554 of file fespace.hpp.
|
protected |
Search row of a DOF table for a DOF set of size 'ndof', return first DOF.
Definition at line 2242 of file fespace.cpp.
|
inlineprotected |
In a variable order space, return edge DOFs associated with a polynomial order that has 'ndof' degrees of freedom.
Definition at line 227 of file fespace.hpp.
|
inlineprotected |
Similar to FindEdgeDof, but used for mixed meshes too.
Definition at line 231 of file fespace.hpp.
|
inlineprotected |
Definition at line 234 of file fespace.hpp.
|
inline |
Definition at line 612 of file fespace.hpp.
|
inline |
Definition at line 614 of file fespace.hpp.
|
virtual |
Returns indices of degrees of freedom for boundary element 'bel'.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 2418 of file fespace.cpp.
|
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 705 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i-th boundary element.
Definition at line 609 of file fespace.hpp.
|
inline |
Returns the type of boundary element i.
Definition at line 596 of file fespace.hpp.
void mfem::FiniteElementSpace::GetBdrElementVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Returns indexes of degrees of freedom for i'th boundary element.
Definition at line 268 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 2686 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 524 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingProlongation | ( | ) | const |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 1168 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingRestriction | ( | ) | const |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 1175 of file fespace.cpp.
|
inline |
Definition at line 549 of file fespace.hpp.
|
protected |
Definition at line 794 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.
Definition at line 2571 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 2749 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 2650 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 298 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 2259 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 280 of file fespace.cpp.
|
virtual |
Returns indices of degrees of freedom of element 'elem'.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 2298 of file fespace.cpp.
|
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 723 of file fespace.hpp.
void mfem::FiniteElementSpace::GetElementInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 2629 of file fespace.cpp.
void mfem::FiniteElementSpace::GetElementInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 292 of file fespace.cpp.
int mfem::FiniteElementSpace::GetElementOrder | ( | int | i | ) | const |
Returns the order of the i'th finite element.
Definition at line 160 of file fespace.cpp.
|
protected |
Return element order: internal version of GetElementOrder without checks.
Definition at line 171 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 1195 of file fespace.cpp.
|
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 700 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i-th element.
Definition at line 600 of file fespace.hpp.
|
inline |
Returns the transformation defining the i-th element in the user-defined variable ElTr.
Definition at line 605 of file fespace.hpp.
|
inline |
Returns the type of element i.
Definition at line 588 of file fespace.hpp.
void 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 262 of file fespace.cpp.
|
inline |
Returns the vertices of element i.
Definition at line 592 of file fespace.hpp.
|
protected |
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition at line 844 of file fespace.cpp.
|
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 506 of file fespace.cpp.
|
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 437 of file fespace.cpp.
|
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.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 2490 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 2719 of file fespace.cpp.
void mfem::FiniteElementSpace::GetFaceInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 2662 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 2270 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 1285 of file fespace.cpp.
|
virtual |
Return an Operator that converts L-vectors to E-vectors on each face.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1228 of file fespace.cpp.
|
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.
Definition at line 713 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 274 of file fespace.cpp.
|
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 2388 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 1182 of file fespace.cpp.
|
inlinevirtual |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 468 of file fespace.hpp.
|
protected |
Definition at line 1685 of file fespace.cpp.
|
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 727 of file fespace.hpp.
|
protected |
Definition at line 1378 of file fespace.cpp.
|
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 1780 of file fespace.cpp.
|
inline |
Return the maximum polynomial order.
Definition at line 428 of file fespace.hpp.
|
inline |
Returns the mesh.
Definition at line 410 of file fespace.hpp.
|
inline |
Returns number of boundary elements in the mesh.
Definition at line 576 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 1189 of file fespace.cpp.
|
inline |
Returns number of degrees of freedom.
Definition at line 537 of file fespace.hpp.
|
inline |
Returns number of elements in the mesh.
Definition at line 567 of file fespace.hpp.
|
inline |
Number of all scalar edge-interior dofs.
Definition at line 559 of file fespace.hpp.
|
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 573 of file fespace.hpp.
|
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 584 of file fespace.hpp.
|
inline |
Number of all scalar face-interior dofs.
Definition at line 561 of file fespace.hpp.
|
protected |
Definition at line 836 of file fespace.cpp.
int mfem::FiniteElementSpace::GetNumElementInteriorDofs | ( | int | i | ) | const |
Definition at line 2644 of file fespace.cpp.
|
inline |
Definition at line 412 of file fespace.hpp.
|
inline |
Definition at line 413 of file fespace.hpp.
|
inline |
Returns number of vertices in the mesh.
Definition at line 564 of file fespace.hpp.
|
protected |
Return number of possible DOF variants for edge/face (var. order spaces).
Definition at line 2286 of file fespace.cpp.
|
inline |
Number of all scalar vertex dofs.
Definition at line 557 of file fespace.hpp.
|
inline |
Returns the polynomial degree of the i'th finite element.
NOTE: it is recommended to use GetElementOrder in new code.
Definition at line 524 of file fespace.hpp.
|
inline |
Return the ordering method.
Definition at line 552 of file fespace.hpp.
|
inlinevirtual |
The returned Operator is owned by the FiniteElementSpace.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 448 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 1256 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 1270 of file fespace.cpp.
|
inlinevirtual |
The returned SparseMatrix is owned by the FiniteElementSpace.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 464 of file fespace.hpp.
|
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 460 of file fespace.hpp.
|
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 454 of file fespace.hpp.
|
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 872 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 2759 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 2817 of file fespace.cpp.
|
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 2842 of file fespace.cpp.
|
inlinevirtual |
Return the number of vector true (conforming) dofs.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 543 of file fespace.hpp.
|
inline |
Get the GridFunction update operator.
Definition at line 847 of file fespace.hpp.
|
inline |
Return the update operator in the given OperatorHandle, T.
Definition at line 850 of file fespace.hpp.
|
inline |
Returns vector dimension.
Definition at line 534 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 177 of file fespace.cpp.
void mfem::FiniteElementSpace::GetVertexDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 2619 of file fespace.cpp.
void mfem::FiniteElementSpace::GetVertexVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 286 of file fespace.cpp.
|
inline |
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition at line 540 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 652 of file fespace.cpp.
|
inline |
Return whether or not the space is discontinuous (L2)
Definition at line 875 of file fespace.hpp.
|
inline |
Returns true if the space contains elements of varying polynomial orders.
Definition at line 432 of file fespace.hpp.
|
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 559 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 3077 of file fespace.cpp.
|
protected |
Build the table var_edge_dofs (or var_face_dofs) in a variable order space; return total edge/face DOFs.
Definition at line 2187 of file fespace.cpp.
|
protected |
Replicate 'mat' in the vector dimension, according to vdim ordering mode.
Definition at line 1140 of file fespace.cpp.
|
static |
Convert a Boolean marker array to a list containing all marked indices.
Definition at line 540 of file fespace.cpp.
|
staticprotected |
Return the minimum order (least significant bit set) in the bit mask.
Definition at line 2060 of file fespace.cpp.
|
inline |
Definition at line 417 of file fespace.hpp.
void mfem::FiniteElementSpace::RebuildElementToDofTable | ( | ) |
(
Definition at line 375 of file fespace.cpp.
|
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 1401 of file fespace.cpp.
|
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.
Definition at line 1316 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 382 of file fespace.cpp.
void mfem::FiniteElementSpace::Save | ( | std::ostream & | out | ) | const |
Save finite element space to output stream out.
Definition at line 3008 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 133 of file fespace.cpp.
|
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 889 of file fespace.hpp.
|
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 857 of file fespace.hpp.
|
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.
Definition at line 865 of file fespace.hpp.
NURBSExtension * mfem::FiniteElementSpace::StealNURBSext | ( | ) |
Definition at line 1854 of file fespace.cpp.
|
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 2905 of file fespace.cpp.
|
protected |
Resize the elem_order array on mesh change.
Definition at line 2883 of file fespace.cpp.
|
protectedvirtual |
Updates the internal mesh pointer.
Definition at line 3003 of file fespace.cpp.
|
protected |
Definition at line 1865 of file fespace.cpp.
|
inlinevirtual |
Free the GridFunction update operator (if any), to save memory.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 868 of file fespace.hpp.
|
inline |
Definition at line 663 of file fespace.hpp.
|
friend |
Definition at line 89 of file fespace.hpp.
|
friend |
Definition at line 92 of file fespace.hpp.
|
friend |
|
friend |
Definition at line 90 of file fespace.hpp.
|
protected |
internal DOFs of elements if mixed/var-order; NULL otherwise
Definition at line 118 of file fespace.hpp.
|
mutableprotected |
Definition at line 131 of file fespace.hpp.
|
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 143 of file fespace.hpp.
|
mutableprotected |
Definition at line 148 of file fespace.hpp.
|
mutableprotected |
Conforming restriction matrix such that cR.cP=I.
Definition at line 145 of file fespace.hpp.
|
mutableprotected |
A version of the conforming restriction matrix for variable-order spaces.
Definition at line 147 of file fespace.hpp.
|
protected |
Definition at line 134 of file fespace.hpp.
|
protected |
Definition at line 134 of file fespace.hpp.
|
mutableprotected |
Definition at line 172 of file fespace.hpp.
|
mutableprotected |
Definition at line 171 of file fespace.hpp.
|
mutableprotected |
Definition at line 170 of file fespace.hpp.
|
mutableprotected |
Definition at line 130 of file fespace.hpp.
|
protected |
Polynomial order for each element. If empty, all elements are assumed to be of the default order (fec->GetOrder()).
Definition at line 114 of file fespace.hpp.
|
mutableprotected |
Definition at line 132 of file fespace.hpp.
|
mutableprotected |
Definition at line 138 of file fespace.hpp.
|
protected |
Associated FE collection (not owned).
Definition at line 99 of file fespace.hpp.
|
mutableprotected |
Definition at line 154 of file fespace.hpp.
|
mutableprotected |
The element restriction operators, see GetElementRestriction().
Definition at line 154 of file fespace.hpp.
|
mutableprotected |
Definition at line 168 of file fespace.hpp.
|
staticprotected |
Definition at line 204 of file fespace.hpp.
|
protected |
The mesh that FE space lives on (not owned).
Definition at line 96 of file fespace.hpp.
|
protected |
Mesh sequence number last seen when constructing the space. The space needs updating if Mesh::GetSequence() is larger than this.
Definition at line 180 of file fespace.hpp.
|
protected |
Definition at line 116 of file fespace.hpp.
|
protected |
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition at line 110 of file fespace.hpp.
|
protected |
Definition at line 116 of file fespace.hpp.
|
protected |
Definition at line 116 of file fespace.hpp.
|
protected |
Definition at line 136 of file fespace.hpp.
|
protected |
Definition at line 116 of file fespace.hpp.
|
protected |
Type of ordering of the vector dofs when vdim > 1.
Definition at line 107 of file fespace.hpp.
|
protected |
True if at least one element order changed (variable-order space only).
Definition at line 183 of file fespace.hpp.
|
protected |
Definition at line 137 of file fespace.hpp.
|
protected |
Definition at line 185 of file fespace.hpp.
|
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 176 of file fespace.hpp.
|
protected |
Transformation to apply to GridFunctions after space Update().
Definition at line 151 of file fespace.hpp.
|
protected |
Definition at line 117 of file fespace.hpp.
|
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 122 of file fespace.hpp.
|
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 127 of file fespace.hpp.
|
protected |
NOTE: also used for spaces with mixed faces.
Definition at line 123 of file fespace.hpp.
|
protected |
Definition at line 127 of file fespace.hpp.
|
protected |
Vector dimension (number of unknowns per degree of freedom).
Definition at line 102 of file fespace.hpp.