MFEM  v4.4.0
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-2022, 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  /** Transpose of R or Rconf. For conforming mesh, this is a matrix-free
86  (Device)ConformingProlongationOperator, for a non-conforming mesh
87  this is a TransposeOperator wrapping R. */
88  mutable Operator *R_transpose;
89 
90  /// Flag indicating the existence of shared triangles with interior ND dofs
91  bool nd_strias;
92 
93  /// Resets nd_strias flag at constuction or after rebalancing
94  void CheckNDSTriaDofs();
95 
96  ParNURBSExtension *pNURBSext() const
97  { return dynamic_cast<ParNURBSExtension *>(NURBSext); }
98 
99  static ParNURBSExtension *MakeLocalNURBSext(
100  const NURBSExtension *globNURBSext, const NURBSExtension *parNURBSext);
101 
102  GroupTopology &GetGroupTopo() const
103  { return (NURBSext) ? pNURBSext()->gtopo : pmesh->gtopo; }
104 
105  // Auxiliary method used in constructors
106  void ParInit(ParMesh *pm);
107 
108  void Construct();
109  void Destroy();
110 
111  // ldof_type = 0 : DOFs communicator, otherwise VDOFs communicator
112  void GetGroupComm(GroupCommunicator &gcomm, int ldof_type,
113  Array<int> *ldof_sign = NULL);
114 
115  /// Construct dof_offsets and tdof_offsets using global communication.
116  void GenerateGlobalOffsets() const;
117 
118  /// Construct ldof_group and ldof_ltdof.
119  void ConstructTrueDofs();
120  void ConstructTrueNURBSDofs();
121 
122  void ApplyLDofSigns(Array<int> &dofs) const;
123  void ApplyLDofSigns(Table &el_dof) const;
124 
125  typedef NCMesh::MeshId MeshId;
126  typedef ParNCMesh::GroupId GroupId;
127 
128  void GetGhostVertexDofs(const MeshId &id, Array<int> &dofs) const;
129  void GetGhostEdgeDofs(const MeshId &edge_id, Array<int> &dofs) const;
130  void GetGhostFaceDofs(const MeshId &face_id, Array<int> &dofs) const;
131 
132  void GetGhostDofs(int entity, const MeshId &id, Array<int> &dofs) const;
133  /// Return the dofs associated with the interior of the given mesh entity.
134  void GetBareDofs(int entity, int index, Array<int> &dofs) const;
135 
136  int PackDof(int entity, int index, int edof) const;
137  void UnpackDof(int dof, int &entity, int &index, int &edof) const;
138 
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;
142 #endif
143 
144  void ScheduleSendRow(const struct PMatrixRow &row, int dof, GroupId group_id,
145  std::map<int, class NeighborRowMessage> &send_msg) const;
146 
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;
150 
151 #ifdef MFEM_DEBUG_PMATRIX
152  void DebugDumpDOFs(std::ostream &os,
153  const SparseMatrix &deps,
154  const Array<GroupId> &dof_group,
155  const Array<GroupId> &dof_owner,
156  const Array<bool> &finalized) const;
157 #endif
158 
159  /// Helper: create a HypreParMatrix from a list of PMatrixRows.
161  MakeVDimHypreMatrix(const std::vector<struct PMatrixRow> &rows,
162  int local_rows, int local_cols,
163  Array<HYPRE_BigInt> &row_starts,
164  Array<HYPRE_BigInt> &col_starts) const;
165 
166  /// Build the P and R matrices.
167  void Build_Dof_TrueDof_Matrix() const;
168 
169  /** Used when the ParMesh is non-conforming, i.e. pmesh->pncmesh != NULL.
170  Constructs the matrices P and R, the DOF and true DOF offset arrays,
171  and the DOF -> true DOF map ('dof_tdof'). Returns the number of
172  vector true DOFs. All pointer arguments are optional and can be NULL. */
173  int BuildParallelConformingInterpolation(HypreParMatrix **P, SparseMatrix **R,
174  Array<HYPRE_BigInt> &dof_offs,
175  Array<HYPRE_BigInt> &tdof_offs,
176  Array<int> *dof_tdof,
177  bool partial = false) const;
178 
179  /** Calculate a GridFunction migration matrix after mesh load balancing.
180  The result is a parallel permutation matrix that can be used to update
181  all grid functions defined on this space. */
182  HypreParMatrix* RebalanceMatrix(int old_ndofs,
183  const Table* old_elem_dof,
184  const Table* old_elem_fos);
185 
186  /** Calculate a GridFunction restriction matrix after mesh derefinement.
187  The matrix is constructed so that the new grid function interpolates
188  the original function, i.e., the original function is evaluated at the
189  nodes of the coarse function. */
190  HypreParMatrix* ParallelDerefinementMatrix(int old_ndofs,
191  const Table *old_elem_dof,
192  const Table *old_elem_fos);
193 
194  /// Updates the internal mesh pointer. @warning @a new_mesh must be
195  /// <b>topologically identical</b> to the existing mesh. Used if the address
196  /// of the Mesh object has changed, e.g. in @a Mesh::Swap.
197  virtual void UpdateMeshPointer(Mesh *new_mesh);
198 
199  /// Copies the prolongation and restriction matrices from @a fes.
200  ///
201  /// Used for low order preconditioning on non-conforming meshes. If the DOFs
202  /// require a permutation, it will be supplied by non-NULL @a perm. NULL @a
203  /// perm indicates that no permutation is required.
204  virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes,
205  const Array<int> *perm);
206 
207 public:
208  // Face-neighbor data
209  // Number of face-neighbor dofs
211  // Face-neighbor-element to face-neighbor dof
213  // Face-neighbor-element face orientations
215  // Face-neighbor to ldof in the face-neighbor numbering
217  // The global ldof indices of the face-neighbor dofs
219  // Local face-neighbor data: face-neighbor to ldof
221 
222  /** @brief Copy constructor: deep copy all data from @a orig except the
223  ParMesh, the FiniteElementCollection, and some derived data. */
224  /** If the @a pmesh or @a fec pointers are NULL (default), then the new
225  ParFiniteElementSpace will reuse the respective pointers from @a orig. If
226  any of these pointers is not NULL, the given pointer will be used instead
227  of the one used by @a orig.
228 
229  @note The objects pointed to by the @a pmesh and @a fec parameters must
230  be either the same objects as the ones used by @a orig, or copies of
231  them. Otherwise, the behavior is undefined.
232 
233  @note Derived data objects, such as the parallel prolongation and
234  restriction operators, the update operator, and any of the face-neighbor
235  data, will not be copied, even if they are created in the @a orig object.
236  */
238  ParMesh *pmesh = NULL,
239  const FiniteElementCollection *fec = NULL);
240 
241  /** @brief Convert/copy the *local* (Par)FiniteElementSpace @a orig to
242  ParFiniteElementSpace: deep copy all data from @a orig except the Mesh,
243  the FiniteElementCollection, and some derived data. */
245  const FiniteElementCollection *fec = NULL);
246 
247  /** @brief Construct the *local* ParFiniteElementSpace corresponding to the
248  global FE space, @a global_fes. */
249  /** The parameter @a pm is the *local* ParMesh obtained by decomposing the
250  global Mesh used by @a global_fes. The array @a partitioning represents
251  the parallel decomposition - it maps global element ids to MPI ranks.
252  If the FiniteElementCollection, @a f, is NULL (default), the FE
253  collection used by @a global_fes will be reused. If @a f is not NULL, it
254  must be the same as, or a copy of, the FE collection used by
255  @a global_fes. */
256  ParFiniteElementSpace(ParMesh *pm, const FiniteElementSpace *global_fes,
257  const int *partitioning,
258  const FiniteElementCollection *f = NULL);
259 
261  int dim = 1, int ordering = Ordering::byNODES);
262 
263  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
264  /** The parameter @a ext will be deleted by this constructor, replaced by a
265  ParNURBSExtension owned by the ParFiniteElementSpace.
266  @note If the pointer @a ext is NULL, this constructor is equivalent to
267  the standard constructor with the same arguments minus the
268  NURBSExtension, @a ext. */
270  const FiniteElementCollection *f,
271  int dim = 1, int ordering = Ordering::byNODES);
272 
273  MPI_Comm GetComm() const { return MyComm; }
274  int GetNRanks() const { return NRanks; }
275  int GetMyRank() const { return MyRank; }
276 
277  inline ParMesh *GetParMesh() const { return pmesh; }
278 
279  int GetDofSign(int i)
280  { return NURBSext || Nonconforming() ? 1 : ldof_sign[VDofToDof(i)]; }
281  HYPRE_BigInt *GetDofOffsets() const { return dof_offsets; }
282  HYPRE_BigInt *GetTrueDofOffsets() const { return tdof_offsets; }
284  { return Dof_TrueDof_Matrix()->GetGlobalNumRows(); }
286  { return Dof_TrueDof_Matrix()->GetGlobalNumCols(); }
287 
288  /// Return the number of local vector true dofs.
289  virtual int GetTrueVSize() const { return ltdof_size; }
290 
291  /// Returns indexes of degrees of freedom in array dofs for i'th element.
292  virtual DofTransformation *GetElementDofs(int i, Array<int> &dofs) const;
293 
294  /// Returns indexes of degrees of freedom for i'th boundary element.
295  virtual DofTransformation *GetBdrElementDofs(int i, Array<int> &dofs) const;
296 
297  /** Returns the indexes of the degrees of freedom for i'th face
298  including the dofs for the edges and the vertices of the face. */
299  virtual int GetFaceDofs(int i, Array<int> &dofs, int variant = 0) const;
300 
301  /** Returns pointer to the FiniteElement in the FiniteElementCollection
302  associated with i'th element in the mesh object. If @a i is greater than
303  or equal to the number of local mesh elements, @a i will be interpreted
304  as a shifted index of a face neighbor element. */
305  virtual const FiniteElement *GetFE(int i) const;
306 
307  /** Returns an Operator that converts L-vectors to E-vectors on each face.
308  The parallel version is different from the serial one because of the
309  presence of shared faces. Shared faces are treated as interior faces,
310  the returned operator handles the communication needed to get the
311  shared face values from other MPI ranks */
312  virtual const FaceRestriction *GetFaceRestriction(
313  ElementDofOrdering e_ordering, FaceType type,
315 
316  void GetSharedEdgeDofs(int group, int ei, Array<int> &dofs) const;
317  void GetSharedTriangleDofs(int group, int fi, Array<int> &dofs) const;
318  void GetSharedQuadrilateralDofs(int group, int fi, Array<int> &dofs) const;
319 
320  /// The true dof-to-dof interpolation matrix
322  { if (!P) { Build_Dof_TrueDof_Matrix(); } return P; }
323 
324  /** @brief For a non-conforming mesh, construct and return the interpolation
325  matrix from the partially conforming true dofs to the local dofs. */
326  /** @note The returned pointer must be deleted by the caller. */
328 
329  /** Create and return a new HypreParVector on the true dofs, which is
330  owned by (i.e. it must be destroyed by) the calling function. */
332  { return (new HypreParVector(MyComm,GlobalTrueVSize(),GetTrueDofOffsets()));}
333 
334  /// Scale a vector of true dofs
335  void DivideByGroupSize(double *vec);
336 
337  /// Return a reference to the internal GroupCommunicator (on VDofs)
338  GroupCommunicator &GroupComm() { return *gcomm; }
339 
340  /// Return a const reference to the internal GroupCommunicator (on VDofs)
341  const GroupCommunicator &GroupComm() const { return *gcomm; }
342 
343  /// Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
344  /** @note The returned pointer must be deleted by the caller. */
346 
347  /** @brief Given an integer array on the local degrees of freedom, perform
348  a bitwise OR between the shared dofs. */
349  /** For non-conforming mesh, synchronization is performed on the cut (aka
350  "partially conforming") space. */
351  void Synchronize(Array<int> &ldof_marker) const;
352 
353  /// Determine the boundary degrees of freedom
354  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
355  Array<int> &ess_dofs,
356  int component = -1) const;
357 
358  /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
359  boundary attributes marked in the array bdr_attr_is_ess. */
360  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
361  Array<int> &ess_tdof_list,
362  int component = -1);
363 
364  /** If the given ldof is owned by the current processor, return its local
365  tdof number, otherwise return -1 */
366  int GetLocalTDofNumber(int ldof) const;
367  /// Returns the global tdof number of the given local degree of freedom
368  HYPRE_BigInt GetGlobalTDofNumber(int ldof) const;
369  /** Returns the global tdof number of the given local degree of freedom in
370  the scalar version of the current finite element space. The input should
371  be a scalar local dof. */
373 
376 
377  virtual const Operator *GetProlongationMatrix() const;
378  /** @brief Return logical transpose of restriction matrix, but in
379  non-assembled optimized matrix-free form.
380 
381  The implementation is like GetProlongationMatrix, but it sets local
382  DOFs to the true DOF values if owned locally, otherwise zero. */
383  virtual const Operator *GetRestrictionTransposeOperator() const;
384  /** Get an Operator that performs the action of GetRestrictionMatrix(),
385  but potentially with a non-assembled optimized matrix-free
386  implementation. */
387  virtual const Operator *GetRestrictionOperator() const;
388  /// Get the R matrix which restricts a local dof vector to true dof vector.
389  virtual const SparseMatrix *GetRestrictionMatrix() const
390  { Dof_TrueDof_Matrix(); return R; }
391 
392  // Face-neighbor functions
393  void ExchangeFaceNbrData();
394  int GetFaceNbrVSize() const { return num_face_nbr_dofs; }
396  void GetFaceNbrFaceVDofs(int i, Array<int> &vdofs) const;
397  const FiniteElement *GetFaceNbrFE(int i) const;
398  const FiniteElement *GetFaceNbrFaceFE(int i) const;
401  { return pmesh->GetFaceNbrElementTransformation(i); }
402 
404  void LoseDofOffsets() { dof_offsets.LoseData(); }
405  void LoseTrueDofOffsets() { tdof_offsets.LoseData(); }
406 
407  bool Conforming() const { return pmesh->pncmesh == NULL && !nonconf_P; }
408  bool Nonconforming() const { return pmesh->pncmesh != NULL || nonconf_P; }
409 
410  bool SharedNDTriangleDofs() const { return nd_strias; }
411 
412  // Transfer parallel true-dof data from coarse_fes, defined on a coarse mesh,
413  // to this FE space, defined on a refined mesh. See full documentation in the
414  // base class, FiniteElementSpace::GetTrueTransferOperator.
415  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
416  OperatorHandle &T) const;
417 
418  /** Reflect changes in the mesh. Calculate one of the refinement/derefinement
419  /rebalance matrices, unless want_transform is false. */
420  virtual void Update(bool want_transform = true);
421 
422  /// Free ParGridFunction transformation matrix (if any), to save memory.
423  virtual void UpdatesFinished()
424  {
426  old_dof_offsets.DeleteAll();
427  }
428 
429  virtual ~ParFiniteElementSpace() { Destroy(); }
430 
431  void PrintPartitionStats();
432 
433  /// Obsolete, kept for backward compatibility
434  int TrueVSize() const { return ltdof_size; }
435 };
436 
437 
438 /// Auxiliary class used by ParFiniteElementSpace.
440 {
441 protected:
444  bool local;
445 
446 public:
448  bool local_=false);
449 
451  bool local_=false);
452 
454 
455  virtual void Mult(const Vector &x, Vector &y) const;
456 
457  virtual void MultTranspose(const Vector &x, Vector &y) const;
458 };
459 
460 /// Auxiliary device class used by ParFiniteElementSpace.
463 {
464 protected:
471  MPI_Request *requests;
472 
473  // Kernel: copy ltdofs from 'src' to 'shr_buf' - prepare for send.
474  // shr_buf[i] = src[shr_ltdof[i]]
475  void BcastBeginCopy(const Vector &src) const;
476 
477  // Kernel: copy ltdofs from 'src' to ldofs in 'dst'.
478  // dst[ltdof_ldof[i]] = src[i]
479  void BcastLocalCopy(const Vector &src, Vector &dst) const;
480 
481  // Kernel: copy ext. dofs from 'ext_buf' to 'dst' - after recv.
482  // dst[ext_ldof[i]] = ext_buf[i]
483  void BcastEndCopy(Vector &dst) const;
484 
485  // Kernel: copy ext. dofs from 'src' to 'ext_buf' - prepare for send.
486  // ext_buf[i] = src[ext_ldof[i]]
487  void ReduceBeginCopy(const Vector &src) const;
488 
489  // Kernel: copy owned ldofs from 'src' to ltdofs in 'dst'.
490  // dst[i] = src[ltdof_ldof[i]]
491  void ReduceLocalCopy(const Vector &src, Vector &dst) const;
492 
493  // Kernel: assemble dofs from 'shr_buf' into to 'dst' - after recv.
494  // dst[shr_ltdof[i]] += shr_buf[i]
495  void ReduceEndAssemble(Vector &dst) const;
496 
497 public:
499  const GroupCommunicator &gc_, const SparseMatrix *R, bool local_=false);
500 
502  bool local_=false);
503 
505 
506  virtual void Mult(const Vector &x, Vector &y) const;
507 
508  virtual void MultTranspose(const Vector &x, Vector &y) const;
509 };
510 
511 }
512 
513 #endif // MFEM_USE_MPI
514 
515 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:235
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
HYPRE_BigInt GetMyTDofOffset() const
Definition: pfespace.cpp:1148
HypreParVector * NewTrueDofVector()
Definition: pfespace.hpp:331
HYPRE_BigInt GetGlobalScalarTDofNumber(int sldof)
Definition: pfespace.cpp:1106
ParMesh * GetParMesh() const
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:1021
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:1153
void ReduceEndAssemble(Vector &dst) const
Definition: pfespace.cpp:3754
Auxiliary class used by ParFiniteElementSpace.
Definition: pfespace.hpp:439
Ordering::Type ordering
Definition: fespace.hpp:108
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:434
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:3284
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:3535
int VDofToDof(int vdof) const
Definition: fespace.hpp:687
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:341
virtual DofTransformation * GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: pfespace.cpp:495
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:3257
ConformingProlongationOperator(int lsize, const GroupCommunicator &gc_, bool local_=false)
Definition: pfespace.cpp:3389
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:975
void BcastEndCopy(Vector &dst) const
Definition: pfespace.cpp:3651
void BcastBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3621
Abstract parallel finite element space.
Definition: pfespace.hpp:28
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:1009
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:1062
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
Definition: pfespace.cpp:964
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:100
const GroupCommunicator & GetGroupCommunicator() const
Definition: pfespace.cpp:3419
void DeleteAll()
Delete the whole array.
Definition: array.hpp:846
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
Definition: pfespace.cpp:584
ParNCMesh * pncmesh
Definition: pmesh.hpp:325
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:289
HYPRE_BigInt GlobalVSize() const
Definition: pfespace.hpp:283
FaceType
Definition: mesh.hpp:45
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1517
MPI_Comm GetComm() const
Definition: pfespace.hpp:273
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:535
ElementTransformation * GetFaceNbrElementTransformation(int i) const
Definition: pfespace.hpp:400
Data type sparse matrix.
Definition: sparsemat.hpp:46
const GroupCommunicator & gc
Definition: pfespace.hpp:443
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:3761
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
Definition: pfespace.hpp:399
void LoseData()
NULL-ifies the data.
Definition: array.hpp:132
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:321
void ReduceLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3725
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:631
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:187
HYPRE_BigInt GetMyDofOffset() const
Definition: pfespace.cpp:1143
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:607
HYPRE_BigInt * GetDofOffsets() const
Definition: pfespace.hpp:281
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:471
Auxiliary device class used by ParFiniteElementSpace.
Definition: pfespace.hpp:461
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:633
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:148
HYPRE_BigInt * GetTrueDofOffsets() const
Definition: pfespace.hpp:282
HYPRE_BigInt GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
Definition: pfespace.cpp:1083
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
bool SharedNDTriangleDofs() const
Definition: pfespace.hpp:410
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:3504
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:338
HYPRE_Int HYPRE_BigInt
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
void BcastLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3643
bool Nonconforming() const
Definition: pfespace.hpp:408
Array< HYPRE_BigInt > face_nbr_glob_dof_map
Definition: pfespace.hpp:218
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: pmesh.cpp:1893
virtual const Operator * GetRestrictionOperator() const
Definition: pfespace.cpp:1184
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:629
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:389
virtual int GetFaceDofs(int i, Array< int > &dofs, int variant=0) const
Definition: pfespace.cpp:519
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:66
int dim
Definition: ex24.cpp:53
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const
Definition: pfespace.cpp:542
int index(int i, int j, int nx, int ny)
Definition: life.cpp:237
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3469
const FiniteElement * GetFaceNbrFaceFE(int i) const
Definition: pfespace.cpp:1531
int GetFaceNbrVSize() const
Definition: pfespace.hpp:394
virtual const Operator * GetRestrictionTransposeOperator() const
Return logical transpose of restriction matrix, but in non-assembled optimized matrix-free form...
Definition: pfespace.cpp:1217
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1489
Vector data type.
Definition: vector.hpp:60
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:894
void ReduceBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3714
NURBSExtension * NURBSext
Definition: fespace.hpp:139
Base class for operators that extracts Face degrees of freedom.
GroupTopology gtopo
Definition: pmesh.hpp:312
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
virtual void UpdatesFinished()
Free ParGridFunction transformation matrix (if any), to save memory.
Definition: pfespace.hpp:423
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
Definition: pfespace.cpp:995
Class for parallel meshes.
Definition: pmesh.hpp:32
GroupTopology gtopo
Definition: nurbs.hpp:423
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3659
DofTransformation * GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1464