15 #include "../config/config.hpp" 19 #include "../linalg/hypre.hpp" 20 #include "../mesh/pmesh.hpp" 21 #include "../mesh/nurbs.hpp" 45 mutable int ltdof_size;
48 int ngvdofs, ngedofs, ngfdofs, ngdofs;
90 void CheckNDSTriaDofs();
112 void GenerateGlobalOffsets()
const;
115 void ConstructTrueDofs();
116 void ConstructTrueNURBSDofs();
119 void ApplyLDofSigns(
Table &el_dof)
const;
132 int PackDof(
int entity,
int index,
int edof)
const;
133 void UnpackDof(
int dof,
int &entity,
int &
index,
int &edof)
const;
135 #ifdef MFEM_PMATRIX_STATS 136 mutable int n_msgs_sent, n_msgs_recv;
137 mutable int n_rows_sent, n_rows_recv, n_rows_fwd;
140 void ScheduleSendRow(
const struct PMatrixRow &row,
int dof, GroupId group_id,
141 std::map<int, class NeighborRowMessage> &send_msg)
const;
143 void ForwardRow(
const struct PMatrixRow &row,
int dof,
144 GroupId group_sent_id, GroupId group_id,
145 std::map<int, class NeighborRowMessage> &send_msg)
const;
147 #ifdef MFEM_DEBUG_PMATRIX 148 void DebugDumpDOFs(std::ostream &os,
157 MakeVDimHypreMatrix(
const std::vector<struct PMatrixRow> &rows,
158 int local_rows,
int local_cols,
163 void Build_Dof_TrueDof_Matrix()
const;
173 bool partial =
false)
const;
179 const Table* old_elem_dof,
180 const Table* old_elem_fos);
187 const Table *old_elem_dof,
188 const Table *old_elem_fos);
193 virtual void UpdateMeshPointer(
Mesh *new_mesh);
253 const int *partitioning,
318 {
if (!P) { Build_Dof_TrueDof_Matrix(); }
return P; }
352 int component = -1)
const;
410 virtual void Update(
bool want_transform =
true);
503 #endif // MFEM_USE_MPI Abstract class for all finite elements.
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
HYPRE_BigInt GetMyTDofOffset() const
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
HypreParVector * NewTrueDofVector()
const FiniteElement * GetFaceNbrFaceFE(int i) const
HYPRE_BigInt GetGlobalScalarTDofNumber(int sldof)
const FiniteElement * GetFaceNbrFE(int i) const
int TrueVSize() const
Obsolete, kept for backward compatibility.
virtual const Operator * GetRestrictionOperator() const
virtual void Update(bool want_transform=true)
Pointer to an Operator of a specified type.
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
HYPRE_BigInt GetMyDofOffset() const
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Abstract parallel finite element space.
int GetFaceNbrVSize() const
void Lose_Dof_TrueDof_Matrix()
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
const FiniteElementCollection * fec
Associated FE collection (not owned).
void DeleteAll()
Delete the whole array.
A parallel extension of the NCMesh class.
std::function< double(const Vector &)> f(double mass_coeff)
DofTransformation * GetFaceNbrElementVDofs(int i, Array< int > &vdofs) 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.
Table face_nbr_element_fos
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
void ExchangeFaceNbrData()
void LoseTrueDofOffsets()
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
void LoseData()
NULL-ifies the data.
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...
ParMesh * GetParMesh() const
Identifies a vertex/edge/face in both Mesh and NCMesh.
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.
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
ElementTransformation * GetFaceNbrElementTransformation(int i) const
HYPRE_BigInt GlobalTrueVSize() const
virtual DofTransformation * GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
Wrapper for hypre's parallel vector class.
virtual int GetFaceDofs(int i, Array< int > &dofs, int variant=0) const
virtual const FiniteElement * GetFE(int i) const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
HYPRE_BigInt GlobalVSize() const
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
int GetLocalTDofNumber(int ldof) const
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.
Array< HYPRE_BigInt > face_nbr_glob_dof_map
HYPRE_BigInt * GetDofOffsets() const
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Table face_nbr_element_dof
void PrintPartitionStats()
virtual ~ParFiniteElementSpace()
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const
int index(int i, int j, int nx, int ny)
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
bool Nonconforming() const
bool SharedNDTriangleDofs() const
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
NURBSExtension * NURBSext
Base class for operators that extracts Face degrees of freedom.
Wrapper for hypre's ParCSR matrix class.
virtual void UpdatesFinished()
Free ParGridFunction transformation matrix (if any), to save memory.
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
HYPRE_BigInt GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
Class for parallel meshes.
HYPRE_BigInt * GetTrueDofOffsets() const
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
virtual DofTransformation * GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i'th boundary element.
int VDofToDof(int vdof) const
Compute the inverse of the Dof to VDof mapping for a single index vdof.