MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
pfespace.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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
24namespace mfem
25{
26
27/// Abstract parallel finite element space.
29{
30private:
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 void GetGhostDofs(int entity, const MeshId &id, Array<int> &dofs) const;
128
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 void UpdateMeshPointer(Mesh *new_mesh) override;
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 void CopyProlongationAndRestriction(const FiniteElementSpace &fes,
201 const Array<int> *perm) override;
202
203public:
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
253 @note Currently the @a partitioning array is not used by this
254 constructor, it is required for general parallel variable-order support.
255 */
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. */
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; }
287
288 /// Return the number of local vector true dofs.
289 int GetTrueVSize() const override { return ltdof_size; }
290
291 /// Returns indexes of degrees of freedom in array dofs for i'th element and
292 /// returns the DofTransformation data in a user-provided object.
294 void GetElementDofs(int i, Array<int> &dofs,
295 DofTransformation &doftrans) const override;
296
297 /// Returns indexes of degrees of freedom for i'th boundary element and
298 /// returns the DofTransformation data in a user-provided object.
300 void GetBdrElementDofs(int i, Array<int> &dofs,
301 DofTransformation &doftrans) const override;
302
303 /** Returns the indexes of the degrees of freedom for i'th face
304 including the dofs for the edges and the vertices of the face. */
305 int GetFaceDofs(int i, Array<int> &dofs, int variant = 0) const override;
306
307 /** Returns pointer to the FiniteElement in the FiniteElementCollection
308 associated with i'th element in the mesh object. If @a i is greater than
309 or equal to the number of local mesh elements, @a i will be interpreted
310 as a shifted index of a face neighbor element. */
311 const FiniteElement *GetFE(int i) const override;
312
313 /** Returns an Operator that converts L-vectors to E-vectors on each face.
314 The parallel version is different from the serial one because of the
315 presence of shared faces. Shared faces are treated as interior faces,
316 the returned operator handles the communication needed to get the
317 shared face values from other MPI ranks */
319 ElementDofOrdering f_ordering, FaceType type,
320 L2FaceValues mul = L2FaceValues::DoubleValued) const override;
321
322 void GetSharedEdgeDofs(int group, int ei, Array<int> &dofs) const;
323 void GetSharedTriangleDofs(int group, int fi, Array<int> &dofs) const;
324 void GetSharedQuadrilateralDofs(int group, int fi, Array<int> &dofs) const;
325
326 /// The true dof-to-dof interpolation matrix
328 { if (!P) { Build_Dof_TrueDof_Matrix(); } return P; }
329
330 /** @brief For a non-conforming mesh, construct and return the interpolation
331 matrix from the partially conforming true dofs to the local dofs. */
332 /** @note The returned pointer must be deleted by the caller. */
334
335 /** Create and return a new HypreParVector on the true dofs, which is
336 owned by (i.e. it must be destroyed by) the calling function. */
339
340 /// Scale a vector of true dofs
341 void DivideByGroupSize(real_t *vec);
342
343 /// Return a reference to the internal GroupCommunicator (on VDofs)
344 GroupCommunicator &GroupComm() { return *gcomm; }
345
346 /// Return a const reference to the internal GroupCommunicator (on VDofs)
347 const GroupCommunicator &GroupComm() const { return *gcomm; }
348
349 /// Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
350 /** @note The returned pointer must be deleted by the caller. */
352
353 /** @brief Given an integer array on the local degrees of freedom, perform
354 a bitwise OR between the shared dofs. */
355 /** For non-conforming mesh, synchronization is performed on the cut (aka
356 "partially conforming") space. */
357 void Synchronize(Array<int> &ldof_marker) const;
358
359 /// Determine the boundary degrees of freedom
360 void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
361 Array<int> &ess_dofs,
362 int component = -1) const override;
363
364 /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
365 boundary attributes marked in the array bdr_attr_is_ess. */
366 void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
367 Array<int> &ess_tdof_list,
368 int component = -1) const override;
369
370 /** If the given ldof is owned by the current processor, return its local
371 tdof number, otherwise return -1 */
372 int GetLocalTDofNumber(int ldof) const;
373 /// Returns the global tdof number of the given local degree of freedom
374 HYPRE_BigInt GetGlobalTDofNumber(int ldof) const;
375 /** Returns the global tdof number of the given local degree of freedom in
376 the scalar version of the current finite element space. The input should
377 be a scalar local dof. */
379
382
383 const Operator *GetProlongationMatrix() const override;
384 /** Get an Operator that performs the action of GetRestrictionMatrix(),
385 but potentially with a non-assembled optimized matrix-free
386 implementation. */
387 const Operator *GetRestrictionOperator() const override;
388 /// Get the R matrix which restricts a local dof vector to true dof vector.
389 const SparseMatrix *GetRestrictionMatrix() const override
390 { Dof_TrueDof_Matrix(); return R; }
391
392 // Face-neighbor functions
393 void ExchangeFaceNbrData();
394 int GetFaceNbrVSize() const { return num_face_nbr_dofs; }
395 void GetFaceNbrElementVDofs(int i, Array<int> &vdofs,
396 DofTransformation &doftrans) const;
398 void GetFaceNbrFaceVDofs(int i, Array<int> &vdofs) const;
399 const FiniteElement *GetFaceNbrFE(int i) const;
400 const FiniteElement *GetFaceNbrFaceFE(int i) const;
404
406 void LoseDofOffsets() { dof_offsets.LoseData(); }
407 void LoseTrueDofOffsets() { tdof_offsets.LoseData(); }
408
409 bool Conforming() const { return pmesh->pncmesh == NULL && !nonconf_P; }
410 bool Nonconforming() const { return pmesh->pncmesh != NULL || nonconf_P; }
411
412 bool SharedNDTriangleDofs() const { return nd_strias; }
413
414 // Transfer parallel true-dof data from coarse_fes, defined on a coarse mesh,
415 // to this FE space, defined on a refined mesh. See full documentation in the
416 // base class, FiniteElementSpace::GetTrueTransferOperator.
417 void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
418 OperatorHandle &T) const override;
419
420 /** Reflect changes in the mesh. Calculate one of the refinement/derefinement
421 /rebalance matrices, unless want_transform is false. */
422 void Update(bool want_transform = true) override;
423
424 /// Free ParGridFunction transformation matrix (if any), to save memory.
425 void UpdatesFinished() override
426 {
428 old_dof_offsets.DeleteAll();
429 }
430
431 virtual ~ParFiniteElementSpace() { Destroy(); }
432
433 void PrintPartitionStats();
434
435 /// Obsolete, kept for backward compatibility
436 int TrueVSize() const { return ltdof_size; }
437};
438
439
440/// Auxiliary class used by ParFiniteElementSpace.
442{
443protected:
446 bool local;
447
448public:
450 bool local_=false);
451
453 bool local_=false);
454
456
457 void Mult(const Vector &x, Vector &y) const override;
458
459 void MultTranspose(const Vector &x, Vector &y) const override;
460};
461
462/// Auxiliary device class used by ParFiniteElementSpace.
465{
466protected:
473 MPI_Request *requests;
474
475 // Kernel: copy ltdofs from 'src' to 'shr_buf' - prepare for send.
476 // shr_buf[i] = src[shr_ltdof[i]]
477 void BcastBeginCopy(const Vector &src) const;
478
479 // Kernel: copy ltdofs from 'src' to ldofs in 'dst'.
480 // dst[ltdof_ldof[i]] = src[i]
481 void BcastLocalCopy(const Vector &src, Vector &dst) const;
482
483 // Kernel: copy ext. dofs from 'ext_buf' to 'dst' - after recv.
484 // dst[ext_ldof[i]] = ext_buf[i]
485 void BcastEndCopy(Vector &dst) const;
486
487 // Kernel: copy ext. dofs from 'src' to 'ext_buf' - prepare for send.
488 // ext_buf[i] = src[ext_ldof[i]]
489 void ReduceBeginCopy(const Vector &src) const;
490
491 // Kernel: copy owned ldofs from 'src' to ltdofs in 'dst'.
492 // dst[i] = src[ltdof_ldof[i]]
493 void ReduceLocalCopy(const Vector &src, Vector &dst) const;
494
495 // Kernel: assemble dofs from 'shr_buf' into to 'dst' - after recv.
496 // dst[shr_ltdof[i]] += shr_buf[i]
497 void ReduceEndAssemble(Vector &dst) const;
498
499public:
501 const GroupCommunicator &gc_, const SparseMatrix *R, bool local_=false);
502
504 bool local_=false);
505
507
508 void Mult(const Vector &x, Vector &y) const override;
509
510 void MultTranspose(const Vector &x, Vector &y) const override;
511};
512
513}
514
515#endif // MFEM_USE_MPI
516
517#endif
void LoseData()
NULL-ifies the data.
Definition array.hpp:138
void DeleteAll()
Delete the whole array.
Definition array.hpp:864
Auxiliary class used by ParFiniteElementSpace.
Definition pfespace.hpp:442
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
const GroupCommunicator & GetGroupCommunicator() const
ConformingProlongationOperator(int lsize, const GroupCommunicator &gc_, bool local_=false)
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
const GroupCommunicator & gc
Definition pfespace.hpp:445
Auxiliary device class used by ParFiniteElementSpace.
Definition pfespace.hpp:465
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void ReduceEndAssemble(Vector &dst) const
void BcastBeginCopy(const Vector &src) const
void ReduceLocalCopy(const Vector &src, Vector &dst) const
void ReduceBeginCopy(const Vector &src) const
DeviceConformingProlongationOperator(const GroupCommunicator &gc_, const SparseMatrix *R, bool local_=false)
void BcastLocalCopy(const Vector &src, Vector &dst) const
Base class for operators that extracts Face degrees of freedom.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
Definition fespace.cpp:2846
NURBSExtension * NURBSext
Definition fespace.hpp:270
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition fespace.hpp:231
int VDofToDof(int vdof) const
Compute the inverse of the Dof to VDof mapping for a single index vdof.
Definition fespace.hpp:995
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition fespace.hpp:1306
Ordering::Type ordering
Definition fespace.hpp:239
DofTransformation * GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element 'bel'. The returned indices are offsets in...
Definition fespace.cpp:2950
Abstract class for all finite elements.
Definition fe_base.hpp:239
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition hypre.hpp:679
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
Definition hypre.hpp:683
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:206
Class used by MFEM to store pointers to host and/or device memory.
Mesh data type.
Definition mesh.hpp:56
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Abstract operator.
Definition operator.hpp:25
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const override
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes,...
HYPRE_BigInt GetGlobalScalarTDofNumber(int sldof)
const GroupCommunicator & GroupComm() const
Return a const reference to the internal GroupCommunicator (on VDofs)
Definition pfespace.hpp:347
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
Definition pfespace.cpp:606
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
Definition pfespace.cpp:583
MPI_Comm GetComm() const
Definition pfespace.hpp:273
const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const override
Definition pfespace.cpp:541
HYPRE_BigInt * GetTrueDofOffsets() const
Definition pfespace.hpp:282
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalVSize() const
Definition pfespace.hpp:283
const Operator * GetRestrictionOperator() const override
int GetLocalTDofNumber(int ldof) const
void UpdatesFinished() override
Free ParGridFunction transformation matrix (if any), to save memory.
Definition pfespace.hpp:425
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.
bool SharedNDTriangleDofs() const
Definition pfespace.hpp:412
ParFiniteElementSpace(const ParFiniteElementSpace &orig, ParMesh *pmesh=NULL, const FiniteElementCollection *fec=NULL)
Copy constructor: deep copy all data from orig except the ParMesh, the FiniteElementCollection,...
Definition pfespace.cpp:29
void DivideByGroupSize(real_t *vec)
Scale a vector of true dofs.
Definition pfespace.cpp:974
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
HypreParVector * NewTrueDofVector()
Definition pfespace.hpp:337
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:289
HYPRE_BigInt GetMyDofOffset() const
HYPRE_BigInt * GetDofOffsets() const
Definition pfespace.hpp:281
Array< HYPRE_BigInt > face_nbr_glob_dof_map
Definition pfespace.hpp:214
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition pfespace.hpp:344
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
const Operator * GetProlongationMatrix() const override
The returned Operator is owned by the FiniteElementSpace.
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
Definition pfespace.cpp:994
void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const override
Determine the boundary degrees of freedom.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:327
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs, DofTransformation &doftrans) const
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
Definition pfespace.hpp:401
int GetFaceDofs(int i, Array< int > &dofs, int variant=0) const override
Definition pfespace.cpp:518
HYPRE_BigInt GetMyTDofOffset() const
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
Definition pfespace.cpp:630
void GetElementDofs(int i, Array< int > &dofs, DofTransformation &doftrans) const override
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be al...
Definition pfespace.cpp:468
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
Definition pfespace.hpp:389
HYPRE_BigInt GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
void Update(bool want_transform=true) override
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
Definition pfespace.cpp:963
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition pfespace.hpp:436
const FiniteElement * GetFaceNbrFaceFE(int i) const
const FiniteElement * GetFaceNbrFE(int i) const
const FiniteElement * GetFE(int i) const override
Definition pfespace.cpp:534
ElementTransformation * GetFaceNbrElementTransformation(int i) const
Definition pfespace.hpp:402
void GetBdrElementDofs(int i, Array< int > &dofs, DofTransformation &doftrans) const override
The same as GetBdrElementDofs(), but with a user-allocated DofTransformation object....
Definition pfespace.cpp:493
Class for parallel meshes.
Definition pmesh.hpp:34
ElementTransformation * GetFaceNbrElementTransformation(int FaceNo)
Returns a pointer to the transformation defining the i-th face neighbor.
Definition pmesh.cpp:3088
GroupTopology gtopo
Definition pmesh.hpp:426
ParNCMesh * pncmesh
Definition pmesh.hpp:439
A parallel extension of the NCMesh class.
Definition pncmesh.hpp:65
GroupTopology gtopo
Definition nurbs.hpp:621
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
int dim
Definition ex24.cpp:53
HYPRE_Int HYPRE_BigInt
int index(int i, int j, int nx, int ny)
Definition life.cpp:235
float real_t
Definition config.hpp:43
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:75
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
FaceType
Definition mesh.hpp:47
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition ncmesh.hpp:189