MFEM
v3.4
Finite element discretization library
|
Abstract parallel finite element space. More...
#include <pfespace.hpp>
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 corresponing 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 |
ParMesh * | GetParMesh () |
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 | GetSharedFaceDofs (int group, int fi, Array< int > &dofs) const |
HypreParMatrix * | Dof_TrueDof_Matrix () const |
The true dof-to-dof interpolation matrix. More... | |
HypreParMatrix * | GetPartialConformingInterpolation () |
For a non-conforming mesh, construct and return the interpolation matrix from the partially conforming true dofs to the local dofs. More... | |
HypreParVector * | NewTrueDofVector () |
void | DivideByGroupSize (double *vec) |
Scale a vector of true dofs. More... | |
GroupCommunicator & | GroupComm () |
Return a reference to the internal GroupCommunicator (on VDofs) More... | |
const GroupCommunicator & | GroupComm () const |
Return a const reference to the internal GroupCommunicator (on VDofs) More... | |
GroupCommunicator * | ScalarGroupComm () |
Return a new GroupCommunicator on Dofs. More... | |
void | Synchronize (Array< int > &ldof_marker) const |
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 Operator * | GetProlongationMatrix () const |
virtual const SparseMatrix * | GetRestrictionMatrix () 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 FiniteElement * | GetFaceNbrFE (int i) const |
const FiniteElement * | GetFaceNbrFaceFE (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 () |
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... | |
Mesh * | GetMesh () const |
Returns the mesh. More... | |
const NURBSExtension * | GetNURBSext () const |
NURBSExtension * | GetNURBSext () |
NURBSExtension * | StealNURBSext () |
bool | Conforming () const |
bool | Nonconforming () const |
const SparseMatrix * | GetConformingProlongation () const |
const SparseMatrix * | GetConformingRestriction () const |
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 FiniteElementCollection * | FEColl () const |
int | GetNVDofs () const |
int | GetNEDofs () const |
int | GetNFDofs () const |
int | GetNV () const |
Returns number of vertices in the mesh. More... | |
int | GetNE () const |
Returns number of elements in the mesh. More... | |
int | GetNF () const |
Returns number of faces (i.e. co-dimension 1 entities) in the mesh. More... | |
int | GetNBE () const |
Returns number of boundary elements in the mesh. More... | |
int | GetElementType (int i) const |
Returns the type of element i. More... | |
void | GetElementVertices (int i, Array< int > &vertices) const |
Returns the vertices of element i. More... | |
int | GetBdrElementType (int i) const |
Returns the type of boundary element i. More... | |
ElementTransformation * | GetElementTransformation (int i) const |
Returns ElementTransformation for the i-th element. More... | |
void | GetElementTransformation (int i, IsoparametricTransformation *ElTr) |
Returns the transformation defining the i-th element in the user-defined variable ElTr. More... | |
ElementTransformation * | GetBdrElementTransformation (int i) const |
Returns ElementTransformation for the i-th boundary element. More... | |
int | GetAttribute (int i) const |
int | GetBdrAttribute (int i) const |
void | GetEdgeDofs (int i, Array< int > &dofs) const |
void | GetVertexDofs (int i, Array< int > &dofs) const |
void | GetElementInteriorDofs (int i, Array< int > &dofs) const |
void | GetFaceInteriorDofs (int i, Array< int > &dofs) const |
int | GetNumElementInteriorDofs (int i) const |
void | GetEdgeInteriorDofs (int i, Array< int > &dofs) const |
void | DofsToVDofs (Array< int > &dofs, int ndofs=-1) const |
void | DofsToVDofs (int vd, Array< int > &dofs, int ndofs=-1) const |
int | DofToVDof (int dof, int vd, int ndofs=-1) const |
int | VDofToDof (int vdof) const |
void | GetElementVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom in array dofs for i'th element. More... | |
void | GetBdrElementVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom for i'th boundary element. More... | |
void | GetFaceVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom for i'th face element (2D and 3D). More... | |
void | GetEdgeVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom for i'th edge. More... | |
void | GetVertexVDofs (int i, Array< int > &vdofs) const |
void | GetElementInteriorVDofs (int i, Array< int > &vdofs) const |
void | GetEdgeInteriorVDofs (int i, Array< int > &vdofs) const |
void | RebuildElementToDofTable () |
void | ReorderElementToDofTable () |
Reorder the scalar DOFs based on the element ordering. More... | |
void | BuildDofToArrays () |
const Table & | GetElementToDofTable () const |
const Table & | GetBdrElementToDofTable () const |
int | GetElementForDof (int i) const |
int | GetLocalDofForDof (int i) const |
const FiniteElement * | GetFE (int i) const |
Returns pointer to the FiniteElement associated with i'th element. More... | |
const FiniteElement * | GetBE (int i) const |
Returns pointer to the FiniteElement for the i'th boundary element. More... | |
const FiniteElement * | GetFaceElement (int i) const |
const FiniteElement * | GetEdgeElement (int i) const |
const FiniteElement * | GetTraceElement (int i, int 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) |
SparseMatrix * | D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | H2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes) |
void | GetTransferOperator (const FiniteElementSpace &coarse_fes, OperatorHandle &T) const |
Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh. More... | |
const Operator * | GetUpdateOperator () |
Get the GridFunction update operator. More... | |
void | GetUpdateOperator (OperatorHandle &T) |
Return the update operator in the given OperatorHandle, T. More... | |
void | SetUpdateOperatorOwner (bool own) |
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator() must be deleted outside the FiniteElementSpace. More... | |
void | SetUpdateOperatorType (Operator::Type tid) |
Specify the Operator::Type to be used by the update operators. More... | |
long | GetSequence () const |
Return update counter (see Mesh::sequence) More... | |
void | Save (std::ostream &out) const |
FiniteElementCollection * | Load (Mesh *m, std::istream &input) |
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller. More... | |
virtual | ~FiniteElementSpace () |
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 |
SparseMatrix * | RefinementMatrix_main (const int coarse_ndofs, const Table &coarse_elem_dof, const DenseTensor &localP) const |
void | GetLocalRefinementMatrices (DenseTensor &localP) const |
void | GetLocalDerefinementMatrices (DenseTensor &localR) const |
SparseMatrix * | RefinementMatrix (int old_ndofs, const Table *old_elem_dof) |
SparseMatrix * | DerefinementMatrix (int old_ndofs, const Table *old_elem_dof) |
Calculate GridFunction restriction matrix after mesh derefinement. More... | |
void | GetLocalRefinementMatrices (const FiniteElementSpace &coarse_fes, 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 | |
Mesh * | mesh |
The mesh that FE space lives on (not owned). More... | |
const FiniteElementCollection * | fec |
Associated FE collection (not owned). More... | |
int | vdim |
Vector dimension (number of unknowns per degree of freedom). More... | |
Ordering::Type | ordering |
int | ndofs |
Number of degrees of freedom. Number of unknowns is ndofs * vdim. More... | |
int | nvdofs |
int | nedofs |
int | nfdofs |
int | nbdofs |
int * | fdofs |
int * | bdofs |
Table * | elem_dof |
Table * | bdrElem_dof |
Array< int > | dof_elem_array |
Array< int > | dof_ldof_array |
NURBSExtension * | NURBSext |
int | own_ext |
SparseMatrix * | cP |
SparseMatrix * | cR |
Conforming restriction matrix such that cR.cP=I. More... | |
bool | cP_is_set |
OperatorHandle | Th |
Transformation to apply to GridFunctions after space Update(). More... | |
long | sequence |
Abstract parallel finite element space.
Definition at line 28 of file pfespace.hpp.
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 poiters 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.
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 corresponing 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.
Definition at line 67 of file pfespace.cpp.
|
inlinevirtual |
Definition at line 358 of file pfespace.hpp.
|
inline |
Definition at line 338 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::DivideByGroupSize | ( | double * | vec | ) |
Scale a vector of true dofs.
Definition at line 557 of file pfespace.cpp.
|
inline |
The true dof-to-dof interpolation matrix.
Definition at line 267 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::ExchangeFaceNbrData | ( | ) |
Definition at line 755 of file pfespace.cpp.
|
virtual |
Returns indexes of degrees of freedom for i'th boundary element.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 375 of file pfespace.cpp.
|
inline |
Definition at line 235 of file pfespace.hpp.
|
inline |
Definition at line 243 of file pfespace.hpp.
|
inline |
Definition at line 241 of file pfespace.hpp.
|
virtual |
Returns indexes of degrees of freedom in array dofs for i'th element.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 361 of file pfespace.cpp.
|
virtual |
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in the array bdr_attr_is_ess.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 626 of file pfespace.cpp.
|
virtual |
Determine the boundary degrees of freedom.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 612 of file pfespace.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 from mfem::FiniteElementSpace.
Definition at line 389 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetFaceNbrElementVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 974 of file pfespace.cpp.
const FiniteElement * mfem::ParFiniteElementSpace::GetFaceNbrFaceFE | ( | int | i | ) | const |
Definition at line 1022 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetFaceNbrFaceVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 980 of file pfespace.cpp.
const FiniteElement * mfem::ParFiniteElementSpace::GetFaceNbrFE | ( | int | i | ) | const |
Definition at line 1008 of file pfespace.cpp.
|
inline |
Definition at line 332 of file pfespace.hpp.
|
inline |
Definition at line 327 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 695 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 672 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 651 of file pfespace.cpp.
HYPRE_Int mfem::ParFiniteElementSpace::GetMyDofOffset | ( | ) | const |
Definition at line 732 of file pfespace.cpp.
|
inline |
Definition at line 237 of file pfespace.hpp.
HYPRE_Int mfem::ParFiniteElementSpace::GetMyTDofOffset | ( | ) | const |
Definition at line 737 of file pfespace.cpp.
|
inline |
Definition at line 236 of file pfespace.hpp.
|
inline |
Definition at line 239 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.
Definition at line 546 of file pfespace.cpp.
|
virtual |
Reimplemented from mfem::FiniteElementSpace.
Definition at line 742 of file pfespace.cpp.
|
inlinevirtual |
Get the R matrix which restricts a local dof vector to true dof vector.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 322 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::GetSharedEdgeDofs | ( | int | group, |
int | ei, | ||
Array< int > & | dofs | ||
) | const |
Definition at line 398 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetSharedFaceDofs | ( | int | group, |
int | fi, | ||
Array< int > & | dofs | ||
) | const |
Definition at line 421 of file pfespace.cpp.
|
inline |
Definition at line 244 of file pfespace.hpp.
|
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 2541 of file pfespace.cpp.
|
inlinevirtual |
Return the number of local vector true dofs.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 251 of file pfespace.hpp.
|
inline |
Definition at line 247 of file pfespace.hpp.
|
inline |
Definition at line 245 of file pfespace.hpp.
|
inline |
Return a reference to the internal GroupCommunicator (on VDofs)
Definition at line 284 of file pfespace.hpp.
|
inline |
Return a const reference to the internal GroupCommunicator (on VDofs)
Definition at line 287 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::Lose_Dof_TrueDof_Matrix | ( | ) |
Definition at line 1031 of file pfespace.cpp.
|
inline |
Definition at line 335 of file pfespace.hpp.
|
inline |
Definition at line 336 of file pfespace.hpp.
|
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 277 of file pfespace.hpp.
|
inline |
Definition at line 339 of file pfespace.hpp.
GroupCommunicator * mfem::ParFiniteElementSpace::ScalarGroupComm | ( | ) |
Return a new GroupCommunicator on Dofs.
Definition at line 575 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.
Definition at line 595 of file pfespace.cpp.
|
inline |
Definition at line 361 of file pfespace.hpp.
|
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 2568 of file pfespace.cpp.
|
inlinevirtual |
Free ParGridFunction transformation matrix (if any), to save memory.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 352 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::face_nbr_element_dof |
Definition at line 176 of file pfespace.hpp.
Array<HYPRE_Int> mfem::ParFiniteElementSpace::face_nbr_glob_dof_map |
Definition at line 180 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::face_nbr_ldof |
Definition at line 178 of file pfespace.hpp.
int mfem::ParFiniteElementSpace::num_face_nbr_dofs |
Definition at line 174 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::send_face_nbr_ldof |
Definition at line 182 of file pfespace.hpp.