MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Public Attributes | List of all members
mfem::ParFiniteElementSpace Class Reference

Abstract parallel finite element space. More...

#include <pfespace.hpp>

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

Public Member Functions

 ParFiniteElementSpace (const ParFiniteElementSpace &orig, ParMesh *pmesh=NULL, const FiniteElementCollection *fec=NULL)
 Copy constructor: deep copy all data from orig except the ParMesh, the FiniteElementCollection, and some derived data. More...
 
 ParFiniteElementSpace (const FiniteElementSpace &orig, ParMesh &pmesh, const FiniteElementCollection *fec=NULL)
 Convert/copy the local (Par)FiniteElementSpace orig to ParFiniteElementSpace: deep copy all data from orig except the Mesh, the FiniteElementCollection, and some derived data. More...
 
 ParFiniteElementSpace (ParMesh *pm, const FiniteElementSpace *global_fes, const int *partitioning, const FiniteElementCollection *f=NULL)
 Construct the local ParFiniteElementSpace corresponding to the global FE space, global_fes. More...
 
 ParFiniteElementSpace (ParMesh *pm, const FiniteElementCollection *f, int dim=1, int ordering=Ordering::byNODES)
 
 ParFiniteElementSpace (ParMesh *pm, NURBSExtension *ext, const FiniteElementCollection *f, int dim=1, int ordering=Ordering::byNODES)
 Construct a NURBS FE space based on the given NURBSExtension, ext. More...
 
MPI_Comm GetComm () const
 
int GetNRanks () const
 
int GetMyRank () const
 
ParMeshGetParMesh ()
 
int GetDofSign (int i)
 
HYPRE_Int * GetDofOffsets () const
 
HYPRE_Int * GetTrueDofOffsets () const
 
HYPRE_Int GlobalVSize () const
 
HYPRE_Int GlobalTrueVSize () const
 
virtual int GetTrueVSize () const
 Return the number of local vector true dofs. More...
 
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 GetSharedEdgeDofs (int group, int ei, Array< int > &dofs) const
 
void GetSharedTriangleDofs (int group, int fi, Array< int > &dofs) const
 
void GetSharedQuadrilateralDofs (int group, int fi, Array< int > &dofs) const
 
HypreParMatrixDof_TrueDof_Matrix () const
 The true dof-to-dof interpolation matrix. More...
 
HypreParMatrixGetPartialConformingInterpolation ()
 For a non-conforming mesh, construct and return the interpolation matrix from the partially conforming true dofs to the local dofs. More...
 
HypreParVectorNewTrueDofVector ()
 
void DivideByGroupSize (double *vec)
 Scale a vector of true dofs. More...
 
GroupCommunicatorGroupComm ()
 Return a reference to the internal GroupCommunicator (on VDofs) More...
 
const GroupCommunicatorGroupComm () const
 Return a const reference to the internal GroupCommunicator (on VDofs) More...
 
GroupCommunicatorScalarGroupComm ()
 Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1. More...
 
void Synchronize (Array< int > &ldof_marker) const
 Given an integer array on the local degrees of freedom, perform a bitwise OR between the shared dofs. More...
 
virtual void GetEssentialVDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const
 Determine the boundary degrees of freedom. More...
 
virtual void GetEssentialTrueDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
 
int GetLocalTDofNumber (int ldof) const
 
HYPRE_Int GetGlobalTDofNumber (int ldof) const
 Returns the global tdof number of the given local degree of freedom. More...
 
HYPRE_Int GetGlobalScalarTDofNumber (int sldof)
 
HYPRE_Int GetMyDofOffset () const
 
HYPRE_Int GetMyTDofOffset () const
 
virtual const OperatorGetProlongationMatrix () const
 The returned Operator is owned by the FiniteElementSpace. More...
 
virtual const SparseMatrixGetRestrictionMatrix () const
 Get the R matrix which restricts a local dof vector to true dof vector. More...
 
void ExchangeFaceNbrData ()
 
int GetFaceNbrVSize () const
 
void GetFaceNbrElementVDofs (int i, Array< int > &vdofs) const
 
void GetFaceNbrFaceVDofs (int i, Array< int > &vdofs) const
 
const FiniteElementGetFaceNbrFE (int i) const
 
const FiniteElementGetFaceNbrFaceFE (int i) const
 
const HYPRE_Int * GetFaceNbrGlobalDofMap ()
 
void Lose_Dof_TrueDof_Matrix ()
 
void LoseDofOffsets ()
 
void LoseTrueDofOffsets ()
 
bool Conforming () const
 
bool Nonconforming () const
 
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)
 
virtual void UpdatesFinished ()
 Free ParGridFunction transformation matrix (if any), to save memory. More...
 
virtual ~ParFiniteElementSpace ()
 
void PrintPartitionStats ()
 
int TrueVSize () const
 
- Public Member Functions inherited from mfem::FiniteElementSpace
 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...
 
MeshGetMesh () const
 Returns the mesh. More...
 
const NURBSExtensionGetNURBSext () const
 
NURBSExtensionGetNURBSext ()
 
NURBSExtensionStealNURBSext ()
 
bool Conforming () const
 
bool Nonconforming () const
 
const SparseMatrixGetConformingProlongation () const
 The returned SparseMatrix is owned by the FiniteElementSpace. More...
 
const SparseMatrixGetConformingRestriction () const
 The returned SparseMatrix is owned by the FiniteElementSpace. More...
 
const OperatorGetElementRestriction (ElementDofOrdering e_ordering) const
 Return an Operator that converts L-vectors to E-vectors. More...
 
const QuadratureInterpolatorGetQuadratureInterpolator (const IntegrationRule &ir) const
 Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). More...
 
const QuadratureInterpolatorGetQuadratureInterpolator (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...
 
int GetNConformingDofs () const
 
int GetConformingVSize () const
 
Ordering::Type GetOrdering () const
 Return the ordering method. More...
 
const FiniteElementCollectionFEColl () 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...
 
ElementTransformationGetElementTransformation (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...
 
ElementTransformationGetBdrElementTransformation (int i) const
 Returns ElementTransformation for the i-th boundary element. More...
 
int GetAttribute (int i) const
 
int GetBdrAttribute (int i) 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, Geometry::Type geom_type) const
 Return the trace element from element 'i' to the given 'geom_type'. More...
 
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)
 
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...
 
const OperatorGetUpdateOperator ()
 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...
 
long GetSequence () const
 Return update counter (see Mesh::sequence) More...
 
void Save (std::ostream &out) const
 
FiniteElementCollectionLoad (Mesh *m, std::istream &input)
 Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller. More...
 
virtual ~FiniteElementSpace ()
 

Public Attributes

int num_face_nbr_dofs
 
Table face_nbr_element_dof
 
Table face_nbr_ldof
 
Array< HYPRE_Int > face_nbr_glob_dof_map
 
Table send_face_nbr_ldof
 

Additional Inherited Members

- Static Public Member Functions inherited from mfem::FiniteElementSpace
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 inherited from mfem::FiniteElementSpace
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
 
SparseMatrixRefinementMatrix_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
 
SparseMatrixRefinementMatrix (int old_ndofs, const Table *old_elem_dof)
 
SparseMatrixDerefinementMatrix (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 inherited from mfem::FiniteElementSpace
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 inherited from mfem::FiniteElementSpace
Meshmesh
 The mesh that FE space lives on (not owned). More...
 
const FiniteElementCollectionfec
 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
 
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
 
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
 

Detailed Description

Abstract parallel finite element space.

Definition at line 28 of file pfespace.hpp.

Constructor & Destructor Documentation

mfem::ParFiniteElementSpace::ParFiniteElementSpace ( const ParFiniteElementSpace orig,
ParMesh pmesh = NULL,
const FiniteElementCollection fec = NULL 
)

Copy constructor: deep copy all data from orig except the ParMesh, the FiniteElementCollection, and some derived data.

If the pmesh or fec pointers are NULL (default), then the new ParFiniteElementSpace 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.

Note
The objects pointed to by the pmesh and fec parameters must be either the same objects as the ones used by orig, or copies of them. Otherwise, the behavior is undefined.
Derived data objects, such as the parallel prolongation and restriction operators, the update operator, and any of the face-neighbor data, will not be copied, even if they are created in the orig object.

Definition at line 28 of file pfespace.cpp.

mfem::ParFiniteElementSpace::ParFiniteElementSpace ( const FiniteElementSpace orig,
ParMesh pmesh,
const FiniteElementCollection fec = NULL 
)

Convert/copy the local (Par)FiniteElementSpace orig to ParFiniteElementSpace: deep copy all data from orig except the Mesh, the FiniteElementCollection, and some derived data.

Definition at line 36 of file pfespace.cpp.

mfem::ParFiniteElementSpace::ParFiniteElementSpace ( ParMesh pm,
const FiniteElementSpace global_fes,
const int *  partitioning,
const FiniteElementCollection f = NULL 
)

Construct the local ParFiniteElementSpace corresponding to the global FE space, global_fes.

The parameter pm is the local ParMesh obtained by decomposing the global Mesh used by global_fes. The array partitioning represents the parallel decomposition - it maps global element ids to MPI ranks. If the FiniteElementCollection, f, is NULL (default), the FE collection used by global_fes will be reused. If f is not NULL, it must be the same as, or a copy of, the FE collection used by global_fes.

Definition at line 44 of file pfespace.cpp.

mfem::ParFiniteElementSpace::ParFiniteElementSpace ( ParMesh pm,
const FiniteElementCollection f,
int  dim = 1,
int  ordering = Ordering::byNODES 
)

Definition at line 60 of file pfespace.cpp.

mfem::ParFiniteElementSpace::ParFiniteElementSpace ( ParMesh pm,
NURBSExtension ext,
const FiniteElementCollection f,
int  dim = 1,
int  ordering = Ordering::byNODES 
)

Construct a NURBS FE space based on the given NURBSExtension, ext.

The parameter ext will be deleted by this constructor, replaced by a ParNURBSExtension owned by the ParFiniteElementSpace.

Note
If the pointer ext is NULL, this constructor is equivalent to the standard constructor with the same arguments minus the NURBSExtension, ext.

Definition at line 67 of file pfespace.cpp.

virtual mfem::ParFiniteElementSpace::~ParFiniteElementSpace ( )
inlinevirtual

Definition at line 366 of file pfespace.hpp.

Member Function Documentation

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

Definition at line 346 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::DivideByGroupSize ( double *  vec)

Scale a vector of true dofs.

Definition at line 677 of file pfespace.cpp.

HypreParMatrix* mfem::ParFiniteElementSpace::Dof_TrueDof_Matrix ( ) const
inline

The true dof-to-dof interpolation matrix.

Definition at line 272 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::ExchangeFaceNbrData ( )

Definition at line 868 of file pfespace.cpp.

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

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

Reimplemented from mfem::FiniteElementSpace.

Definition at line 470 of file pfespace.cpp.

MPI_Comm mfem::ParFiniteElementSpace::GetComm ( ) const
inline

Definition at line 239 of file pfespace.hpp.

HYPRE_Int* mfem::ParFiniteElementSpace::GetDofOffsets ( ) const
inline

Definition at line 247 of file pfespace.hpp.

int mfem::ParFiniteElementSpace::GetDofSign ( int  i)
inline

Definition at line 245 of file pfespace.hpp.

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

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

Reimplemented from mfem::FiniteElementSpace.

Definition at line 456 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::GetEssentialTrueDofs ( const Array< int > &  bdr_attr_is_ess,
Array< int > &  ess_tdof_list,
int  component = -1 
)
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 from mfem::FiniteElementSpace.

Definition at line 737 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::GetEssentialVDofs ( const Array< int > &  bdr_attr_is_ess,
Array< int > &  ess_dofs,
int  component = -1 
) const
virtual

Determine the boundary degrees of freedom.

Reimplemented from mfem::FiniteElementSpace.

Definition at line 723 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::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 from mfem::FiniteElementSpace.

Definition at line 484 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::GetFaceNbrElementVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 1109 of file pfespace.cpp.

const FiniteElement * mfem::ParFiniteElementSpace::GetFaceNbrFaceFE ( int  i) const

Definition at line 1157 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::GetFaceNbrFaceVDofs ( int  i,
Array< int > &  vdofs 
) const

Definition at line 1115 of file pfespace.cpp.

const FiniteElement * mfem::ParFiniteElementSpace::GetFaceNbrFE ( int  i) const

Definition at line 1143 of file pfespace.cpp.

const HYPRE_Int* mfem::ParFiniteElementSpace::GetFaceNbrGlobalDofMap ( )
inline

Definition at line 340 of file pfespace.hpp.

int mfem::ParFiniteElementSpace::GetFaceNbrVSize ( ) const
inline

Definition at line 335 of file pfespace.hpp.

HYPRE_Int mfem::ParFiniteElementSpace::GetGlobalScalarTDofNumber ( int  sldof)

Returns the global tdof number of the given local degree of freedom in the scalar version of the current finite element space. The input should be a scalar local dof.

Definition at line 808 of file pfespace.cpp.

HYPRE_Int mfem::ParFiniteElementSpace::GetGlobalTDofNumber ( int  ldof) const

Returns the global tdof number of the given local degree of freedom.

Definition at line 785 of file pfespace.cpp.

int mfem::ParFiniteElementSpace::GetLocalTDofNumber ( int  ldof) const

If the given ldof is owned by the current processor, return its local tdof number, otherwise return -1

Definition at line 764 of file pfespace.cpp.

HYPRE_Int mfem::ParFiniteElementSpace::GetMyDofOffset ( ) const

Definition at line 845 of file pfespace.cpp.

int mfem::ParFiniteElementSpace::GetMyRank ( ) const
inline

Definition at line 241 of file pfespace.hpp.

HYPRE_Int mfem::ParFiniteElementSpace::GetMyTDofOffset ( ) const

Definition at line 850 of file pfespace.cpp.

int mfem::ParFiniteElementSpace::GetNRanks ( ) const
inline

Definition at line 240 of file pfespace.hpp.

ParMesh* mfem::ParFiniteElementSpace::GetParMesh ( )
inline

Definition at line 243 of file pfespace.hpp.

HypreParMatrix * mfem::ParFiniteElementSpace::GetPartialConformingInterpolation ( )

For a non-conforming mesh, construct and return the interpolation matrix from the partially conforming true dofs to the local dofs.

Note
The returned pointer must be deleted by the caller.

Definition at line 666 of file pfespace.cpp.

const Operator * mfem::ParFiniteElementSpace::GetProlongationMatrix ( ) const
virtual

The returned Operator is owned by the FiniteElementSpace.

Reimplemented from mfem::FiniteElementSpace.

Definition at line 855 of file pfespace.cpp.

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

Get the R matrix which restricts a local dof vector to true dof vector.

Reimplemented from mfem::FiniteElementSpace.

Definition at line 330 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::GetSharedEdgeDofs ( int  group,
int  ei,
Array< int > &  dofs 
) const

Definition at line 493 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::GetSharedQuadrilateralDofs ( int  group,
int  fi,
Array< int > &  dofs 
) const

Definition at line 540 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::GetSharedTriangleDofs ( int  group,
int  fi,
Array< int > &  dofs 
) const

Definition at line 516 of file pfespace.cpp.

HYPRE_Int* mfem::ParFiniteElementSpace::GetTrueDofOffsets ( ) const
inline

Definition at line 248 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::GetTrueTransferOperator ( const FiniteElementSpace coarse_fes,
OperatorHandle T 
) const
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 from mfem::FiniteElementSpace.

Definition at line 2727 of file pfespace.cpp.

virtual int mfem::ParFiniteElementSpace::GetTrueVSize ( ) const
inlinevirtual

Return the number of local vector true dofs.

Reimplemented from mfem::FiniteElementSpace.

Definition at line 255 of file pfespace.hpp.

HYPRE_Int mfem::ParFiniteElementSpace::GlobalTrueVSize ( ) const
inline

Definition at line 251 of file pfespace.hpp.

HYPRE_Int mfem::ParFiniteElementSpace::GlobalVSize ( ) const
inline

Definition at line 249 of file pfespace.hpp.

GroupCommunicator& mfem::ParFiniteElementSpace::GroupComm ( )
inline

Return a reference to the internal GroupCommunicator (on VDofs)

Definition at line 289 of file pfespace.hpp.

const GroupCommunicator& mfem::ParFiniteElementSpace::GroupComm ( ) const
inline

Return a const reference to the internal GroupCommunicator (on VDofs)

Definition at line 292 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::Lose_Dof_TrueDof_Matrix ( )

Definition at line 1166 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::LoseDofOffsets ( )
inline

Definition at line 343 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::LoseTrueDofOffsets ( )
inline

Definition at line 344 of file pfespace.hpp.

HypreParVector* mfem::ParFiniteElementSpace::NewTrueDofVector ( )
inline

Create and return a new HypreParVector on the true dofs, which is owned by (i.e. it must be destroyed by) the calling function.

Definition at line 282 of file pfespace.hpp.

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

Definition at line 347 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::PrintPartitionStats ( )

Definition at line 181 of file pfespace.cpp.

GroupCommunicator * mfem::ParFiniteElementSpace::ScalarGroupComm ( )

Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.

Note
The returned pointer must be deleted by the caller.

Definition at line 697 of file pfespace.cpp.

void mfem::ParFiniteElementSpace::Synchronize ( Array< int > &  ldof_marker) const

Given an integer array on the local degrees of freedom, perform a bitwise OR between the shared dofs.

For non-conforming mesh, synchronization is performed on the cut (aka "partially conforming") space.

Definition at line 711 of file pfespace.cpp.

int mfem::ParFiniteElementSpace::TrueVSize ( ) const
inline

Definition at line 371 of file pfespace.hpp.

void mfem::ParFiniteElementSpace::Update ( bool  want_transform = true)
virtual

Reflect changes in the mesh. Calculate one of the refinement/derefinement /rebalance matrices, unless want_transform is false.

Reimplemented from mfem::FiniteElementSpace.

Definition at line 2754 of file pfespace.cpp.

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

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

Reimplemented from mfem::FiniteElementSpace.

Definition at line 360 of file pfespace.hpp.

Member Data Documentation

Table mfem::ParFiniteElementSpace::face_nbr_element_dof

Definition at line 180 of file pfespace.hpp.

Array<HYPRE_Int> mfem::ParFiniteElementSpace::face_nbr_glob_dof_map

Definition at line 184 of file pfespace.hpp.

Table mfem::ParFiniteElementSpace::face_nbr_ldof

Definition at line 182 of file pfespace.hpp.

int mfem::ParFiniteElementSpace::num_face_nbr_dofs

Definition at line 178 of file pfespace.hpp.

Table mfem::ParFiniteElementSpace::send_face_nbr_ldof

Definition at line 186 of file pfespace.hpp.


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