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 void UpdateMeshPointer(
Mesh *new_mesh)
override;
257 const int *partitioning,
328 {
if (!P) { Build_Dof_TrueDof_Matrix(); }
return P; }
362 int component = -1)
const override;
368 int component = -1)
const override;
422 void Update(
bool want_transform =
true)
override;
void LoseData()
NULL-ifies the data.
void DeleteAll()
Delete the whole array.
Base class for operators that extracts Face degrees of freedom.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
NURBSExtension * NURBSext
const FiniteElementCollection * fec
Associated FE collection (not owned).
int VDofToDof(int vdof) const
Compute the inverse of the Dof to VDof mapping for a single index vdof.
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
DofTransformation * GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element 'bel'. The returned indices are offsets in...
Abstract class for all finite elements.
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
Wrapper for hypre's ParCSR matrix class.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
Wrapper for hypre's parallel vector class.
Class used by MFEM to store pointers to host and/or device memory.
Pointer to an Operator of a specified type.
Abstract parallel finite element space.
void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const override
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes,...
HYPRE_BigInt GetGlobalScalarTDofNumber(int sldof)
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const override
HYPRE_BigInt * GetTrueDofOffsets() const
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
void PrintPartitionStats()
HYPRE_BigInt GlobalVSize() const
const Operator * GetRestrictionOperator() const override
int GetLocalTDofNumber(int ldof) const
void UpdatesFinished() override
Free ParGridFunction transformation matrix (if any), to save memory.
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.
void LoseTrueDofOffsets()
bool SharedNDTriangleDofs() 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,...
void DivideByGroupSize(real_t *vec)
Scale a vector of true dofs.
void ExchangeFaceNbrData()
HYPRE_BigInt GlobalTrueVSize() const
HypreParVector * NewTrueDofVector()
int GetTrueVSize() const override
Return the number of local vector true dofs.
HYPRE_BigInt GetMyDofOffset() const
HYPRE_BigInt * GetDofOffsets() const
Array< HYPRE_BigInt > face_nbr_glob_dof_map
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
bool Nonconforming() const
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
const Operator * GetProlongationMatrix() const override
The returned Operator is owned by the FiniteElementSpace.
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const override
Determine the boundary degrees of freedom.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs, DofTransformation &doftrans) const
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
int GetFaceDofs(int i, Array< int > &dofs, int variant=0) const override
int GetFaceNbrVSize() const
Table face_nbr_element_fos
HYPRE_BigInt GetMyTDofOffset() const
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
void GetElementDofs(int i, Array< int > &dofs, DofTransformation &doftrans) const override
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be al...
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
HYPRE_BigInt GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
virtual ~ParFiniteElementSpace()
ParMesh * GetParMesh() const
void Update(bool want_transform=true) override
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
int TrueVSize() const
Obsolete, kept for backward compatibility.
void Lose_Dof_TrueDof_Matrix()
const FiniteElement * GetFaceNbrFaceFE(int i) const
const FiniteElement * GetFaceNbrFE(int i) const
const FiniteElement * GetFE(int i) const override
Table face_nbr_element_dof
ElementTransformation * GetFaceNbrElementTransformation(int i) const
void GetBdrElementDofs(int i, Array< int > &dofs, DofTransformation &doftrans) const override
The same as GetBdrElementDofs(), but with a user-allocated DofTransformation object....
Class for parallel meshes.
ElementTransformation * GetFaceNbrElementTransformation(int FaceNo)
Returns a pointer to the transformation defining the i-th face neighbor.
A parallel extension of the NCMesh class.
int index(int i, int j, int nx, int ny)
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Identifies a vertex/edge/face in both Mesh and NCMesh.