MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
pmesh.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_PMESH
13#define MFEM_PMESH
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_MPI
18
21#include "mesh.hpp"
22#include "pncmesh.hpp"
23#include <iostream>
24
25namespace mfem
26{
27
28#ifdef MFEM_USE_PUMI
29class ParPumiMesh;
30#endif
31
32/// Class for parallel meshes
33class ParMesh : public Mesh
34{
35 friend class ParNCMesh;
36 friend class ParSubMesh;
37#ifdef MFEM_USE_PUMI
38 friend class ParPumiMesh;
39#endif
40#ifdef MFEM_USE_ADIOS2
41 friend class adios2stream;
42#endif
43
44protected:
45 MPI_Comm MyComm;
47
48 struct Vert3
49 {
50 int v[3];
51 Vert3() = default;
52 Vert3(int v0, int v1, int v2) { v[0] = v0; v[1] = v1; v[2] = v2; }
53 void Set(int v0, int v1, int v2) { v[0] = v0; v[1] = v1; v[2] = v2; }
54 void Set(const int *w) { v[0] = w[0]; v[1] = w[1]; v[2] = w[2]; }
55 };
56
57 struct Vert4
58 {
59 int v[4];
60 Vert4() = default;
61 Vert4(int v0, int v1, int v2, int v3)
62 { v[0] = v0; v[1] = v1; v[2] = v2; v[3] = v3; }
63 void Set(int v0, int v1, int v2, int v3)
64 { v[0] = v0; v[1] = v1; v[2] = v2; v[3] = v3; }
65 void Set(const int *w)
66 { v[0] = w[0]; v[1] = w[1]; v[2] = w[2]; v[3] = w[3]; }
67 };
68
70 // shared face id 'i' is:
71 // * triangle id 'i', if i < shared_trias.Size()
72 // * quad id 'i-shared_trias.Size()', otherwise
75
76 /// Shared objects in each group.
79 Table group_stria; // contains shared triangle indices
80 Table group_squad; // contains shared quadrilateral indices
81
82 /// Shared to local index mapping.
85 // sface ids: all triangles first, then all quads
87
88 /// Table that maps from face neighbor element number, to the face numbers of
89 /// that element.
90 std::unique_ptr<Table> face_nbr_el_to_face;
91 /// orientations for each face (from nbr processor)
92 std::unique_ptr<Table> face_nbr_el_ori;
93
94 // glob_elem_offset + local element number defines a global element numbering
95 mutable long long glob_elem_offset;
97 void ComputeGlobalElementOffset() const;
98
99 // Enable Print() to add the parallel interface as boundary (typically used
100 // for visualization purposes)
101 bool print_shared = true;
102
103 /// Create from a nonconforming mesh.
104 ParMesh(const ParNCMesh &pncmesh);
105
106 // Convert the local 'meshgen' to a global one.
107 void ReduceMeshGen();
108
109 // Determine sedge_ledge and sface_lface.
110 void FinalizeParTopo();
111
112 // Mark all tets to ensure consistency across MPI tasks; also mark the
113 // shared and boundary triangle faces using the consistently marked tets.
114 void MarkTetMeshForRefinement(const DSTable &v_to_v) override;
115
116 /// Return a number(0-1) identifying how the given edge has been split
117 int GetEdgeSplittings(Element *edge, const DSTable &v_to_v, int *middle);
118 /// Append codes identifying how the given face has been split to @a codes
119 void GetFaceSplittings(const int *fv, const HashTable<Hashed2> &v_to_v,
120 Array<unsigned> &codes);
121
122 bool DecodeFaceSplittings(HashTable<Hashed2> &v_to_v, const int *v,
123 const Array<unsigned> &codes, int &pos);
124
125 // Given a completed FacesTable and SharedFacesTable, construct a table that
126 // maps from face neighbor element number, to the set of faces of that
127 // element. Store the resulting data in the member variable
128 // face_nbr_el_to_face. If the mesh is nonconforming, this also builds the
129 // the face_nbr_el_ori variable from the faces_info.
131
132 /**
133 * @brief Helper function for adding triangle face neighbor element to face
134 * table entries. Have to use a template here rather than lambda capture
135 * because the FaceVert entries in Geometry have inner size of 3 for tets and
136 * 4 for everything else.
137 *
138 * @tparam N Inner dimension on the fvert variable, 3 for tet, 4 otherwise
139 * @param[in] v Set of vertices for this element
140 * @param[in] faces Table of faces interior to this rank
141 * @param[in] shared_faces Table of faces shared by this rank and another
142 * @param[in] elem The face neighbor element
143 * @param[in] start Starting index into fverts
144 * @param[in] end End index into fverts
145 * @param[in] fverts Array of face vertices for this particular geometry.
146 */
147 template <int N>
148 void AddTriFaces(const Array<int> &v, const std::unique_ptr<STable3D> &faces,
149 const std::unique_ptr<STable3D> &shared_faces,
150 int elem, int start, int end, const int fverts[][N]);
151
154 Geometry::Type face_geom) const;
157 Geometry::Type face_geom) const
158 {
159 MFEM_ASSERT(FElTr, "Missing FaceElementTransformations object!");
160 GetGhostFaceTransformation(*FElTr, face_type, face_geom);
161 }
162
163 /// Update the groups after triangle refinement
164 void RefineGroups(const DSTable &v_to_v, int *middle);
165
166 /// Update the groups after tetrahedron refinement
167 void RefineGroups(int old_nv, const HashTable<Hashed2> &v_to_v);
168
169 void UniformRefineGroups2D(int old_nv);
170
171 // f2qf can be NULL if all faces are quads or there are no quad faces
172 void UniformRefineGroups3D(int old_nv, int old_nedges,
173 const DSTable &old_v_to_v,
174 const STable3D &old_faces,
175 Array<int> *f2qf);
176
177 void ExchangeFaceNbrData(Table *gr_sface, int *s2l_face);
178
179 /// Refine a mixed 2D mesh uniformly.
180 void UniformRefinement2D() override;
181
182 /// Refine a mixed 3D mesh uniformly.
183 void UniformRefinement3D() override;
184
185 /** @brief Refine NURBS mesh, with an optional refinement factor.
186
187 @param[in] rf Optional refinement factor. If scalar, the factor is used
188 for all dimensions. If an array, factors can be specified
189 for each dimension.
190 @param[in] tol NURBS geometry deviation tolerance. */
191 void NURBSUniformRefinement(int rf = 2, real_t tol=1.0e-12) override;
192 void NURBSUniformRefinement(const Array<int> &rf, real_t tol=1.e-12) override;
193
194 /// This function is not public anymore. Use GeneralRefinement instead.
195 void LocalRefinement(const Array<int> &marked_el, int type = 3) override;
196
197 /// This function is not public anymore. Use GeneralRefinement instead.
198 void NonconformingRefinement(const Array<Refinement> &refinements,
199 int nc_limit = 0) override;
200
202 real_t threshold, int nc_limit = 0,
203 int op = 1) override;
204
205 void RebalanceImpl(const Array<int> *partition);
206
207 void DeleteFaceNbrData();
208
209 bool WantSkipSharedMaster(const NCMesh::Master &master) const;
210
211 /// Fills out partitioned Mesh::vertices
212 int BuildLocalVertices(const Mesh& global_mesh, const int *partitioning,
213 Array<int> &vert_global_local);
214
215 /// Fills out partitioned Mesh::elements
216 int BuildLocalElements(const Mesh& global_mesh, const int *partitioning,
217 const Array<int> &vert_global_local);
218
219 /// Fills out partitioned Mesh::boundary
220 int BuildLocalBoundary(const Mesh& global_mesh, const int *partitioning,
221 const Array<int> &vert_global_local,
222 Array<bool>& activeBdrElem,
223 Table* &edge_element);
224
225 void FindSharedFaces(const Mesh &mesh, const int* partition,
226 Array<int>& face_group,
227 ListOfIntegerSets& groups);
228
229 int FindSharedEdges(const Mesh &mesh, const int* partition,
230 Table* &edge_element, ListOfIntegerSets& groups);
231
232 int FindSharedVertices(const int *partition, Table* vertex_element,
233 ListOfIntegerSets& groups);
234
235 void BuildFaceGroup(int ngroups, const Mesh &mesh,
236 const Array<int>& face_group,
237 int &nstria, int &nsquad);
238
239 void BuildEdgeGroup(int ngroups, const Table& edge_element);
240
241 void BuildVertexGroup(int ngroups, const Table& vert_element);
242
243 void BuildSharedFaceElems(int ntri_faces, int nquad_faces,
244 const Mesh &mesh, int *partitioning,
245 const STable3D *faces_tbl,
246 const Array<int> &face_group,
247 const Array<int> &vert_global_local);
248
249 void BuildSharedEdgeElems(int nedges, Mesh &mesh,
250 const Array<int> &vert_global_local,
251 const Table *edge_element);
252
253 void BuildSharedVertMapping(int nvert, const Table* vert_element,
254 const Array<int> &vert_global_local);
255
256 /**
257 * @brief Get the shared edges GroupCommunicator.
258 *
259 * The output of the shared edges is chosen by the @a ordering parameter with
260 * the following options
261 * 0: Internal ordering. Not exposed to public interfaces.
262 * 1: Contiguous ordering.
263 *
264 * @param[in] ordering Ordering for the shared edges.
265 * @param[out] sedge_comm
266 */
267 void GetSharedEdgeCommunicator(int ordering,
268 GroupCommunicator& sedge_comm) const;
269
270 /**
271 * @brief Get the shared vertices GroupCommunicator.
272 *
273 * The output of the shared vertices is chosen by the @a ordering parameter
274 * with the following options
275 * 0: Internal ordering. Not exposed to public interfaces.
276 * 1: Contiguous ordering.
277 *
278 * @param[in] ordering
279 * @param[out] svert_comm
280 */
281 void GetSharedVertexCommunicator(int ordering,
282 GroupCommunicator& svert_comm) const;
283
284 /**
285 * @brief Get the shared face quadrilaterals GroupCommunicator.
286 *
287 * The output of the shared face quadrilaterals is chosen by the @a ordering
288 * parameter with the following options
289 * 0: Internal ordering. Not exposed to public interfaces.
290 * 1: Contiguous ordering.
291 *
292 * @param[in] ordering
293 * @param[out] squad_comm
294 */
295 void GetSharedQuadCommunicator(int ordering,
296 GroupCommunicator& squad_comm) const;
297
298 /**
299 * @brief Get the shared face triangles GroupCommunicator.
300 *
301 * The output of the shared face triangles is chosen by the @a ordering
302 * parameter with the following options
303 * 0: Internal ordering. Not exposed to public interfaces.
304 * 1: Contiguous ordering.
305 *
306 * @param[in] ordering
307 * @param[out] stria_comm
308 */
309 void GetSharedTriCommunicator(int ordering,
310 GroupCommunicator& stria_comm) const;
311
312 // Similar to Mesh::GetFacesTable()
314
315 /// Ensure that bdr_attributes and attributes agree across processors
317
318 void LoadSharedEntities(std::istream &input);
319
320 /// If the mesh is curved, make sure 'Nodes' is ParGridFunction.
321 /** Note that this method is not related to the public 'Mesh::EnsureNodes`.*/
322 void EnsureParNodes();
323
324 /// Internal function used in ParMesh::MakeRefined (and related constructor)
325 void MakeRefined_(ParMesh &orig_mesh, int ref_factor, int ref_type);
326
327 // Mark Mesh::Swap as protected, should use ParMesh::Swap to swap @a ParMesh
328 // objects.
329 using Mesh::Swap;
330
331 void Destroy();
332
333public:
334 /// Default constructor. Create an empty @a ParMesh.
335 ParMesh() : MyComm(0), NRanks(0), MyRank(-1),
337 have_face_nbr_data(false), pncmesh(NULL) { }
338
339 /// Create a parallel mesh by partitioning a serial Mesh.
340 /** The mesh is partitioned automatically or using external partitioning
341 data (the optional parameter 'partitioning_[i]' contains the desired MPI
342 rank for element 'i'). Automatic partitioning uses METIS for conforming
343 meshes and quick space-filling curve equipartitioning for nonconforming
344 meshes (elements of nonconforming meshes should ideally be ordered as a
345 sequence of face-neighbors). */
346 ParMesh(MPI_Comm comm, Mesh &mesh, int *partitioning_ = NULL,
347 int part_method = 1);
348
349 /** Copy constructor. Performs a deep copy of (almost) all data, so that the
350 source mesh can be modified (e.g. deleted, refined) without affecting the
351 new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the
352 nodes, if present. */
353 explicit ParMesh(const ParMesh &pmesh, bool copy_nodes = true);
354
355 /// Read a parallel mesh, each MPI rank from its own file/stream.
356 /** The @a generate_edges parameter is passed to Mesh::Loader. The @a refine
357 and @a fix_orientation parameters are passed to the method
358 Mesh::Finalize().
359
360 @note The order of arguments and their default values are different than
361 for the Mesh class. */
362 ParMesh(MPI_Comm comm, std::istream &input, bool refine = true,
363 int generate_edges = 1, bool fix_orientation = true);
364
365 /// Deprecated: see @a ParMesh::MakeRefined
366 MFEM_DEPRECATED
367 ParMesh(ParMesh *orig_mesh, int ref_factor, int ref_type);
368
369 /// Move constructor. Used for named constructors.
370 ParMesh(ParMesh &&mesh);
371
372 /// Move assignment operator.
373 ParMesh& operator=(ParMesh &&mesh);
374
375 /// Explicitly delete the copy assignment operator.
376 ParMesh& operator=(const ParMesh &mesh) = delete;
377
378 /// Create a uniformly refined (by any factor) version of @a orig_mesh.
379 /** @param[in] orig_mesh The starting coarse mesh.
380 @param[in] ref_factor The refinement factor, an integer > 1.
381 @param[in] ref_type Specify the positions of the new vertices. The
382 options are BasisType::ClosedUniform or
383 BasisType::GaussLobatto.
384
385 The refinement data which can be accessed with GetRefinementTransforms()
386 is set to reflect the performed refinements.
387
388 @note The constructed ParMesh is linear, i.e. it does not have nodes. */
389 static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type);
390
391 /** Create a mesh by splitting each element of @a orig_mesh into simplices.
392 See @a Mesh::MakeSimplicial for more details. */
393 static ParMesh MakeSimplicial(ParMesh &orig_mesh);
394
395 void Finalize(bool refine = false, bool fix_orientation = false) override;
396
397 void SetAttributes() override;
398
399 /// Checks if any rank in the mesh has boundary elements
400 bool HasBoundaryElements() const override;
401
402 MPI_Comm GetComm() const { return MyComm; }
403 int GetNRanks() const { return NRanks; }
404 int GetMyRank() const { return MyRank; }
405
406 /** Map a global element number to a local element number. If the global
407 element is not on this processor, return -1. */
408 int GetLocalElementNum(long long global_element_num) const;
409
410 /// Map a local element number to a global element number.
411 long long GetGlobalElementNum(int local_element_num) const;
412
413 /** The following functions define global indices for all local vertices,
414 edges, faces, or elements. The global indices have no meaning or
415 significance for ParMesh, but can be used for purposes beyond this class.
416 */
417 /// AMR meshes are not supported.
419 /// AMR meshes are not supported.
421 /// AMR meshes are not supported.
423 /// AMR meshes are supported.
425
427
428 // Face-neighbor elements and vertices
435 // Local face-neighbor elements and vertices ordered by face-neighbor
438
440
441 int *partitioning_cache = nullptr;
442
443 int GetNGroups() const { return gtopo.NGroups(); }
444
445 ///@{ @name These methods require group > 0
446 int GroupNVertices(int group) const { return group_svert.RowSize(group-1); }
447 int GroupNEdges(int group) const { return group_sedge.RowSize(group-1); }
448 int GroupNTriangles(int group) const { return group_stria.RowSize(group-1); }
449 int GroupNQuadrilaterals(int group) const { return group_squad.RowSize(group-1); }
450
451 int GroupVertex(int group, int i) const
452 { return svert_lvert[group_svert.GetRow(group-1)[i]]; }
453 void GroupEdge(int group, int i, int &edge, int &o) const;
454 void GroupTriangle(int group, int i, int &face, int &o) const;
455 void GroupQuadrilateral(int group, int i, int &face, int &o) const;
456 ///@}
457
458 /**
459 * @brief Get the shared edges GroupCommunicator.
460 *
461 * @param[out] sedge_comm
462 */
464 {
465 GetSharedEdgeCommunicator(1, sedge_comm);
466 }
467
468 /**
469 * @brief Get the shared vertices GroupCommunicator.
470 *
471 * @param[out] svert_comm
472 */
474 {
475 GetSharedVertexCommunicator(1, svert_comm);
476 }
477
478 /**
479 * @brief Get the shared face quadrilaterals GroupCommunicator.
480 *
481 * @param[out] squad_comm
482 */
484 {
485 GetSharedQuadCommunicator(1, squad_comm);
486 }
487
488 /**
489 * @brief Get the shared face triangles GroupCommunicator.
490 *
491 * @param[out] stria_comm
492 */
494 {
495 GetSharedTriCommunicator(1, stria_comm);
496 }
497
498 void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[],
499 Array<HYPRE_BigInt> *offsets[]) const;
500
501 /** Return true if the face is interior or shared. In parallel, this
502 method only works if the face neighbor data is exchanged. */
503 inline bool FaceIsTrueInterior(int FaceNo) const { return Mesh::FaceIsTrueInterior(FaceNo); }
504
505 void ExchangeFaceNbrData();
507
508 void SetCurvature(int order, bool discont = false, int space_dim = -1,
509 int ordering = 1) override;
510
511 /** Replace the internal node GridFunction with a new GridFunction defined
512 on the given FiniteElementSpace. The new node coordinates are projected
513 (derived) from the current nodes/vertices. */
514 void SetNodalFESpace(FiniteElementSpace *nfes) override;
516
517 int GetNFaceNeighbors() const { return face_nbr_group.Size(); }
518 int GetNFaceNeighborElements() const { return face_nbr_elements.Size(); }
519 int GetFaceNbrGroup(int fn) const { return face_nbr_group[fn]; }
520 int GetFaceNbrRank(int fn) const;
521
522 /** Similar to Mesh::GetElementFaces */
524 Array<int> &orientation) const;
525
526 /** Similar to Mesh::GetFaceToElementTable with added face-neighbor elements
527 with indices offset by the local number of elements. */
528 /// @note The returned Table should be deleted by the caller
530
531 /// Returns (a pointer to an object containing) the following data:
532 ///
533 /// 1) Elem1No - the index of the first element that contains this face this
534 /// is the element that has the same outward unit normal vector as the
535 /// face;
536 ///
537 /// 2) Elem2No - the index of the second element that contains this face this
538 /// element has outward unit normal vector as the face multiplied with -1;
539 ///
540 /// 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first
541 /// and the second element respectively;
542 ///
543 /// 4) Face - pointer to the ElementTransformation of the face;
544 ///
545 /// 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face
546 /// coordinate system to the element coordinate system (both in their
547 /// reference elements). Used to transform IntegrationPoints from face to
548 /// element. More formally, let:
549 /// TL1, TL2 be the transformations represented by Loc1, Loc2,
550 /// TE1, TE2 - the transformations represented by Elem1, Elem2,
551 /// TF - the transformation represented by Face, then
552 /// TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
553 ///
554 /// 6) FaceGeom - the base geometry for the face.
555 ///
556 /// The mask specifies which fields in the structure to return:
557 /// mask & 1 - Elem1, mask & 2 - Elem2
558 /// mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.
559 /// These mask values are defined in the ConfigMasks enum type as part of the
560 /// FaceElementTransformations class in fem/eltrans.hpp.
561 ///
562 /// @note The returned object is owned by the class and is shared, i.e.,
563 /// calling this function resets pointers obtained from previous calls.
564 /// Also, the returned object should NOT be deleted by the caller.
566 GetFaceElementTransformations(int FaceNo, int mask = 31) override;
567
568 /// @brief Variant of GetFaceElementTransformations using a user allocated
569 /// FaceElementTransformations object.
570 void GetFaceElementTransformations(int FaceNo,
574 int mask = 31) const override;
575
576 /// @brief Get the FaceElementTransformations for the given shared face
577 /// (edge 2D) using the shared face index @a sf. @a fill2 specify if the
578 /// information for elem2 of the face should be computed or not.
579 /// In the returned object, 1 and 2 refer to the local and the neighbor
580 /// elements, respectively.
581 ///
582 /// @note The returned object is owned by the class and is shared, i.e.,
583 /// calling this function resets pointers obtained from previous calls.
584 /// Also, the returned object should NOT be deleted by the caller.
586 GetSharedFaceTransformations(int sf, bool fill2 = true);
587
588 /// @brief Variant of GetSharedFaceTransformations using a user allocated
589 /// FaceElementTransformations object.
594 bool fill2 = true) const;
595
596 /// @brief Get the FaceElementTransformations for the given shared face
597 /// (edge 2D) using the face index @a FaceNo. @a fill2 specify if the
598 /// information for elem2 of the face should be computed or not.
599 /// In the returned object, 1 and 2 refer to the local and the neighbor
600 /// elements, respectively.
601 ///
602 /// @note The returned object is owned by the class and is shared, i.e.,
603 /// calling this function resets pointers obtained from previous calls.
604 /// Also, the returned object should NOT be deleted by the caller.
606 GetSharedFaceTransformationsByLocalIndex(int FaceNo, bool fill2 = true);
607
608 /// @brief Variant of GetSharedFaceTransformationsByLocalIndex using a user
609 /// allocated FaceElementTransformations object.
614 bool fill2 = true) const;
615
616 /// @brief Returns a pointer to the transformation defining the i-th face
617 /// neighbor.
618 ///
619 /// @note The returned object is owned by the class and is shared, i.e.,
620 /// calling this function resets pointers obtained from previous calls.
621 /// Also, the returned object should NOT be deleted by the caller.
623
624 /// @brief Variant of GetFaceNbrElementTransformation using a user allocated
625 /// IsoparametricTransformation object.
626 void GetFaceNbrElementTransformation(int FaceNo,
627 IsoparametricTransformation &ElTr) const;
628
629 /// Get the size of the i-th face neighbor element relative to the reference
630 /// element.
631 real_t GetFaceNbrElementSize(int i, int type = 0);
632
633 /// Return the number of shared faces (3D), edges (2D), vertices (1D)
634 int GetNSharedFaces() const;
635
636 /// Return the local face index for the given shared face.
637 int GetSharedFace(int sface) const;
638
639 /** @brief Returns the number of local faces according to the requested type,
640 does not count master non-conforming faces.
641
642 If type==Boundary returns only the number of true boundary faces
643 contrary to GetNBE() that returns all "boundary" elements which may
644 include actual interior faces.
645 Similarly, if type==Interior, only the true interior faces (including
646 shared faces) are counted excluding all master non-conforming faces. */
647 int GetNFbyType(FaceType type) const override;
648
650 { MFEM_ABORT("Generation of boundary elements works properly only on serial meshes."); }
651
652 /// See the remarks for the serial version in mesh.hpp
653 MFEM_DEPRECATED void ReorientTetMesh() override;
654
655 /// Utility function: sum integers from all processors (Allreduce).
656 long long ReduceInt(int value) const override;
657
658 /** Load balance the mesh by equipartitioning the global space-filling
659 sequence of elements. Works for nonconforming meshes only. */
660 void Rebalance();
661
662 /** Load balance a nonconforming mesh using a user-defined partition.
663 Each local element 'i' is migrated to processor rank 'partition[i]',
664 for 0 <= i < GetNE(). */
665 void Rebalance(const Array<int> &partition);
666
667 /** Save the mesh in a parallel mesh format. If @a comments is non-empty, it
668 will be printed after the first line of the file, and each line should
669 begin with '#'. */
670 void ParPrint(std::ostream &out, const std::string &comments = "") const;
671
672 // Enable Print() to add the parallel interface as boundary (typically used
673 // for visualization purposes)
674 void SetPrintShared(bool print) { print_shared = print; }
675
676 /** Print the part of the mesh in the calling processor using the mfem v1.0
677 format. Depending on SetPrintShared(), the parallel interface can be
678 added as boundary for visualization (true by default). If @a comments is
679 non-empty, it will be printed after the first line of the file, and each
680 line should begin with '#'. */
681 void Print(std::ostream &out = mfem::out,
682 const std::string &comments = "") const override;
683
684 /// Save the ParMesh to files (one for each MPI rank). The files will be
685 /// given suffixes according to the MPI rank. The mesh will be written to the
686 /// files using ParMesh::Print. The given @a precision will be used for ASCII
687 /// output.
688 void Save(const std::string &fname, int precision=16) const override;
689
690#ifdef MFEM_USE_ADIOS2
691 /** Print the part of the mesh in the calling processor using adios2 bp
692 format. */
693 void Print(adios2stream &out) const override;
694#endif
695
696 /** Print the part of the mesh in the calling processor adding the interface
697 as boundary (for visualization purposes) using Netgen/Truegrid format .*/
698 void PrintXG(std::ostream &out = mfem::out) const override;
699
700 /** Write the mesh to the stream 'out' on Process 0 in a form suitable for
701 visualization: the mesh is written as a disjoint mesh and the shared
702 boundary is added to the actual boundary; both the element and boundary
703 attributes are set to the processor number. If @a comments is non-empty,
704 it will be printed after the first line of the file, and each line should
705 begin with '#'. */
706 void PrintAsOne(std::ostream &out = mfem::out,
707 const std::string &comments = "") const;
708
709 /** Write the mesh to the stream 'out' on Process 0 as a serial mesh. The
710 output mesh does not have any duplication of vertices/nodes at processor
711 boundaries. If @a comments is non-empty, it will be printed after the
712 first line of the file, and each line should begin with '#'. */
713 void PrintAsSerial(std::ostream &out = mfem::out,
714 const std::string &comments = "") const;
715
716 /** Returns a Serial mesh on MPI rank @a save_rank that does not have any
717 duplication of vertices/nodes at processor boundaries. */
718 Mesh GetSerialMesh(int save_rank) const;
719
720 /// Save the mesh as a single file (using ParMesh::PrintAsOne). The given
721 /// @a precision is used for ASCII output.
722 void SaveAsOne(const std::string &fname, int precision=16) const;
723
724 /// Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'
725 void PrintAsOneXG(std::ostream &out = mfem::out);
726
727 /** Print the mesh in parallel PVTU format. The PVTU and VTU files will be
728 stored in the directory specified by @a pathname. If the directory does
729 not exist, it will be created. */
730 void PrintVTU(std::string pathname,
732 bool high_order_output=false,
733 int compression_level=0,
734 bool bdr=false) override;
735
736 /// Parallel version of Mesh::Load().
737 void Load(std::istream &input, int generate_edges = 0,
738 int refine = 1, bool fix_orientation = true) override;
739
740 /// Returns the minimum and maximum corners of the mesh bounding box. For
741 /// high-order meshes, the geometry is refined first "ref" times.
742 void GetBoundingBox(Vector &p_min, Vector &p_max, int ref = 2);
743
744 void GetCharacteristics(real_t &h_min, real_t &h_max,
745 real_t &kappa_min, real_t &kappa_max);
746
747 /// Swaps internal data with another ParMesh, including non-geometry members.
748 /// See @a Mesh::Swap
749 void Swap(ParMesh &other);
750
751 /// Print various parallel mesh stats
752 void PrintInfo(std::ostream &out = mfem::out) override;
753
754 int FindPoints(DenseMatrix& point_mat, Array<int>& elem_ids,
755 Array<IntegrationPoint>& ips, bool warn = true,
756 InverseElementTransformation *inv_trans = NULL) override;
757
758 /// Debugging method
759 void PrintSharedEntities(const std::string &fname_prefix) const;
760
761 virtual ~ParMesh();
762};
763
764}
765
766#endif // MFEM_USE_MPI
767
768#endif
int Size() const
Return the logical size of the array.
Definition array.hpp:144
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Abstract data type element.
Definition element.hpp:29
Type
Constants for the classes derived from Element.
Definition element.hpp:41
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:484
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
int NGroups() const
Return the number of groups.
The inverse transformation of a given ElementTransformation.
Definition eltrans.hpp:187
A standard isoparametric element transformation.
Definition eltrans.hpp:363
List of integer sets.
Definition sets.hpp:63
Mesh data type.
Definition mesh.hpp:56
Array< Element * > faces
Definition mesh.hpp:99
bool FaceIsTrueInterior(int FaceNo) const
Definition mesh.hpp:540
void Swap(Mesh &other, bool non_geometry)
Definition mesh.cpp:10378
Abstract parallel finite element space.
Definition pfespace.hpp:29
Class for parallel meshes.
Definition pmesh.hpp:34
Mesh GetSerialMesh(int save_rank) const
Definition pmesh.cpp:5292
void GetCharacteristics(real_t &h_min, real_t &h_max, real_t &kappa_min, real_t &kappa_max)
Definition pmesh.cpp:6174
void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0) override
This function is not public anymore. Use GeneralRefinement instead.
Definition pmesh.cpp:3889
int GroupNQuadrilaterals(int group) const
Definition pmesh.hpp:449
ElementTransformation * GetFaceNbrElementTransformation(int FaceNo)
Returns a pointer to the transformation defining the i-th face neighbor.
Definition pmesh.cpp:3088
void GetFaceSplittings(const int *fv, const HashTable< Hashed2 > &v_to_v, Array< unsigned > &codes)
Append codes identifying how the given face has been split to codes.
Definition pmesh.cpp:1872
void GetGhostFaceTransformation(FaceElementTransformations *FElTr, Element::Type face_type, Geometry::Type face_geom) const
Definition pmesh.hpp:155
Table send_face_nbr_elements
Definition pmesh.hpp:436
Array< int > face_nbr_vertices_offset
Definition pmesh.hpp:432
void BuildVertexGroup(int ngroups, const Table &vert_element)
Definition pmesh.cpp:706
void PrintVTU(std::string pathname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr=false) override
Definition pmesh.cpp:6398
void NURBSUniformRefinement(int rf=2, real_t tol=1.0e-12) override
Refine NURBS mesh, with an optional refinement factor.
Definition pmesh.cpp:4570
MPI_Comm GetComm() const
Definition pmesh.hpp:402
void PrintXG(std::ostream &out=mfem::out) const override
Definition pmesh.cpp:4586
Array< Element * > shared_edges
Definition pmesh.hpp:69
int GetMyRank() const
Definition pmesh.hpp:404
int GetEdgeSplittings(Element *edge, const DSTable &v_to_v, int *middle)
Return a number(0-1) identifying how the given edge has been split.
Definition pmesh.cpp:1857
void RefineGroups(const DSTable &v_to_v, int *middle)
Update the groups after triangle refinement.
Definition pmesh.cpp:4061
void ParPrint(std::ostream &out, const std::string &comments="") const
Definition pmesh.cpp:6314
void GetSharedTriCommunicator(int ordering, GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
Definition pmesh.cpp:1718
bool NonconformingDerefinement(Array< real_t > &elem_error, real_t threshold, int nc_limit=0, int op=1) override
NC version of GeneralDerefinement.
Definition pmesh.cpp:3947
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition pmesh.cpp:3153
Array< int > sface_lface
Definition pmesh.hpp:86
void SaveAsOne(const std::string &fname, int precision=16) const
Definition pmesh.cpp:5609
void ExchangeFaceNbrData()
Definition pmesh.cpp:2069
Table group_sedge
Definition pmesh.hpp:78
int GetNRanks() const
Definition pmesh.hpp:403
void BuildEdgeGroup(int ngroups, const Table &edge_element)
Definition pmesh.cpp:680
Table group_svert
Shared objects in each group.
Definition pmesh.hpp:77
MFEM_DEPRECATED void ReorientTetMesh() override
See the remarks for the serial version in mesh.hpp.
Definition pmesh.cpp:3226
int GroupVertex(int group, int i) const
Definition pmesh.hpp:451
void UniformRefinement2D() override
Refine a mixed 2D mesh uniformly.
Definition pmesh.cpp:4518
void BuildSharedFaceElems(int ntri_faces, int nquad_faces, const Mesh &mesh, int *partitioning, const STable3D *faces_tbl, const Array< int > &face_group, const Array< int > &vert_global_local)
Definition pmesh.cpp:732
bool WantSkipSharedMaster(const NCMesh::Master &master) const
Definition pmesh.cpp:4788
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31) override
Definition pmesh.cpp:2897
void BuildSharedVertMapping(int nvert, const Table *vert_element, const Array< int > &vert_global_local)
Definition pmesh.cpp:846
long long GetGlobalElementNum(int local_element_num) const
Map a local element number to a global element number.
Definition pmesh.cpp:1547
Table * GetFaceToAllElementTable() const
Definition pmesh.cpp:2839
virtual ~ParMesh()
Definition pmesh.cpp:6741
void ReduceMeshGen()
Definition pmesh.cpp:892
void GetGlobalElementIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are supported.
Definition pmesh.cpp:6668
void AddTriFaces(const Array< int > &v, const std::unique_ptr< STable3D > &faces, const std::unique_ptr< STable3D > &shared_faces, int elem, int start, int end, const int fverts[][N])
Helper function for adding triangle face neighbor element to face table entries. Have to use a templa...
Definition pmesh.cpp:2677
void FindSharedFaces(const Mesh &mesh, const int *partition, Array< int > &face_group, ListOfIntegerSets &groups)
Definition pmesh.cpp:516
long long glob_elem_offset
Definition pmesh.hpp:95
void GetSharedEdgeCommunicator(int ordering, GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
Definition pmesh.cpp:1646
void PrintInfo(std::ostream &out=mfem::out) override
Print various parallel mesh stats.
Definition pmesh.cpp:6191
int GetNFbyType(FaceType type) const override
Returns the number of local faces according to the requested type, does not count master non-conformi...
Definition pmesh.cpp:3194
void GetGlobalVertexIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Definition pmesh.cpp:6601
int GetNFaceNeighborElements() const
Definition pmesh.hpp:518
bool HasBoundaryElements() const override
Checks if any rank in the mesh has boundary elements.
Definition pmesh.cpp:1607
bool print_shared
Definition pmesh.hpp:101
void LocalRefinement(const Array< int > &marked_el, int type=3) override
This function is not public anymore. Use GeneralRefinement instead.
Definition pmesh.cpp:3385
void GetFaceNbrElementFaces(int i, Array< int > &faces, Array< int > &orientation) const
Definition pmesh.cpp:2814
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1) override
Set the curvature of the mesh nodes using the given polynomial degree.
Definition pmesh.cpp:2007
void SetPrintShared(bool print)
Definition pmesh.hpp:674
void GetSharedTriCommunicator(GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
Definition pmesh.hpp:493
long glob_offset_sequence
Definition pmesh.hpp:96
int GetLocalElementNum(long long global_element_num) const
Definition pmesh.cpp:1539
void SetAttributes() override
Determine the sets of unique attribute values in domain and boundary elements.
Definition pmesh.cpp:1589
int FindPoints(DenseMatrix &point_mat, Array< int > &elem_ids, Array< IntegrationPoint > &ips, bool warn=true, InverseElementTransformation *inv_trans=NULL) override
Find the ids of the elements that contain the given points, and their corresponding reference coordin...
Definition pmesh.cpp:6465
void GetSharedQuadCommunicator(int ordering, GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
Definition pmesh.cpp:1694
void BuildFaceNbrElementToFaceTable()
Definition pmesh.cpp:2710
void MakeRefined_(ParMesh &orig_mesh, int ref_factor, int ref_type)
Internal function used in ParMesh::MakeRefined (and related constructor)
Definition pmesh.cpp:1140
void GetGlobalFaceIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Definition pmesh.cpp:6640
Array< Vertex > face_nbr_vertices
Definition pmesh.hpp:434
void UniformRefineGroups2D(int old_nv)
Definition pmesh.cpp:4318
void GetGlobalEdgeIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Definition pmesh.cpp:6617
void ExchangeFaceNbrNodes()
Definition pmesh.cpp:2552
void SetNodalFESpace(FiniteElementSpace *nfes) override
Definition pmesh.cpp:2028
void UniformRefinement3D() override
Refine a mixed 3D mesh uniformly.
Definition pmesh.cpp:4542
void Rebalance()
Definition pmesh.cpp:4009
bool have_face_nbr_data
Definition pmesh.hpp:429
long long ReduceInt(int value) const override
Utility function: sum integers from all processors (Allreduce).
Definition pmesh.cpp:6307
Table send_face_nbr_vertices
Definition pmesh.hpp:437
ParMesh()
Default constructor. Create an empty ParMesh.
Definition pmesh.hpp:335
int BuildLocalElements(const Mesh &global_mesh, const int *partitioning, const Array< int > &vert_global_local)
Fills out partitioned Mesh::elements.
Definition pmesh.cpp:366
MPI_Comm MyComm
Definition pmesh.hpp:45
bool DecodeFaceSplittings(HashTable< Hashed2 > &v_to_v, const int *v, const Array< unsigned > &codes, int &pos)
Definition pmesh.cpp:1907
static ParMesh MakeSimplicial(ParMesh &orig_mesh)
Definition pmesh.cpp:1382
Array< Vert4 > shared_quads
Definition pmesh.hpp:74
void GenerateBoundaryElements() override
Definition pmesh.hpp:649
void LoadSharedEntities(std::istream &input)
Definition pmesh.cpp:986
void Finalize(bool refine=false, bool fix_orientation=false) override
Finalize the construction of a general Mesh.
Definition pmesh.cpp:1524
int GetFaceNbrGroup(int fn) const
Definition pmesh.hpp:519
Table group_squad
Definition pmesh.hpp:80
void GetSharedEdgeCommunicator(GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
Definition pmesh.hpp:463
void BuildFaceGroup(int ngroups, const Mesh &mesh, const Array< int > &face_group, int &nstria, int &nsquad)
Definition pmesh.cpp:634
void BuildSharedEdgeElems(int nedges, Mesh &mesh, const Array< int > &vert_global_local, const Table *edge_element)
Definition pmesh.cpp:809
Array< int > svert_lvert
Shared to local index mapping.
Definition pmesh.hpp:83
Array< int > sedge_ledge
Definition pmesh.hpp:84
Array< Element * > face_nbr_elements
Definition pmesh.hpp:433
GroupTopology gtopo
Definition pmesh.hpp:426
void Destroy()
Definition pmesh.cpp:6725
FaceElementTransformations * GetSharedFaceTransformationsByLocalIndex(int FaceNo, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the face index FaceNo....
Definition pmesh.cpp:2942
void PrintSharedEntities(const std::string &fname_prefix) const
Debugging method.
Definition pmesh.cpp:6514
int GroupNTriangles(int group) const
Definition pmesh.hpp:448
void GetSharedVertexCommunicator(int ordering, GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
Definition pmesh.cpp:1670
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
Definition pmesh.cpp:3172
void EnsureParNodes()
If the mesh is curved, make sure 'Nodes' is ParGridFunction.
Definition pmesh.cpp:2049
int GroupNEdges(int group) const
Definition pmesh.hpp:447
void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[], Array< HYPRE_BigInt > *offsets[]) const
Definition pmesh.cpp:1943
void MarkTetMeshForRefinement(const DSTable &v_to_v) override
Definition pmesh.cpp:1742
Array< int > face_nbr_group
Definition pmesh.hpp:430
Array< int > face_nbr_elements_offset
Definition pmesh.hpp:431
ParMesh & operator=(ParMesh &&mesh)
Move assignment operator.
Definition pmesh.cpp:100
void DeleteFaceNbrData()
Definition pmesh.cpp:1986
void Load(std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true) override
Parallel version of Mesh::Load().
Definition pmesh.cpp:950
int GetNFaceNeighbors() const
Definition pmesh.hpp:517
void UniformRefineGroups3D(int old_nv, int old_nedges, const DSTable &old_v_to_v, const STable3D &old_faces, Array< int > *f2qf)
Definition pmesh.cpp:4369
void GetGhostFaceTransformation(FaceElementTransformations &FElTr, Element::Type face_type, Geometry::Type face_geom) const
Definition pmesh.cpp:3056
void RebalanceImpl(const Array< int > *partition)
Definition pmesh.cpp:4019
int FindSharedEdges(const Mesh &mesh, const int *partition, Table *&edge_element, ListOfIntegerSets &groups)
Definition pmesh.cpp:543
void FinalizeParTopo()
Definition pmesh.cpp:898
void PrintAsOneXG(std::ostream &out=mfem::out)
Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'.
Definition pmesh.cpp:5620
std::unique_ptr< Table > face_nbr_el_to_face
Definition pmesh.hpp:90
Array< Vert3 > shared_trias
Definition pmesh.hpp:73
void GroupQuadrilateral(int group, int i, int &face, int &o) const
Definition pmesh.cpp:1635
void Save(const std::string &fname, int precision=16) const override
Definition pmesh.cpp:4940
void DistributeAttributes(Array< int > &attr)
Ensure that bdr_attributes and attributes agree across processors.
Definition pmesh.cpp:1553
int * partitioning_cache
Definition pmesh.hpp:441
ParMesh & operator=(const ParMesh &mesh)=delete
Explicitly delete the copy assignment operator.
STable3D * GetSharedFacesTable()
Definition pmesh.cpp:2612
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
Definition pmesh.cpp:1375
int GetNGroups() const
Definition pmesh.hpp:443
real_t GetFaceNbrElementSize(int i, int type=0)
Definition pmesh.cpp:3148
ParNCMesh * pncmesh
Definition pmesh.hpp:439
std::unique_ptr< Table > face_nbr_el_ori
orientations for each face (from nbr processor)
Definition pmesh.hpp:92
bool FaceIsTrueInterior(int FaceNo) const
Definition pmesh.hpp:503
int GetFaceNbrRank(int fn) const
Definition pmesh.cpp:2796
int FindSharedVertices(const int *partition, Table *vertex_element, ListOfIntegerSets &groups)
Definition pmesh.cpp:597
void GroupTriangle(int group, int i, int &face, int &o) const
Definition pmesh.cpp:1624
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the shared face index sf...
Definition pmesh.cpp:2923
void PrintAsSerial(std::ostream &out=mfem::out, const std::string &comments="") const
Definition pmesh.cpp:5281
void GetSharedVertexCommunicator(GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
Definition pmesh.hpp:473
void GetSharedQuadCommunicator(GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
Definition pmesh.hpp:483
int BuildLocalBoundary(const Mesh &global_mesh, const int *partitioning, const Array< int > &vert_global_local, Array< bool > &activeBdrElem, Table *&edge_element)
Fills out partitioned Mesh::boundary.
Definition pmesh.cpp:396
Table group_stria
Definition pmesh.hpp:79
int BuildLocalVertices(const Mesh &global_mesh, const int *partitioning, Array< int > &vert_global_local)
Fills out partitioned Mesh::vertices.
Definition pmesh.cpp:318
void ComputeGlobalElementOffset() const
Definition pmesh.cpp:879
void GetBoundingBox(Vector &p_min, Vector &p_max, int ref=2)
Definition pmesh.cpp:6154
void PrintAsOne(std::ostream &out=mfem::out, const std::string &comments="") const
Definition pmesh.cpp:4968
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
void Swap(ParMesh &other)
Definition pmesh.cpp:6682
void GroupEdge(int group, int i, int &edge, int &o) const
Definition pmesh.cpp:1616
int GroupNVertices(int group) const
Definition pmesh.hpp:446
A parallel extension of the NCMesh class.
Definition pncmesh.hpp:65
Class for PUMI parallel meshes.
Definition pumi.hpp:70
Subdomain representation of a topological parent in another ParMesh.
Definition psubmesh.hpp:52
Symmetric 3D Table stored as an array of rows each of which has a stack of column,...
Definition stable3d.hpp:35
int RowSize(int i) const
Definition table.hpp:108
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:187
Vector data type.
Definition vector.hpp:80
HYPRE_Int HYPRE_BigInt
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
VTKFormat
Data array format for VTK and VTU files.
Definition vtk.hpp:99
@ ASCII
Data arrays will be written in ASCII format.
float real_t
Definition config.hpp:43
FaceType
Definition mesh.hpp:47
void Set(const int *w)
Definition pmesh.hpp:54
void Set(int v0, int v1, int v2)
Definition pmesh.hpp:53
Vert3(int v0, int v1, int v2)
Definition pmesh.hpp:52
void Set(int v0, int v1, int v2, int v3)
Definition pmesh.hpp:63
void Set(const int *w)
Definition pmesh.hpp:65
Vert4(int v0, int v1, int v2, int v3)
Definition pmesh.hpp:61