MFEM
v3.1
Finite element discretization library
|
Abstract finite element space. More...
#include <fespace.hpp>
Public Member Functions | |
FiniteElementSpace (Mesh *m, const FiniteElementCollection *f, int vdim=1, int ordering=Ordering::byNODES) | |
Mesh * | GetMesh () const |
Returns the mesh. More... | |
NURBSExtension * | GetNURBSext () |
NURBSExtension * | StealNURBSext () |
bool | Conforming () const |
bool | Nonconforming () const |
const SparseMatrix * | GetConformingProlongation () |
const SparseMatrix * | GetConformingRestriction () |
virtual const SparseMatrix * | GetRestrictionMatrix () |
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 |
virtual int | GetTrueVSize () |
Return the number of vector true (conforming) dofs. More... | |
int | GetNConformingDofs () |
int | GetConformingVSize () |
int | GetOrdering () const |
Return the ordering method. More... | |
const FiniteElementCollection * | FEColl () const |
int | GetNVDofs () const |
int | GetNEDofs () const |
int | GetNFDofs () const |
int | GetNE () const |
Returns number of elements in the mesh. More... | |
int | GetNV () const |
Returns number of nodes 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) |
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 |
int | GetNumElementInteriorDofs (int i) const |
void | GetEdgeInteriorDofs (int i, Array< int > &dofs) const |
void | DofsToVDofs (Array< int > &dofs) 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 | BuildElementToDofTable () |
void | BuildDofToArrays () |
const Table & | GetElementToDofTable () const |
const Table & | GetBdrElementToDofTable () const |
int | GetElementForDof (int i) |
int | GetLocalDofForDof (int i) |
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, int geom_type) const |
Return the trace element from element 'i' to the given 'geom_type'. More... | |
SparseMatrix * | GlobalRestrictionMatrix (FiniteElementSpace *cfes, int one_vdim=-1) |
virtual void | GetEssentialVDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs) const |
virtual void | GetEssentialTrueDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list) |
void | ConvertToConformingVDofs (const Array< int > &dofs, Array< int > &cdofs) |
void | ConvertFromConformingVDofs (const Array< int > &cdofs, Array< int > &dofs) |
void | EliminateEssentialBCFromGRM (FiniteElementSpace *cfes, Array< int > &bdr_attr_is_ess, SparseMatrix *R) |
SparseMatrix * | GlobalRestrictionMatrix (FiniteElementSpace *cfes, Array< int > &bdr_attr_is_ess, int one_vdim=-1) |
Generate the global restriction matrix with eliminated essential bc. More... | |
SparseMatrix * | D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | H2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes) |
virtual void | Update () |
virtual void | UpdateAndInterpolate (int num_grid_fns,...) |
void | UpdateAndInterpolate (GridFunction *gf) |
A shortcut for passing only one GridFunction to UndateAndInterpolate. More... | |
virtual FiniteElementSpace * | SaveUpdate () |
Return a copy of the current FE space and update. More... | |
void | Save (std::ostream &out) const |
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 | Constructor () |
void | Destructor () |
FiniteElementSpace (FiniteElementSpace &) | |
void | ConstructRefinementData (int k, int cdofs, RefinementType type) |
Constructs new refinement data using coarse element k as a template. More... | |
DenseMatrix * | LocalInterpolation (int k, int cdofs, RefinementType type, Array< int > &rows) |
Generates the local interpolation matrix for coarse element k. More... | |
SparseMatrix * | NC_GlobalRestrictionMatrix (FiniteElementSpace *cfes, NCMesh *ncmesh) |
void | GetEdgeFaceDofs (int type, int index, Array< int > &dofs) |
void | GetConformingInterpolation () |
void | MakeVDimMatrix (SparseMatrix &mat) const |
Protected Attributes | |
Mesh * | mesh |
The mesh that FE space lives on. More... | |
int | vdim |
Vector dimension (number of unknowns per degree of freedom). More... | |
int | ndofs |
Number of degrees of freedom. Number of unknowns are ndofs*vdim. More... | |
int | ordering |
const FiniteElementCollection * | fec |
int | nvdofs |
int | nedofs |
int | nfdofs |
int | nbdofs |
int * | fdofs |
int * | bdofs |
Array< RefinementData * > | RefData |
Collection of currently known refinement data. More... | |
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 |
Abstract finite element space.
Definition at line 62 of file fespace.hpp.
|
protected |
Definition at line 932 of file fespace.cpp.
mfem::FiniteElementSpace::FiniteElementSpace | ( | Mesh * | m, |
const FiniteElementCollection * | f, | ||
int | vdim = 1 , |
||
int | ordering = Ordering::byNODES |
||
) |
Definition at line 965 of file fespace.cpp.
|
virtual |
Definition at line 1479 of file fespace.cpp.
|
static |
Definition at line 148 of file fespace.cpp.
void mfem::FiniteElementSpace::BuildDofToArrays | ( | ) |
Definition at line 225 of file fespace.cpp.
void mfem::FiniteElementSpace::BuildElementToDofTable | ( | ) |
Definition at line 203 of file fespace.cpp.
|
inline |
Definition at line 145 of file fespace.hpp.
|
protected |
Definition at line 1037 of file fespace.cpp.
|
protected |
Constructs new refinement data using coarse element k as a template.
Definition at line 1568 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 523 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 515 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 583 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 615 of file fespace.cpp.
|
protected |
Definition at line 1489 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | Array< int > & | dofs | ) | const |
Definition at line 35 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | int | vd, |
Array< int > & | dofs, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 74 of file fespace.cpp.
int mfem::FiniteElementSpace::DofToVDof | ( | int | dof, |
int | vd, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 116 of file fespace.cpp.
void mfem::FiniteElementSpace::EliminateEssentialBCFromGRM | ( | FiniteElementSpace * | cfes, |
Array< int > & | bdr_attr_is_ess, | ||
SparseMatrix * | R | ||
) |
Definition at line 532 of file fespace.cpp.
|
inline |
Definition at line 178 of file fespace.hpp.
|
inline |
Definition at line 218 of file fespace.hpp.
|
inline |
Definition at line 220 of file fespace.hpp.
|
virtual |
Returns indexes of degrees of freedom for i'th boundary element.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1212 of file fespace.cpp.
|
inline |
Definition at line 279 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i'th boundary element.
Definition at line 215 of file fespace.hpp.
|
inline |
Returns the type of boundary element i.
Definition at line 202 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 167 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 1419 of file fespace.cpp.
|
protected |
Calculate the cP and cR matrices for a nonconforming mesh.
Definition at line 733 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingProlongation | ( | ) |
Definition at line 912 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingRestriction | ( | ) |
Definition at line 919 of file fespace.cpp.
|
inline |
Definition at line 173 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 1354 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetEdgeElement | ( | int | i | ) | const |
Definition at line 1468 of file fespace.cpp.
|
protected |
This is a helper function to get edge (type == 0) or face (type == 1) DOFs. The function is aware of ghost edges/faces in parallel, for which an empty DOF list is returned.
Definition at line 720 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1407 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 197 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 179 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 1100 of file fespace.cpp.
|
inline |
Definition at line 281 of file fespace.hpp.
void mfem::FiniteElementSpace::GetElementInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1395 of file fespace.cpp.
void mfem::FiniteElementSpace::GetElementInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 191 of file fespace.cpp.
|
inline |
Definition at line 278 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i'th element.
Definition at line 206 of file fespace.hpp.
|
inline |
Returns the transformation defining the i-th element in the user-defined variable.
Definition at line 211 of file fespace.hpp.
|
inline |
Returns the type of element i.
Definition at line 194 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 161 of file fespace.cpp.
|
inline |
Returns the vertices of element i.
Definition at line 198 of file fespace.hpp.
|
virtual |
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in the array bdr_attr_is_ess.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 469 of file fespace.cpp.
|
virtual |
Mark degrees of freedom associated with boundary elements with the specified boundary attributes (marked in 'bdr_attr_is_ess').
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 432 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 1295 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetFaceElement | ( | int | i | ) | const |
Definition at line 1445 of file fespace.cpp.
int mfem::FiniteElementSpace::GetFaceOrder | ( | int | i | ) | const |
Returns the order of the i'th face finite element.
Definition at line 29 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 173 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 1199 of file fespace.cpp.
|
inline |
Definition at line 282 of file fespace.hpp.
|
inline |
Returns the mesh.
Definition at line 140 of file fespace.hpp.
|
inline |
Returns number of boundary elements in the mesh.
Definition at line 191 of file fespace.hpp.
int mfem::FiniteElementSpace::GetNConformingDofs | ( | ) |
Returns the number of conforming ("true") degrees of freedom (if the space is on a nonconforming mesh with hanging nodes).
Definition at line 926 of file fespace.cpp.
|
inline |
Returns number of degrees of freedom.
Definition at line 162 of file fespace.hpp.
|
inline |
Returns number of elements in the mesh.
Definition at line 185 of file fespace.hpp.
|
inline |
Definition at line 181 of file fespace.hpp.
|
inline |
Definition at line 182 of file fespace.hpp.
|
inline |
Definition at line 240 of file fespace.hpp.
|
inline |
Definition at line 142 of file fespace.hpp.
|
inline |
Returns number of nodes in the mesh.
Definition at line 188 of file fespace.hpp.
|
inline |
Definition at line 180 of file fespace.hpp.
int mfem::FiniteElementSpace::GetOrder | ( | int | i | ) | const |
Returns the order of the i'th finite element.
Definition at line 23 of file fespace.cpp.
|
inline |
Return the ordering method.
Definition at line 176 of file fespace.hpp.
|
inlinevirtual |
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 150 of file fespace.hpp.
const FiniteElement * mfem::FiniteElementSpace::GetTraceElement | ( | int | i, |
int | geom_type | ||
) | const |
Return the trace element from element 'i' to the given 'geom_type'.
Definition at line 1473 of file fespace.cpp.
|
inlinevirtual |
Return the number of vector true (conforming) dofs.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 167 of file fespace.hpp.
|
inline |
Returns vector dimension.
Definition at line 154 of file fespace.hpp.
void mfem::FiniteElementSpace::GetVertexDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1383 of file fespace.cpp.
void mfem::FiniteElementSpace::GetVertexVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 185 of file fespace.cpp.
|
inline |
Definition at line 164 of file fespace.hpp.
SparseMatrix * mfem::FiniteElementSpace::GlobalRestrictionMatrix | ( | FiniteElementSpace * | cfes, |
int | one_vdim = -1 |
||
) |
Return the restriction matrix from this FE space to the coarse FE space 'cfes'. Both FE spaces must use the same FE collection and be defined on the same Mesh which must be in TWO_LEVEL_* state. When vdim > 1, 'one_vdim' specifies whether the restriction matrix built should be the scalar restriction (one_vdim=1) or the full vector restriction (one_vdim=0); if one_vdim=-1 then the behavior depends on the ordering of this FE space: if ordering=byNodes then the scalar restriction matrix is built and if ordering=byVDim – the full vector restriction matrix.
Definition at line 302 of file fespace.cpp.
SparseMatrix * mfem::FiniteElementSpace::GlobalRestrictionMatrix | ( | FiniteElementSpace * | cfes, |
Array< int > & | bdr_attr_is_ess, | ||
int | one_vdim = -1 |
||
) |
Generate the global restriction matrix with eliminated essential bc.
Definition at line 572 of file fespace.cpp.
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 646 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 504 of file fespace.cpp.
|
protected |
Generates the local interpolation matrix for coarse element k.
Definition at line 250 of file fespace.cpp.
|
protected |
Definition at line 885 of file fespace.cpp.
|
static |
Convert a Boolean marker array to a list containing all marked indices.
Definition at line 487 of file fespace.cpp.
|
protected |
Construct the restriction matrix from the coarse FE space 'cfes' to (*this) space, where both spaces use the same FE collection and their meshes are obtained from different levels of a single NCMesh. (Also, the coarse level must have been marked in 'ncmesh' before refinement).
Definition at line 346 of file fespace.cpp.
|
inline |
Definition at line 146 of file fespace.hpp.
void mfem::FiniteElementSpace::Save | ( | std::ostream & | out | ) | const |
Definition at line 1635 of file fespace.cpp.
|
virtual |
Return a copy of the current FE space and update.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1524 of file fespace.cpp.
NURBSExtension * mfem::FiniteElementSpace::StealNURBSext | ( | ) |
Definition at line 1008 of file fespace.cpp.
|
virtual |
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1511 of file fespace.cpp.
|
virtual |
Updates the space after the underlying mesh has been refined and interpolates one or more GridFunctions so that they represent the same functions on the new mesh. The grid functions are passed as pointers after 'num_grid_fns'.
Definition at line 1531 of file fespace.cpp.
|
inline |
A shortcut for passing only one GridFunction to UndateAndInterpolate.
Definition at line 371 of file fespace.hpp.
|
protected |
Definition at line 1019 of file fespace.cpp.
|
inline |
Definition at line 251 of file fespace.hpp.
|
protected |
Definition at line 81 of file fespace.hpp.
|
protected |
Definition at line 87 of file fespace.hpp.
|
protected |
Definition at line 96 of file fespace.hpp.
|
protected |
Definition at line 98 of file fespace.hpp.
|
protected |
Definition at line 88 of file fespace.hpp.
|
protected |
Definition at line 88 of file fespace.hpp.
|
protected |
Definition at line 86 of file fespace.hpp.
|
protected |
Definition at line 81 of file fespace.hpp.
|
protected |
Definition at line 79 of file fespace.hpp.
|
protected |
The mesh that FE space lives on.
Definition at line 66 of file fespace.hpp.
|
protected |
Definition at line 80 of file fespace.hpp.
|
protected |
Number of degrees of freedom. Number of unknowns are ndofs*vdim.
Definition at line 72 of file fespace.hpp.
|
protected |
Definition at line 80 of file fespace.hpp.
|
protected |
Definition at line 80 of file fespace.hpp.
|
protected |
Definition at line 90 of file fespace.hpp.
|
protected |
Definition at line 80 of file fespace.hpp.
|
protected |
Type of ordering of dofs. Ordering::byNODES - first nodes, then vector dimension, Ordering::byVDIM - first vector dimension, then nodes
Definition at line 77 of file fespace.hpp.
|
protected |
Definition at line 91 of file fespace.hpp.
|
protected |
Collection of currently known refinement data.
Definition at line 84 of file fespace.hpp.
|
protected |
Vector dimension (number of unknowns per degree of freedom).
Definition at line 69 of file fespace.hpp.