MFEM  v4.6.0
Finite element discretization library
pfespace.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
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_BigInt> dof_offsets;
58 
59  /// Offsets for the true dofs in each processor in global numbering.
60  mutable Array<HYPRE_BigInt> tdof_offsets;
61 
62  /// Offsets for the true dofs in neighbor processor in global numbering.
63  mutable Array<HYPRE_BigInt> tdof_nb_offsets;
64 
65  /// Previous 'dof_offsets' (before Update()), column partition of T.
66  Array<HYPRE_BigInt> 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 Operator *Pconf;
75 
76  /** Used to indicate that the space is nonconforming (even if the underlying
77  mesh has NULL @a ncmesh). This occurs in low-order preconditioning on
78  nonconforming meshes. */
79  bool nonconf_P;
80 
81  /// The (block-diagonal) matrix R (restriction of dof to true dof). Owned.
82  mutable SparseMatrix *R;
83  /// Optimized action-only restriction operator for conforming meshes. Owned.
84  mutable Operator *Rconf;
85 
86  /// Flag indicating the existence of shared triangles with interior ND dofs
87  bool nd_strias;
88 
89  /// Resets nd_strias flag at construction or after rebalancing
90  void CheckNDSTriaDofs();
91 
92  ParNURBSExtension *pNURBSext() const
93  { return dynamic_cast<ParNURBSExtension *>(NURBSext); }
94 
95  static ParNURBSExtension *MakeLocalNURBSext(
96  const NURBSExtension *globNURBSext, const NURBSExtension *parNURBSext);
97 
98  GroupTopology &GetGroupTopo() const
99  { return (NURBSext) ? pNURBSext()->gtopo : pmesh->gtopo; }
100 
101  // Auxiliary method used in constructors
102  void ParInit(ParMesh *pm);
103 
104  void Construct();
105  void Destroy();
106 
107  // ldof_type = 0 : DOFs communicator, otherwise VDOFs communicator
108  void GetGroupComm(GroupCommunicator &gcomm, int ldof_type,
109  Array<int> *ldof_sign = NULL);
110 
111  /// Construct dof_offsets and tdof_offsets using global communication.
112  void GenerateGlobalOffsets() const;
113 
114  /// Construct ldof_group and ldof_ltdof.
115  void ConstructTrueDofs();
116  void ConstructTrueNURBSDofs();
117 
118  void ApplyLDofSigns(Array<int> &dofs) const;
119  void ApplyLDofSigns(Table &el_dof) const;
120 
121  typedef NCMesh::MeshId MeshId;
122  typedef ParNCMesh::GroupId GroupId;
123 
124  void GetGhostVertexDofs(const MeshId &id, Array<int> &dofs) const;
125  void GetGhostEdgeDofs(const MeshId &edge_id, Array<int> &dofs) const;
126  void GetGhostFaceDofs(const MeshId &face_id, Array<int> &dofs) const;
127 
128  void GetGhostDofs(int entity, const MeshId &id, Array<int> &dofs) const;
129  /// Return the dofs associated with the interior of the given mesh entity.
130  void GetBareDofs(int entity, int index, Array<int> &dofs) const;
131 
132  int PackDof(int entity, int index, int edof) const;
133  void UnpackDof(int dof, int &entity, int &index, int &edof) const;
134 
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;
138 #endif
139 
140  void ScheduleSendRow(const struct PMatrixRow &row, int dof, GroupId group_id,
141  std::map<int, class NeighborRowMessage> &send_msg) const;
142 
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;
146 
147 #ifdef MFEM_DEBUG_PMATRIX
148  void DebugDumpDOFs(std::ostream &os,
149  const SparseMatrix &deps,
150  const Array<GroupId> &dof_group,
151  const Array<GroupId> &dof_owner,
152  const Array<bool> &finalized) const;
153 #endif
154 
155  /// Helper: create a HypreParMatrix from a list of PMatrixRows.
157  MakeVDimHypreMatrix(const std::vector<struct PMatrixRow> &rows,
158  int local_rows, int local_cols,
159  Array<HYPRE_BigInt> &row_starts,
160  Array<HYPRE_BigInt> &col_starts) const;
161 
162  /// Build the P and R matrices.
163  void Build_Dof_TrueDof_Matrix() const;
164 
165  /** Used when the ParMesh is non-conforming, i.e. pmesh->pncmesh != NULL.
166  Constructs the matrices P and R, the DOF and true DOF offset arrays,
167  and the DOF -> true DOF map ('dof_tdof'). Returns the number of
168  vector true DOFs. All pointer arguments are optional and can be NULL. */
169  int BuildParallelConformingInterpolation(HypreParMatrix **P, SparseMatrix **R,
170  Array<HYPRE_BigInt> &dof_offs,
171  Array<HYPRE_BigInt> &tdof_offs,
172  Array<int> *dof_tdof,
173  bool partial = false) const;
174 
175  /** Calculate a GridFunction migration matrix after mesh load balancing.
176  The result is a parallel permutation matrix that can be used to update
177  all grid functions defined on this space. */
178  HypreParMatrix* RebalanceMatrix(int old_ndofs,
179  const Table* old_elem_dof,
180  const Table* old_elem_fos);
181 
182  /** Calculate a GridFunction restriction matrix after mesh derefinement.
183  The matrix is constructed so that the new grid function interpolates
184  the original function, i.e., the original function is evaluated at the
185  nodes of the coarse function. */
186  HypreParMatrix* ParallelDerefinementMatrix(int old_ndofs,
187  const Table *old_elem_dof,
188  const Table *old_elem_fos);
189 
190  /// Updates the internal mesh pointer. @warning @a new_mesh must be
191  /// <b>topologically identical</b> to the existing mesh. Used if the address
192  /// of the Mesh object has changed, e.g. in @a Mesh::Swap.
193  virtual void UpdateMeshPointer(Mesh *new_mesh);
194 
195  /// Copies the prolongation and restriction matrices from @a fes.
196  ///
197  /// Used for low order preconditioning on non-conforming meshes. If the DOFs
198  /// require a permutation, it will be supplied by non-NULL @a perm. NULL @a
199  /// perm indicates that no permutation is required.
200  virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes,
201  const Array<int> *perm);
202 
203 public:
204  // Face-neighbor data
205  // Number of face-neighbor dofs
207  // Face-neighbor-element to face-neighbor dof
209  // Face-neighbor-element face orientations
211  // Face-neighbor to ldof in the face-neighbor numbering
213  // The global ldof indices of the face-neighbor dofs
215  // Local face-neighbor data: face-neighbor to ldof
217 
218  /** @brief Copy constructor: deep copy all data from @a orig except the
219  ParMesh, the FiniteElementCollection, and some derived data. */
220  /** If the @a pmesh or @a fec pointers are NULL (default), then the new
221  ParFiniteElementSpace will reuse the respective pointers from @a orig. If
222  any of these pointers is not NULL, the given pointer will be used instead
223  of the one used by @a orig.
224 
225  @note The objects pointed to by the @a pmesh and @a fec parameters must
226  be either the same objects as the ones used by @a orig, or copies of
227  them. Otherwise, the behavior is undefined.
228 
229  @note Derived data objects, such as the parallel prolongation and
230  restriction operators, the update operator, and any of the face-neighbor
231  data, will not be copied, even if they are created in the @a orig object.
232  */
234  ParMesh *pmesh = NULL,
235  const FiniteElementCollection *fec = NULL);
236 
237  /** @brief Convert/copy the *local* (Par)FiniteElementSpace @a orig to
238  ParFiniteElementSpace: deep copy all data from @a orig except the Mesh,
239  the FiniteElementCollection, and some derived data. */
241  const FiniteElementCollection *fec = NULL);
242 
243  /** @brief Construct the *local* ParFiniteElementSpace corresponding to the
244  global FE space, @a global_fes. */
245  /** The parameter @a pm is the *local* ParMesh obtained by decomposing the
246  global Mesh used by @a global_fes. The array @a partitioning represents
247  the parallel decomposition - it maps global element ids to MPI ranks.
248  If the FiniteElementCollection, @a f, is NULL (default), the FE
249  collection used by @a global_fes will be reused. If @a f is not NULL, it
250  must be the same as, or a copy of, the FE collection used by
251  @a global_fes. */
252  ParFiniteElementSpace(ParMesh *pm, const FiniteElementSpace *global_fes,
253  const int *partitioning,
254  const FiniteElementCollection *f = NULL);
255 
257  int dim = 1, int ordering = Ordering::byNODES);
258 
259  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
260  /** The parameter @a ext will be deleted by this constructor, replaced by a
261  ParNURBSExtension owned by the ParFiniteElementSpace.
262  @note If the pointer @a ext is NULL, this constructor is equivalent to
263  the standard constructor with the same arguments minus the
264  NURBSExtension, @a ext. */
266  const FiniteElementCollection *f,
267  int dim = 1, int ordering = Ordering::byNODES);
268 
269  MPI_Comm GetComm() const { return MyComm; }
270  int GetNRanks() const { return NRanks; }
271  int GetMyRank() const { return MyRank; }
272 
273  inline ParMesh *GetParMesh() const { return pmesh; }
274 
275  int GetDofSign(int i)
276  { return NURBSext || Nonconforming() ? 1 : ldof_sign[VDofToDof(i)]; }
277  HYPRE_BigInt *GetDofOffsets() const { return dof_offsets; }
278  HYPRE_BigInt *GetTrueDofOffsets() const { return tdof_offsets; }
280  { return Dof_TrueDof_Matrix()->GetGlobalNumRows(); }
282  { return Dof_TrueDof_Matrix()->GetGlobalNumCols(); }
283 
284  /// Return the number of local vector true dofs.
285  virtual int GetTrueVSize() const { return ltdof_size; }
286 
287  /// Returns indexes of degrees of freedom in array dofs for i'th element.
288  virtual DofTransformation *GetElementDofs(int i, Array<int> &dofs) const;
289 
290  /// Returns indexes of degrees of freedom for i'th boundary element.
291  virtual DofTransformation *GetBdrElementDofs(int i, Array<int> &dofs) const;
292 
293  /** Returns the indexes of the degrees of freedom for i'th face
294  including the dofs for the edges and the vertices of the face. */
295  virtual int GetFaceDofs(int i, Array<int> &dofs, int variant = 0) const;
296 
297  /** Returns pointer to the FiniteElement in the FiniteElementCollection
298  associated with i'th element in the mesh object. If @a i is greater than
299  or equal to the number of local mesh elements, @a i will be interpreted
300  as a shifted index of a face neighbor element. */
301  virtual const FiniteElement *GetFE(int i) const;
302 
303  /** Returns an Operator that converts L-vectors to E-vectors on each face.
304  The parallel version is different from the serial one because of the
305  presence of shared faces. Shared faces are treated as interior faces,
306  the returned operator handles the communication needed to get the
307  shared face values from other MPI ranks */
308  virtual const FaceRestriction *GetFaceRestriction(
309  ElementDofOrdering f_ordering, FaceType type,
311 
312  void GetSharedEdgeDofs(int group, int ei, Array<int> &dofs) const;
313  void GetSharedTriangleDofs(int group, int fi, Array<int> &dofs) const;
314  void GetSharedQuadrilateralDofs(int group, int fi, Array<int> &dofs) const;
315 
316  /// The true dof-to-dof interpolation matrix
318  { if (!P) { Build_Dof_TrueDof_Matrix(); } return P; }
319 
320  /** @brief For a non-conforming mesh, construct and return the interpolation
321  matrix from the partially conforming true dofs to the local dofs. */
322  /** @note The returned pointer must be deleted by the caller. */
324 
325  /** Create and return a new HypreParVector on the true dofs, which is
326  owned by (i.e. it must be destroyed by) the calling function. */
328  { return (new HypreParVector(MyComm,GlobalTrueVSize(),GetTrueDofOffsets()));}
329 
330  /// Scale a vector of true dofs
331  void DivideByGroupSize(double *vec);
332 
333  /// Return a reference to the internal GroupCommunicator (on VDofs)
334  GroupCommunicator &GroupComm() { return *gcomm; }
335 
336  /// Return a const reference to the internal GroupCommunicator (on VDofs)
337  const GroupCommunicator &GroupComm() const { return *gcomm; }
338 
339  /// Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
340  /** @note The returned pointer must be deleted by the caller. */
342 
343  /** @brief Given an integer array on the local degrees of freedom, perform
344  a bitwise OR between the shared dofs. */
345  /** For non-conforming mesh, synchronization is performed on the cut (aka
346  "partially conforming") space. */
347  void Synchronize(Array<int> &ldof_marker) const;
348 
349  /// Determine the boundary degrees of freedom
350  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
351  Array<int> &ess_dofs,
352  int component = -1) const;
353 
354  /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
355  boundary attributes marked in the array bdr_attr_is_ess. */
356  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
357  Array<int> &ess_tdof_list,
358  int component = -1);
359 
360  /** If the given ldof is owned by the current processor, return its local
361  tdof number, otherwise return -1 */
362  int GetLocalTDofNumber(int ldof) const;
363  /// Returns the global tdof number of the given local degree of freedom
364  HYPRE_BigInt GetGlobalTDofNumber(int ldof) const;
365  /** Returns the global tdof number of the given local degree of freedom in
366  the scalar version of the current finite element space. The input should
367  be a scalar local dof. */
369 
372 
373  virtual const Operator *GetProlongationMatrix() const;
374  /** Get an Operator that performs the action of GetRestrictionMatrix(),
375  but potentially with a non-assembled optimized matrix-free
376  implementation. */
377  virtual const Operator *GetRestrictionOperator() const;
378  /// Get the R matrix which restricts a local dof vector to true dof vector.
379  virtual const SparseMatrix *GetRestrictionMatrix() const
380  { Dof_TrueDof_Matrix(); return R; }
381 
382  // Face-neighbor functions
383  void ExchangeFaceNbrData();
384  int GetFaceNbrVSize() const { return num_face_nbr_dofs; }
386  void GetFaceNbrFaceVDofs(int i, Array<int> &vdofs) const;
387  const FiniteElement *GetFaceNbrFE(int i) const;
388  const FiniteElement *GetFaceNbrFaceFE(int i) const;
391  { return pmesh->GetFaceNbrElementTransformation(i); }
392 
394  void LoseDofOffsets() { dof_offsets.LoseData(); }
395  void LoseTrueDofOffsets() { tdof_offsets.LoseData(); }
396 
397  bool Conforming() const { return pmesh->pncmesh == NULL && !nonconf_P; }
398  bool Nonconforming() const { return pmesh->pncmesh != NULL || nonconf_P; }
399 
400  bool SharedNDTriangleDofs() const { return nd_strias; }
401 
402  // Transfer parallel true-dof data from coarse_fes, defined on a coarse mesh,
403  // to this FE space, defined on a refined mesh. See full documentation in the
404  // base class, FiniteElementSpace::GetTrueTransferOperator.
405  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
406  OperatorHandle &T) const;
407 
408  /** Reflect changes in the mesh. Calculate one of the refinement/derefinement
409  /rebalance matrices, unless want_transform is false. */
410  virtual void Update(bool want_transform = true);
411 
412  /// Free ParGridFunction transformation matrix (if any), to save memory.
413  virtual void UpdatesFinished()
414  {
416  old_dof_offsets.DeleteAll();
417  }
418 
419  virtual ~ParFiniteElementSpace() { Destroy(); }
420 
421  void PrintPartitionStats();
422 
423  /// Obsolete, kept for backward compatibility
424  int TrueVSize() const { return ltdof_size; }
425 };
426 
427 
428 /// Auxiliary class used by ParFiniteElementSpace.
430 {
431 protected:
434  bool local;
435 
436 public:
438  bool local_=false);
439 
441  bool local_=false);
442 
444 
445  virtual void Mult(const Vector &x, Vector &y) const;
446 
447  virtual void MultTranspose(const Vector &x, Vector &y) const;
448 };
449 
450 /// Auxiliary device class used by ParFiniteElementSpace.
453 {
454 protected:
461  MPI_Request *requests;
462 
463  // Kernel: copy ltdofs from 'src' to 'shr_buf' - prepare for send.
464  // shr_buf[i] = src[shr_ltdof[i]]
465  void BcastBeginCopy(const Vector &src) const;
466 
467  // Kernel: copy ltdofs from 'src' to ldofs in 'dst'.
468  // dst[ltdof_ldof[i]] = src[i]
469  void BcastLocalCopy(const Vector &src, Vector &dst) const;
470 
471  // Kernel: copy ext. dofs from 'ext_buf' to 'dst' - after recv.
472  // dst[ext_ldof[i]] = ext_buf[i]
473  void BcastEndCopy(Vector &dst) const;
474 
475  // Kernel: copy ext. dofs from 'src' to 'ext_buf' - prepare for send.
476  // ext_buf[i] = src[ext_ldof[i]]
477  void ReduceBeginCopy(const Vector &src) const;
478 
479  // Kernel: copy owned ldofs from 'src' to ltdofs in 'dst'.
480  // dst[i] = src[ltdof_ldof[i]]
481  void ReduceLocalCopy(const Vector &src, Vector &dst) const;
482 
483  // Kernel: assemble dofs from 'shr_buf' into to 'dst' - after recv.
484  // dst[shr_ltdof[i]] += shr_buf[i]
485  void ReduceEndAssemble(Vector &dst) const;
486 
487 public:
489  const GroupCommunicator &gc_, const SparseMatrix *R, bool local_=false);
490 
492  bool local_=false);
493 
495 
496  virtual void Mult(const Vector &x, Vector &y) const;
497 
498  virtual void MultTranspose(const Vector &x, Vector &y) const;
499 };
500 
501 }
502 
503 #endif // MFEM_USE_MPI
504 
505 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:233
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:630
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1031
HYPRE_BigInt GetMyTDofOffset() const
Definition: pfespace.cpp:1147
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:1152
HypreParVector * NewTrueDofVector()
Definition: pfespace.hpp:327
const FiniteElement * GetFaceNbrFaceFE(int i) const
Definition: pfespace.cpp:1524
HYPRE_BigInt GetGlobalScalarTDofNumber(int sldof)
Definition: pfespace.cpp:1105
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1510
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:424
Auxiliary class used by ParFiniteElementSpace.
Definition: pfespace.hpp:429
virtual const Operator * GetRestrictionOperator() const
Definition: pfespace.cpp:1183
Ordering::Type ordering
Definition: fespace.hpp:239
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:3323
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
DeviceConformingProlongationOperator(const GroupCommunicator &gc_, const SparseMatrix *R, bool local_=false)
Definition: pfespace.cpp:3574
void BcastLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3688
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1482
void ReduceLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3770
ConformingProlongationOperator(int lsize, const GroupCommunicator &gc_, bool local_=false)
Definition: pfespace.cpp:3428
HYPRE_BigInt GetMyDofOffset() const
Definition: pfespace.cpp:1142
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:3543
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:974
Abstract parallel finite element space.
Definition: pfespace.hpp:28
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
Definition: pfespace.cpp:963
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:379
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:231
void DeleteAll()
Delete the whole array.
Definition: array.hpp:854
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
ParNCMesh * pncmesh
Definition: pmesh.hpp:388
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3508
DofTransformation * GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1457
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:1020
FaceType
Definition: mesh.hpp:45
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
Data type sparse matrix.
Definition: sparsemat.hpp:50
void ReduceBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3759
const GroupCommunicator & gc
Definition: pfespace.hpp:433
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
Definition: pfespace.hpp:389
void LoseData()
NULL-ifies the data.
Definition: array.hpp:135
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...
Definition: pfespace.cpp:1008
ParMesh * GetParMesh() const
Definition: pfespace.hpp:273
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3704
const GroupCommunicator & GetGroupCommunicator() const
Definition: pfespace.cpp:3458
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:187
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:3296
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:606
void BcastBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3663
ElementTransformation * GetFaceNbrElementTransformation(int i) const
Definition: pfespace.hpp:390
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
Auxiliary device class used by ParFiniteElementSpace.
Definition: pfespace.hpp:451
virtual DofTransformation * GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: pfespace.cpp:470
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:161
virtual int GetFaceDofs(int i, Array< int > &dofs, int variant=0) const
Definition: pfespace.cpp:518
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:534
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
HYPRE_BigInt GlobalVSize() const
Definition: pfespace.hpp:279
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:334
HYPRE_Int HYPRE_BigInt
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:1061
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:29
Array< HYPRE_BigInt > face_nbr_glob_dof_map
Definition: pfespace.hpp:214
HYPRE_BigInt * GetDofOffsets() const
Definition: pfespace.hpp:277
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: pmesh.cpp:1976
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:285
void ReduceEndAssemble(Vector &dst) const
Definition: pfespace.cpp:3799
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:337
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
int dim
Definition: ex24.cpp:53
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const
Definition: pfespace.cpp:541
int index(int i, int j, int nx, int ny)
Definition: life.cpp:235
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
Definition: pfespace.cpp:583
MPI_Comm GetComm() const
Definition: pfespace.hpp:269
Vector data type.
Definition: vector.hpp:58
bool SharedNDTriangleDofs() const
Definition: pfespace.hpp:400
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:1267
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:639
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:317
NURBSExtension * NURBSext
Definition: fespace.hpp:270
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:3806
Base class for operators that extracts Face degrees of freedom.
GroupTopology gtopo
Definition: pmesh.hpp:375
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
virtual void UpdatesFinished()
Free ParGridFunction transformation matrix (if any), to save memory.
Definition: pfespace.hpp:413
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
Definition: pfespace.cpp:994
HYPRE_BigInt GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
Definition: pfespace.cpp:1082
Class for parallel meshes.
Definition: pmesh.hpp:32
HYPRE_BigInt * GetTrueDofOffsets() const
Definition: pfespace.hpp:278
GroupTopology gtopo
Definition: nurbs.hpp:521
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:635
virtual DofTransformation * GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: pfespace.cpp:494
int VDofToDof(int vdof) const
Compute the inverse of the Dof to VDof mapping for a single index vdof.
Definition: fespace.hpp:974