MFEM
v3.2
Finite element discretization library
|
#include <fespace.hpp>
Public Member Functions | |
FiniteElementSpace (Mesh *mesh, const FiniteElementCollection *fec, 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 () |
Ordering::Type | 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, 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 | 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... | |
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) |
SparseMatrix * | D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | H2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes) |
virtual void | Update (bool want_transform=true) |
const Operator * | GetUpdateOperator () |
Get the GridFunction update matrix. 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... | |
virtual void | UpdatesFinished () |
Free GridFunction transformation matrix (if any), to save memory. More... | |
long | GetSequence () const |
Return update counter (see Mesh::sequence) 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 | Construct () |
void | Destroy () |
void | BuildElementToDofTable () const |
void | GetEdgeFaceDofs (int type, int index, Array< int > &dofs) |
void | GetConformingInterpolation () |
Calculate the cP and cR matrices for a nonconforming mesh. More... | |
void | MakeVDimMatrix (SparseMatrix &mat) const |
SparseMatrix * | RefinementMatrix (int old_ndofs, const Table *old_elem_dof) |
Calculate GridFunction interpolation matrix after mesh refinement. More... | |
void | GetLocalDerefinementMatrices (int geom, const CoarseFineTransformations &dt, DenseTensor &localR) |
SparseMatrix * | DerefinementMatrix (int old_ndofs, const Table *old_elem_dof) |
Calculate GridFunction restriction matrix after mesh derefinement. More... | |
Protected Attributes | |
Mesh * | mesh |
The mesh that FE space lives on. More... | |
const FiniteElementCollection * | fec |
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 are 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 |
Operator * | T |
Transformation to apply to GridFunctions after space Update(). More... | |
bool | own_T |
long | sequence |
Class FiniteElementSpace - responsible for providing FEM view of the mesh (mainly managing the set of degrees of freedom).
Definition at line 60 of file fespace.hpp.
mfem::FiniteElementSpace::FiniteElementSpace | ( | Mesh * | mesh, |
const FiniteElementCollection * | fec, | ||
int | vdim = 1 , |
||
int | ordering = Ordering::byNODES |
||
) |
Definition at line 872 of file fespace.cpp.
|
virtual |
Definition at line 1393 of file fespace.cpp.
|
static |
Definition at line 119 of file fespace.cpp.
void mfem::FiniteElementSpace::BuildDofToArrays | ( | ) |
Definition at line 203 of file fespace.cpp.
|
protected |
Definition at line 174 of file fespace.cpp.
|
inline |
Definition at line 141 of file fespace.hpp.
|
protected |
Definition at line 950 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 328 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 320 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 337 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 369 of file fespace.cpp.
|
protected |
Calculate GridFunction restriction matrix after mesh derefinement.
Definition at line 815 of file fespace.cpp.
|
protected |
Definition at line 1398 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | Array< int > & | dofs, |
int | ndofs = -1 |
||
) | const |
Definition at line 67 of file fespace.cpp.
void mfem::FiniteElementSpace::DofsToVDofs | ( | int | vd, |
Array< int > & | dofs, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 82 of file fespace.cpp.
int mfem::FiniteElementSpace::DofToVDof | ( | int | dof, |
int | vd, | ||
int | ndofs = -1 |
||
) | const |
Definition at line 103 of file fespace.cpp.
|
inline |
Definition at line 175 of file fespace.hpp.
|
inline |
Definition at line 215 of file fespace.hpp.
|
inline |
Definition at line 217 of file fespace.hpp.
|
virtual |
Returns indexes of degrees of freedom for i'th boundary element.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 1126 of file fespace.cpp.
|
inline |
Definition at line 276 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i'th boundary element.
Definition at line 212 of file fespace.hpp.
|
inline |
Returns the type of boundary element i.
Definition at line 199 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 138 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 1333 of file fespace.cpp.
|
protected |
Calculate the cP and cR matrices for a nonconforming mesh.
Definition at line 487 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingProlongation | ( | ) |
Definition at line 669 of file fespace.cpp.
const SparseMatrix * mfem::FiniteElementSpace::GetConformingRestriction | ( | ) |
Definition at line 676 of file fespace.cpp.
|
inline |
Definition at line 170 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 1268 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetEdgeElement | ( | int | i | ) | const |
Definition at line 1382 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 474 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1321 of file fespace.cpp.
void mfem::FiniteElementSpace::GetEdgeInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 168 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 150 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 1014 of file fespace.cpp.
|
inline |
Definition at line 278 of file fespace.hpp.
void mfem::FiniteElementSpace::GetElementInteriorDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1309 of file fespace.cpp.
void mfem::FiniteElementSpace::GetElementInteriorVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 162 of file fespace.cpp.
|
inline |
Definition at line 275 of file fespace.hpp.
|
inline |
Returns ElementTransformation for the i'th element.
Definition at line 203 of file fespace.hpp.
|
inline |
Returns the transformation defining the i-th element in the user-defined variable.
Definition at line 208 of file fespace.hpp.
|
inline |
Returns the type of element i.
Definition at line 191 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 132 of file fespace.cpp.
|
inline |
Returns the vertices of element i.
Definition at line 195 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 274 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 237 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 1209 of file fespace.cpp.
const FiniteElement * mfem::FiniteElementSpace::GetFaceElement | ( | int | i | ) | const |
Definition at line 1359 of file fespace.cpp.
int mfem::FiniteElementSpace::GetFaceOrder | ( | int | i | ) | const |
Returns the order of the i'th face finite element.
Definition at line 61 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 144 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 1113 of file fespace.cpp.
|
protected |
Definition at line 770 of file fespace.cpp.
|
inline |
Definition at line 279 of file fespace.hpp.
|
inline |
Returns the mesh.
Definition at line 136 of file fespace.hpp.
|
inline |
Returns number of boundary elements in the mesh.
Definition at line 188 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 683 of file fespace.cpp.
|
inline |
Returns number of degrees of freedom.
Definition at line 159 of file fespace.hpp.
|
inline |
Returns number of elements in the mesh.
Definition at line 182 of file fespace.hpp.
|
inline |
Definition at line 178 of file fespace.hpp.
|
inline |
Definition at line 179 of file fespace.hpp.
|
inline |
Definition at line 237 of file fespace.hpp.
|
inline |
Definition at line 138 of file fespace.hpp.
|
inline |
Returns number of nodes in the mesh.
Definition at line 185 of file fespace.hpp.
|
inline |
Definition at line 177 of file fespace.hpp.
int mfem::FiniteElementSpace::GetOrder | ( | int | i | ) | const |
Returns the order of the i'th finite element.
Definition at line 55 of file fespace.cpp.
|
inline |
Return the ordering method.
Definition at line 173 of file fespace.hpp.
|
inlinevirtual |
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 147 of file fespace.hpp.
|
inline |
Return update counter (see Mesh::sequence)
Definition at line 357 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 1387 of file fespace.cpp.
|
inlinevirtual |
Return the number of vector true (conforming) dofs.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 164 of file fespace.hpp.
|
inline |
Get the GridFunction update matrix.
Definition at line 346 of file fespace.hpp.
|
inline |
Returns vector dimension.
Definition at line 151 of file fespace.hpp.
void mfem::FiniteElementSpace::GetVertexDofs | ( | int | i, |
Array< int > & | dofs | ||
) | const |
Definition at line 1297 of file fespace.cpp.
void mfem::FiniteElementSpace::GetVertexVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 156 of file fespace.cpp.
|
inline |
Definition at line 161 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 400 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 309 of file fespace.cpp.
|
protected |
Definition at line 642 of file fespace.cpp.
|
static |
Convert a Boolean marker array to a list containing all marked indices.
Definition at line 292 of file fespace.cpp.
|
inline |
Definition at line 142 of file fespace.hpp.
void mfem::FiniteElementSpace::RebuildElementToDofTable | ( | ) |
Definition at line 196 of file fespace.cpp.
|
protected |
Calculate GridFunction interpolation matrix after mesh refinement.
Definition at line 689 of file fespace.cpp.
void mfem::FiniteElementSpace::Save | ( | std::ostream & | out | ) | const |
Definition at line 1485 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.
Definition at line 351 of file fespace.hpp.
NURBSExtension * mfem::FiniteElementSpace::StealNURBSext | ( | ) |
Definition at line 923 of file fespace.cpp.
|
virtual |
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation matrix (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 1421 of file fespace.cpp.
|
protected |
Definition at line 934 of file fespace.cpp.
|
inlinevirtual |
Free GridFunction transformation matrix (if any), to save memory.
Reimplemented in mfem::ParFiniteElementSpace.
Definition at line 354 of file fespace.hpp.
|
inline |
Definition at line 248 of file fespace.hpp.
|
protected |
Definition at line 80 of file fespace.hpp.
|
protected |
Definition at line 83 of file fespace.hpp.
|
protected |
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 93 of file fespace.hpp.
|
protected |
Definition at line 96 of file fespace.hpp.
|
protected |
Conforming restriction matrix such that cR.cP=I.
Definition at line 95 of file fespace.hpp.
|
protected |
Definition at line 85 of file fespace.hpp.
|
protected |
Definition at line 85 of file fespace.hpp.
|
mutableprotected |
Definition at line 82 of file fespace.hpp.
|
protected |
Definition at line 80 of file fespace.hpp.
|
protected |
Definition at line 66 of file fespace.hpp.
|
protected |
The mesh that FE space lives on.
Definition at line 64 of file fespace.hpp.
|
protected |
Definition at line 79 of file fespace.hpp.
|
protected |
Number of degrees of freedom. Number of unknowns are ndofs*vdim.
Definition at line 77 of file fespace.hpp.
|
protected |
Definition at line 79 of file fespace.hpp.
|
protected |
Definition at line 79 of file fespace.hpp.
|
protected |
Definition at line 87 of file fespace.hpp.
|
protected |
Definition at line 79 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 74 of file fespace.hpp.
|
protected |
Definition at line 88 of file fespace.hpp.
|
protected |
Definition at line 100 of file fespace.hpp.
|
protected |
Definition at line 102 of file fespace.hpp.
|
protected |
Transformation to apply to GridFunctions after space Update().
Definition at line 99 of file fespace.hpp.
|
protected |
Vector dimension (number of unknowns per degree of freedom).
Definition at line 69 of file fespace.hpp.