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;
99 void GenerateGlobalOffsets()
const;
102 void ConstructTrueDofs();
103 void ConstructTrueNURBSDofs();
106 void ApplyLDofSigns(
Table &el_dof)
const;
119 int PackDof(
int entity,
int index,
int edof)
const;
120 void UnpackDof(
int dof,
int &entity,
int &index,
int &edof)
const;
122 #ifdef MFEM_PMATRIX_STATS
123 mutable int n_msgs_sent, n_msgs_recv;
124 mutable int n_rows_sent, n_rows_recv, n_rows_fwd;
127 void ScheduleSendRow(
const struct PMatrixRow &row,
int dof, GroupId group_id,
128 std::map<int, class NeighborRowMessage> &send_msg)
const;
130 void ForwardRow(
const struct PMatrixRow &row,
int dof,
131 GroupId group_sent_id, GroupId group_id,
132 std::map<int, class NeighborRowMessage> &send_msg)
const;
134 #ifdef MFEM_DEBUG_PMATRIX
135 void DebugDumpDOFs(std::ostream &os,
144 MakeVDimHypreMatrix(
const std::vector<struct PMatrixRow> &rows,
145 int local_rows,
int local_cols,
150 void Build_Dof_TrueDof_Matrix()
const;
160 bool partial =
false)
const;
166 const Table* old_elem_dof);
173 const Table *old_elem_dof);
223 const int *partitioning,
282 {
if (!P) { Build_Dof_TrueDof_Matrix(); }
return P; }
316 int component = -1)
const;
368 virtual void Update(
bool want_transform =
true);
449 #endif // MFEM_USE_MPI
Abstract class for all finite elements.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
HYPRE_Int GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
HypreParVector * NewTrueDofVector()
ParMesh * GetParMesh() 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.
HYPRE_Int * GetDofOffsets() const
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
int TrueVSize() const
Obsolete, kept for backward compatibility.
virtual void Update(bool want_transform=true)
Pointer to an Operator of a specified type.
int VDofToDof(int vdof) const
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
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 DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Abstract parallel finite element space.
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 Lose_Dof_TrueDof_Matrix()
int GetLocalTDofNumber(int ldof) const
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
const FiniteElementCollection * fec
Associated FE collection (not owned).
const HYPRE_Int * GetFaceNbrGlobalDofMap()
HYPRE_Int GetMyDofOffset() const
void DeleteAll()
Delete the whole array.
A parallel extension of the NCMesh class.
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i'th boundary element.
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const FiniteElement * GetFaceNbrFE(int i) const
void ExchangeFaceNbrData()
HYPRE_Int GetGlobalNumCols() const
Return the global number of columns.
void LoseTrueDofOffsets()
ElementTransformation * GetFaceNbrElementTransformation(int i) const
void LoseData()
NULL-ifies the data.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
HYPRE_Int GlobalTrueVSize() const
HYPRE_Int GetMyTDofOffset() const
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
Identifies a vertex/edge/face in both Mesh and NCMesh.
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
virtual void GetFaceDofs(int i, Array< int > &dofs) const
HYPRE_Int GlobalVSize() const
Wrapper for hypre's parallel vector class.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
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)
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.
HYPRE_Int GetGlobalScalarTDofNumber(int sldof)
virtual const Operator * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const
bool Nonconforming() const
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
Table face_nbr_element_dof
void PrintPartitionStats()
virtual ~ParFiniteElementSpace()
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
int index(int i, int j, int nx, int ny)
const FiniteElement * GetFaceNbrFaceFE(int i) const
int GetFaceNbrVSize() const
HYPRE_Int * GetTrueDofOffsets() const
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
NURBSExtension * NURBSext
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.
Class for parallel meshes.
double f(const Vector &p)
Array< HYPRE_Int > face_nbr_glob_dof_map