MFEM v4.7.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(). | |
FiniteElementSpace (const FiniteElementSpace &orig, Mesh *mesh=NULL, const FiniteElementCollection *fec=NULL) | |
Copy constructor: deep copy all data from orig except the Mesh, the FiniteElementCollection, and some derived data. | |
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. | |
FiniteElementSpace & | operator= (const FiniteElementSpace &)=delete |
Copy assignment not supported. | |
Mesh * | GetMesh () const |
Returns the mesh. | |
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. | |
int | GetElementOrder (int i) const |
Returns the order of the i'th finite element. | |
int | GetMaxElementOrder () const |
Return the maximum polynomial order. | |
bool | IsVariableOrder () const |
Returns true if the space contains elements of varying polynomial orders. | |
const SparseMatrix * | GetConformingProlongation () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
const SparseMatrix * | GetConformingRestriction () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
const SparseMatrix * | GetHpConformingRestriction () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
virtual const Operator * | GetProlongationMatrix () const |
The returned Operator is owned by the FiniteElementSpace. | |
const Operator * | GetRestrictionTransposeOperator () const |
Return an operator that performs the transpose of GetRestrictionOperator. | |
virtual const Operator * | GetRestrictionOperator () const |
An abstract operator that performs the same action as GetRestrictionMatrix. | |
virtual const SparseMatrix * | GetRestrictionMatrix () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
virtual const SparseMatrix * | GetHpRestrictionMatrix () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
const ElementRestrictionOperator * | GetElementRestriction (ElementDofOrdering e_ordering) const |
Return an Operator that converts L-vectors to E-vectors. | |
virtual const FaceRestriction * | GetFaceRestriction (ElementDofOrdering f_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const |
Return an Operator that converts L-vectors to E-vectors on each face. | |
const QuadratureInterpolator * | GetQuadratureInterpolator (const IntegrationRule &ir) const |
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). | |
const QuadratureInterpolator * | GetQuadratureInterpolator (const QuadratureSpace &qs) const |
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). | |
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). | |
int | GetOrder (int i) const |
Returns the polynomial degree of the i'th finite element. | |
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. | |
int | GetVDim () const |
Returns vector dimension. | |
int | GetNDofs () const |
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom. | |
int | GetVSize () const |
Return the number of vector dofs, i.e. GetNDofs() x GetVDim(). | |
virtual int | GetTrueVSize () const |
Return the number of vector true (conforming) dofs. | |
int | GetNConformingDofs () const |
int | GetConformingVSize () const |
Ordering::Type | GetOrdering () const |
Return the ordering method. | |
const FiniteElementCollection * | FEColl () const |
int | GetNVDofs () const |
Number of all scalar vertex dofs. | |
int | GetNEDofs () const |
Number of all scalar edge-interior dofs. | |
int | GetNFDofs () const |
Number of all scalar face-interior dofs. | |
int | GetNV () const |
Returns number of vertices in the mesh. | |
int | GetNE () const |
Returns number of elements in the mesh. | |
int | GetNF () const |
Returns number of faces (i.e. co-dimension 1 entities) in the mesh. | |
int | GetNBE () const |
Returns number of boundary elements in the mesh. | |
int | GetNFbyType (FaceType type) const |
Returns the number of faces according to the requested type. | |
int | GetElementType (int i) const |
Returns the type of element i. | |
void | GetElementVertices (int i, Array< int > &vertices) const |
Returns the vertices of element i. | |
int | GetBdrElementType (int i) const |
Returns the type of boundary element i. | |
ElementTransformation * | GetElementTransformation (int i) const |
Returns ElementTransformation for the i-th element. | |
void | GetElementTransformation (int i, IsoparametricTransformation *ElTr) |
Returns the transformation defining the i-th element in the user-defined variable ElTr. | |
ElementTransformation * | GetBdrElementTransformation (int i) const |
Returns ElementTransformation for the i-th boundary element. | |
int | GetAttribute (int i) const |
int | GetBdrAttribute (int i) const |
void | GetPatchDofs (int patch, Array< int > &dofs) const |
Returns indices of degrees of freedom for NURBS patch index patch. Cartesian ordering is used, for the tensor-product degrees of freedom. | |
MFEM_DEPRECATED void | RebuildElementToDofTable () |
( | |
void | ReorderElementToDofTable () |
Reorder the scalar DOFs based on the element ordering. | |
const Table * | GetElementToFaceOrientationTable () const |
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(). | |
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(). | |
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. | |
void | BuildDofToArrays () |
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof(). | |
int | GetElementForDof (int i) const |
Return the index of the first element that contains dof i. | |
int | GetLocalDofForDof (int i) const |
Return the local dof index in the first element that contains dof i. | |
virtual const FiniteElement * | GetFE (int i) const |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in the mesh object. Note: The method has been updated to abort instead of returning NULL for an empty partition. | |
const FiniteElement * | GetBE (int i) const |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary face in the mesh object. | |
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. | |
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. | |
const FiniteElement * | GetTraceElement (int i, Geometry::Type geom_type) const |
Return the trace element from element 'i' to the given 'geom_type'. | |
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. | |
virtual void | GetEssentialTrueDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const |
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. | |
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. | |
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. | |
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. | |
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. | |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE space. | |
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. | |
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. | |
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. | |
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. | |
const Operator * | GetUpdateOperator () |
Get the GridFunction update operator. | |
void | GetUpdateOperator (OperatorHandle &T) |
Return the update operator in the given OperatorHandle, T. | |
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. | |
void | SetUpdateOperatorType (Operator::Type tid) |
Specify the Operator::Type to be used by the update operators. | |
virtual void | UpdatesFinished () |
Free the GridFunction update operator (if any), to save memory. | |
long | GetSequence () const |
bool | IsDGSpace () const |
Return whether or not the space is discontinuous (L2) | |
void | SetRelaxedHpConformity (bool relaxed=true) |
void | Save (std::ostream &out) const |
Save finite element space to output stream out. | |
FiniteElementCollection * | Load (Mesh *m, std::istream &input) |
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller. | |
virtual | ~FiniteElementSpace () |
Local DoF Access Members | |
These member functions produce arrays of local degree of freedom indices, see ldof. If vdim == 1 these indices can be used to access entries in GridFunction, LinearForm, and BilinearForm objects. If vdim != 1 the corresponding Get*VDofs methods should be used instead or one of the DofToVDof methods could be used to produce the appropriate offsets from these local dofs. | |
DofTransformation * | GetElementDofs (int elem, Array< int > &dofs) const |
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldof vector. See also GetElementVDofs(). | |
virtual void | GetElementDofs (int elem, Array< int > &dofs, DofTransformation &doftrans) const |
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
DofTransformation * | GetBdrElementDofs (int bel, Array< int > &dofs) const |
Returns indices of degrees of freedom for boundary element 'bel'. The returned indices are offsets into an ldof vector. See also GetBdrElementVDofs(). | |
virtual void | GetBdrElementDofs (int bel, Array< int > &dofs, DofTransformation &doftrans) const |
The same as GetBdrElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
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. | |
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. | |
void | GetVertexDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the specified vertices. | |
void | GetElementInteriorDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the interior of the specified element. | |
int | GetNumElementInteriorDofs (int i) const |
Returns the number of degrees of freedom associated with the interior of the specified element. | |
void | GetFaceInteriorDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the interior of the specified face. | |
void | GetEdgeInteriorDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the interior of the specified edge. | |
DoF To VDoF Conversion methods | |
These methods convert between local dof and local vector dof using the appropriate relationship based on the Ordering::Type defined in this FiniteElementSpace object. These methods assume the index set has a range [0, GetNDofs()) which will be mapped to the range [0, GetVSize()). This assumption can be changed in the forward mappings by passing a value for ndofs which differs from that returned by GetNDofs().
| |
void | GetVDofs (int vd, Array< int > &dofs, int ndofs=-1) const |
Returns the indices of all of the VDofs for the specified dimension 'vd'. | |
void | DofsToVDofs (Array< int > &dofs, int ndofs=-1) const |
Compute the full set of vdofs corresponding to each entry in dofs. | |
void | DofsToVDofs (int vd, Array< int > &dofs, int ndofs=-1) const |
Compute the set of vdofs corresponding to each entry in dofs for the given vector index vd. | |
int | DofToVDof (int dof, int vd, int ndofs=-1) const |
Compute a single vdof corresponding to the index dof and the vector index vd. | |
int | VDofToDof (int vdof) const |
Compute the inverse of the Dof to VDof mapping for a single index vdof. | |
Local Vector DoF Access Members | |
These member functions produce arrays of local vector degree of freedom indices, see ldof and vdof. These indices can be used to access entries in GridFunction, LinearForm, and BilinearForm objects regardless of the value of vdim. | |
DofTransformation * | GetElementVDofs (int i, Array< int > &vdofs) const |
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. The returned indices are always ordered byNODES, irrespective of whether the space is byNODES or byVDIM. See also GetElementDofs(). | |
void | GetElementVDofs (int i, Array< int > &vdofs, DofTransformation &doftrans) const |
The same as GetElementVDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
DofTransformation * | GetBdrElementVDofs (int i, Array< int > &vdofs) const |
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See also GetBdrElementDofs(). | |
void | GetBdrElementVDofs (int i, Array< int > &vdofs, DofTransformation &doftrans) const |
The same as GetBdrElementVDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
void | GetPatchVDofs (int i, Array< int > &vdofs) const |
Returns indices of degrees of freedom in vdofs for NURBS patch i. | |
void | GetFaceVDofs (int i, Array< int > &vdofs) 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. | |
void | GetEdgeVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vertices of the edge. | |
void | GetVertexVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the specified vertices. | |
void | GetElementInteriorVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the interior of the specified element. | |
void | GetEdgeInteriorVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the interior of the specified edge. | |
Static Public Member Functions | |
static void | AdjustVDofs (Array< int > &vdofs) |
Remove the orientation information encoded into an array of dofs Some basis function types have a relative orientation associated with degrees of freedom shared between neighboring elements, see ldof for more information. An orientation mismatch is indicated in the dof indices by a negative index value. This method replaces such negative indices with the corresponding positive offsets. | |
static int | EncodeDof (int entity_base, int idx) |
Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes). | |
static int | DecodeDof (int dof) |
Helper to return the DOF associated with a sign encoded DOF. | |
static int | DecodeDof (int dof, real_t &sign) |
Helper to determine the DOF and sign of a sign encoded DOF. | |
static void | MarkerToList (const Array< int > &marker, Array< int > &list) |
Convert a Boolean marker array to a list containing all marked indices. | |
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. | |
Protected Types | |
using | key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues> |
The face restriction operators, see GetFaceRestriction(). | |
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. | |
Protected Member Functions | |
void | UpdateNURBS () |
void | Construct () |
void | Destroy () |
void | ConstructDoFTransArray () |
void | DestroyDoFTransArray () |
void | BuildElementToDofTable () const |
void | BuildBdrElementToDofTable () const |
void | BuildFaceToDofTable () const |
void | BuildNURBSFaceToDofTable () const |
Generates partial face_dof table for a NURBS space. | |
int | GetElementOrderImpl (int i) const |
Return element order: internal version of GetElementOrder without checks. | |
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. | |
int | FindEdgeDof (int edge, int ndof) const |
int | FindFaceDof (int face, int ndof) const |
Similar to FindEdgeDof, but used for mixed meshes too. | |
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). | |
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.). | |
int | GetEntityVDofs (int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID, int variant=0) const |
Helper to get vertex, edge or face VDOFs (entity=0,1,2 resp.). | |
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. | |
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. | |
SparseMatrix * | RefinementMatrix_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 |
SparseMatrix * | RefinementMatrix (int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos) |
SparseMatrix * | DerefinementMatrix (int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos) |
Calculate GridFunction restriction matrix after mesh derefinement. | |
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. | |
void | Constructor (Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES) |
Help function for constructors + Load(). | |
virtual void | UpdateMeshPointer (Mesh *new_mesh) |
void | UpdateElementOrders () |
Resize the elem_order array on mesh change. | |
virtual void | CopyProlongationAndRestriction (const FiniteElementSpace &fes, const Array< int > *perm) |
Copies the prolongation and restriction matrices from fes. | |
Static Protected Member Functions | |
static int | MinOrder (VarOrderBits bits) |
Return the minimum order (least significant bit set) in the bit mask. | |
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). | |
const FiniteElementCollection * | fec |
Associated FE collection (not owned). | |
int | vdim |
Vector dimension (number of unknowns per degree of freedom). | |
Ordering::Type | ordering |
int | ndofs |
Number of degrees of freedom. Number of unknowns is ndofs * vdim. | |
Array< char > | elem_order |
int | nvdofs |
int | nedofs |
int | nfdofs |
int | nbdofs |
int | uni_fdof |
int * | bdofs |
internal DOFs of elements if mixed/var-order; NULL otherwise | |
Table | var_edge_dofs |
Table | var_face_dofs |
NOTE: also used for spaces with mixed faces. | |
Array< char > | var_edge_orders |
Array< char > | var_face_orders |
Table * | elem_dof |
Table * | elem_fos |
Table * | bdr_elem_dof |
Table * | bdr_elem_fos |
Table * | face_dof |
Array< int > | dof_elem_array |
Array< int > | dof_ldof_array |
NURBSExtension * | NURBSext |
int | own_ext |
Array< int > | face_to_be |
Array< StatelessDofTransformation * > | DoFTransArray |
DofTransformation | DoFTrans |
std::unique_ptr< SparseMatrix > | cP |
std::unique_ptr< SparseMatrix > | cR |
Conforming restriction matrix such that cR.cP=I. | |
std::unique_ptr< SparseMatrix > | cR_hp |
A version of the conforming restriction matrix for variable-order spaces. | |
bool | cP_is_set |
std::unique_ptr< Operator > | R_transpose |
Operator computing the action of the transpose of the restriction. | |
OperatorHandle | Th |
Transformation to apply to GridFunctions after space Update(). | |
OperatorHandle | L2E_nat |
The element restriction operators, see GetElementRestriction(). | |
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). | |
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.
The term "degree of freedom", or "dof" for short, can mean different things in different contexts. In MFEM we use "dof" to refer to four closely related types of data; edofs, ldofs, tdofs, and vdofs.
Definition at line 219 of file fespace.hpp.
|
protected |
The face restriction operators, see GetFaceRestriction().
Definition at line 295 of file fespace.hpp.
|
protected |
Definition at line 306 of file fespace.hpp.
|
protected |
Bit-mask representing a set of orders needed by an edge/face.
Definition at line 345 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, and 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 71 of file fespace.cpp.
|
inline |
Definition at line 541 of file fespace.hpp.
|
inline |
Construct a NURBS FE space based on the given NURBSExtension, ext.
Definition at line 550 of file fespace.hpp.
|
virtual |
Definition at line 3282 of file fespace.cpp.
|
staticprotected |
Definition at line 807 of file fespace.cpp.
|
protected |
Definition at line 832 of file fespace.cpp.
|
static |
Remove the orientation information encoded into an array of dofs Some basis function types have a relative orientation associated with degrees of freedom shared between neighboring elements, see ldof for more information. An orientation mismatch is indicated in the dof indices by a negative index value. This method replaces such negative indices with the corresponding positive offsets.
Definition at line 265 of file fespace.cpp.
|
protected |
Definition at line 387 of file fespace.cpp.
|
protected |
Calculate the cP and cR matrices for a nonconforming mesh.
Definition at line 985 of file fespace.cpp.
void mfem::FiniteElementSpace::BuildDofToArrays | ( | ) |
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof().
Definition at line 481 of file fespace.cpp.
|
protected |
Definition at line 346 of file fespace.cpp.
|
protected |
Definition at line 426 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 2325 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 2513 of file fespace.cpp.
|
inline |
Definition at line 565 of file fespace.hpp.
|
protected |
Definition at line 2374 of file fespace.cpp.
|
protected |
Definition at line 2252 of file fespace.cpp.
|
protected |
Help function for constructors + Load().
Definition at line 2198 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 688 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 680 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 98 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 697 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 729 of file fespace.cpp.
|
inlinestatic |
Helper to return the DOF associated with a sign encoded DOF.
Definition at line 1017 of file fespace.hpp.
|
inlinestatic |
Helper to determine the DOF and sign of a sign encoded DOF.
Definition at line 1021 of file fespace.hpp.
|
protected |
Calculate GridFunction restriction matrix after mesh derefinement.
TODO: Implement DofTransformation support
Definition at line 2095 of file fespace.cpp.
|
protected |
Definition at line 3287 of file fespace.cpp.
|
protected |
Definition at line 3340 of file fespace.cpp.
|
staticprotected |
Definition at line 886 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | Array< int > & | dofs, |
int | ndofs = -1 ) const |
Compute the full set of vdofs corresponding to each entry in dofs.
Produces a set of vdofs of length vdim * dofs.Size() corresponding to the entries contained in the dofs array.
The ndofs parameter can be used to indicate the total number of Dofs associated with each component of vdim. If ndofs is -1 (the default value), then the number of Dofs is <determined by the FiniteElementSpace::GetNDofs().
Definition at line 213 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | int | vd, |
Array< int > & | dofs, | ||
int | ndofs = -1 ) const |
Compute the set of vdofs corresponding to each entry in dofs for the given vector index vd.
The ndofs parameter can be used to indicate the total number of Dofs associated with each component of vdim. If ndofs is -1 (the default value), then the number of Dofs is <determined by the FiniteElementSpace::GetNDofs().
Definition at line 228 of file fespace.cpp.
int mfem::FiniteElementSpace::DofToVDof | ( | int | dof, |
int | vd, | ||
int | ndofs = -1 ) const |
Compute a single vdof corresponding to the index dof and the vector index vd.
The ndofs parameter can be used to indicate the total number of Dofs associated with each component of vdim. If ndofs is -1 (the default value), then the number of Dofs is <determined by the FiniteElementSpace::GetNDofs().
Definition at line 249 of file fespace.cpp.
|
inlinestatic |
Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
Definition at line 1013 of file fespace.hpp.
|
inline |
Definition at line 727 of file fespace.hpp.
|
protected |
Search row of a DOF table for a DOF set of size 'ndof', return first DOF.
Definition at line 2682 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 369 of file fespace.hpp.
|
inlineprotected |
Similar to FindEdgeDof, but used for mixed meshes too.
Definition at line 373 of file fespace.hpp.
|
inlineprotected |
Definition at line 376 of file fespace.hpp.
|
inline |
Definition at line 785 of file fespace.hpp.
|
inline |
Definition at line 787 of file fespace.hpp.
DofTransformation * mfem::FiniteElementSpace::GetBdrElementDofs | ( | int | bel, |
Array< int > & | dofs ) const |
Returns indices of degrees of freedom for boundary element 'bel'. The returned indices are offsets into an ldof vector. See also GetBdrElementVDofs().
Definition at line 2950 of file fespace.cpp.
|
virtual |
The same as GetBdrElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 2854 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 1141 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i-th boundary element.
Definition at line 782 of file fespace.hpp.
|
inline |
Returns the type of boundary element i.
Definition at line 769 of file fespace.hpp.
DofTransformation * mfem::FiniteElementSpace::GetBdrElementVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See also GetBdrElementDofs().
Definition at line 303 of file fespace.cpp.
void mfem::FiniteElementSpace::GetBdrElementVDofs | ( | int | i, |
Array< int > & | vdofs, | ||
DofTransformation & | doftrans ) const |
The same as GetBdrElementVDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used.
Definition at line 294 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 3204 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 632 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingProlongation | ( | ) | const |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 1302 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingRestriction | ( | ) | const |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 1309 of file fespace.cpp.
|
inline |
Definition at line 722 of file fespace.hpp.
|
protected |
Definition at line 900 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.
The returned indices are offsets into an ldof vector. See also GetEdgeVDofs().
Definition at line 3046 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 3267 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorDofs | ( | int | i, |
Array< int > & | dofs ) const |
Returns the indices of the degrees of freedom for the interior of the specified edge.
The returned indices are offsets into an ldof vector. See also GetEdgeInteriorVDofs().
Definition at line 3149 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Returns the indices of the degrees of freedom for the interior of the specified edge.
The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See also GetEdgeInteriorDofs().
Definition at line 340 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 2699 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vertices of the edge.
The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See GetEdgeDofs() for more information.
Definition at line 322 of file fespace.cpp.
DofTransformation * mfem::FiniteElementSpace::GetElementDofs | ( | int | elem, |
Array< int > & | dofs ) const |
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldof vector. See also GetElementVDofs().
Definition at line 2846 of file fespace.cpp.
|
virtual |
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 2738 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 1159 of file fespace.hpp.
void mfem::FiniteElementSpace::GetElementInteriorDofs | ( | int | i, |
Array< int > & | dofs ) const |
Returns the indices of the degrees of freedom for the interior of the specified element.
Specifically this refers to degrees of freedom which are not associated with the vertices, edges, or faces of the mesh. This method may be useful in conjunction with schemes which process shared and non-shared degrees of freedom differently such as static condensation.
The returned indices are offsets into an ldof vector. See also GetElementInteriorVDofs().
Definition at line 3104 of file fespace.cpp.
void mfem::FiniteElementSpace::GetElementInteriorVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Returns the indices of the degrees of freedom for the interior of the specified element.
The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See GetElementInteriorDofs() for more information.
Definition at line 334 of file fespace.cpp.
int mfem::FiniteElementSpace::GetElementOrder | ( | int | i | ) | const |
Returns the order of the i'th finite element.
Definition at line 176 of file fespace.cpp.
|
protected |
Return element order: internal version of GetElementOrder without checks.
Definition at line 187 of file fespace.cpp.
const ElementRestrictionOperator * 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 in the E-vector, 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 1336 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 1136 of file fespace.hpp.
|
inline |
Definition at line 1132 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i-th element.
Definition at line 773 of file fespace.hpp.
|
inline |
Returns the transformation defining the i-th element in the user-defined variable ElTr.
Definition at line 778 of file fespace.hpp.
|
inline |
Returns the type of element i.
Definition at line 761 of file fespace.hpp.
DofTransformation * mfem::FiniteElementSpace::GetElementVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. The returned indices are always ordered byNODES, irrespective of whether the space is byNODES or byVDIM. See also GetElementDofs().
Definition at line 287 of file fespace.cpp.
void mfem::FiniteElementSpace::GetElementVDofs | ( | int | i, |
Array< int > & | vdofs, | ||
DofTransformation & | doftrans ) const |
The same as GetElementVDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used.
Definition at line 278 of file fespace.cpp.
|
inline |
Returns the vertices of element i.
Definition at line 765 of file fespace.hpp.
|
protected |
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition at line 951 of file fespace.cpp.
|
protected |
Helper to get vertex, edge or face VDOFs (entity=0,1,2 resp.).
Definition at line 976 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 588 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 514 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.
The returned indices are offsets into an ldof vector. See also GetFaceVDofs().
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 2958 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 3237 of file fespace.cpp.
void mfem::FiniteElementSpace::GetFaceInteriorDofs | ( | int | i, |
Array< int > & | dofs ) const |
Returns the indices of the degrees of freedom for the interior of the specified face.
Specifically this refers to degrees of freedom which are not associated with the vertices, edges, or cell interiors of the mesh. This method may be useful in conjunction with schemes which process shared and non-shared degrees of freedom differently such as static condensation.
The returned indices are offsets into an ldof vector. See also GetFaceInteriorVDofs().
Definition at line 3125 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 2710 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 1433 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 1369 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 1149 of file fespace.hpp.
void mfem::FiniteElementSpace::GetFaceVDofs | ( | int | i, |
Array< int > & | vdofs ) 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.
The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See GetFaceDofs() for more information.
Definition at line 316 of file fespace.cpp.
|
virtual |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in the mesh object. Note: The method has been updated to abort instead of returning NULL for an empty partition.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 3168 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 1317 of file fespace.cpp.
|
inlinevirtual |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 624 of file fespace.hpp.
|
protected |
Definition at line 2071 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 1163 of file fespace.hpp.
|
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 2171 of file fespace.cpp.
|
protected |
Definition at line 1528 of file fespace.cpp.
|
inline |
Return the maximum polynomial order.
Definition at line 577 of file fespace.hpp.
|
inline |
Returns the mesh.
Definition at line 559 of file fespace.hpp.
|
inline |
Returns number of boundary elements in the mesh.
Definition at line 749 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 1330 of file fespace.cpp.
|
inline |
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition at line 710 of file fespace.hpp.
|
inline |
Returns number of elements in the mesh.
Definition at line 740 of file fespace.hpp.
|
inline |
Number of all scalar edge-interior dofs.
Definition at line 732 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 746 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 757 of file fespace.hpp.
|
inline |
Number of all scalar face-interior dofs.
Definition at line 734 of file fespace.hpp.
|
protected |
Definition at line 942 of file fespace.cpp.
int mfem::FiniteElementSpace::GetNumElementInteriorDofs | ( | int | i | ) | const |
Returns the number of degrees of freedom associated with the interior of the specified element.
See GetElementInteriorDofs() for more information or to obtain the relevant indices.
Definition at line 3119 of file fespace.cpp.
|
inline |
Definition at line 562 of file fespace.hpp.
|
inline |
Definition at line 561 of file fespace.hpp.
|
inline |
Returns number of vertices in the mesh.
Definition at line 737 of file fespace.hpp.
|
protected |
Return number of possible DOF variants for edge/face (var. order spaces).
Definition at line 2726 of file fespace.cpp.
|
inline |
Number of all scalar vertex dofs.
Definition at line 730 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 696 of file fespace.hpp.
|
inline |
Return the ordering method.
Definition at line 725 of file fespace.hpp.
void mfem::FiniteElementSpace::GetPatchDofs | ( | int | patch, |
Array< int > & | dofs ) const |
Returns indices of degrees of freedom for NURBS patch index patch. Cartesian ordering is used, for the tensor-product degrees of freedom.
Definition at line 3161 of file fespace.cpp.
void mfem::FiniteElementSpace::GetPatchVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Returns indices of degrees of freedom in vdofs for NURBS patch i.
Definition at line 310 of file fespace.cpp.
|
inlinevirtual |
The returned Operator is owned by the FiniteElementSpace.
Reimplemented in mfem::ceed::AlgebraicCoarseSpace, mfem::ceed::ParAlgebraicCoarseSpace, and mfem::ParFiniteElementSpace.
Definition at line 597 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 1404 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 1418 of file fespace.cpp.
|
inlinevirtual |
The returned SparseMatrix is owned by the FiniteElementSpace.
Reimplemented in mfem::ceed::AlgebraicCoarseSpace, mfem::ceed::ParAlgebraicCoarseSpace, and mfem::ParFiniteElementSpace.
Definition at line 620 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 616 of file fespace.hpp.
const Operator * mfem::FiniteElementSpace::GetRestrictionTransposeOperator | ( | ) | const |
Return an operator that performs the transpose of GetRestrictionOperator.
The returned operator is owned by the FiniteElementSpace.
For a serial conforming space, this returns NULL, indicating the identity operator.
For a parallel conforming space, this will return a matrix-free (Device)ConformingProlongationOperator.
For a non-conforming mesh this will return a TransposeOperator wrapping the restriction matrix.
Definition at line 1324 of file fespace.cpp.
|
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 1310 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 3276 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 3349 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 3376 of file fespace.cpp.
|
inlinevirtual |
Return the number of vector true (conforming) dofs.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 716 of file fespace.hpp.
|
inline |
Get the GridFunction update operator.
Definition at line 1285 of file fespace.hpp.
|
inline |
Return the update operator in the given OperatorHandle, T.
Definition at line 1288 of file fespace.hpp.
|
inline |
Returns vector dimension.
Definition at line 706 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 can be used to indicate the total number of Dofs associated with each component of vdim. If ndofs is -1 (the default value), then the number of Dofs is determined by the FiniteElementSpace::GetNDofs().
Definition at line 193 of file fespace.cpp.
void mfem::FiniteElementSpace::GetVertexDofs | ( | int | i, |
Array< int > & | dofs ) const |
Returns the indices of the degrees of freedom for the specified vertices.
The returned indices are offsets into an ldof vector. See also GetVertexVDofs().
Definition at line 3094 of file fespace.cpp.
void mfem::FiniteElementSpace::GetVertexVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Returns the indices of the degrees of freedom for the specified vertices.
The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See also GetVertexDofs().
Definition at line 328 of file fespace.cpp.
|
inline |
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition at line 713 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 760 of file fespace.cpp.
|
inline |
Return whether or not the space is discontinuous (L2)
Definition at line 1313 of file fespace.hpp.
|
inline |
Returns true if the space contains elements of varying polynomial orders.
Definition at line 581 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 667 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 3618 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 2629 of file fespace.cpp.
|
protected |
Replicate 'mat' in the vector dimension, according to vdim ordering mode.
Definition at line 1274 of file fespace.cpp.
|
static |
Convert a Boolean marker array to a list containing all marked indices.
Definition at line 648 of file fespace.cpp.
|
staticprotected |
Return the minimum order (least significant bit set) in the bit mask.
Definition at line 2503 of file fespace.cpp.
|
inline |
Definition at line 566 of file fespace.hpp.
|
delete |
Copy assignment not supported.
void mfem::FiniteElementSpace::RebuildElementToDofTable | ( | ) |
(
Definition at line 452 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 1551 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.
TODO: Implement DofTransformation support
Definition at line 1464 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 461 of file fespace.cpp.
void mfem::FiniteElementSpace::Save | ( | std::ostream & | out | ) | const |
Save finite element space to output stream out.
Definition at line 3549 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 149 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 1327 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 1295 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 1303 of file fespace.hpp.
NURBSExtension * mfem::FiniteElementSpace::StealNURBSext | ( | ) |
Definition at line 2290 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 3439 of file fespace.cpp.
|
protected |
Resize the elem_order array on mesh change.
Definition at line 3417 of file fespace.cpp.
|
protectedvirtual |
Updates the internal mesh pointer.
Definition at line 3544 of file fespace.cpp.
|
protected |
Definition at line 2301 of file fespace.cpp.
|
inlinevirtual |
Free the GridFunction update operator (if any), to save memory.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1306 of file fespace.hpp.
|
inline |
Compute the inverse of the Dof to VDof mapping for a single index vdof.
Definition at line 995 of file fespace.hpp.
|
friend |
Definition at line 221 of file fespace.hpp.
|
friend |
Definition at line 224 of file fespace.hpp.
|
friend |
|
friend |
Definition at line 222 of file fespace.hpp.
|
protected |
internal DOFs of elements if mixed/var-order; NULL otherwise
Definition at line 250 of file fespace.hpp.
|
mutableprotected |
Definition at line 264 of file fespace.hpp.
|
mutableprotected |
Definition at line 265 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 280 of file fespace.hpp.
|
mutableprotected |
Definition at line 285 of file fespace.hpp.
|
mutableprotected |
Conforming restriction matrix such that cR.cP=I.
Definition at line 282 of file fespace.hpp.
|
mutableprotected |
A version of the conforming restriction matrix for variable-order spaces.
Definition at line 284 of file fespace.hpp.
|
protected |
Definition at line 268 of file fespace.hpp.
|
protected |
Definition at line 268 of file fespace.hpp.
|
mutableprotected |
Definition at line 275 of file fespace.hpp.
|
protected |
Definition at line 274 of file fespace.hpp.
|
mutableprotected |
Definition at line 311 of file fespace.hpp.
|
mutableprotected |
Definition at line 310 of file fespace.hpp.
|
mutableprotected |
Definition at line 309 of file fespace.hpp.
|
mutableprotected |
Definition at line 262 of file fespace.hpp.
|
mutableprotected |
Definition at line 263 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 246 of file fespace.hpp.
|
mutableprotected |
Definition at line 266 of file fespace.hpp.
|
mutableprotected |
Definition at line 272 of file fespace.hpp.
|
protected |
Associated FE collection (not owned).
Definition at line 231 of file fespace.hpp.
|
protected |
Definition at line 293 of file fespace.hpp.
|
mutableprotected |
The element restriction operators, see GetElementRestriction().
Definition at line 293 of file fespace.hpp.
|
mutableprotected |
Definition at line 307 of file fespace.hpp.
|
staticconstexprprotected |
Definition at line 346 of file fespace.hpp.
|
protected |
The mesh that FE space lives on (not owned).
Definition at line 228 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 319 of file fespace.hpp.
|
protected |
Definition at line 248 of file fespace.hpp.
|
protected |
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition at line 242 of file fespace.hpp.
|
protected |
Definition at line 248 of file fespace.hpp.
|
protected |
Definition at line 248 of file fespace.hpp.
|
protected |
Definition at line 270 of file fespace.hpp.
|
protected |
Definition at line 248 of file fespace.hpp.
|
protected |
Type of ordering of the vector dofs when vdim > 1.
Definition at line 239 of file fespace.hpp.
|
protected |
True if at least one element order changed (variable-order space only).
Definition at line 322 of file fespace.hpp.
|
protected |
Definition at line 271 of file fespace.hpp.
|
mutableprotected |
Operator computing the action of the transpose of the restriction.
Definition at line 287 of file fespace.hpp.
|
protected |
Definition at line 324 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 315 of file fespace.hpp.
|
protected |
Transformation to apply to GridFunctions after space Update().
Definition at line 290 of file fespace.hpp.
|
protected |
Definition at line 249 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 254 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 259 of file fespace.hpp.
|
protected |
NOTE: also used for spaces with mixed faces.
Definition at line 255 of file fespace.hpp.
|
protected |
Definition at line 259 of file fespace.hpp.
|
protected |
Vector dimension (number of unknowns per degree of freedom).
Definition at line 234 of file fespace.hpp.