MFEM  v4.5.2
Finite element discretization library
pmesh.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_PMESH
13 #define MFEM_PMESH
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include "../general/communication.hpp"
20 #include "../general/globals.hpp"
21 #include "mesh.hpp"
22 #include "pncmesh.hpp"
23 #include <iostream>
24 
25 namespace mfem
26 {
27 #ifdef MFEM_USE_PUMI
28 class ParPumiMesh;
29 #endif
30 
31 /// Class for parallel meshes
32 class ParMesh : public Mesh
33 {
34 protected:
35  friend class ParSubMesh;
36 
37  MPI_Comm MyComm;
38  int NRanks, MyRank;
39 
40  struct Vert3
41  {
42  int v[3];
43  Vert3() = default;
44  Vert3(int v0, int v1, int v2) { v[0] = v0; v[1] = v1; v[2] = v2; }
45  void Set(int v0, int v1, int v2) { v[0] = v0; v[1] = v1; v[2] = v2; }
46  void Set(const int *w) { v[0] = w[0]; v[1] = w[1]; v[2] = w[2]; }
47  };
48 
49  struct Vert4
50  {
51  int v[4];
52  Vert4() = default;
53  Vert4(int v0, int v1, int v2, int v3)
54  { v[0] = v0; v[1] = v1; v[2] = v2; v[3] = v3; }
55  void Set(int v0, int v1, int v2, int v3)
56  { v[0] = v0; v[1] = v1; v[2] = v2; v[3] = v3; }
57  void Set(const int *w)
58  { v[0] = w[0]; v[1] = w[1]; v[2] = w[2]; v[3] = w[3]; }
59  };
60 
62  // shared face id 'i' is:
63  // * triangle id 'i', if i < shared_trias.Size()
64  // * quad id 'i-shared_trias.Size()', otherwise
67 
68  /// Shared objects in each group.
71  Table group_stria; // contains shared triangle indices
72  Table group_squad; // contains shared quadrilateral indices
73 
74  /// Shared to local index mapping.
77  // sface ids: all triangles first, then all quads
79 
81  Table face_nbr_el_ori; // orientations for each face (from nbr processor)
82 
84 
85  // glob_elem_offset + local element number defines a global element numbering
86  mutable long long glob_elem_offset;
87  mutable long glob_offset_sequence;
88  void ComputeGlobalElementOffset() const;
89 
90  // Enable Print() to add the parallel interface as boundary (typically used
91  // for visualization purposes)
92  bool print_shared = true;
93 
94  /// Create from a nonconforming mesh.
95  ParMesh(const ParNCMesh &pncmesh);
96 
97  // Convert the local 'meshgen' to a global one.
98  void ReduceMeshGen();
99 
100  // Determine sedge_ledge and sface_lface.
101  void FinalizeParTopo();
102 
103  // Mark all tets to ensure consistency across MPI tasks; also mark the
104  // shared and boundary triangle faces using the consistently marked tets.
105  void MarkTetMeshForRefinement(DSTable &v_to_v) override;
106 
107  /// Return a number(0-1) identifying how the given edge has been split
108  int GetEdgeSplittings(Element *edge, const DSTable &v_to_v, int *middle);
109  /// Append codes identifying how the given face has been split to @a codes
110  void GetFaceSplittings(const int *fv, const HashTable<Hashed2> &v_to_v,
111  Array<unsigned> &codes);
112 
113  bool DecodeFaceSplittings(HashTable<Hashed2> &v_to_v, const int *v,
114  const Array<unsigned> &codes, int &pos);
115 
116  STable3D *GetFaceNbrElementToFaceTable(int ret_ftbl = 0);
117 
119  int i, IsoparametricTransformation *ElTr);
120 
122  FaceElementTransformations* FETr, Element::Type face_type,
123  Geometry::Type face_geom);
124 
125  /// Update the groups after triangle refinement
126  void RefineGroups(const DSTable &v_to_v, int *middle);
127 
128  /// Update the groups after tetrahedron refinement
129  void RefineGroups(int old_nv, const HashTable<Hashed2> &v_to_v);
130 
131  void UniformRefineGroups2D(int old_nv);
132 
133  // f2qf can be NULL if all faces are quads or there are no quad faces
134  void UniformRefineGroups3D(int old_nv, int old_nedges,
135  const DSTable &old_v_to_v,
136  const STable3D &old_faces,
137  Array<int> *f2qf);
138 
139  void ExchangeFaceNbrData(Table *gr_sface, int *s2l_face);
140 
141  /// Refine a mixed 2D mesh uniformly.
142  void UniformRefinement2D() override;
143 
144  /// Refine a mixed 3D mesh uniformly.
145  void UniformRefinement3D() override;
146 
147  void NURBSUniformRefinement() override;
148 
149  /// This function is not public anymore. Use GeneralRefinement instead.
150  void LocalRefinement(const Array<int> &marked_el, int type = 3) override;
151 
152  /// This function is not public anymore. Use GeneralRefinement instead.
153  void NonconformingRefinement(const Array<Refinement> &refinements,
154  int nc_limit = 0) override;
155 
156  bool NonconformingDerefinement(Array<double> &elem_error,
157  double threshold, int nc_limit = 0,
158  int op = 1) override;
159 
160  void RebalanceImpl(const Array<int> *partition);
161 
162  void DeleteFaceNbrData();
163 
164  bool WantSkipSharedMaster(const NCMesh::Master &master) const;
165 
166  /// Fills out partitioned Mesh::vertices
167  int BuildLocalVertices(const Mesh& global_mesh, const int *partitioning,
168  Array<int> &vert_global_local);
169 
170  /// Fills out partitioned Mesh::elements
171  int BuildLocalElements(const Mesh& global_mesh, const int *partitioning,
172  const Array<int> &vert_global_local);
173 
174  /// Fills out partitioned Mesh::boundary
175  int BuildLocalBoundary(const Mesh& global_mesh, const int *partitioning,
176  const Array<int> &vert_global_local,
177  Array<bool>& activeBdrElem,
178  Table* &edge_element);
179 
180  void FindSharedFaces(const Mesh &mesh, const int* partition,
181  Array<int>& face_group,
182  ListOfIntegerSets& groups);
183 
184  int FindSharedEdges(const Mesh &mesh, const int* partition,
185  Table* &edge_element, ListOfIntegerSets& groups);
186 
187  int FindSharedVertices(const int *partition, Table* vertex_element,
188  ListOfIntegerSets& groups);
189 
190  void BuildFaceGroup(int ngroups, const Mesh &mesh,
191  const Array<int>& face_group,
192  int &nstria, int &nsquad);
193 
194  void BuildEdgeGroup(int ngroups, const Table& edge_element);
195 
196  void BuildVertexGroup(int ngroups, const Table& vert_element);
197 
198  void BuildSharedFaceElems(int ntri_faces, int nquad_faces,
199  const Mesh &mesh, int *partitioning,
200  const STable3D *faces_tbl,
201  const Array<int> &face_group,
202  const Array<int> &vert_global_local);
203 
204  void BuildSharedEdgeElems(int nedges, Mesh &mesh,
205  const Array<int> &vert_global_local,
206  const Table *edge_element);
207 
208  void BuildSharedVertMapping(int nvert, const Table* vert_element,
209  const Array<int> &vert_global_local);
210 
211  /**
212  * @brief Get the shared edges GroupCommunicator.
213  *
214  * The output of the shared edges is chosen by the @a ordering parameter with
215  * the following options
216  * 0: Internal ordering. Not exposed to public interfaces.
217  * 1: Contiguous ordering.
218  *
219  * @param[in] ordering Ordering for the shared edges.
220  * @param[out] sedge_comm
221  */
222  void GetSharedEdgeCommunicator(int ordering,
223  GroupCommunicator& sedge_comm) const;
224 
225  /**
226  * @brief Get the shared vertices GroupCommunicator.
227  *
228  * The output of the shared vertices is chosen by the @a ordering parameter
229  * with the following options
230  * 0: Internal ordering. Not exposed to public interfaces.
231  * 1: Contiguous ordering.
232  *
233  * @param[in] ordering
234  * @param[out] svert_comm
235  */
236  void GetSharedVertexCommunicator(int ordering,
237  GroupCommunicator& svert_comm) const;
238 
239  /**
240  * @brief Get the shared face quadrilaterals GroupCommunicator.
241  *
242  * The output of the shared face quadrilaterals is chosen by the @a ordering
243  * parameter with the following options
244  * 0: Internal ordering. Not exposed to public interfaces.
245  * 1: Contiguous ordering.
246  *
247  * @param[in] ordering
248  * @param[out] squad_comm
249  */
250  void GetSharedQuadCommunicator(int ordering,
251  GroupCommunicator& squad_comm) const;
252 
253  /**
254  * @brief Get the shared face triangles GroupCommunicator.
255  *
256  * The output of the shared face triangles is chosen by the @a ordering
257  * parameter with the following options
258  * 0: Internal ordering. Not exposed to public interfaces.
259  * 1: Contiguous ordering.
260  *
261  * @param[in] ordering
262  * @param[out] stria_comm
263  */
264  void GetSharedTriCommunicator(int ordering,
265  GroupCommunicator& stria_comm) const;
266 
267  // Similar to Mesh::GetFacesTable()
269 
270  /// Ensure that bdr_attributes and attributes agree across processors
271  void DistributeAttributes(Array<int> &attr);
272 
273  void LoadSharedEntities(std::istream &input);
274 
275  /// If the mesh is curved, make sure 'Nodes' is ParGridFunction.
276  /** Note that this method is not related to the public 'Mesh::EnsureNodes`.*/
277  void EnsureParNodes();
278 
279  /// Internal function used in ParMesh::MakeRefined (and related constructor)
280  void MakeRefined_(ParMesh &orig_mesh, int ref_factor, int ref_type);
281 
282  // Mark Mesh::Swap as protected, should use ParMesh::Swap to swap @a ParMesh
283  // objects.
284  using Mesh::Swap;
285 
286  void Destroy();
287 
288 public:
289  /// Default constructor. Create an empty @a ParMesh.
292  have_face_nbr_data(false), pncmesh(NULL) { }
293 
294  /// Create a parallel mesh by partitioning a serial Mesh.
295  /** The mesh is partitioned automatically or using external partitioning
296  data (the optional parameter 'partitioning_[i]' contains the desired MPI
297  rank for element 'i'). Automatic partitioning uses METIS for conforming
298  meshes and quick space-filling curve equipartitioning for nonconforming
299  meshes (elements of nonconforming meshes should ideally be ordered as a
300  sequence of face-neighbors). */
301  ParMesh(MPI_Comm comm, Mesh &mesh, int *partitioning_ = NULL,
302  int part_method = 1);
303 
304  /** Copy constructor. Performs a deep copy of (almost) all data, so that the
305  source mesh can be modified (e.g. deleted, refined) without affecting the
306  new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the
307  nodes, if present. */
308  explicit ParMesh(const ParMesh &pmesh, bool copy_nodes = true);
309 
310  /// Read a parallel mesh, each MPI rank from its own file/stream.
311  /** The @a refine parameter is passed to the method Mesh::Finalize(). */
312  ParMesh(MPI_Comm comm, std::istream &input, bool refine = true);
313 
314  /// Deprecated: see @a ParMesh::MakeRefined
315  MFEM_DEPRECATED
316  ParMesh(ParMesh *orig_mesh, int ref_factor, int ref_type);
317 
318  /// Move constructor. Used for named constructors.
319  ParMesh(ParMesh &&mesh);
320 
321  /// Move assignment operator.
322  ParMesh& operator=(ParMesh &&mesh);
323 
324  /// Explicitly delete the copy assignment operator.
325  ParMesh& operator=(const ParMesh &mesh) = delete;
326 
327  /// Create a uniformly refined (by any factor) version of @a orig_mesh.
328  /** @param[in] orig_mesh The starting coarse mesh.
329  @param[in] ref_factor The refinement factor, an integer > 1.
330  @param[in] ref_type Specify the positions of the new vertices. The
331  options are BasisType::ClosedUniform or
332  BasisType::GaussLobatto.
333 
334  The refinement data which can be accessed with GetRefinementTransforms()
335  is set to reflect the performed refinements.
336 
337  @note The constructed ParMesh is linear, i.e. it does not have nodes. */
338  static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type);
339 
340  /** Create a mesh by splitting each element of @a orig_mesh into simplices.
341  See @a Mesh::MakeSimplicial for more details. */
342  static ParMesh MakeSimplicial(ParMesh &orig_mesh);
343 
344  void Finalize(bool refine = false, bool fix_orientation = false) override;
345 
346  void SetAttributes() override;
347 
348  /// Checks if any rank in the mesh has boundary elements
349  bool HasBoundaryElements() const override;
350 
351  MPI_Comm GetComm() const { return MyComm; }
352  int GetNRanks() const { return NRanks; }
353  int GetMyRank() const { return MyRank; }
354 
355  /** Map a global element number to a local element number. If the global
356  element is not on this processor, return -1. */
357  int GetLocalElementNum(long long global_element_num) const;
358 
359  /// Map a local element number to a global element number.
360  long long GetGlobalElementNum(int local_element_num) const;
361 
362  /** The following functions define global indices for all local vertices,
363  edges, faces, or elements. The global indices have no meaning or
364  significance for ParMesh, but can be used for purposes beyond this class.
365  */
366  /// AMR meshes are not supported.
368  /// AMR meshes are not supported.
370  /// AMR meshes are not supported.
372  /// AMR meshes are supported.
374 
376 
377  // Face-neighbor elements and vertices
384  // Local face-neighbor elements and vertices ordered by face-neighbor
387 
389 
390  int *partitioning_cache = nullptr;
391 
392  int GetNGroups() const { return gtopo.NGroups(); }
393 
394  ///@{ @name These methods require group > 0
395  int GroupNVertices(int group) const { return group_svert.RowSize(group-1); }
396  int GroupNEdges(int group) const { return group_sedge.RowSize(group-1); }
397  int GroupNTriangles(int group) const { return group_stria.RowSize(group-1); }
398  int GroupNQuadrilaterals(int group) const { return group_squad.RowSize(group-1); }
399 
400  int GroupVertex(int group, int i) const
401  { return svert_lvert[group_svert.GetRow(group-1)[i]]; }
402  void GroupEdge(int group, int i, int &edge, int &o) const;
403  void GroupTriangle(int group, int i, int &face, int &o) const;
404  void GroupQuadrilateral(int group, int i, int &face, int &o) const;
405  ///@}
406 
407  /**
408  * @brief Get the shared edges GroupCommunicator.
409  *
410  * @param[out] sedge_comm
411  */
413  {
414  GetSharedEdgeCommunicator(1, sedge_comm);
415  }
416 
417  /**
418  * @brief Get the shared vertices GroupCommunicator.
419  *
420  * @param[out] svert_comm
421  */
423  {
424  GetSharedVertexCommunicator(1, svert_comm);
425  }
426 
427  /**
428  * @brief Get the shared face quadrilaterals GroupCommunicator.
429  *
430  * @param[out] squad_comm
431  */
433  {
434  GetSharedQuadCommunicator(1, squad_comm);
435  }
436 
437  /**
438  * @brief Get the shared face triangles GroupCommunicator.
439  *
440  * @param[out] stria_comm
441  */
443  {
444  GetSharedTriCommunicator(1, stria_comm);
445  }
446 
447  void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[],
448  Array<HYPRE_BigInt> *offsets[]) const;
449 
450  void ExchangeFaceNbrData();
451  void ExchangeFaceNbrNodes();
452 
453  void SetCurvature(int order, bool discont = false, int space_dim = -1,
454  int ordering = 1) override;
455 
456  /** Replace the internal node GridFunction with a new GridFunction defined
457  on the given FiniteElementSpace. The new node coordinates are projected
458  (derived) from the current nodes/vertices. */
459  void SetNodalFESpace(FiniteElementSpace *nfes) override;
461 
462  int GetNFaceNeighbors() const { return face_nbr_group.Size(); }
463  int GetNFaceNeighborElements() const { return face_nbr_elements.Size(); }
464  int GetFaceNbrGroup(int fn) const { return face_nbr_group[fn]; }
465  int GetFaceNbrRank(int fn) const;
466 
467  /** Similar to Mesh::GetElementFaces */
468  void GetFaceNbrElementFaces(int i, Array<int> &fcs, Array<int> &cor) const;
469 
470  /** Similar to Mesh::GetFaceToElementTable with added face-neighbor elements
471  with indices offset by the local number of elements. */
473 
474  /// Returns (a pointer to an object containing) the following data:
475  ///
476  /// 1) Elem1No - the index of the first element that contains this face this
477  /// is the element that has the same outward unit normal vector as the
478  /// face;
479  ///
480  /// 2) Elem2No - the index of the second element that contains this face this
481  /// element has outward unit normal vector as the face multiplied with -1;
482  ///
483  /// 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first
484  /// and the second element respectively;
485  ///
486  /// 4) Face - pointer to the ElementTransformation of the face;
487  ///
488  /// 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face
489  /// coordinate system to the element coordinate system (both in their
490  /// reference elements). Used to transform IntegrationPoints from face to
491  /// element. More formally, let:
492  /// TL1, TL2 be the transformations represented by Loc1, Loc2,
493  /// TE1, TE2 - the transformations represented by Elem1, Elem2,
494  /// TF - the transformation represented by Face, then
495  /// TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
496  ///
497  /// 6) FaceGeom - the base geometry for the face.
498  ///
499  /// The mask specifies which fields in the structure to return:
500  /// mask & 1 - Elem1, mask & 2 - Elem2
501  /// mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.
502  /// These mask values are defined in the ConfigMasks enum type as part of the
503  /// FaceElementTransformations class in fem/eltrans.hpp.
504  ///
505  /// Note that the pointer is owned by the class and is shared, i.e., calling
506  /// this function resets pointers obtained from previous calls.
508  int FaceNo,
509  int mask = 31) override;
510 
511  /** Get the FaceElementTransformations for the given shared face (edge 2D)
512  using the shared face index @a sf. @a fill2 specify if the information
513  for elem2 of the face should be computed or not.
514  In the returned object, 1 and 2 refer to the local and the neighbor
515  elements, respectively.
516  Note that the pointer is owned by the class and is shared, i.e., calling
517  this function resets pointers obtained from previous calls. */
519  GetSharedFaceTransformations(int sf, bool fill2 = true);
520 
521  /** Get the FaceElementTransformations for the given shared face (edge 2D)
522  using the face index @a FaceNo. @a fill2 specify if the information
523  for elem2 of the face should be computed or not.
524  In the returned object, 1 and 2 refer to the local and the neighbor
525  elements, respectively.
526  Note that the pointer is owned by the class and is shared, i.e., calling
527  this function resets pointers obtained from previous calls. */
529  GetSharedFaceTransformationsByLocalIndex(int FaceNo, bool fill2 = true);
530 
531  /// Returns a pointer to the transformation defining the i-th face neighbor.
532  /// Note that the pointer is owned by the class and is shared, i.e., calling
533  /// this function resets pointers obtained from previous calls.
535  {
537 
538  return &FaceNbrTransformation;
539  }
540 
541  /// Get the size of the i-th face neighbor element relative to the reference
542  /// element.
543  double GetFaceNbrElementSize(int i, int type=0);
544 
545  /// Return the number of shared faces (3D), edges (2D), vertices (1D)
546  int GetNSharedFaces() const;
547 
548  /// Return the local face index for the given shared face.
549  int GetSharedFace(int sface) const;
550 
551  /** @brief Returns the number of local faces according to the requested type,
552  does not count master non-conforming faces.
553 
554  If type==Boundary returns only the number of true boundary faces
555  contrary to GetNBE() that returns all "boundary" elements which may
556  include actual interior faces.
557  Similarly, if type==Interior, only the true interior faces (including
558  shared faces) are counted excluding all master non-conforming faces. */
559  int GetNFbyType(FaceType type) const override;
560 
561  /// See the remarks for the serial version in mesh.hpp
562  MFEM_DEPRECATED void ReorientTetMesh() override;
563 
564  /// Utility function: sum integers from all processors (Allreduce).
565  long long ReduceInt(int value) const override;
566 
567  /** Load balance the mesh by equipartitioning the global space-filling
568  sequence of elements. Works for nonconforming meshes only. */
569  void Rebalance();
570 
571  /** Load balance a nonconforming mesh using a user-defined partition.
572  Each local element 'i' is migrated to processor rank 'partition[i]',
573  for 0 <= i < GetNE(). */
574  void Rebalance(const Array<int> &partition);
575 
576  /// Save the mesh in a parallel mesh format.
577  void ParPrint(std::ostream &out) const;
578 
579  // Enable Print() to add the parallel interface as boundary (typically used
580  // for visualization purposes)
581  void SetPrintShared(bool print) { print_shared = print; }
582 
583  /** Print the part of the mesh in the calling processor using the mfem v1.0
584  format. Depending on SetPrintShared(), the parallel interface can be
585  added as boundary for visualization (true by default) . */
586  void Print(std::ostream &out = mfem::out) const override;
587 
588  /// Save the ParMesh to files (one for each MPI rank). The files will be
589  /// given suffixes according to the MPI rank. The mesh will be written to the
590  /// files using ParMesh::Print. The given @a precision will be used for ASCII
591  /// output.
592  void Save(const char *fname, int precision=16) const override;
593 
594 #ifdef MFEM_USE_ADIOS2
595  /** Print the part of the mesh in the calling processor using adios2 bp
596  format. */
597  void Print(adios2stream &out) const override;
598 #endif
599 
600  /** Print the part of the mesh in the calling processor adding the interface
601  as boundary (for visualization purposes) using Netgen/Truegrid format .*/
602  void PrintXG(std::ostream &out = mfem::out) const override;
603 
604  /** Write the mesh to the stream 'out' on Process 0 in a form suitable for
605  visualization: the mesh is written as a disjoint mesh and the shared
606  boundary is added to the actual boundary; both the element and boundary
607  attributes are set to the processor number. */
608  void PrintAsOne(std::ostream &out = mfem::out) const;
609 
610  /** Write the mesh to the stream 'out' on Process 0 as a serial mesh. The
611  output mesh does not have any duplication of vertices/nodes at
612  processor boundaries. */
613  void PrintAsSerial(std::ostream &out = mfem::out) const;
614 
615  /** Returns a Serial mesh on MPI rank @a save_rank that does not have any
616  duplication of vertices/nodes at processor boundaries. */
617  Mesh GetSerialMesh(int save_rank) const;
618 
619  /// Save the mesh as a single file (using ParMesh::PrintAsOne). The given
620  /// @a precision is used for ASCII output.
621  void SaveAsOne(const char *fname, int precision=16) const;
622 
623  /// Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'
624  void PrintAsOneXG(std::ostream &out = mfem::out);
625 
626  /** Print the mesh in parallel PVTU format. The PVTU and VTU files will be
627  stored in the directory specified by @a pathname. If the directory does
628  not exist, it will be created. */
629  void PrintVTU(std::string pathname,
631  bool high_order_output=false,
632  int compression_level=0,
633  bool bdr=false) override;
634 
635  /// Parallel version of Mesh::Load().
636  void Load(std::istream &input, int generate_edges = 0,
637  int refine = 1, bool fix_orientation = true) override;
638 
639  /// Returns the minimum and maximum corners of the mesh bounding box. For
640  /// high-order meshes, the geometry is refined first "ref" times.
641  void GetBoundingBox(Vector &p_min, Vector &p_max, int ref = 2);
642 
643  void GetCharacteristics(double &h_min, double &h_max,
644  double &kappa_min, double &kappa_max);
645 
646  /// Swaps internal data with another ParMesh, including non-geometry members.
647  /// See @a Mesh::Swap
648  void Swap(ParMesh &other);
649 
650  /// Print various parallel mesh stats
651  void PrintInfo(std::ostream &out = mfem::out) override;
652 
653  int FindPoints(DenseMatrix& point_mat, Array<int>& elem_ids,
654  Array<IntegrationPoint>& ips, bool warn = true,
655  InverseElementTransformation *inv_trans = NULL) override;
656 
657  /// Debugging method
658  void PrintSharedEntities(const char *fname_prefix) const;
659 
660  virtual ~ParMesh();
661 
662  friend class ParNCMesh;
663 #ifdef MFEM_USE_PUMI
664  friend class ParPumiMesh;
665 #endif
666 #ifdef MFEM_USE_ADIOS2
667  friend class adios2stream;
668 #endif
669 };
670 
671 }
672 
673 #endif // MFEM_USE_MPI
674 
675 #endif
void UniformRefinement3D() override
Refine a mixed 3D mesh uniformly.
Definition: pmesh.cpp:4587
void PrintAsOneXG(std::ostream &out=mfem::out)
Old mesh format (Netgen/Truegrid) version of &#39;PrintAsOne&#39;.
Definition: pmesh.cpp:5629
void GetGhostFaceTransformation(FaceElementTransformations *FETr, Element::Type face_type, Geometry::Type face_geom)
Definition: pmesh.cpp:3032
void GetSharedVertexCommunicator(int ordering, GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
Definition: pmesh.cpp:1661
virtual ~ParMesh()
Definition: pmesh.cpp:6734
void LoadSharedEntities(std::istream &input)
Definition: pmesh.cpp:978
long glob_offset_sequence
Definition: pmesh.hpp:87
void UniformRefineGroups2D(int old_nv)
Definition: pmesh.cpp:4363
int * partitioning_cache
Definition: pmesh.hpp:390
FaceElementTransformations * GetSharedFaceTransformationsByLocalIndex(int FaceNo, bool fill2=true)
Definition: pmesh.cpp:3088
int NRanks
Definition: pmesh.hpp:38
ParMesh & operator=(ParMesh &&mesh)
Move assignment operator.
Definition: pmesh.cpp:101
void GetSharedQuadCommunicator(int ordering, GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
Definition: pmesh.cpp:1685
void MakeRefined_(ParMesh &orig_mesh, int ref_factor, int ref_type)
Internal function used in ParMesh::MakeRefined (and related constructor)
Definition: pmesh.cpp:1132
int FindSharedVertices(const int *partition, Table *vertex_element, ListOfIntegerSets &groups)
Definition: pmesh.cpp:585
void PrintSharedEntities(const char *fname_prefix) const
Debugging method.
Definition: pmesh.cpp:6508
ElementTransformation * GetFaceNbrElementTransformation(int i)
Definition: pmesh.hpp:534
void PrintInfo(std::ostream &out=mfem::out) override
Print various parallel mesh stats.
Definition: pmesh.cpp:6190
void Set(const int *w)
Definition: pmesh.hpp:57
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:1367
bool print_shared
Definition: pmesh.hpp:92
void GetSharedTriCommunicator(int ordering, GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
Definition: pmesh.cpp:1709
int BuildLocalVertices(const Mesh &global_mesh, const int *partitioning, Array< int > &vert_global_local)
Fills out partitioned Mesh::vertices.
Definition: pmesh.cpp:302
void UniformRefinement2D() override
Refine a mixed 2D mesh uniformly.
Definition: pmesh.cpp:4563
Array< Element * > face_nbr_elements
Definition: pmesh.hpp:382
Array< int > sface_lface
Definition: pmesh.hpp:78
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:720
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int NGroups() const
Return the number of groups.
void GetCharacteristics(double &h_min, double &h_max, double &kappa_min, double &kappa_max)
Definition: pmesh.cpp:6177
void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0) override
This function is not public anymore. Use GeneralRefinement instead.
Definition: pmesh.cpp:3948
void EnsureParNodes()
If the mesh is curved, make sure &#39;Nodes&#39; is ParGridFunction.
Definition: pmesh.cpp:2095
void UniformRefineGroups3D(int old_nv, int old_nedges, const DSTable &old_v_to_v, const STable3D &old_faces, Array< int > *f2qf)
Definition: pmesh.cpp:4414
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:6459
Abstract parallel finite element space.
Definition: pfespace.hpp:28
bool have_face_nbr_data
Definition: pmesh.hpp:378
Array< Vert3 > shared_trias
Definition: pmesh.hpp:65
Array< int > face_nbr_vertices_offset
Definition: pmesh.hpp:381
Array< int > face_nbr_group
Definition: pmesh.hpp:379
bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1) override
NC version of GeneralDerefinement.
Definition: pmesh.cpp:4000
Data arrays will be written in ASCII format.
The inverse transformation of a given ElementTransformation.
Definition: eltrans.hpp:185
void RebalanceImpl(const Array< int > *partition)
Definition: pmesh.cpp:4068
void GroupEdge(int group, int i, int &edge, int &o) const
Definition: pmesh.cpp:1607
void MarkTetMeshForRefinement(DSTable &v_to_v) override
Definition: pmesh.cpp:1733
void GetSharedEdgeCommunicator(GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
Definition: pmesh.hpp:412
void BuildFaceGroup(int ngroups, const Mesh &mesh, const Array< int > &face_group, int &nstria, int &nsquad)
Definition: pmesh.cpp:622
void Set(int v0, int v1, int v2)
Definition: pmesh.hpp:45
int RowSize(int i) const
Definition: table.hpp:108
void Set(const int *w)
Definition: pmesh.hpp:46
STable3D * GetSharedFacesTable()
Definition: pmesh.cpp:2661
Array< int > sedge_ledge
Definition: pmesh.hpp:76
bool HasBoundaryElements() const override
Checks if any rank in the mesh has boundary elements.
Definition: pmesh.cpp:1598
Table group_stria
Definition: pmesh.hpp:71
void GetSharedTriCommunicator(GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
Definition: pmesh.hpp:442
Subdomain representation of a topological parent in another ParMesh.
Definition: psubmesh.hpp:51
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
ParNCMesh * pncmesh
Definition: pmesh.hpp:388
Class for PUMI parallel meshes.
Definition: pumi.hpp:69
int FindSharedEdges(const Mesh &mesh, const int *partition, Table *&edge_element, ListOfIntegerSets &groups)
Definition: pmesh.cpp:531
void ExchangeFaceNbrData()
Definition: pmesh.cpp:2117
void BuildSharedVertMapping(int nvert, const Table *vert_element, const Array< int > &vert_global_local)
Definition: pmesh.cpp:834
Vert4(int v0, int v1, int v2, int v3)
Definition: pmesh.hpp:53
void ReduceMeshGen()
Definition: pmesh.cpp:881
void NURBSUniformRefinement() override
Refine NURBS mesh.
Definition: pmesh.cpp:4616
void Rebalance()
Definition: pmesh.cpp:4058
void GetGlobalEdgeIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Definition: pmesh.cpp:6611
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:383
FaceType
Definition: mesh.hpp:45
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
void ParPrint(std::ostream &out) const
Save the mesh in a parallel mesh format.
Definition: pmesh.cpp:6309
Array< int > face_nbr_elements_offset
Definition: pmesh.hpp:380
MPI_Comm MyComm
Definition: pmesh.hpp:37
void Finalize(bool refine=false, bool fix_orientation=false) override
Finalize the construction of a general Mesh.
Definition: pmesh.cpp:1516
Symmetric 3D Table stored as an array of rows each of which has a stack of column, floor, number nodes. The number of the node is assigned by counting the nodes from zero as they are pushed into the table. Diagonals of any kind are not allowed so the row, column and floor must all be different for each node. Only one node is stored for all 6 symmetric entries that are indexable by unique triplets of row, column, and floor.
Definition: stable3d.hpp:34
void GetBoundingBox(Vector &p_min, Vector &p_max, int ref=2)
Definition: pmesh.cpp:6159
Array< Element * > shared_edges
Definition: pmesh.hpp:61
void Destroy()
Definition: pmesh.cpp:6717
void BuildSharedEdgeElems(int nedges, Mesh &mesh, const Array< int > &vert_global_local, const Table *edge_element)
Definition: pmesh.cpp:797
int GroupNVertices(int group) const
Definition: pmesh.hpp:395
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31) override
Definition: pmesh.cpp:3064
int GroupVertex(int group, int i) const
Definition: pmesh.hpp:400
void GetSharedQuadCommunicator(GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
Definition: pmesh.hpp:432
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:3080
Table send_face_nbr_vertices
Definition: pmesh.hpp:386
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
int GetNGroups() const
Definition: pmesh.hpp:392
VTKFormat
Data array format for VTK and VTU files.
Definition: vtk.hpp:96
void GetGlobalVertexIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Definition: pmesh.cpp:6595
void GetSharedVertexCommunicator(GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
Definition: pmesh.hpp:422
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:6392
int GetNFaceNeighbors() const
Definition: pmesh.hpp:462
void Set(int v0, int v1, int v2, int v3)
Definition: pmesh.hpp:55
void GetSharedEdgeCommunicator(int ordering, GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
Definition: pmesh.cpp:1637
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:3237
Table * face_nbr_el_to_face
Definition: pmesh.hpp:80
long long ReduceInt(int value) const override
Utility function: sum integers from all processors (Allreduce).
Definition: pmesh.cpp:6302
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:187
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
Definition: pmesh.cpp:3215
STable3D * GetFaceNbrElementToFaceTable(int ret_ftbl=0)
Definition: pmesh.cpp:2724
void GetGlobalElementIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are supported.
Definition: pmesh.cpp:6662
MPI_Comm GetComm() const
Definition: pmesh.hpp:351
Array< Vert4 > shared_quads
Definition: pmesh.hpp:66
int GetNFaceNeighborElements() const
Definition: pmesh.hpp:463
int GroupNTriangles(int group) const
Definition: pmesh.hpp:397
void PrintXG(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4624
Table face_nbr_el_ori
Definition: pmesh.hpp:81
int GetMyRank() const
Definition: pmesh.hpp:353
void FindSharedFaces(const Mesh &mesh, const int *partition, Array< int > &face_group, ListOfIntegerSets &groups)
Definition: pmesh.cpp:504
void LocalRefinement(const Array< int > &marked_el, int type=3) override
This function is not public anymore. Use GeneralRefinement instead.
Definition: pmesh.cpp:3444
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:3196
void RefineGroups(const DSTable &v_to_v, int *middle)
Update the groups after triangle refinement.
Definition: pmesh.cpp:4106
Array< Vertex > face_nbr_vertices
Definition: pmesh.hpp:383
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
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
void Swap(ParMesh &other)
Definition: pmesh.cpp:6676
void SaveAsOne(const char *fname, int precision=16) const
Definition: pmesh.cpp:5618
HYPRE_Int HYPRE_BigInt
long long glob_elem_offset
Definition: pmesh.hpp:86
Table * GetFaceToAllElementTable() const
Definition: pmesh.cpp:2974
void PrintAsOne(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4983
bool WantSkipSharedMaster(const NCMesh::Master &master) const
Definition: pmesh.cpp:4826
ParMesh()
Default constructor. Create an empty ParMesh.
Definition: pmesh.hpp:290
void SetAttributes() override
Definition: pmesh.cpp:1580
void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[], Array< HYPRE_BigInt > *offsets[]) const
Definition: pmesh.cpp:1933
void BuildVertexGroup(int ngroups, const Table &vert_element)
Definition: pmesh.cpp:694
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: pmesh.cpp:1976
int GetLocalElementNum(long long global_element_num) const
Definition: pmesh.cpp:1530
void GetFaceNbrElementFaces(int i, Array< int > &fcs, Array< int > &cor) const
Definition: pmesh.cpp:2944
void DistributeAttributes(Array< int > &attr)
Ensure that bdr_attributes and attributes agree across processors.
Definition: pmesh.cpp:1544
int GetNRanks() const
Definition: pmesh.hpp:352
int GetFaceNbrGroup(int fn) const
Definition: pmesh.hpp:464
void Swap(Mesh &other, bool non_geometry)
Definition: mesh.cpp:9346
void BuildEdgeGroup(int ngroups, const Table &edge_element)
Definition: pmesh.cpp:668
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:1862
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
Table group_sedge
Definition: pmesh.hpp:70
void ExchangeFaceNbrNodes()
Definition: pmesh.cpp:2601
int BuildLocalElements(const Mesh &global_mesh, const int *partitioning, const Array< int > &vert_global_local)
Fills out partitioned Mesh::elements.
Definition: pmesh.cpp:350
Table group_svert
Shared objects in each group.
Definition: pmesh.hpp:69
void PrintAsSerial(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:5292
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:2056
bool DecodeFaceSplittings(HashTable< Hashed2 > &v_to_v, const int *v, const Array< unsigned > &codes, int &pos)
Definition: pmesh.cpp:1897
double GetFaceNbrElementSize(int i, int type=0)
Definition: pmesh.cpp:2030
long long GetGlobalElementNum(int local_element_num) const
Map a local element number to a global element number.
Definition: pmesh.cpp:1538
void SetPrintShared(bool print)
Definition: pmesh.hpp:581
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
void GetGlobalFaceIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Definition: pmesh.cpp:6634
IsoparametricTransformation FaceNbrTransformation
Definition: pmesh.hpp:83
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:1847
Vector data type.
Definition: vector.hpp:60
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:2926
int GroupNEdges(int group) const
Definition: pmesh.hpp:396
static ParMesh MakeSimplicial(ParMesh &orig_mesh)
Definition: pmesh.cpp:1374
void SetNodalFESpace(FiniteElementSpace *nfes) override
Definition: pmesh.cpp:2076
Array< int > svert_lvert
Shared to local index mapping.
Definition: pmesh.hpp:75
int GroupNQuadrilaterals(int group) const
Definition: pmesh.hpp:398
int MyRank
Definition: pmesh.hpp:38
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4839
void GroupQuadrilateral(int group, int i, int &face, int &o) const
Definition: pmesh.cpp:1626
Table send_face_nbr_elements
Definition: pmesh.hpp:385
List of integer sets.
Definition: sets.hpp:62
void FinalizeParTopo()
Definition: pmesh.cpp:887
void DeleteFaceNbrData()
Definition: pmesh.cpp:2035
GroupTopology gtopo
Definition: pmesh.hpp:375
Mesh GetSerialMesh(int save_rank) const
Definition: pmesh.cpp:5303
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:942
void GroupTriangle(int group, int i, int &face, int &o) const
Definition: pmesh.cpp:1615
Class for parallel meshes.
Definition: pmesh.hpp:32
Abstract data type element.
Definition: element.hpp:28
void ComputeGlobalElementOffset() const
Definition: pmesh.cpp:868
Table group_squad
Definition: pmesh.hpp:72
void Save(const char *fname, int precision=16) const override
Definition: pmesh.cpp:4955
Vert3(int v0, int v1, int v2)
Definition: pmesh.hpp:44
MFEM_DEPRECATED void ReorientTetMesh() override
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:3269