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;
94 void CheckNDSTriaDofs();
116 void GenerateGlobalOffsets()
const;
119 void ConstructTrueDofs();
120 void ConstructTrueNURBSDofs();
123 void ApplyLDofSigns(
Table &el_dof)
const;
136 int PackDof(
int entity,
int index,
int edof)
const;
137 void UnpackDof(
int dof,
int &entity,
int &index,
int &edof)
const;
139 #ifdef MFEM_PMATRIX_STATS
140 mutable int n_msgs_sent, n_msgs_recv;
141 mutable int n_rows_sent, n_rows_recv, n_rows_fwd;
144 void ScheduleSendRow(
const struct PMatrixRow &row,
int dof, GroupId group_id,
145 std::map<int, class NeighborRowMessage> &send_msg)
const;
147 void ForwardRow(
const struct PMatrixRow &row,
int dof,
148 GroupId group_sent_id, GroupId group_id,
149 std::map<int, class NeighborRowMessage> &send_msg)
const;
151 #ifdef MFEM_DEBUG_PMATRIX
152 void DebugDumpDOFs(std::ostream &os,
161 MakeVDimHypreMatrix(
const std::vector<struct PMatrixRow> &rows,
162 int local_rows,
int local_cols,
167 void Build_Dof_TrueDof_Matrix()
const;
177 bool partial =
false)
const;
183 const Table* old_elem_dof,
184 const Table* old_elem_fos);
191 const Table *old_elem_dof,
192 const Table *old_elem_fos);
197 virtual void UpdateMeshPointer(
Mesh *new_mesh);
257 const int *partitioning,
322 {
if (!P) { Build_Dof_TrueDof_Matrix(); }
return P; }
356 int component = -1)
const;
420 virtual void Update(
bool want_transform =
true);
513 #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_BigInt GetMyTDofOffset() const
HypreParVector * NewTrueDofVector()
HYPRE_BigInt GetGlobalScalarTDofNumber(int sldof)
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.
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
HYPRE_BigInt GlobalTrueVSize() const
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
virtual DofTransformation * GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i'th boundary element.
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).
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.
Table face_nbr_element_fos
HYPRE_BigInt GlobalVSize() const
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const FiniteElement * GetFaceNbrFE(int i) const
void ExchangeFaceNbrData()
double f(const Vector &xvec)
virtual const FiniteElement * GetFE(int i) const
void LoseTrueDofOffsets()
ElementTransformation * GetFaceNbrElementTransformation(int i) const
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
void LoseData()
NULL-ifies the data.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
Identifies a vertex/edge/face in both Mesh and NCMesh.
HYPRE_BigInt GetMyDofOffset() const
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
HYPRE_BigInt * GetDofOffsets() const
virtual DofTransformation * GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
Wrapper for hypre's parallel vector class.
HYPRE_BigInt * GetTrueDofOffsets() const
HYPRE_BigInt GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
bool SharedNDTriangleDofs() 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)
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.
bool Nonconforming() const
Array< HYPRE_BigInt > face_nbr_glob_dof_map
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
virtual const Operator * GetRestrictionOperator() const
Table face_nbr_element_dof
void PrintPartitionStats()
virtual ~ParFiniteElementSpace()
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
virtual int GetFaceDofs(int i, Array< int > &dofs, int variant=0) const
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const
int index(int i, int j, int nx, int ny)
const FiniteElement * GetFaceNbrFaceFE(int i) const
int GetFaceNbrVSize() const
virtual const Operator * GetRestrictionTransposeOperator() const
Return logical transpose of restriction matrix, but in non-assembled optimized matrix-free form...
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
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.
Class for parallel meshes.
DofTransformation * GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const