MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pfespace.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_PFESPACE
13 #define MFEM_PFESPACE
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include "../linalg/hypre.hpp"
20 #include "../mesh/pmesh.hpp"
21 #include "../mesh/nurbs.hpp"
22 #include "fespace.hpp"
23 
24 namespace mfem
25 {
26 
27 /// Abstract parallel finite element space.
29 {
30 private:
31  /// MPI data.
32  MPI_Comm MyComm;
33  int NRanks, MyRank;
34 
35  /// Parallel mesh; #mesh points to this object as well. Not owned.
36  ParMesh *pmesh;
37  /** Parallel non-conforming mesh extension object; same as pmesh->pncmesh.
38  Not owned. */
39  ParNCMesh *pncmesh;
40 
41  /// GroupCommunicator on the local VDofs. Owned.
42  GroupCommunicator *gcomm;
43 
44  /// Number of true dofs in this processor (local true dofs).
45  mutable int ltdof_size;
46 
47  /// Number of vertex/edge/face/total ghost DOFs (nonconforming case).
48  int ngvdofs, ngedofs, ngfdofs, ngdofs;
49 
50  /// The group of each local dof.
51  Array<int> ldof_group;
52 
53  /// For a local dof: the local true dof number in the master of its group.
54  mutable Array<int> ldof_ltdof;
55 
56  /// Offsets for the dofs in each processor in global numbering.
57  mutable Array<HYPRE_Int> dof_offsets;
58 
59  /// Offsets for the true dofs in each processor in global numbering.
60  mutable Array<HYPRE_Int> tdof_offsets;
61 
62  /// Offsets for the true dofs in neighbor processor in global numbering.
63  mutable Array<HYPRE_Int> tdof_nb_offsets;
64 
65  /// Previous 'dof_offsets' (before Update()), column partition of T.
66  Array<HYPRE_Int> old_dof_offsets;
67 
68  /// The sign of the basis functions at the scalar local dofs.
69  Array<int> ldof_sign;
70 
71  /// The matrix P (interpolation from true dof to dof). Owned.
72  mutable HypreParMatrix *P;
73  /// Optimized action-only prolongation operator for conforming meshes. Owned.
74  mutable class ConformingProlongationOperator *Pconf;
75 
76  /// The (block-diagonal) matrix R (restriction of dof to true dof). Owned.
77  mutable SparseMatrix *R;
78 
79  ParNURBSExtension *pNURBSext() const
80  { return dynamic_cast<ParNURBSExtension *>(NURBSext); }
81 
82  static ParNURBSExtension *MakeLocalNURBSext(
83  const NURBSExtension *globNURBSext, const NURBSExtension *parNURBSext);
84 
85  GroupTopology &GetGroupTopo() const
86  { return (NURBSext) ? pNURBSext()->gtopo : pmesh->gtopo; }
87 
88  // Auxiliary method used in constructors
89  void ParInit(ParMesh *pm);
90 
91  void Construct();
92  void Destroy();
93 
94  // ldof_type = 0 : DOFs communicator, otherwise VDOFs communicator
95  void GetGroupComm(GroupCommunicator &gcomm, int ldof_type,
96  Array<int> *ldof_sign = NULL);
97 
98  /// Construct dof_offsets and tdof_offsets using global communication.
99  void GenerateGlobalOffsets() const;
100 
101  /// Construct ldof_group and ldof_ltdof.
102  void ConstructTrueDofs();
103  void ConstructTrueNURBSDofs();
104 
105  void ApplyLDofSigns(Array<int> &dofs) const;
106  void ApplyLDofSigns(Table &el_dof) const;
107 
108  typedef NCMesh::MeshId MeshId;
109  typedef ParNCMesh::GroupId GroupId;
110 
111  void GetGhostVertexDofs(const MeshId &id, Array<int> &dofs) const;
112  void GetGhostEdgeDofs(const MeshId &edge_id, Array<int> &dofs) const;
113  void GetGhostFaceDofs(const MeshId &face_id, Array<int> &dofs) const;
114 
115  void GetGhostDofs(int entity, const MeshId &id, Array<int> &dofs) const;
116  // Return the dofs associated with the interior of the given mesh entity.
117  // The MeshId may be the id of a regular or a ghost mesh entity.
118  void GetBareDofs(int entity, const MeshId &id, Array<int> &dofs) const;
119 
120  int PackDof(int entity, int index, int edof) const;
121  void UnpackDof(int dof, int &entity, int &index, int &edof) const;
122 
123  void ScheduleSendRow(const struct PMatrixRow &row, int dof, GroupId group_id,
124  std::map<int, class NeighborRowMessage> &send_msg) const;
125 
126  void ForwardRow(const struct PMatrixRow &row, int dof,
127  GroupId group_sent_id, GroupId group_id,
128  std::map<int, class NeighborRowMessage> &send_msg) const;
129 
130 #ifdef MFEM_DEBUG_PMATRIX
131  void DebugDumpDOFs(std::ofstream &os,
132  const SparseMatrix &deps,
133  const Array<GroupId> &dof_group,
134  const Array<GroupId> &dof_owner,
135  const Array<bool> &finalized) const;
136 #endif
137 
138  /// Helper: create a HypreParMatrix from a list of PMatrixRows.
140  MakeVDimHypreMatrix(const std::vector<struct PMatrixRow> &rows,
141  int local_rows, int local_cols,
142  Array<HYPRE_Int> &row_starts,
143  Array<HYPRE_Int> &col_starts) const;
144 
145  /// Build the P and R matrices.
146  void Build_Dof_TrueDof_Matrix() const;
147 
148  /** Used when the ParMesh is non-conforming, i.e. pmesh->pncmesh != NULL.
149  Constructs the matrices P and R, the DOF and true DOF offset arrays,
150  and the DOF -> true DOF map ('dof_tdof'). Returns the number of
151  vector true DOFs. All pointer arguments are optional and can be NULL. */
152  int BuildParallelConformingInterpolation(HypreParMatrix **P, SparseMatrix **R,
153  Array<HYPRE_Int> &dof_offs,
154  Array<HYPRE_Int> &tdof_offs,
155  Array<int> *dof_tdof,
156  bool partial = false) const;
157 
158  /** Calculate a GridFunction migration matrix after mesh load balancing.
159  The result is a parallel permutation matrix that can be used to update
160  all grid functions defined on this space. */
161  HypreParMatrix* RebalanceMatrix(int old_ndofs,
162  const Table* old_elem_dof);
163 
164  /** Calculate a GridFunction restriction matrix after mesh derefinement.
165  The matrix is constructed so that the new grid function interpolates
166  the original function, i.e., the original function is evaluated at the
167  nodes of the coarse function. */
168  HypreParMatrix* ParallelDerefinementMatrix(int old_ndofs,
169  const Table *old_elem_dof);
170 
171 public:
172  // Face-neighbor data
173  // Number of face-neighbor dofs
175  // Face-neighbor-element to face-neighbor dof
177  // Face-neighbor to ldof in the face-neighbor numbering
179  // The global ldof indices of the face-neighbor dofs
181  // Local face-neighbor data: face-neighbor to ldof
183 
184  /** @brief Copy constructor: deep copy all data from @a orig except the
185  ParMesh, the FiniteElementCollection, and some derived data. */
186  /** If the @a pmesh or @a fec poiters are NULL (default), then the new
187  ParFiniteElementSpace will reuse the respective pointers from @a orig. If
188  any of these pointers is not NULL, the given pointer will be used instead
189  of the one used by @a orig.
190 
191  @note The objects pointed to by the @a pmesh and @a fec parameters must
192  be either the same objects as the ones used by @a orig, or copies of
193  them. Otherwise, the behavior is undefined.
194 
195  @note Derived data objects, such as the parallel prolongation and
196  restriction operators, the update operator, and any of the face-neighbor
197  data, will not be copied, even if they are created in the @a orig object.
198  */
200  ParMesh *pmesh = NULL,
201  const FiniteElementCollection *fec = NULL);
202 
203  /** @brief Convert/copy the *local* (Par)FiniteElementSpace @a orig to
204  ParFiniteElementSpace: deep copy all data from @a orig except the Mesh,
205  the FiniteElementCollection, and some derived data. */
207  const FiniteElementCollection *fec = NULL);
208 
209  /** @brief Construct the *local* ParFiniteElementSpace corresponing to the
210  global FE space, @a global_fes. */
211  /** The parameter @a pm is the *local* ParMesh obtained by decomposing the
212  global Mesh used by @a global_fes. The array @a partitioning represents
213  the parallel decomposition - it maps global element ids to MPI ranks.
214  If the FiniteElementCollection, @a f, is NULL (default), the FE
215  collection used by @a global_fes will be reused. If @a f is not NULL, it
216  must be the same as, or a copy of, the FE collection used by
217  @a global_fes. */
218  ParFiniteElementSpace(ParMesh *pm, const FiniteElementSpace *global_fes,
219  const int *partitioning,
220  const FiniteElementCollection *f = NULL);
221 
223  int dim = 1, int ordering = Ordering::byNODES);
224 
225  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
226  /** The parameter @a ext will be deleted by this constructor, replaced by a
227  ParNURBSExtension owned by the ParFiniteElementSpace.
228  @note If the pointer @a ext is NULL, this constructor is equivalent to
229  the standard constructor with the same arguments minus the
230  NURBSExtension, @a ext. */
232  const FiniteElementCollection *f,
233  int dim = 1, int ordering = Ordering::byNODES);
234 
235  MPI_Comm GetComm() const { return MyComm; }
236  int GetNRanks() const { return NRanks; }
237  int GetMyRank() const { return MyRank; }
238 
239  inline ParMesh *GetParMesh() { return pmesh; }
240 
241  int GetDofSign(int i)
242  { return NURBSext || Nonconforming() ? 1 : ldof_sign[VDofToDof(i)]; }
243  HYPRE_Int *GetDofOffsets() const { return dof_offsets; }
244  HYPRE_Int *GetTrueDofOffsets() const { return tdof_offsets; }
245  HYPRE_Int GlobalVSize() const
246  { return Dof_TrueDof_Matrix()->GetGlobalNumRows(); }
247  HYPRE_Int GlobalTrueVSize() const
248  { return Dof_TrueDof_Matrix()->GetGlobalNumCols(); }
249 
250  /// Return the number of local vector true dofs.
251  virtual int GetTrueVSize() const { return ltdof_size; }
252 
253  /// Returns indexes of degrees of freedom in array dofs for i'th element.
254  virtual void GetElementDofs(int i, Array<int> &dofs) const;
255 
256  /// Returns indexes of degrees of freedom for i'th boundary element.
257  virtual void GetBdrElementDofs(int i, Array<int> &dofs) const;
258 
259  /** Returns the indexes of the degrees of freedom for i'th face
260  including the dofs for the edges and the vertices of the face. */
261  virtual void GetFaceDofs(int i, Array<int> &dofs) const;
262 
263  void GetSharedEdgeDofs(int group, int ei, Array<int> &dofs) const;
264  void GetSharedFaceDofs(int group, int fi, Array<int> &dofs) const;
265 
266  /// The true dof-to-dof interpolation matrix
268  { if (!P) { Build_Dof_TrueDof_Matrix(); } return P; }
269 
270  /** @brief For a non-conforming mesh, construct and return the interpolation
271  matrix from the partially conforming true dofs to the local dofs. */
272  /** @note The returned pointer must be deleted by the caller. */
274 
275  /** Create and return a new HypreParVector on the true dofs, which is
276  owned by (i.e. it must be destroyed by) the calling function. */
278  { return (new HypreParVector(MyComm,GlobalTrueVSize(),GetTrueDofOffsets()));}
279 
280  /// Scale a vector of true dofs
281  void DivideByGroupSize(double *vec);
282 
283  /// Return a reference to the internal GroupCommunicator (on VDofs)
284  GroupCommunicator &GroupComm() { return *gcomm; }
285 
286  /// Return a const reference to the internal GroupCommunicator (on VDofs)
287  const GroupCommunicator &GroupComm() const { return *gcomm; }
288 
289  /// Return a new GroupCommunicator on Dofs
291 
292  /** Given an integer array on the local degrees of freedom, perform
293  a bitwise OR between the shared dofs. */
294  void Synchronize(Array<int> &ldof_marker) const;
295 
296  /// Determine the boundary degrees of freedom
297  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
298  Array<int> &ess_dofs,
299  int component = -1) const;
300 
301  /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
302  boundary attributes marked in the array bdr_attr_is_ess. */
303  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
304  Array<int> &ess_tdof_list,
305  int component = -1);
306 
307  /** If the given ldof is owned by the current processor, return its local
308  tdof number, otherwise return -1 */
309  int GetLocalTDofNumber(int ldof) const;
310  /// Returns the global tdof number of the given local degree of freedom
311  HYPRE_Int GetGlobalTDofNumber(int ldof) const;
312  /** Returns the global tdof number of the given local degree of freedom in
313  the scalar version of the current finite element space. The input should
314  be a scalar local dof. */
315  HYPRE_Int GetGlobalScalarTDofNumber(int sldof);
316 
317  HYPRE_Int GetMyDofOffset() const;
318  HYPRE_Int GetMyTDofOffset() const;
319 
320  virtual const Operator *GetProlongationMatrix() const;
321  /// Get the R matrix which restricts a local dof vector to true dof vector.
322  virtual const SparseMatrix *GetRestrictionMatrix() const
323  { Dof_TrueDof_Matrix(); return R; }
324 
325  // Face-neighbor functions
326  void ExchangeFaceNbrData();
327  int GetFaceNbrVSize() const { return num_face_nbr_dofs; }
328  void GetFaceNbrElementVDofs(int i, Array<int> &vdofs) const;
329  void GetFaceNbrFaceVDofs(int i, Array<int> &vdofs) const;
330  const FiniteElement *GetFaceNbrFE(int i) const;
331  const FiniteElement *GetFaceNbrFaceFE(int i) const;
332  const HYPRE_Int *GetFaceNbrGlobalDofMap() { return face_nbr_glob_dof_map; }
333 
335  void LoseDofOffsets() { dof_offsets.LoseData(); }
336  void LoseTrueDofOffsets() { tdof_offsets.LoseData(); }
337 
338  bool Conforming() const { return pmesh->pncmesh == NULL; }
339  bool Nonconforming() const { return pmesh->pncmesh != NULL; }
340 
341  // Transfer parallel true-dof data from coarse_fes, defined on a coarse mesh,
342  // to this FE space, defined on a refined mesh. See full documentation in the
343  // base class, FiniteElementSpace::GetTrueTransferOperator.
344  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
345  OperatorHandle &T) const;
346 
347  /** Reflect changes in the mesh. Calculate one of the refinement/derefinement
348  /rebalance matrices, unless want_transform is false. */
349  virtual void Update(bool want_transform = true);
350 
351  /// Free ParGridFunction transformation matrix (if any), to save memory.
352  virtual void UpdatesFinished()
353  {
355  old_dof_offsets.DeleteAll();
356  }
357 
358  virtual ~ParFiniteElementSpace() { Destroy(); }
359 
360  // Obsolete, kept for backward compatibility
361  int TrueVSize() const { return ltdof_size; }
362 };
363 
364 
365 /// Auxiliary class used by ParFiniteElementSpace.
367 {
368 protected:
371 
372 public:
374 
375  virtual void Mult(const Vector &x, Vector &y) const;
376 
377  virtual void MultTranspose(const Vector &x, Vector &y) const;
378 };
379 
380 }
381 
382 #endif // MFEM_USE_MPI
383 
384 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:140
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:626
HYPRE_Int GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
Definition: pfespace.cpp:672
HypreParVector * NewTrueDofVector()
Definition: pfespace.hpp:277
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const
Determine the boundary degrees of freedom.
Definition: pfespace.cpp:612
ConformingProlongationOperator(const ParFiniteElementSpace &pfes)
Definition: pfespace.cpp:2656
HYPRE_Int * GetDofOffsets() const
Definition: pfespace.hpp:243
virtual const Operator * GetProlongationMatrix() const
Definition: pfespace.cpp:742
Auxiliary class used by ParFiniteElementSpace.
Definition: pfespace.hpp:366
Ordering::Type ordering
Definition: fespace.hpp:81
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:2568
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int VDofToDof(int vdof) const
Definition: fespace.hpp:348
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:287
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.
Definition: pfespace.cpp:2541
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:557
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void Synchronize(Array< int > &ldof_marker) const
Definition: pfespace.cpp:595
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:651
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
Definition: pfespace.cpp:546
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:73
const HYPRE_Int * GetFaceNbrGlobalDofMap()
Definition: pfespace.hpp:332
HYPRE_Int GetMyDofOffset() const
Definition: pfespace.cpp:732
void DeleteAll()
Delete whole array.
Definition: array.hpp:709
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:61
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
Definition: pfespace.cpp:398
ParNCMesh * pncmesh
Definition: pmesh.hpp:140
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:251
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: pfespace.cpp:375
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:974
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:379
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1008
void GetSharedFaceDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:421
int dim
Definition: ex3.cpp:47
MPI_Comm GetComm() const
Definition: pfespace.hpp:235
HYPRE_Int GetGlobalNumCols() const
Definition: hypre.hpp:382
Data type sparse matrix.
Definition: sparsemat.hpp:38
const GroupCommunicator & gc
Definition: pfespace.hpp:370
void LoseData()
NULL-ifies the data.
Definition: array.hpp:127
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:267
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:247
HYPRE_Int GetMyTDofOffset() const
Definition: pfespace.cpp:737
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:134
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: pfespace.cpp:389
HYPRE_Int GlobalVSize() const
Definition: pfespace.hpp:245
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: pfespace.cpp:2724
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:284
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.
Definition: pfespace.cpp:28
HYPRE_Int GetGlobalScalarTDofNumber(int sldof)
Definition: pfespace.cpp:695
bool Nonconforming() const
Definition: pfespace.hpp:339
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: pfespace.cpp:361
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:322
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:2699
const FiniteElement * GetFaceNbrFaceFE(int i) const
Definition: pfespace.cpp:1022
int GetFaceNbrVSize() const
Definition: pfespace.hpp:327
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:244
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:980
Vector data type.
Definition: vector.hpp:48
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:514
NURBSExtension * NURBSext
Definition: fespace.hpp:94
GroupTopology gtopo
Definition: pmesh.hpp:127
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
virtual void UpdatesFinished()
Free ParGridFunction transformation matrix (if any), to save memory.
Definition: pfespace.hpp:352
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on Dofs.
Definition: pfespace.cpp:575
Class for parallel meshes.
Definition: pmesh.hpp:32
GroupTopology gtopo
Definition: nurbs.hpp:386
Array< HYPRE_Int > face_nbr_glob_dof_map
Definition: pfespace.hpp:180