MFEM
v4.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 |
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 |
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... | |
virtual const Operator * | GetProlongationMatrix () const |
The returned Operator is owned by the FiniteElementSpace. More... | |
virtual const SparseMatrix * | GetRestrictionMatrix () 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... | |
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... | |
int | GetVDim () const |
Returns vector dimension. More... | |
int | GetOrder (int i) const |
Returns the order of the i'th finite element. More... | |
int | GetFaceOrder (int i) const |
Returns the order of the i'th face finite element. 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 | 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 i, Array< int > &dofs) const |
Returns indexes of degrees of freedom in array dofs for i'th element. More... | |
virtual void | GetBdrElementDofs (int i, Array< int > &dofs) const |
Returns indexes of degrees of freedom for i'th boundary element. More... | |
virtual void | GetFaceDofs (int i, Array< int > &dofs) const |
void | GetEdgeDofs (int i, Array< int > &dofs) const |
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 | 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 |
void | RebuildElementToDofTable () |
void | ReorderElementToDofTable () |
Reorder the scalar DOFs based on the element ordering. More... | |
void | BuildDofToArrays () |
const Table & | GetElementToDofTable () const |
const Table & | GetBdrElementToDofTable () const |
int | GetElementForDof (int i) const |
int | GetLocalDofForDof (int i) const |
const FiniteElement * | GetFE (int i) const |
Returns pointer to the FiniteElement associated with i'th element. More... | |
const FiniteElement * | GetBE (int i) const |
Returns pointer to the FiniteElement for the i'th boundary element. More... | |
const FiniteElement * | GetFaceElement (int i) const |
const FiniteElement * | GetEdgeElement (int i) const |
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 |
virtual void | GetEssentialTrueDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) |
void | ConvertToConformingVDofs (const Array< int > &dofs, Array< int > &cdofs) |
void | ConvertFromConformingVDofs (const Array< int > &cdofs, Array< int > &dofs) |
SparseMatrix * | D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | H2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes) |
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) |
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 |
Return update counter (see Mesh::sequence) More... | |
void | Save (std::ostream &out) const |
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) |
Protected Member Functions | |
void | UpdateNURBS () |
void | Construct () |
void | Destroy () |
void | BuildElementToDofTable () const |
void | GetEntityDofs (int entity, int index, Array< int > &dofs) const |
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.). More... | |
void | BuildConformingInterpolation () const |
Calculate the cP and cR matrices for a nonconforming mesh. More... | |
void | MakeVDimMatrix (SparseMatrix &mat) const |
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 |
void | Constructor (Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES) |
Help function for constructors + Load(). More... | |
Static Protected Member Functions | |
static int | DecodeDof (int dof, double &sign) |
Helper to remove encoded sign from a DOF. More... | |
static void | AddDependencies (SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I) |
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... | |
int | nvdofs |
int | nedofs |
int | nfdofs |
int | nbdofs |
int * | fdofs |
int * | bdofs |
Table * | elem_dof |
Table * | bdrElem_dof |
Array< int > | dof_elem_array |
Array< int > | dof_ldof_array |
NURBSExtension * | NURBSext |
int | own_ext |
SparseMatrix * | cP |
SparseMatrix * | cR |
Conforming restriction matrix such that cR.cP=I. 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 |
Array< QuadratureInterpolator * > | E2Q_array |
long | sequence |
Friends | |
class | InterpolationGridTransfer |
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of degrees of freedom.
Definition at line 85 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 70 of file fespace.cpp.
|
inline |
Definition at line 257 of file fespace.hpp.
|
inline |
Construct a NURBS FE space based on the given NURBSExtension, ext.
Definition at line 266 of file fespace.hpp.
|
virtual |
Definition at line 1841 of file fespace.cpp.
|
staticprotected |
Definition at line 532 of file fespace.cpp.
|
static |
Definition at line 159 of file fespace.cpp.
|
protected |
Calculate the cP and cR matrices for a nonconforming mesh.
Definition at line 582 of file fespace.cpp.
void mfem::FiniteElementSpace::BuildDofToArrays | ( | ) |
Definition at line 263 of file fespace.cpp.
|
protected |
Definition at line 214 of file fespace.cpp.
|
inline |
Definition at line 278 of file fespace.hpp.
|
protected |
Definition at line 1368 of file fespace.cpp.
|
protected |
Help function for constructors + Load().
Definition at line 1295 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 422 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 414 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 431 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 463 of file fespace.cpp.
|
inlinestaticprotected |
Helper to remove encoded sign from a DOF.
Definition at line 144 of file fespace.hpp.
|
protected |
Calculate GridFunction restriction matrix after mesh derefinement.
Definition at line 1192 of file fespace.cpp.
|
protected |
Definition at line 1846 of file fespace.cpp.
|
staticprotected |
Definition at line 556 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | Array< int > & | dofs, |
int | ndofs = -1 |
||
) | const |
Definition at line 107 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | int | vd, |
Array< int > & | dofs, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 122 of file fespace.cpp.
int mfem::FiniteElementSpace::DofToVDof | ( | int | dof, |
int | vd, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 143 of file fespace.cpp.
|
inline |
Definition at line 361 of file fespace.hpp.
|
inline |
Definition at line 410 of file fespace.hpp.
|
inline |
Definition at line 412 of file fespace.hpp.
|
virtual |
Returns indexes of degrees of freedom for i'th boundary element.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1558 of file fespace.cpp.
|
inline |
Definition at line 482 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i-th boundary element.
Definition at line 407 of file fespace.hpp.
|
inline |
Returns the type of boundary element i.
Definition at line 394 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 178 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetBE | ( | int | i | ) | const |
Returns pointer to the FiniteElement for the i'th boundary element.
Definition at line 1781 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingProlongation | ( | ) | const |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 769 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingRestriction | ( | ) | const |
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition at line 776 of file fespace.cpp.
|
inline |
Definition at line 356 of file fespace.hpp.
void mfem::FiniteElementSpace::GetEdgeDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Returns the indexes of the degrees of freedom for i'th edge including the dofs for the vertices of the edge.
Definition at line 1700 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetEdgeElement | ( | int | i | ) | const |
Definition at line 1830 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1754 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 208 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 190 of file fespace.cpp.
|
virtual |
Returns indexes of degrees of freedom in array dofs for i'th element.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1439 of file fespace.cpp.
|
inline |
Definition at line 484 of file fespace.hpp.
void mfem::FiniteElementSpace::GetElementInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1741 of file fespace.cpp.
void mfem::FiniteElementSpace::GetElementInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 202 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, where the element-restriction is the identity, this method will return NULL.
The returned Operator is owned by the FiniteElementSpace.
Definition at line 789 of file fespace.cpp.
|
inline |
Definition at line 481 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i-th element.
Definition at line 398 of file fespace.hpp.
|
inline |
Returns the transformation defining the i-th element in the user-defined variable ElTr.
Definition at line 403 of file fespace.hpp.
|
inline |
Returns the type of element i.
Definition at line 386 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 172 of file fespace.cpp.
|
inline |
Returns the vertices of element i.
Definition at line 390 of file fespace.hpp.
|
protected |
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition at line 571 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 366 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 297 of file fespace.cpp.
|
virtual |
Returns the indexes of the degrees of freedom for i'th face including the dofs for the edges and the vertices of the face.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1641 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetFaceElement | ( | int | i | ) | const |
Definition at line 1807 of file fespace.cpp.
void mfem::FiniteElementSpace::GetFaceInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1766 of file fespace.cpp.
int mfem::FiniteElementSpace::GetFaceOrder | ( | int | i | ) | const |
Returns the order of the i'th face finite element.
Definition at line 101 of file fespace.cpp.
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 184 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetFE | ( | int | i | ) | const |
Returns pointer to the FiniteElement associated with i'th element.
Definition at line 1541 of file fespace.cpp.
|
protected |
Definition at line 1166 of file fespace.cpp.
|
inline |
Definition at line 485 of file fespace.hpp.
|
protected |
Definition at line 902 of file fespace.cpp.
|
protected |
Definition at line 1267 of file fespace.cpp.
|
inline |
Returns the mesh.
Definition at line 272 of file fespace.hpp.
|
inline |
Returns number of boundary elements in the mesh.
Definition at line 383 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 783 of file fespace.cpp.
|
inline |
Returns number of degrees of freedom.
Definition at line 344 of file fespace.hpp.
|
inline |
Returns number of elements in the mesh.
Definition at line 374 of file fespace.hpp.
|
inline |
Number of all scalar edge-interior dofs.
Definition at line 366 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 380 of file fespace.hpp.
|
inline |
Number of all scalar face-interior dofs.
Definition at line 368 of file fespace.hpp.
|
inline |
Definition at line 434 of file fespace.hpp.
|
inline |
Definition at line 274 of file fespace.hpp.
|
inline |
Definition at line 275 of file fespace.hpp.
|
inline |
Returns number of vertices in the mesh.
Definition at line 371 of file fespace.hpp.
|
inline |
Number of all scalar vertex dofs.
Definition at line 364 of file fespace.hpp.
int mfem::FiniteElementSpace::GetOrder | ( | int | i | ) | const |
Returns the order of the i'th finite element.
Definition at line 95 of file fespace.cpp.
|
inline |
Return the ordering method.
Definition at line 359 of file fespace.hpp.
|
inlinevirtual |
The returned Operator is owned by the FiniteElementSpace.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 288 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 812 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 826 of file fespace.cpp.
|
inlinevirtual |
The returned SparseMatrix is owned by the FiniteElementSpace.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 292 of file fespace.hpp.
|
inline |
Return update counter (see Mesh::sequence)
Definition at line 614 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 1835 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 1876 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 1901 of file fespace.cpp.
|
inlinevirtual |
Return the number of vector true (conforming) dofs.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 350 of file fespace.hpp.
|
inline |
Get the GridFunction update operator.
Definition at line 590 of file fespace.hpp.
|
inline |
Return the update operator in the given OperatorHandle, T.
Definition at line 593 of file fespace.hpp.
|
inline |
Returns vector dimension.
Definition at line 336 of file fespace.hpp.
void mfem::FiniteElementSpace::GetVertexDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1729 of file fespace.cpp.
void mfem::FiniteElementSpace::GetVertexVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 196 of file fespace.cpp.
|
inline |
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition at line 347 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 494 of file fespace.cpp.
|
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 403 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 2088 of file fespace.cpp.
|
protected |
Definition at line 741 of file fespace.cpp.
|
static |
Convert a Boolean marker array to a list containing all marked indices.
Definition at line 385 of file fespace.cpp.
|
inline |
Definition at line 279 of file fespace.hpp.
void mfem::FiniteElementSpace::RebuildElementToDofTable | ( | ) |
Definition at line 236 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 926 of file fespace.cpp.
|
protected |
Definition at line 840 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 243 of file fespace.cpp.
void mfem::FiniteElementSpace::Save | ( | std::ostream & | out | ) | const |
Definition at line 2019 of file fespace.cpp.
|
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 600 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 608 of file fespace.hpp.
NURBSExtension * mfem::FiniteElementSpace::StealNURBSext | ( | ) |
Definition at line 1341 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 1942 of file fespace.cpp.
|
protected |
Definition at line 1352 of file fespace.cpp.
|
inlinevirtual |
Free the GridFunction update operator (if any), to save memory.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 611 of file fespace.hpp.
|
inline |
Definition at line 445 of file fespace.hpp.
|
friend |
Definition at line 87 of file fespace.hpp.
|
protected |
Definition at line 108 of file fespace.hpp.
|
protected |
Definition at line 111 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 121 of file fespace.hpp.
|
mutableprotected |
Definition at line 124 of file fespace.hpp.
|
mutableprotected |
Conforming restriction matrix such that cR.cP=I.
Definition at line 123 of file fespace.hpp.
|
protected |
Definition at line 113 of file fespace.hpp.
|
protected |
Definition at line 113 of file fespace.hpp.
|
mutableprotected |
Definition at line 132 of file fespace.hpp.
|
mutableprotected |
Definition at line 110 of file fespace.hpp.
|
protected |
Definition at line 108 of file fespace.hpp.
|
protected |
Associated FE collection (not owned).
Definition at line 94 of file fespace.hpp.
|
mutableprotected |
Definition at line 130 of file fespace.hpp.
|
mutableprotected |
The element restriction operators, see GetElementRestriction().
Definition at line 130 of file fespace.hpp.
|
protected |
The mesh that FE space lives on (not owned).
Definition at line 91 of file fespace.hpp.
|
protected |
Definition at line 107 of file fespace.hpp.
|
protected |
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition at line 105 of file fespace.hpp.
|
protected |
Definition at line 107 of file fespace.hpp.
|
protected |
Definition at line 107 of file fespace.hpp.
|
protected |
Definition at line 115 of file fespace.hpp.
|
protected |
Definition at line 107 of file fespace.hpp.
|
protected |
Type of ordering of the vector dofs when vdim > 1.
Definition at line 102 of file fespace.hpp.
|
protected |
Definition at line 116 of file fespace.hpp.
|
protected |
Definition at line 134 of file fespace.hpp.
|
protected |
Transformation to apply to GridFunctions after space Update().
Definition at line 127 of file fespace.hpp.
|
protected |
Vector dimension (number of unknowns per degree of freedom).
Definition at line 97 of file fespace.hpp.