MFEM  v4.2.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-2020, 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_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 Operator *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  void GetBareDofs(int entity, int index, Array<int> &dofs) const;
118 
119  int PackDof(int entity, int index, int edof) const;
120  void UnpackDof(int dof, int &entity, int &index, int &edof) const;
121 
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;
125 #endif
126 
127  void ScheduleSendRow(const struct PMatrixRow &row, int dof, GroupId group_id,
128  std::map<int, class NeighborRowMessage> &send_msg) const;
129 
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;
133 
134 #ifdef MFEM_DEBUG_PMATRIX
135  void DebugDumpDOFs(std::ostream &os,
136  const SparseMatrix &deps,
137  const Array<GroupId> &dof_group,
138  const Array<GroupId> &dof_owner,
139  const Array<bool> &finalized) const;
140 #endif
141 
142  /// Helper: create a HypreParMatrix from a list of PMatrixRows.
144  MakeVDimHypreMatrix(const std::vector<struct PMatrixRow> &rows,
145  int local_rows, int local_cols,
146  Array<HYPRE_Int> &row_starts,
147  Array<HYPRE_Int> &col_starts) const;
148 
149  /// Build the P and R matrices.
150  void Build_Dof_TrueDof_Matrix() const;
151 
152  /** Used when the ParMesh is non-conforming, i.e. pmesh->pncmesh != NULL.
153  Constructs the matrices P and R, the DOF and true DOF offset arrays,
154  and the DOF -> true DOF map ('dof_tdof'). Returns the number of
155  vector true DOFs. All pointer arguments are optional and can be NULL. */
156  int BuildParallelConformingInterpolation(HypreParMatrix **P, SparseMatrix **R,
157  Array<HYPRE_Int> &dof_offs,
158  Array<HYPRE_Int> &tdof_offs,
159  Array<int> *dof_tdof,
160  bool partial = false) const;
161 
162  /** Calculate a GridFunction migration matrix after mesh load balancing.
163  The result is a parallel permutation matrix that can be used to update
164  all grid functions defined on this space. */
165  HypreParMatrix* RebalanceMatrix(int old_ndofs,
166  const Table* old_elem_dof);
167 
168  /** Calculate a GridFunction restriction matrix after mesh derefinement.
169  The matrix is constructed so that the new grid function interpolates
170  the original function, i.e., the original function is evaluated at the
171  nodes of the coarse function. */
172  HypreParMatrix* ParallelDerefinementMatrix(int old_ndofs,
173  const Table *old_elem_dof);
174 
175 public:
176  // Face-neighbor data
177  // Number of face-neighbor dofs
179  // Face-neighbor-element to face-neighbor dof
181  // Face-neighbor to ldof in the face-neighbor numbering
183  // The global ldof indices of the face-neighbor dofs
185  // Local face-neighbor data: face-neighbor to ldof
187 
188  /** @brief Copy constructor: deep copy all data from @a orig except the
189  ParMesh, the FiniteElementCollection, and some derived data. */
190  /** If the @a pmesh or @a fec pointers are NULL (default), then the new
191  ParFiniteElementSpace will reuse the respective pointers from @a orig. If
192  any of these pointers is not NULL, the given pointer will be used instead
193  of the one used by @a orig.
194 
195  @note The objects pointed to by the @a pmesh and @a fec parameters must
196  be either the same objects as the ones used by @a orig, or copies of
197  them. Otherwise, the behavior is undefined.
198 
199  @note Derived data objects, such as the parallel prolongation and
200  restriction operators, the update operator, and any of the face-neighbor
201  data, will not be copied, even if they are created in the @a orig object.
202  */
204  ParMesh *pmesh = NULL,
205  const FiniteElementCollection *fec = NULL);
206 
207  /** @brief Convert/copy the *local* (Par)FiniteElementSpace @a orig to
208  ParFiniteElementSpace: deep copy all data from @a orig except the Mesh,
209  the FiniteElementCollection, and some derived data. */
211  const FiniteElementCollection *fec = NULL);
212 
213  /** @brief Construct the *local* ParFiniteElementSpace corresponding to the
214  global FE space, @a global_fes. */
215  /** The parameter @a pm is the *local* ParMesh obtained by decomposing the
216  global Mesh used by @a global_fes. The array @a partitioning represents
217  the parallel decomposition - it maps global element ids to MPI ranks.
218  If the FiniteElementCollection, @a f, is NULL (default), the FE
219  collection used by @a global_fes will be reused. If @a f is not NULL, it
220  must be the same as, or a copy of, the FE collection used by
221  @a global_fes. */
222  ParFiniteElementSpace(ParMesh *pm, const FiniteElementSpace *global_fes,
223  const int *partitioning,
224  const FiniteElementCollection *f = NULL);
225 
227  int dim = 1, int ordering = Ordering::byNODES);
228 
229  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
230  /** The parameter @a ext will be deleted by this constructor, replaced by a
231  ParNURBSExtension owned by the ParFiniteElementSpace.
232  @note If the pointer @a ext is NULL, this constructor is equivalent to
233  the standard constructor with the same arguments minus the
234  NURBSExtension, @a ext. */
236  const FiniteElementCollection *f,
237  int dim = 1, int ordering = Ordering::byNODES);
238 
239  MPI_Comm GetComm() const { return MyComm; }
240  int GetNRanks() const { return NRanks; }
241  int GetMyRank() const { return MyRank; }
242 
243  inline ParMesh *GetParMesh() const { return pmesh; }
244 
245  int GetDofSign(int i)
246  { return NURBSext || Nonconforming() ? 1 : ldof_sign[VDofToDof(i)]; }
247  HYPRE_Int *GetDofOffsets() const { return dof_offsets; }
248  HYPRE_Int *GetTrueDofOffsets() const { return tdof_offsets; }
249  HYPRE_Int GlobalVSize() const
250  { return Dof_TrueDof_Matrix()->GetGlobalNumRows(); }
251  HYPRE_Int GlobalTrueVSize() const
252  { return Dof_TrueDof_Matrix()->GetGlobalNumCols(); }
253 
254  /// Return the number of local vector true dofs.
255  virtual int GetTrueVSize() const { return ltdof_size; }
256 
257  /// Returns indexes of degrees of freedom in array dofs for i'th element.
258  virtual void GetElementDofs(int i, Array<int> &dofs) const;
259 
260  /// Returns indexes of degrees of freedom for i'th boundary element.
261  virtual void GetBdrElementDofs(int i, Array<int> &dofs) const;
262 
263  /** Returns the indexes of the degrees of freedom for i'th face
264  including the dofs for the edges and the vertices of the face. */
265  virtual void GetFaceDofs(int i, Array<int> &dofs) const;
266 
267  /** Returns an Operator that converts L-vectors to E-vectors on each face.
268  The parallel version is different from the serial one because of the
269  presence of shared faces. Shared faces are treated as interior faces,
270  the returned operator handles the communication needed to get the
271  shared face values from other MPI ranks */
272  virtual const Operator *GetFaceRestriction(
273  ElementDofOrdering e_ordering, FaceType type,
275 
276  void GetSharedEdgeDofs(int group, int ei, Array<int> &dofs) const;
277  void GetSharedTriangleDofs(int group, int fi, Array<int> &dofs) const;
278  void GetSharedQuadrilateralDofs(int group, int fi, Array<int> &dofs) const;
279 
280  /// The true dof-to-dof interpolation matrix
282  { if (!P) { Build_Dof_TrueDof_Matrix(); } return P; }
283 
284  /** @brief For a non-conforming mesh, construct and return the interpolation
285  matrix from the partially conforming true dofs to the local dofs. */
286  /** @note The returned pointer must be deleted by the caller. */
288 
289  /** Create and return a new HypreParVector on the true dofs, which is
290  owned by (i.e. it must be destroyed by) the calling function. */
292  { return (new HypreParVector(MyComm,GlobalTrueVSize(),GetTrueDofOffsets()));}
293 
294  /// Scale a vector of true dofs
295  void DivideByGroupSize(double *vec);
296 
297  /// Return a reference to the internal GroupCommunicator (on VDofs)
298  GroupCommunicator &GroupComm() { return *gcomm; }
299 
300  /// Return a const reference to the internal GroupCommunicator (on VDofs)
301  const GroupCommunicator &GroupComm() const { return *gcomm; }
302 
303  /// Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
304  /** @note The returned pointer must be deleted by the caller. */
306 
307  /** @brief Given an integer array on the local degrees of freedom, perform
308  a bitwise OR between the shared dofs. */
309  /** For non-conforming mesh, synchronization is performed on the cut (aka
310  "partially conforming") space. */
311  void Synchronize(Array<int> &ldof_marker) const;
312 
313  /// Determine the boundary degrees of freedom
314  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
315  Array<int> &ess_dofs,
316  int component = -1) const;
317 
318  /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
319  boundary attributes marked in the array bdr_attr_is_ess. */
320  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
321  Array<int> &ess_tdof_list,
322  int component = -1);
323 
324  /** If the given ldof is owned by the current processor, return its local
325  tdof number, otherwise return -1 */
326  int GetLocalTDofNumber(int ldof) const;
327  /// Returns the global tdof number of the given local degree of freedom
328  HYPRE_Int GetGlobalTDofNumber(int ldof) const;
329  /** Returns the global tdof number of the given local degree of freedom in
330  the scalar version of the current finite element space. The input should
331  be a scalar local dof. */
332  HYPRE_Int GetGlobalScalarTDofNumber(int sldof);
333 
334  HYPRE_Int GetMyDofOffset() const;
335  HYPRE_Int GetMyTDofOffset() const;
336 
337  virtual const Operator *GetProlongationMatrix() const;
338  /// Get the R matrix which restricts a local dof vector to true dof vector.
339  virtual const SparseMatrix *GetRestrictionMatrix() const
340  { Dof_TrueDof_Matrix(); return R; }
341 
342  // Face-neighbor functions
343  void ExchangeFaceNbrData();
344  int GetFaceNbrVSize() const { return num_face_nbr_dofs; }
345  void GetFaceNbrElementVDofs(int i, Array<int> &vdofs) const;
346  void GetFaceNbrFaceVDofs(int i, Array<int> &vdofs) const;
347  const FiniteElement *GetFaceNbrFE(int i) const;
348  const FiniteElement *GetFaceNbrFaceFE(int i) const;
349  const HYPRE_Int *GetFaceNbrGlobalDofMap() { return face_nbr_glob_dof_map; }
351  { return pmesh->GetFaceNbrElementTransformation(i); }
352 
354  void LoseDofOffsets() { dof_offsets.LoseData(); }
355  void LoseTrueDofOffsets() { tdof_offsets.LoseData(); }
356 
357  bool Conforming() const { return pmesh->pncmesh == NULL; }
358  bool Nonconforming() const { return pmesh->pncmesh != NULL; }
359 
360  // Transfer parallel true-dof data from coarse_fes, defined on a coarse mesh,
361  // to this FE space, defined on a refined mesh. See full documentation in the
362  // base class, FiniteElementSpace::GetTrueTransferOperator.
363  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
364  OperatorHandle &T) const;
365 
366  /** Reflect changes in the mesh. Calculate one of the refinement/derefinement
367  /rebalance matrices, unless want_transform is false. */
368  virtual void Update(bool want_transform = true);
369 
370  /// Free ParGridFunction transformation matrix (if any), to save memory.
371  virtual void UpdatesFinished()
372  {
374  old_dof_offsets.DeleteAll();
375  }
376 
377  virtual ~ParFiniteElementSpace() { Destroy(); }
378 
379  void PrintPartitionStats();
380 
381  /// Obsolete, kept for backward compatibility
382  int TrueVSize() const { return ltdof_size; }
383 };
384 
385 
386 /// Auxiliary class used by ParFiniteElementSpace.
388 {
389 protected:
392 
393 public:
395 
396  virtual void Mult(const Vector &x, Vector &y) const;
397 
398  virtual void MultTranspose(const Vector &x, Vector &y) const;
399 };
400 
401 /// Auxiliary device class used by ParFiniteElementSpace.
404 {
405 protected:
412  MPI_Request *requests;
413  // Kernel: copy ltdofs from 'src' to 'shr_buf' - prepare for send.
414  // shr_buf[i] = src[shr_ltdof[i]]
415  void BcastBeginCopy(const Vector &src) const;
416 
417  // Kernel: copy ltdofs from 'src' to ldofs in 'dst'.
418  // dst[ltdof_ldof[i]] = src[i]
419  void BcastLocalCopy(const Vector &src, Vector &dst) const;
420 
421  // Kernel: copy ext. dofs from 'ext_buf' to 'dst' - after recv.
422  // dst[ext_ldof[i]] = ext_buf[i]
423  void BcastEndCopy(Vector &dst) const;
424 
425  // Kernel: copy ext. dofs from 'src' to 'ext_buf' - prepare for send.
426  // ext_buf[i] = src[ext_ldof[i]]
427  void ReduceBeginCopy(const Vector &src) const;
428 
429  // Kernel: copy owned ldofs from 'src' to ltdofs in 'dst'.
430  // dst[i] = src[ltdof_ldof[i]]
431  void ReduceLocalCopy(const Vector &src, Vector &dst) const;
432 
433  // Kernel: assemble dofs from 'shr_buf' into to 'dst' - after recv.
434  // dst[shr_ltdof[i]] += shr_buf[i]
435  void ReduceEndAssemble(Vector &dst) const;
436 
437 public:
439 
441 
442  virtual void Mult(const Vector &x, Vector &y) const;
443 
444  virtual void MultTranspose(const Vector &x, Vector &y) const;
445 };
446 
447 }
448 
449 #endif // MFEM_USE_MPI
450 
451 #endif
Abstract class for all finite elements.
Definition: fe.hpp:235
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:775
HYPRE_Int GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
Definition: pfespace.cpp:825
L2FaceValues
Definition: restriction.hpp:26
HypreParVector * NewTrueDofVector()
Definition: pfespace.hpp:291
ParMesh * GetParMesh() const
Definition: pfespace.hpp:243
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:761
ConformingProlongationOperator(const ParFiniteElementSpace &pfes)
Definition: pfespace.cpp:2963
HYPRE_Int * GetDofOffsets() const
Definition: pfespace.hpp:247
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:895
void ReduceEndAssemble(Vector &dst) const
Definition: pfespace.cpp:3252
Auxiliary class used by ParFiniteElementSpace.
Definition: pfespace.hpp:387
Ordering::Type ordering
Definition: fespace.hpp:105
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:382
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:2875
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int VDofToDof(int vdof) const
Definition: fespace.hpp:507
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:301
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:2848
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:715
void BcastEndCopy(Vector &dst) const
Definition: pfespace.cpp:3164
void BcastBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3135
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:749
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:804
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
Definition: pfespace.cpp:704
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:97
const HYPRE_Int * GetFaceNbrGlobalDofMap()
Definition: pfespace.hpp:349
HYPRE_Int GetMyDofOffset() const
Definition: pfespace.cpp:885
void DeleteAll()
Delete the whole array.
Definition: array.hpp:821
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
Definition: pfespace.cpp:531
ParNCMesh * pncmesh
Definition: pmesh.hpp:260
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:255
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: pfespace.cpp:474
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1165
DeviceConformingProlongationOperator(const ParFiniteElementSpace &pfes)
Definition: pfespace.cpp:3056
FaceType
Definition: mesh.hpp:45
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:415
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1199
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
HYPRE_Int GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:419
ElementTransformation * GetFaceNbrElementTransformation(int i) const
Definition: pfespace.hpp:350
Data type sparse matrix.
Definition: sparsemat.hpp:46
const GroupCommunicator & gc
Definition: pfespace.hpp:391
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:3260
void LoseData()
NULL-ifies the data.
Definition: array.hpp:118
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:281
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
HYPRE_Int GetMyTDofOffset() const
Definition: pfespace.cpp:890
void ReduceLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3222
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:578
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:153
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:554
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: pfespace.cpp:488
Auxiliary device class used by ParFiniteElementSpace.
Definition: pfespace.hpp:402
HYPRE_Int GlobalVSize() const
Definition: pfespace.hpp:249
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:70
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
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:3031
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:298
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:30
HYPRE_Int GetGlobalScalarTDofNumber(int sldof)
Definition: pfespace.cpp:848
virtual const Operator * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const
Definition: pfespace.cpp:503
void BcastLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3156
bool Nonconforming() const
Definition: pfespace.hpp:358
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: pmesh.cpp:1685
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:460
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:339
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
int dim
Definition: ex24.cpp:53
int index(int i, int j, int nx, int ny)
Definition: life.cpp:241
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3006
const FiniteElement * GetFaceNbrFaceFE(int i) const
Definition: pfespace.cpp:1213
int GetFaceNbrVSize() const
Definition: pfespace.hpp:344
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:248
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1171
Vector data type.
Definition: vector.hpp:51
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:704
void ReduceBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3211
NURBSExtension * NURBSext
Definition: fespace.hpp:120
GroupTopology gtopo
Definition: pmesh.hpp:247
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
virtual void UpdatesFinished()
Free ParGridFunction transformation matrix (if any), to save memory.
Definition: pfespace.hpp:371
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
Definition: pfespace.cpp:735
Class for parallel meshes.
Definition: pmesh.hpp:32
GroupTopology gtopo
Definition: nurbs.hpp:418
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3172
double f(const Vector &p)
Array< HYPRE_Int > face_nbr_glob_dof_map
Definition: pfespace.hpp:184