MFEM
v3.1
Finite element discretization library
|
Abstract parallel finite element space. More...
#include <pfespace.hpp>
Public Member Functions | |
ParFiniteElementSpace (ParMesh *pm, const FiniteElementCollection *f, int dim=1, int ordering=Ordering::byNODES) | |
MPI_Comm | GetComm () |
int | GetNRanks () |
int | GetMyRank () |
ParMesh * | GetParMesh () |
int | GetDofSign (int i) |
HYPRE_Int * | GetDofOffsets () |
HYPRE_Int * | GetTrueDofOffsets () |
HYPRE_Int | GlobalVSize () |
HYPRE_Int | GlobalTrueVSize () |
virtual int | GetTrueVSize () |
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 |
HypreParMatrix * | Dof_TrueDof_Matrix () |
The true dof-to-dof interpolation matrix. 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... | |
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) const |
Determine the boundary degrees of freedom. More... | |
virtual void | GetEssentialTrueDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list) |
int | GetLocalTDofNumber (int ldof) |
HYPRE_Int | GetGlobalTDofNumber (int ldof) |
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 SparseMatrix * | GetRestrictionMatrix () |
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 |
const FiniteElement * | GetFaceNbrFE (int i) const |
const HYPRE_Int * | GetFaceNbrGlobalDofMap () |
void | Lose_Dof_TrueDof_Matrix () |
void | LoseDofOffsets () |
void | LoseTrueDofOffsets () |
bool | Conforming () const |
bool | Nonconforming () const |
virtual void | Update () |
virtual FiniteElementSpace * | SaveUpdate () |
Return a copy of the current FE space and update. More... | |
virtual | ~ParFiniteElementSpace () |
int | TrueVSize () |
Public Member Functions inherited from mfem::FiniteElementSpace | |
FiniteElementSpace (Mesh *m, const FiniteElementCollection *f, int vdim=1, int ordering=Ordering::byNODES) | |
Mesh * | GetMesh () const |
Returns the mesh. More... | |
NURBSExtension * | GetNURBSext () |
NURBSExtension * | StealNURBSext () |
bool | Conforming () const |
bool | Nonconforming () const |
const SparseMatrix * | GetConformingProlongation () |
const SparseMatrix * | GetConformingRestriction () |
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 |
int | GetNConformingDofs () |
int | GetConformingVSize () |
int | GetOrdering () const |
Return the ordering method. More... | |
const FiniteElementCollection * | FEColl () const |
int | GetNVDofs () const |
int | GetNEDofs () const |
int | GetNFDofs () const |
int | GetNE () const |
Returns number of elements in the mesh. More... | |
int | GetNV () const |
Returns number of nodes in the mesh. More... | |
int | GetNBE () const |
Returns number of boundary elements in the mesh. More... | |
int | GetElementType (int i) const |
Returns the type of element i. More... | |
void | GetElementVertices (int i, Array< int > &vertices) const |
Returns the vertices of element i. More... | |
int | GetBdrElementType (int i) const |
Returns the type of boundary element i. More... | |
ElementTransformation * | GetElementTransformation (int i) const |
Returns ElementTransformation for the i'th element. More... | |
void | GetElementTransformation (int i, IsoparametricTransformation *ElTr) |
ElementTransformation * | GetBdrElementTransformation (int i) const |
Returns ElementTransformation for the i'th boundary element. More... | |
int | GetAttribute (int i) const |
int | GetBdrAttribute (int i) const |
void | GetEdgeDofs (int i, Array< int > &dofs) const |
void | GetVertexDofs (int i, Array< int > &dofs) const |
void | GetElementInteriorDofs (int i, Array< int > &dofs) const |
int | GetNumElementInteriorDofs (int i) const |
void | GetEdgeInteriorDofs (int i, Array< int > &dofs) const |
void | DofsToVDofs (Array< int > &dofs) const |
void | DofsToVDofs (int vd, Array< int > &dofs, int ndofs=-1) const |
int | DofToVDof (int dof, int vd, int ndofs=-1) const |
int | VDofToDof (int vdof) const |
void | GetElementVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom in array dofs for i'th element. More... | |
void | GetBdrElementVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom for i'th boundary element. More... | |
void | GetFaceVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom for i'th face element (2D and 3D). More... | |
void | GetEdgeVDofs (int i, Array< int > &vdofs) const |
Returns indexes of degrees of freedom for i'th edge. More... | |
void | GetVertexVDofs (int i, Array< int > &vdofs) const |
void | GetElementInteriorVDofs (int i, Array< int > &vdofs) const |
void | GetEdgeInteriorVDofs (int i, Array< int > &vdofs) const |
void | BuildElementToDofTable () |
void | BuildDofToArrays () |
const Table & | GetElementToDofTable () const |
const Table & | GetBdrElementToDofTable () const |
int | GetElementForDof (int i) |
int | GetLocalDofForDof (int i) |
const FiniteElement * | GetFE (int i) const |
Returns pointer to the FiniteElement associated with i'th element. More... | |
const FiniteElement * | GetBE (int i) const |
Returns pointer to the FiniteElement for the i'th boundary element. More... | |
const FiniteElement * | GetFaceElement (int i) const |
const FiniteElement * | GetEdgeElement (int i) const |
const FiniteElement * | GetTraceElement (int i, int geom_type) const |
Return the trace element from element 'i' to the given 'geom_type'. More... | |
SparseMatrix * | GlobalRestrictionMatrix (FiniteElementSpace *cfes, int one_vdim=-1) |
void | ConvertToConformingVDofs (const Array< int > &dofs, Array< int > &cdofs) |
void | ConvertFromConformingVDofs (const Array< int > &cdofs, Array< int > &dofs) |
void | EliminateEssentialBCFromGRM (FiniteElementSpace *cfes, Array< int > &bdr_attr_is_ess, SparseMatrix *R) |
SparseMatrix * | GlobalRestrictionMatrix (FiniteElementSpace *cfes, Array< int > &bdr_attr_is_ess, int one_vdim=-1) |
Generate the global restriction matrix with eliminated essential bc. More... | |
SparseMatrix * | D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
SparseMatrix * | H2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes) |
virtual void | UpdateAndInterpolate (int num_grid_fns,...) |
void | UpdateAndInterpolate (GridFunction *gf) |
A shortcut for passing only one GridFunction to UndateAndInterpolate. More... | |
void | Save (std::ostream &out) const |
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 | Constructor () |
void | Destructor () |
FiniteElementSpace (FiniteElementSpace &) | |
void | ConstructRefinementData (int k, int cdofs, RefinementType type) |
Constructs new refinement data using coarse element k as a template. More... | |
DenseMatrix * | LocalInterpolation (int k, int cdofs, RefinementType type, Array< int > &rows) |
Generates the local interpolation matrix for coarse element k. More... | |
SparseMatrix * | NC_GlobalRestrictionMatrix (FiniteElementSpace *cfes, NCMesh *ncmesh) |
void | GetEdgeFaceDofs (int type, int index, Array< int > &dofs) |
void | GetConformingInterpolation () |
void | MakeVDimMatrix (SparseMatrix &mat) const |
Protected Attributes inherited from mfem::FiniteElementSpace | |
Mesh * | mesh |
The mesh that FE space lives on. More... | |
int | vdim |
Vector dimension (number of unknowns per degree of freedom). More... | |
int | ndofs |
Number of degrees of freedom. Number of unknowns are ndofs*vdim. More... | |
int | ordering |
const FiniteElementCollection * | fec |
int | nvdofs |
int | nedofs |
int | nfdofs |
int | nbdofs |
int * | fdofs |
int * | bdofs |
Array< RefinementData * > | RefData |
Collection of currently known refinement data. More... | |
Table * | elem_dof |
Table * | bdrElem_dof |
Array< int > | dof_elem_array |
Array< int > | dof_ldof_array |
NURBSExtension * | NURBSext |
int | own_ext |
SparseMatrix * | cP |
SparseMatrix * | cR |
Abstract parallel finite element space.
Definition at line 28 of file pfespace.hpp.
mfem::ParFiniteElementSpace::ParFiniteElementSpace | ( | ParMesh * | pm, |
const FiniteElementCollection * | f, | ||
int | dim = 1 , |
||
int | ordering = Ordering::byNODES |
||
) |
Definition at line 53 of file pfespace.cpp.
|
inlinevirtual |
Definition at line 236 of file pfespace.hpp.
|
inline |
Definition at line 229 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::DivideByGroupSize | ( | double * | vec | ) |
Scale a vector of true dofs.
Definition at line 479 of file pfespace.cpp.
HypreParMatrix * mfem::ParFiniteElementSpace::Dof_TrueDof_Matrix | ( | ) |
The true dof-to-dof interpolation matrix.
Definition at line 411 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::ExchangeFaceNbrData | ( | ) |
Definition at line 647 of file pfespace.cpp.
|
virtual |
Returns indexes of degrees of freedom for i'th boundary element.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 288 of file pfespace.cpp.
|
inline |
Definition at line 143 of file pfespace.hpp.
|
inline |
Definition at line 151 of file pfespace.hpp.
|
inline |
Definition at line 149 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 274 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 545 of file pfespace.cpp.
|
virtual |
Determine the boundary degrees of freedom.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 532 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 302 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetFaceNbrElementVDofs | ( | int | i, |
Array< int > & | vdofs | ||
) | const |
Definition at line 867 of file pfespace.cpp.
const FiniteElement * mfem::ParFiniteElementSpace::GetFaceNbrFE | ( | int | i | ) | const |
Definition at line 873 of file pfespace.cpp.
|
inline |
Definition at line 223 of file pfespace.hpp.
|
inline |
Definition at line 220 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 600 of file pfespace.cpp.
HYPRE_Int mfem::ParFiniteElementSpace::GetGlobalTDofNumber | ( | int | ldof | ) |
Returns the global tdof number of the given local degree of freedom.
Definition at line 577 of file pfespace.cpp.
int mfem::ParFiniteElementSpace::GetLocalTDofNumber | ( | int | ldof | ) |
If the given ldof is owned by the current processor, return its local tdof number, otherwise return -1
Definition at line 555 of file pfespace.cpp.
HYPRE_Int mfem::ParFiniteElementSpace::GetMyDofOffset | ( | ) | const |
Definition at line 637 of file pfespace.cpp.
|
inline |
Definition at line 145 of file pfespace.hpp.
HYPRE_Int mfem::ParFiniteElementSpace::GetMyTDofOffset | ( | ) | const |
Definition at line 642 of file pfespace.cpp.
|
inline |
Definition at line 144 of file pfespace.hpp.
|
inline |
Definition at line 147 of file pfespace.hpp.
|
inlinevirtual |
Get the R matrix which restricts a local dof vector to true dof vector.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 215 of file pfespace.hpp.
|
inline |
Definition at line 152 of file pfespace.hpp.
|
inlinevirtual |
Return the number of local vector true dofs.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 159 of file pfespace.hpp.
|
inline |
Definition at line 155 of file pfespace.hpp.
|
inline |
Definition at line 153 of file pfespace.hpp.
|
inline |
Return a reference to the internal GroupCommunicator (on VDofs)
Definition at line 183 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::Lose_Dof_TrueDof_Matrix | ( | ) |
Definition at line 887 of file pfespace.cpp.
|
inline |
Definition at line 226 of file pfespace.hpp.
|
inline |
Definition at line 227 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 176 of file pfespace.hpp.
|
inline |
Definition at line 230 of file pfespace.hpp.
|
virtual |
Return a copy of the current FE space and update.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 1554 of file pfespace.cpp.
GroupCommunicator * mfem::ParFiniteElementSpace::ScalarGroupComm | ( | ) |
Return a new GroupCommunicator on Dofs.
Definition at line 495 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 515 of file pfespace.cpp.
|
inline |
Definition at line 239 of file pfespace.hpp.
|
virtual |
Reimplemented from mfem::FiniteElementSpace.
Definition at line 1522 of file pfespace.cpp.
Table mfem::ParFiniteElementSpace::face_nbr_element_dof |
Definition at line 132 of file pfespace.hpp.
Array<HYPRE_Int> mfem::ParFiniteElementSpace::face_nbr_glob_dof_map |
Definition at line 136 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::face_nbr_ldof |
Definition at line 134 of file pfespace.hpp.
int mfem::ParFiniteElementSpace::num_face_nbr_dofs |
Definition at line 130 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::send_face_nbr_ldof |
Definition at line 138 of file pfespace.hpp.