MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::FiniteElementSpace Class Reference

#include <fespace.hpp>

Inheritance diagram for mfem::FiniteElementSpace:
[legend]
Collaboration diagram for mfem::FiniteElementSpace:
[legend]

Public Member Functions

 FiniteElementSpace (Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
 
MeshGetMesh () const
 Returns the mesh. More...
 
NURBSExtensionGetNURBSext ()
 
NURBSExtensionStealNURBSext ()
 
bool Conforming () const
 
bool Nonconforming () const
 
const SparseMatrixGetConformingProlongation ()
 
const SparseMatrixGetConformingRestriction ()
 
virtual const OperatorGetProlongationMatrix ()
 
virtual const SparseMatrixGetRestrictionMatrix ()
 
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 FiniteElementCollectionFEColl () 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...
 
ElementTransformationGetElementTransformation (int i) const
 Returns ElementTransformation for the i'th element. More...
 
void GetElementTransformation (int i, IsoparametricTransformation *ElTr)
 
ElementTransformationGetBdrElementTransformation (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 TableGetElementToDofTable () const
 
const TableGetBdrElementToDofTable () const
 
int GetElementForDof (int i) const
 
int GetLocalDofForDof (int i) const
 
const FiniteElementGetFE (int i) const
 Returns pointer to the FiniteElement associated with i'th element. More...
 
const FiniteElementGetBE (int i) const
 Returns pointer to the FiniteElement for the i'th boundary element. More...
 
const FiniteElementGetFaceElement (int i) const
 
const FiniteElementGetEdgeElement (int i) const
 
const FiniteElementGetTraceElement (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)
 
SparseMatrixD2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes)
 
SparseMatrixD2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes)
 
SparseMatrixH2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes)
 
virtual void Update (bool want_transform=true)
 
const OperatorGetUpdateOperator ()
 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
 
SparseMatrixRefinementMatrix (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)
 
SparseMatrixDerefinementMatrix (int old_ndofs, const Table *old_elem_dof)
 Calculate GridFunction restriction matrix after mesh derefinement. More...
 

Protected Attributes

Meshmesh
 The mesh that FE space lives on. More...
 
const FiniteElementCollectionfec
 
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
 
Tableelem_dof
 
TablebdrElem_dof
 
Array< int > dof_elem_array
 
Array< int > dof_ldof_array
 
NURBSExtensionNURBSext
 
int own_ext
 
SparseMatrixcP
 
SparseMatrixcR
 Conforming restriction matrix such that cR.cP=I. More...
 
bool cP_is_set
 
OperatorT
 Transformation to apply to GridFunctions after space Update(). More...
 
bool own_T
 
long sequence
 

Detailed Description

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.

Constructor & Destructor Documentation

mfem::FiniteElementSpace::FiniteElementSpace ( Mesh mesh,
const FiniteElementCollection fec,
int  vdim = 1,
int  ordering = Ordering::byNODES 
)

Definition at line 893 of file fespace.cpp.

mfem::FiniteElementSpace::~FiniteElementSpace ( )
virtual

Definition at line 1429 of file fespace.cpp.

Member Function Documentation

void mfem::FiniteElementSpace::AdjustVDofs ( Array< int > &  vdofs)
static

Definition at line 120 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildDofToArrays ( )

Definition at line 224 of file fespace.cpp.

void mfem::FiniteElementSpace::BuildElementToDofTable ( ) const
protected

Definition at line 175 of file fespace.cpp.

bool mfem::FiniteElementSpace::Conforming ( ) const
inline

Definition at line 141 of file fespace.hpp.

void mfem::FiniteElementSpace::Construct ( )
protected

Definition at line 971 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 349 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 341 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 358 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 390 of file fespace.cpp.

SparseMatrix * mfem::FiniteElementSpace::DerefinementMatrix ( int  old_ndofs,
const Table old_elem_dof 
)
protected

Calculate GridFunction restriction matrix after mesh derefinement.

Definition at line 836 of file fespace.cpp.

void mfem::FiniteElementSpace::Destroy ( )
protected

Definition at line 1434 of file fespace.cpp.

void mfem::FiniteElementSpace::DofsToVDofs ( Array< int > &  dofs,
int  ndofs = -1 
) const

Definition at line 68 of file fespace.cpp.

void mfem::FiniteElementSpace::DofsToVDofs ( int  vd,
Array< int > &  dofs,
int  ndofs = -1 
) const

Definition at line 83 of file fespace.cpp.

int mfem::FiniteElementSpace::DofToVDof ( int  dof,
int  vd,
int  ndofs = -1 
) const

Definition at line 104 of file fespace.cpp.

const FiniteElementCollection* mfem::FiniteElementSpace::FEColl ( ) const
inline

Definition at line 177 of file fespace.hpp.

int mfem::FiniteElementSpace::GetAttribute ( int  i) const
inline

Definition at line 217 of file fespace.hpp.

int mfem::FiniteElementSpace::GetBdrAttribute ( int  i) const
inline

Definition at line 219 of file fespace.hpp.

void mfem::FiniteElementSpace::GetBdrElementDofs ( int  i,
Array< int > &  dofs 
) const
virtual

Returns indexes of degrees of freedom for i'th boundary element.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 1147 of file fespace.cpp.

const Table& mfem::FiniteElementSpace::GetBdrElementToDofTable ( ) const
inline

Definition at line 289 of file fespace.hpp.

ElementTransformation* mfem::FiniteElementSpace::GetBdrElementTransformation ( int  i) const
inline

Returns ElementTransformation for the i'th boundary element.

Definition at line 214 of file fespace.hpp.

int mfem::FiniteElementSpace::GetBdrElementType ( int  i) const
inline

Returns the type of boundary element i.

Definition at line 201 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 139 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 1369 of file fespace.cpp.

void mfem::FiniteElementSpace::GetConformingInterpolation ( )
protected

Calculate the cP and cR matrices for a nonconforming mesh.

Definition at line 508 of file fespace.cpp.

const SparseMatrix * mfem::FiniteElementSpace::GetConformingProlongation ( )

Definition at line 690 of file fespace.cpp.

const SparseMatrix * mfem::FiniteElementSpace::GetConformingRestriction ( )

Definition at line 697 of file fespace.cpp.

int mfem::FiniteElementSpace::GetConformingVSize ( )
inline

Definition at line 172 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 1289 of file fespace.cpp.

const FiniteElement * mfem::FiniteElementSpace::GetEdgeElement ( int  i) const

Definition at line 1418 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEdgeFaceDofs ( int  type,
int  index,
Array< int > &  dofs 
)
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 495 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEdgeInteriorDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 1342 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEdgeInteriorVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 169 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 151 of file fespace.cpp.

void mfem::FiniteElementSpace::GetElementDofs ( int  i,
Array< int > &  dofs 
) const
virtual

Returns indexes of degrees of freedom in array dofs for i'th element.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 1035 of file fespace.cpp.

int mfem::FiniteElementSpace::GetElementForDof ( int  i) const
inline

Definition at line 291 of file fespace.hpp.

void mfem::FiniteElementSpace::GetElementInteriorDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 1330 of file fespace.cpp.

void mfem::FiniteElementSpace::GetElementInteriorVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 163 of file fespace.cpp.

const Table& mfem::FiniteElementSpace::GetElementToDofTable ( ) const
inline

Definition at line 288 of file fespace.hpp.

ElementTransformation* mfem::FiniteElementSpace::GetElementTransformation ( int  i) const
inline

Returns ElementTransformation for the i'th element.

Definition at line 205 of file fespace.hpp.

void mfem::FiniteElementSpace::GetElementTransformation ( int  i,
IsoparametricTransformation ElTr 
)
inline

Returns the transformation defining the i-th element in the user-defined variable.

Definition at line 210 of file fespace.hpp.

int mfem::FiniteElementSpace::GetElementType ( int  i) const
inline

Returns the type of element i.

Definition at line 193 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 133 of file fespace.cpp.

void mfem::FiniteElementSpace::GetElementVertices ( int  i,
Array< int > &  vertices 
) const
inline

Returns the vertices of element i.

Definition at line 197 of file fespace.hpp.

void mfem::FiniteElementSpace::GetEssentialTrueDofs ( const Array< int > &  bdr_attr_is_ess,
Array< int > &  ess_tdof_list 
)
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 295 of file fespace.cpp.

void mfem::FiniteElementSpace::GetEssentialVDofs ( const Array< int > &  bdr_attr_is_ess,
Array< int > &  ess_vdofs 
) const
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 258 of file fespace.cpp.

void mfem::FiniteElementSpace::GetFaceDofs ( int  i,
Array< int > &  dofs 
) const
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 1230 of file fespace.cpp.

const FiniteElement * mfem::FiniteElementSpace::GetFaceElement ( int  i) const

Definition at line 1395 of file fespace.cpp.

void mfem::FiniteElementSpace::GetFaceInteriorDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 1354 of file fespace.cpp.

int mfem::FiniteElementSpace::GetFaceOrder ( int  i) const

Returns the order of the i'th face finite element.

Definition at line 62 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 145 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 1134 of file fespace.cpp.

void mfem::FiniteElementSpace::GetLocalDerefinementMatrices ( int  geom,
const CoarseFineTransformations dt,
DenseTensor localR 
)
protected

Definition at line 791 of file fespace.cpp.

int mfem::FiniteElementSpace::GetLocalDofForDof ( int  i) const
inline

Definition at line 292 of file fespace.hpp.

Mesh* mfem::FiniteElementSpace::GetMesh ( ) const
inline

Returns the mesh.

Definition at line 136 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNBE ( ) const
inline

Returns number of boundary elements in the mesh.

Definition at line 190 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 704 of file fespace.cpp.

int mfem::FiniteElementSpace::GetNDofs ( ) const
inline

Returns number of degrees of freedom.

Definition at line 161 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNE ( ) const
inline

Returns number of elements in the mesh.

Definition at line 184 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNEDofs ( ) const
inline

Definition at line 180 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNFDofs ( ) const
inline

Definition at line 181 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNumElementInteriorDofs ( int  i) const
inline

Definition at line 241 of file fespace.hpp.

NURBSExtension* mfem::FiniteElementSpace::GetNURBSext ( )
inline

Definition at line 138 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNV ( ) const
inline

Returns number of nodes in the mesh.

Definition at line 187 of file fespace.hpp.

int mfem::FiniteElementSpace::GetNVDofs ( ) const
inline

Definition at line 179 of file fespace.hpp.

int mfem::FiniteElementSpace::GetOrder ( int  i) const

Returns the order of the i'th finite element.

Definition at line 56 of file fespace.cpp.

Ordering::Type mfem::FiniteElementSpace::GetOrdering ( ) const
inline

Return the ordering method.

Definition at line 175 of file fespace.hpp.

virtual const Operator* mfem::FiniteElementSpace::GetProlongationMatrix ( )
inlinevirtual

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 147 of file fespace.hpp.

virtual const SparseMatrix* mfem::FiniteElementSpace::GetRestrictionMatrix ( )
inlinevirtual

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 149 of file fespace.hpp.

long mfem::FiniteElementSpace::GetSequence ( ) const
inline

Return update counter (see Mesh::sequence)

Definition at line 370 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 1423 of file fespace.cpp.

virtual int mfem::FiniteElementSpace::GetTrueVSize ( )
inlinevirtual

Return the number of vector true (conforming) dofs.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 166 of file fespace.hpp.

const Operator* mfem::FiniteElementSpace::GetUpdateOperator ( )
inline

Get the GridFunction update matrix.

Definition at line 359 of file fespace.hpp.

int mfem::FiniteElementSpace::GetVDim ( ) const
inline

Returns vector dimension.

Definition at line 153 of file fespace.hpp.

void mfem::FiniteElementSpace::GetVertexDofs ( int  i,
Array< int > &  dofs 
) const

Definition at line 1318 of file fespace.cpp.

void mfem::FiniteElementSpace::GetVertexVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 157 of file fespace.cpp.

int mfem::FiniteElementSpace::GetVSize ( ) const
inline

Definition at line 163 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 421 of file fespace.cpp.

void mfem::FiniteElementSpace::ListToMarker ( const Array< int > &  list,
int  marker_size,
Array< int > &  marker,
int  mark_val = -1 
)
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 330 of file fespace.cpp.

void mfem::FiniteElementSpace::MakeVDimMatrix ( SparseMatrix mat) const
protected

Definition at line 663 of file fespace.cpp.

void mfem::FiniteElementSpace::MarkerToList ( const Array< int > &  marker,
Array< int > &  list 
)
static

Convert a Boolean marker array to a list containing all marked indices.

Definition at line 313 of file fespace.cpp.

bool mfem::FiniteElementSpace::Nonconforming ( ) const
inline

Definition at line 142 of file fespace.hpp.

void mfem::FiniteElementSpace::RebuildElementToDofTable ( )

Definition at line 197 of file fespace.cpp.

SparseMatrix * mfem::FiniteElementSpace::RefinementMatrix ( int  old_ndofs,
const Table old_elem_dof 
)
protected

Calculate GridFunction interpolation matrix after mesh refinement.

Definition at line 710 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 204 of file fespace.cpp.

void mfem::FiniteElementSpace::Save ( std::ostream &  out) const

Definition at line 1521 of file fespace.cpp.

void mfem::FiniteElementSpace::SetUpdateOperatorOwner ( bool  own)
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 364 of file fespace.hpp.

NURBSExtension * mfem::FiniteElementSpace::StealNURBSext ( )

Definition at line 944 of file fespace.cpp.

void mfem::FiniteElementSpace::Update ( bool  want_transform = true)
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 1457 of file fespace.cpp.

void mfem::FiniteElementSpace::UpdateNURBS ( )
protected

Definition at line 955 of file fespace.cpp.

virtual void mfem::FiniteElementSpace::UpdatesFinished ( )
inlinevirtual

Free GridFunction transformation matrix (if any), to save memory.

Reimplemented in mfem::ParFiniteElementSpace.

Definition at line 367 of file fespace.hpp.

int mfem::FiniteElementSpace::VDofToDof ( int  vdof) const
inline

Definition at line 252 of file fespace.hpp.

Member Data Documentation

int * mfem::FiniteElementSpace::bdofs
protected

Definition at line 80 of file fespace.hpp.

Table* mfem::FiniteElementSpace::bdrElem_dof
protected

Definition at line 83 of file fespace.hpp.

SparseMatrix* mfem::FiniteElementSpace::cP
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.

bool mfem::FiniteElementSpace::cP_is_set
protected

Definition at line 96 of file fespace.hpp.

SparseMatrix* mfem::FiniteElementSpace::cR
protected

Conforming restriction matrix such that cR.cP=I.

Definition at line 95 of file fespace.hpp.

Array<int> mfem::FiniteElementSpace::dof_elem_array
protected

Definition at line 85 of file fespace.hpp.

Array<int> mfem::FiniteElementSpace::dof_ldof_array
protected

Definition at line 85 of file fespace.hpp.

Table* mfem::FiniteElementSpace::elem_dof
mutableprotected

Definition at line 82 of file fespace.hpp.

int* mfem::FiniteElementSpace::fdofs
protected

Definition at line 80 of file fespace.hpp.

const FiniteElementCollection* mfem::FiniteElementSpace::fec
protected

Definition at line 66 of file fespace.hpp.

Mesh* mfem::FiniteElementSpace::mesh
protected

The mesh that FE space lives on.

Definition at line 64 of file fespace.hpp.

int mfem::FiniteElementSpace::nbdofs
protected

Definition at line 79 of file fespace.hpp.

int mfem::FiniteElementSpace::ndofs
protected

Number of degrees of freedom. Number of unknowns are ndofs*vdim.

Definition at line 77 of file fespace.hpp.

int mfem::FiniteElementSpace::nedofs
protected

Definition at line 79 of file fespace.hpp.

int mfem::FiniteElementSpace::nfdofs
protected

Definition at line 79 of file fespace.hpp.

NURBSExtension* mfem::FiniteElementSpace::NURBSext
protected

Definition at line 87 of file fespace.hpp.

int mfem::FiniteElementSpace::nvdofs
protected

Definition at line 79 of file fespace.hpp.

Ordering::Type mfem::FiniteElementSpace::ordering
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.

int mfem::FiniteElementSpace::own_ext
protected

Definition at line 88 of file fespace.hpp.

bool mfem::FiniteElementSpace::own_T
protected

Definition at line 100 of file fespace.hpp.

long mfem::FiniteElementSpace::sequence
protected

Definition at line 102 of file fespace.hpp.

Operator* mfem::FiniteElementSpace::T
protected

Transformation to apply to GridFunctions after space Update().

Definition at line 99 of file fespace.hpp.

int mfem::FiniteElementSpace::vdim
protected

Vector dimension (number of unknowns per degree of freedom).

Definition at line 69 of file fespace.hpp.


The documentation for this class was generated from the following files: