MFEM  v4.5.2
Finite element discretization library
mesh.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_MESH
13 #define MFEM_MESH
14 
15 #include "../config/config.hpp"
16 #include "../general/stable3d.hpp"
17 #include "../general/globals.hpp"
18 #include "triangle.hpp"
19 #include "tetrahedron.hpp"
20 #include "vertex.hpp"
21 #include "vtk.hpp"
22 #include "ncmesh.hpp"
23 #include "../fem/eltrans.hpp"
24 #include "../fem/coefficient.hpp"
25 #include "../general/zstr.hpp"
26 #ifdef MFEM_USE_ADIOS2
27 #include "../general/adios2stream.hpp"
28 #endif
29 #include <iostream>
30 
31 namespace mfem
32 {
33 
34 // Data type mesh
35 
36 class GeometricFactors;
37 class FaceGeometricFactors;
38 class KnotVector;
39 class NURBSExtension;
40 class FiniteElementSpace;
41 class GridFunction;
42 struct Refinement;
43 
44 /** An enum type to specify if interior or boundary faces are desired. */
45 enum class FaceType : bool {Interior, Boundary};
46 
47 #ifdef MFEM_USE_MPI
48 class ParMesh;
49 class ParNCMesh;
50 #endif
51 
52 class Mesh
53 {
54 #ifdef MFEM_USE_MPI
55  friend class ParMesh;
56  friend class ParNCMesh;
57 #endif
58  friend class NCMesh;
59  friend class NURBSExtension;
60 
61 #ifdef MFEM_USE_ADIOS2
62  friend class adios2stream;
63 #endif
64 
65 protected:
66  int Dim;
67  int spaceDim;
68 
71  /** These variables store the number of Interior and Boundary faces. Calling
72  fes->GetMesh()->GetNBE() doesn't return the expected value in 3D because
73  periodic meshes in 3D have some of their faces marked as boundary for
74  visualization purpose in GLVis. */
76 
77  int meshgen; // see MeshGenerator()
78  int mesh_geoms; // sum of (1 << geom) for all geom of all dimensions
79 
80  // Counter for Mesh transformations: refinement, derefinement, rebalancing.
81  // Used for checking during Update operations on objects depending on the
82  // Mesh, such as FiniteElementSpace, GridFunction, etc.
83  long sequence;
84 
86  // Vertices are only at the corners of elements, where you would expect them
87  // in the lowest-order mesh. In some cases, e.g. in a Mesh that defines the
88  // patch topology for a NURBS mesh (see LoadPatchTopo()) the vertices may be
89  // empty while NumOfVertices is positive.
93 
94  /** @brief This structure stores the low level information necessary to
95  interpret the configuration of elements on a specific face. This
96  information can be accessed using methods like GetFaceElements(),
97  GetFaceInfos(), FaceIsInterior(), etc.
98 
99  For accessing higher level deciphered information look at
100  Mesh::FaceInformation, and its accessor Mesh::GetFaceInformation().
101 
102  Each face contains information on the indices, local reference faces,
103  orientations, and potential nonconformity for the two neighboring
104  elements on a face.
105  Each face can either be an interior, boundary, or shared interior face.
106  Each interior face is shared by two elements referred as Elem1 and Elem2.
107  For boundary faces only the information on Elem1 is relevant.
108  Shared interior faces correspond to faces where Elem1 and Elem2 are
109  distributed on different MPI ranks.
110  Regarding conformity, three cases are distinguished, conforming faces,
111  nonconforming slave faces, and nonconforming master faces. Master and
112  slave referring to the coarse and fine elements respectively on a
113  nonconforming face.
114  Nonconforming slave faces always have the slave element as Elem1 and
115  the master element as Elem2. On the other side, nonconforming master
116  faces always have the master element as Elem1, and one of the slave
117  element as Elem2. Except for ghost nonconforming slave faces, where
118  Elem1 is the master side and Elem2 is the slave side.
119 
120  The indices of Elem1 and Elem2 can be indirectly extracted from
121  FaceInfo::Elem1No and FaceInfo::Elem2No, read the note below for special
122  cases on the index of Elem2.
123 
124  The local face identifiers are deciphered from FaceInfo::Elem1Inf and
125  FaceInfo::Elem2Inf through the formula: LocalFaceIndex = ElemInf/64,
126  the semantic of the computed local face identifier can be found in
127  fem/geom.cpp. The local face identifier corresponds to an index
128  in the Constants<Geometry>::Edges arrays for 2D element geometries, and
129  to an index in the Constants<Geometry>::FaceVert arrays for 3D element
130  geometries.
131 
132  The orientation of each element relative to a face is obtained through
133  the formula: Orientation = ElemInf%64, the semantic of the orientation
134  can also be found in fem/geom.cpp. The orientation corresponds to
135  an index in the Constants<Geometry>::Orient arrays, providing the
136  sequence of vertices identifying the orientation of an edge/face. By
137  convention the orientation of Elem1 is always set to 0, serving as the
138  reference orientation. The orientation of Elem2 relatively to Elem1 is
139  therefore determined just by using the orientation of Elem2. An important
140  special case is the one of nonconforming faces, the orientation should
141  be composed with the PointMatrix, which also contains orientation
142  information. A special treatment should be done for 2D, the orientation
143  in the PointMatrix is not included, therefore when applying the
144  PointMatrix transformation, the PointMatrix should be flipped, except for
145  shared nonconforming slave faces where the transformation can be applied
146  as is.
147 
148  Another special case is the case of shared nonconforming faces. Ghost
149  faces use a different design based on so called "ghost" faces.
150  Ghost faces, as their name suggest are very well hidden, and they
151  usually have a separate interface from "standard" faces.
152  */
153  struct FaceInfo
154  {
155  // Inf = 64 * LocalFaceIndex + FaceOrientation
157  int NCFace; /* -1 if this is a regular conforming/boundary face;
158  index into 'nc_faces_info' if >= 0. */
159  };
160  // NOTE: in NC meshes, master faces have Elem2No == -1. Slave faces on the
161  // other hand have Elem2No and Elem2Inf set to the master face's element and
162  // its local face number.
163  //
164  // A local face is one generated from a local element and has index i in
165  // faces_info such that i < GetNumFaces(). Also, Elem1No always refers to the
166  // element (slave or master, in the nonconforming case) that generated the
167  // face.
168  // Classification of a local (non-ghost) face based on its FaceInfo:
169  // - Elem2No >= 0 --> local interior face; can be either:
170  // - NCFace == -1 --> conforming face, or
171  // - NCFace >= 0 --> nonconforming slave face; Elem2No is the index of
172  // the master volume element; Elem2Inf%64 is 0, see the note in
173  // Mesh::GenerateNCFaceInfo().
174  // - Elem2No < 0 --> local "boundary" face; can be one of:
175  // - NCFace == -1 --> conforming face; can be either:
176  // - Elem2Inf < 0 --> true boundary face (no element on side 2)
177  // - Elem2Inf >= 0 --> shared face where element 2 is a face-neighbor
178  // element with index -1-Elem2No. This state is initialized by
179  // ParMesh::ExchangeFaceNbrData().
180  // - NCFace >= 0 --> nonconforming face; can be one of:
181  // - Elem2Inf < 0 --> master nonconforming face, interior or shared;
182  // In this case, Elem2No is -1; see GenerateNCFaceInfo().
183  // - Elem2Inf >= 0 --> shared slave nonconforming face where element 2
184  // is the master face-neighbor element with index -1-Elem2No; see
185  // ParNCMesh::GetFaceNeighbors().
186  //
187  // A ghost face is a nonconforming face that is generated by a non-local,
188  // i.e. ghost, element. A ghost face has index i in faces_info such that
189  // i >= GetNumFaces().
190  // Classification of a ghost (non-local) face based on its FaceInfo:
191  // - Elem1No == -1 --> master ghost face? These ghost faces also have:
192  // Elem2No == -1, Elem1Inf == Elem2Inf == -1, and NCFace == -1.
193  // - Elem1No >= 0 --> slave ghost face; Elem1No is the index of the local
194  // master side element, i.e. side 1 IS NOT the side that generated the
195  // face. Elem2No is < 0 and -1-Elem2No is the index of the ghost
196  // face-neighbor element that generated this slave ghost face. In this
197  // case, Elem2Inf >= 0 and NCFace >= 0.
198  // Relevant methods: GenerateFaces(), GenerateNCFaceInfo(),
199  // ParNCMesh::GetFaceNeighbors(),
200  // ParMesh::ExchangeFaceNbrData()
201 
202  struct NCFaceInfo
203  {
204  bool Slave; // true if this is a slave face, false if master face
205  int MasterFace; // if Slave, this is the index of the master face
206  // If not Slave, 'MasterFace' is the local face index of this master face
207  // as a face in the unique adjacent element.
208  const DenseMatrix* PointMatrix; // if Slave, position within master face
209  // (NOTE: PointMatrix points to a matrix owned by NCMesh.)
210 
211  NCFaceInfo() = default;
212 
213  NCFaceInfo(bool slave, int master, const DenseMatrix* pm)
214  : Slave(slave), MasterFace(master), PointMatrix(pm) {}
215  };
216 
219 
224  Table *bel_to_edge; // for 3D
226 
227  // Note that the following tables are owned by this class and should not be
228  // deleted by the caller. Of these three tables, only face_edge and
229  // edge_vertex are returned by access functions.
230  mutable Table *face_to_elem; // Used by FindFaceNeighbors, not returned.
231  mutable Table *face_edge; // Returned by GetFaceEdgeTable().
232  mutable Table *edge_vertex; // Returned by GetEdgeVertexTable().
233 
238 
239  // refinement embeddings for forward compatibility with NCMesh
241 
242  // Nodes are only active for higher order meshes, and share locations with
243  // the vertices, plus all the higher- order control points within the
244  // element and along the edges and on the faces.
247 
248  static const int vtk_quadratic_tet[10];
249  static const int vtk_quadratic_pyramid[13];
250  static const int vtk_quadratic_wedge[18];
251  static const int vtk_quadratic_hex[27];
252 
253 #ifdef MFEM_USE_MEMALLOC
254  friend class Tetrahedron;
256 #endif
257 
258  // used during NC mesh initialization only
260 
261 public:
269 
271 
272  /// A list of all unique element attributes used by the Mesh.
274  /// A list of all unique boundary attributes used by the Mesh.
276 
277  NURBSExtension *NURBSext; ///< Optional NURBS mesh extension.
278  NCMesh *ncmesh; ///< Optional nonconforming mesh extension.
279  Array<GeometricFactors*> geom_factors; ///< Optional geometric factors.
280  Array<FaceGeometricFactors*> face_geom_factors; /**< Optional face geometric
281  factors. */
282 
283  // Global parameter that can be used to control the removal of unused
284  // vertices performed when reading a mesh in MFEM format. The default value
285  // (true) is set in mesh_readers.cpp.
287 
288 protected:
290 
291  void Init();
292  void InitTables();
293  void SetEmpty(); // Init all data members with empty values
294  void DestroyTables();
296  void DestroyPointers(); // Delete data specifically allocated by class Mesh.
297  void Destroy(); // Delete all owned data.
298  void ResetLazyData();
299 
300  Element *ReadElementWithoutAttr(std::istream &);
301  static void PrintElementWithoutAttr(const Element *, std::ostream &);
302 
303  Element *ReadElement(std::istream &);
304  static void PrintElement(const Element *, std::ostream &);
305 
306  // Readers for different mesh formats, used in the Load() method.
307  // The implementations of these methods are in mesh_readers.cpp.
308  void ReadMFEMMesh(std::istream &input, int version, int &curved);
309  void ReadLineMesh(std::istream &input);
310  void ReadNetgen2DMesh(std::istream &input, int &curved);
311  void ReadNetgen3DMesh(std::istream &input);
312  void ReadTrueGridMesh(std::istream &input);
313  void CreateVTKMesh(const Vector &points, const Array<int> &cell_data,
314  const Array<int> &cell_offsets,
315  const Array<int> &cell_types,
316  const Array<int> &cell_attributes,
317  int &curved, int &read_gf, bool &finalize_topo);
318  void ReadVTKMesh(std::istream &input, int &curved, int &read_gf,
319  bool &finalize_topo);
320  void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf,
321  bool &finalize_topo, const std::string &xml_prefix="");
322  void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf);
323  void ReadInlineMesh(std::istream &input, bool generate_edges = false);
324  void ReadGmshMesh(std::istream &input, int &curved, int &read_gf);
325  /* Note NetCDF (optional library) is used for reading cubit files */
326 #ifdef MFEM_USE_NETCDF
327  void ReadCubit(const char *filename, int &curved, int &read_gf);
328 #endif
329 
330  /// Determine the mesh generator bitmask #meshgen, see MeshGenerator().
331  /** Also, initializes #mesh_geoms. */
332  void SetMeshGen();
333 
334  /// Return the length of the segment from node i to node j.
335  double GetLength(int i, int j) const;
336 
337  /** Compute the Jacobian of the transformation from the perfect
338  reference element at the center of the element. */
339  void GetElementJacobian(int i, DenseMatrix &J);
340 
341  void MarkForRefinement();
343  void GetEdgeOrdering(DSTable &v_to_v, Array<int> &order);
344  virtual void MarkTetMeshForRefinement(DSTable &v_to_v);
345 
346  // Methods used to prepare and apply permutation of the mesh nodes assuming
347  // that the mesh elements may be rotated (e.g. to mark triangle or tet edges
348  // for refinement) between the two calls - PrepareNodeReorder() and
349  // DoNodeReorder(). The latter method assumes that the 'faces' have not been
350  // updated after the element rotations.
351  void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert);
352  void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert);
353 
355  STable3D *GetElementToFaceTable(int ret_ftbl = 0);
356 
357  /** Red refinement. Element with index i is refined. The default
358  red refinement for now is Uniform. */
359  void RedRefinement(int i, const DSTable &v_to_v,
360  int *edge1, int *edge2, int *middle)
361  { UniformRefinement(i, v_to_v, edge1, edge2, middle); }
362 
363  /** Green refinement. Element with index i is refined. The default
364  refinement for now is Bisection. */
365  void GreenRefinement(int i, const DSTable &v_to_v,
366  int *edge1, int *edge2, int *middle)
367  { Bisection(i, v_to_v, edge1, edge2, middle); }
368 
369  /// Bisect a triangle: element with index @a i is bisected.
370  void Bisection(int i, const DSTable &, int *, int *, int *);
371 
372  /// Bisect a tetrahedron: element with index @a i is bisected.
373  void Bisection(int i, HashTable<Hashed2> &);
374 
375  /// Bisect a boundary triangle: boundary element with index @a i is bisected.
376  void BdrBisection(int i, const HashTable<Hashed2> &);
377 
378  /** Uniform Refinement. Element with index i is refined uniformly. */
379  void UniformRefinement(int i, const DSTable &, int *, int *, int *);
380 
381  /** @brief Averages the vertices with given @a indexes and saves the result
382  in #vertices[result]. */
383  void AverageVertices(const int *indexes, int n, int result);
384 
386  int FindCoarseElement(int i);
387 
388  /** @brief Update the nodes of a curved mesh after the topological part of a
389  Mesh::Operation, such as refinement, has been performed. */
390  /** If Nodes GridFunction is defined, i.e. not NULL, this method calls
391  NodesUpdated().
392 
393  @note Unlike the similarly named public method NodesUpdated() this
394  method modifies the mesh nodes (if they exist) and calls NodesUpdated().
395  */
396  void UpdateNodes();
397 
398  /// Helper to set vertex coordinates given a high-order curvature function.
399  void SetVerticesFromNodes(const GridFunction *nodes);
400 
401  void UniformRefinement2D_base(bool update_nodes = true);
402 
403  /// Refine a mixed 2D mesh uniformly.
405 
406  /* If @a f2qf is not NULL, adds all quadrilateral faces to @a f2qf which
407  represents a "face-to-quad-face" index map. When all faces are quads, the
408  array @a f2qf is kept empty since it is not needed. */
409  void UniformRefinement3D_base(Array<int> *f2qf = NULL,
410  DSTable *v_to_v_p = NULL,
411  bool update_nodes = true);
412 
413  /// Refine a mixed 3D mesh uniformly.
415 
416  /// Refine NURBS mesh.
417  virtual void NURBSUniformRefinement();
418 
419  /// This function is not public anymore. Use GeneralRefinement instead.
420  virtual void LocalRefinement(const Array<int> &marked_el, int type = 3);
421 
422  /// This function is not public anymore. Use GeneralRefinement instead.
423  virtual void NonconformingRefinement(const Array<Refinement> &refinements,
424  int nc_limit = 0);
425 
426  /// NC version of GeneralDerefinement.
427  virtual bool NonconformingDerefinement(Array<double> &elem_error,
428  double threshold, int nc_limit = 0,
429  int op = 1);
430  /// Derefinement helper.
431  double AggregateError(const Array<double> &elem_error,
432  const int *fine, int nfine, int op);
433 
434  /// Read NURBS patch/macro-element mesh
435  void LoadPatchTopo(std::istream &input, Array<int> &edge_to_knot);
436 
437  void UpdateNURBS();
438 
439  void PrintTopo(std::ostream &out, const Array<int> &e_to_k) const;
440 
441  /// Used in GetFaceElementTransformations (...)
444  int i);
446  int i);
447  /// Used in GetFaceElementTransformations (...)
449  int i);
450  /// Used in GetFaceElementTransformations (...)
452  int i);
453  /// Used in GetFaceElementTransformations (...)
455  int i);
456  /// Used in GetFaceElementTransformations (...)
458  int i);
459  /// Used in GetFaceElementTransformations (...)
461  int i);
462  /// Used in GetFaceElementTransformations (...)
464  int i);
465 
466  /** Used in GetFaceElementTransformations to account for the fact that a
467  slave face occupies only a portion of its master face. */
469  const FaceInfo &fi, bool is_ghost);
470 
471  bool IsSlaveFace(const FaceInfo &fi) const;
472 
473  /// Returns the orientation of "test" relative to "base"
474  static int GetTriOrientation (const int * base, const int * test);
475  /// Returns the orientation of "test" relative to "base"
476  static int GetQuadOrientation (const int * base, const int * test);
477  /// Returns the orientation of "test" relative to "base"
478  static int GetTetOrientation (const int * base, const int * test);
479 
480  static void GetElementArrayEdgeTable(const Array<Element*> &elem_array,
481  const DSTable &v_to_v,
482  Table &el_to_edge);
483 
484  /** Return element to edge table and the indices for the boundary edges.
485  The entries in the table are ordered according to the order of the
486  nodes in the elements. For example, if T is the element to edge table
487  T(i, 0) gives the index of edge in element i that connects vertex 0
488  to vertex 1, etc. Returns the number of the edges. */
490 
491  /// Used in GenerateFaces()
492  void AddPointFaceElement(int lf, int gf, int el);
493 
494  void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1);
495 
496  void AddTriangleFaceElement (int lf, int gf, int el,
497  int v0, int v1, int v2);
498 
499  void AddQuadFaceElement (int lf, int gf, int el,
500  int v0, int v1, int v2, int v3);
501  /** For a serial Mesh, return true if the face is interior. For a parallel
502  ParMesh return true if the face is interior or shared. In parallel, this
503  method only works if the face neighbor data is exchanged. */
504  bool FaceIsTrueInterior(int FaceNo) const
505  {
506  return FaceIsInterior(FaceNo) || (faces_info[FaceNo].Elem2Inf >= 0);
507  }
508 
509  void FreeElement(Element *E);
510 
511  void GenerateFaces();
512  void GenerateNCFaceInfo();
513 
514  /// Begin construction of a mesh
515  void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem);
516 
517  // Used in the methods FinalizeXXXMesh() and FinalizeTopology()
518  void FinalizeCheck();
519 
520  void Loader(std::istream &input, int generate_edges = 0,
521  std::string parse_tag = "");
522 
523  // If NURBS mesh, write NURBS format. If NCMesh, write mfem v1.1 format.
524  // If section_delimiter is empty, write mfem v1.0 format. Otherwise, write
525  // mfem v1.2 format with the given section_delimiter at the end.
526  void Printer(std::ostream &out = mfem::out,
527  std::string section_delimiter = "") const;
528 
529  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
530  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
531  type=TETRAHEDRON. The parameter @a sfc_ordering controls how the elements
532  (when type=HEXAHEDRON) are ordered: true - use space-filling curve
533  ordering, or false - use lexicographic ordering. */
534  void Make3D(int nx, int ny, int nz, Element::Type type,
535  double sx, double sy, double sz, bool sfc_ordering);
536 
537  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
538  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
539  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
540  if 1 edges are generated. The parameter @a sfc_ordering controls how the
541  elements (when type=QUADRILATERAL) are ordered: true - use space-filling
542  curve ordering, or false - use lexicographic ordering. */
543  void Make2D(int nx, int ny, Element::Type type, double sx, double sy,
544  bool generate_edges, bool sfc_ordering);
545 
546  /// Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.
547  void Make1D(int n, double sx = 1.0);
548 
549  /// Internal function used in Mesh::MakeRefined
550  void MakeRefined_(Mesh &orig_mesh, const Array<int> ref_factors,
551  int ref_type);
552 
553  /// Initialize vertices/elements/boundary/tables from a nonconforming mesh.
554  void InitFromNCMesh(const NCMesh &ncmesh);
555 
556  /// Create from a nonconforming mesh.
557  explicit Mesh(const NCMesh &ncmesh);
558 
559  // used in GetElementData() and GetBdrElementData()
560  void GetElementData(const Array<Element*> &elem_array, int geom,
561  Array<int> &elem_vtx, Array<int> &attr) const;
562 
563  double GetElementSize(ElementTransformation *T, int type = 0);
564 
565  // Internal helper used in MakeSimplicial (and ParMesh::MakeSimplicial).
566  void MakeSimplicial_(const Mesh &orig_mesh, int *vglobal);
567 
568 public:
569 
570  Mesh() { SetEmpty(); }
571 
572  /** Copy constructor. Performs a deep copy of (almost) all data, so that the
573  source mesh can be modified (e.g. deleted, refined) without affecting the
574  new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the
575  nodes, if present. */
576  explicit Mesh(const Mesh &mesh, bool copy_nodes = true);
577 
578  /// Move constructor, useful for using a Mesh as a function return value.
579  Mesh(Mesh &&mesh);
580 
581  /// Move assignment operstor.
582  Mesh& operator=(Mesh &&mesh);
583 
584  /// Explicitly delete the copy assignment operator.
585  Mesh& operator=(const Mesh &mesh) = delete;
586 
587  /** @name Named mesh constructors.
588 
589  Each of these constructors uses the move constructor, and can be used as
590  the right-hand side of an assignment when creating new meshes. */
591  ///@{
592 
593  /** Creates mesh by reading a file in MFEM, Netgen, or VTK format. If
594  generate_edges = 0 (default) edges are not generated, if 1 edges are
595  generated. */
596  static Mesh LoadFromFile(const char *filename,
597  int generate_edges = 0, int refine = 1,
598  bool fix_orientation = true);
599 
600  /** Creates 1D mesh , divided into n equal intervals. */
601  static Mesh MakeCartesian1D(int n, double sx = 1.0);
602 
603  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
604  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
605  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
606  if 1 edges are generated. If scf_ordering = true (default), elements are
607  ordered along a space-filling curve, instead of row by row. */
608  static Mesh MakeCartesian2D(
609  int nx, int ny, Element::Type type, bool generate_edges = false,
610  double sx = 1.0, double sy = 1.0, bool sfc_ordering = true);
611 
612  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
613  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
614  type=TETRAHEDRON. If sfc_ordering = true (default), elements are ordered
615  along a space-filling curve, instead of row by row and layer by layer. */
616  static Mesh MakeCartesian3D(
617  int nx, int ny, int nz, Element::Type type,
618  double sx = 1.0, double sy = 1.0, double sz = 1.0,
619  bool sfc_ordering = true);
620 
621  /// Create a refined (by any factor) version of @a orig_mesh.
622  /** @param[in] orig_mesh The starting coarse mesh.
623  @param[in] ref_factor The refinement factor, an integer > 1.
624  @param[in] ref_type Specify the positions of the new vertices. The
625  options are BasisType::ClosedUniform or
626  BasisType::GaussLobatto.
627 
628  The refinement data which can be accessed with GetRefinementTransforms()
629  is set to reflect the performed refinements.
630 
631  @note The constructed Mesh is straight-sided. */
632  static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type);
633 
634  /// Create a refined mesh, where each element of the original mesh may be
635  /// refined by a different factor.
636  /** @param[in] orig_mesh The starting coarse mesh.
637  @param[in] ref_factors An array of integers whose size is the number of
638  elements of @a orig_mesh. The @a ith element of
639  @a orig_mesh is refined by refinement factor
640  @a ref_factors[i].
641  @param[in] ref_type Specify the positions of the new vertices. The
642  options are BasisType::ClosedUniform or
643  BasisType::GaussLobatto.
644 
645  The refinement data which can be accessed with GetRefinementTransforms()
646  is set to reflect the performed refinements.
647 
648  @note The constructed Mesh is straight-sided. */
649  /// refined @a ref_factors[i] times in each dimension.
650  static Mesh MakeRefined(Mesh &orig_mesh, const Array<int> &ref_factors,
651  int ref_type);
652 
653  /** Create a mesh by splitting each element of @a orig_mesh into simplices.
654  Quadrilaterals are split into two triangles, prisms are split into
655  3 tetrahedra, and hexahedra are split into either 5 or 6 tetrahedra
656  depending on the configuration.
657  @warning The curvature of the original mesh is not carried over to the
658  new mesh. Periodic meshes are not supported. */
659  static Mesh MakeSimplicial(const Mesh &orig_mesh);
660 
661  /// Create a periodic mesh by identifying vertices of @a orig_mesh.
662  /** Each vertex @a i will be mapped to vertex @a v2v[i], such that all
663  vertices that are coincident under the periodic mapping get mapped to
664  the same index. The mapping @a v2v can be generated from translation
665  vectors using Mesh::CreatePeriodicVertexMapping.
666  @note MFEM requires that each edge of the resulting mesh be uniquely
667  identifiable by a pair of distinct vertices. As a consequence, periodic
668  boundaries must be connected by at least three edges. */
669  static Mesh MakePeriodic(const Mesh &orig_mesh, const std::vector<int> &v2v);
670 
671  ///@}
672 
673  /// @brief Creates a mapping @a v2v from the vertex indices of the mesh such
674  /// that coincident vertices under the given @a translations are identified.
675  /** Each Vector in @a translations should be of size @a sdim (the spatial
676  dimension of the mesh). Two vertices are considered coincident if the
677  translated coordinates of one vertex are within the given tolerance (@a
678  tol, relative to the mesh diameter) of the coordinates of the other
679  vertex.
680  @warning This algorithm does not scale well with the number of boundary
681  vertices in the mesh, and may run slowly on very large meshes. */
682  std::vector<int> CreatePeriodicVertexMapping(
683  const std::vector<Vector> &translations, double tol = 1e-8) const;
684 
685  /// Construct a Mesh from the given primary data.
686  /** The array @a vertices is used as external data, i.e. the Mesh does not
687  copy the data and will not delete the pointer.
688 
689  The data from the other arrays is copied into the internal Mesh data
690  structures.
691 
692  This method calls the method FinalizeTopology(). The method Finalize()
693  may be called after this constructor and after optionally setting the
694  Mesh nodes. */
695  Mesh(double *vertices, int num_vertices,
696  int *element_indices, Geometry::Type element_type,
697  int *element_attributes, int num_elements,
698  int *boundary_indices, Geometry::Type boundary_type,
699  int *boundary_attributes, int num_boundary_elements,
700  int dimension, int space_dimension = -1);
701 
702  /** @anchor mfem_Mesh_init_ctor
703  @brief _Init_ constructor: begin the construction of a Mesh object. */
704  Mesh(int Dim_, int NVert, int NElem, int NBdrElem = 0, int spaceDim_ = -1)
705  {
706  if (spaceDim_ == -1) { spaceDim_ = Dim_; }
707  InitMesh(Dim_, spaceDim_, NVert, NElem, NBdrElem);
708  }
709 
710  /** @name Methods for Mesh construction.
711 
712  These methods are intended to be used with the @ref mfem_Mesh_init_ctor
713  "init constructor". */
714  ///@{
715 
716  Element *NewElement(int geom);
717 
718  int AddVertex(double x, double y = 0.0, double z = 0.0);
719  int AddVertex(const double *coords);
720  int AddVertex(const Vector &coords);
721  /// Mark vertex @a i as nonconforming, with parent vertices @a p1 and @a p2.
722  void AddVertexParents(int i, int p1, int p2);
723 
724  int AddSegment(int v1, int v2, int attr = 1);
725  int AddSegment(const int *vi, int attr = 1);
726 
727  int AddTriangle(int v1, int v2, int v3, int attr = 1);
728  int AddTriangle(const int *vi, int attr = 1);
729  int AddTri(const int *vi, int attr = 1) { return AddTriangle(vi, attr); }
730 
731  int AddQuad(int v1, int v2, int v3, int v4, int attr = 1);
732  int AddQuad(const int *vi, int attr = 1);
733 
734  int AddTet(int v1, int v2, int v3, int v4, int attr = 1);
735  int AddTet(const int *vi, int attr = 1);
736 
737  int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr = 1);
738  int AddWedge(const int *vi, int attr = 1);
739 
740  int AddPyramid(int v1, int v2, int v3, int v4, int v5, int attr = 1);
741  int AddPyramid(const int *vi, int attr = 1);
742 
743  int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8,
744  int attr = 1);
745  int AddHex(const int *vi, int attr = 1);
746  void AddHexAsTets(const int *vi, int attr = 1);
747  void AddHexAsWedges(const int *vi, int attr = 1);
748  void AddHexAsPyramids(const int *vi, int attr = 1);
749 
750  /// The parameter @a elem should be allocated using the NewElement() method
751  int AddElement(Element *elem);
752  int AddBdrElement(Element *elem);
753 
754  int AddBdrSegment(int v1, int v2, int attr = 1);
755  int AddBdrSegment(const int *vi, int attr = 1);
756 
757  int AddBdrTriangle(int v1, int v2, int v3, int attr = 1);
758  int AddBdrTriangle(const int *vi, int attr = 1);
759 
760  int AddBdrQuad(int v1, int v2, int v3, int v4, int attr = 1);
761  int AddBdrQuad(const int *vi, int attr = 1);
762  void AddBdrQuadAsTriangles(const int *vi, int attr = 1);
763 
764  int AddBdrPoint(int v, int attr = 1);
765 
767  /// Finalize the construction of a triangular Mesh.
768  void FinalizeTriMesh(int generate_edges = 0, int refine = 0,
769  bool fix_orientation = true);
770  /// Finalize the construction of a quadrilateral Mesh.
771  void FinalizeQuadMesh(int generate_edges = 0, int refine = 0,
772  bool fix_orientation = true);
773  /// Finalize the construction of a tetrahedral Mesh.
774  void FinalizeTetMesh(int generate_edges = 0, int refine = 0,
775  bool fix_orientation = true);
776  /// Finalize the construction of a wedge Mesh.
777  void FinalizeWedgeMesh(int generate_edges = 0, int refine = 0,
778  bool fix_orientation = true);
779  /// Finalize the construction of a hexahedral Mesh.
780  void FinalizeHexMesh(int generate_edges = 0, int refine = 0,
781  bool fix_orientation = true);
782  /// Finalize the construction of any type of Mesh.
783  /** This method calls FinalizeTopology() and Finalize(). */
784  void FinalizeMesh(int refine = 0, bool fix_orientation = true);
785 
786  ///@}
787 
788  /** @brief Finalize the construction of the secondary topology (connectivity)
789  data of a Mesh. */
790  /** This method does not require any actual coordinate data (either vertex
791  coordinates for linear meshes or node coordinates for meshes with nodes)
792  to be available. However, the data generated by this method is generally
793  required by the FiniteElementSpace class.
794 
795  After calling this method, setting the Mesh vertices or nodes, it may be
796  appropriate to call the method Finalize(). */
797  void FinalizeTopology(bool generate_bdr = true);
798 
799  /// Finalize the construction of a general Mesh.
800  /** This method will:
801  - check and optionally fix the orientation of regular elements
802  - check and fix the orientation of boundary elements
803  - assume that #vertices are defined, if #Nodes == NULL
804  - assume that #Nodes are defined, if #Nodes != NULL.
805  @param[in] refine If true, prepare the Mesh for conforming refinement of
806  triangular or tetrahedral meshes.
807  @param[in] fix_orientation
808  If true, fix the orientation of inverted mesh elements
809  by permuting their vertices.
810 
811  Before calling this method, call FinalizeTopology() and ensure that the
812  Mesh vertices or nodes are set. */
813  virtual void Finalize(bool refine = false, bool fix_orientation = false);
814 
815  virtual void SetAttributes();
816 
817  /** This is our integration with the Gecko library. The method finds an
818  element ordering that will increase memory coherency by putting elements
819  that are in physical proximity closer in memory. It can also be used to
820  obtain a space-filling curve ordering for ParNCMesh partitioning.
821  @param[out] ordering Output element ordering.
822  @param iterations Total number of V cycles. The ordering may improve with
823  more iterations. The best iteration is returned at the end.
824  @param window Initial window size. This determines the number of
825  permutations tested at each multigrid level and strongly influences the
826  quality of the result, but the cost of increasing 'window' is exponential.
827  @param period The window size is incremented every 'period' iterations.
828  @param seed Seed for initial random ordering (0 = skip random reorder).
829  @param verbose Print the progress of the optimization to mfem::out.
830  @param time_limit Optional time limit for the optimization, in seconds.
831  When reached, ordering from the best iteration so far is returned
832  (0 = no limit).
833  @return The final edge product cost of the ordering. The function may be
834  called in an external loop with different seeds, and the best ordering can
835  then be retained. */
836  double GetGeckoElementOrdering(Array<int> &ordering,
837  int iterations = 4, int window = 4,
838  int period = 2, int seed = 0,
839  bool verbose = false, double time_limit = 0);
840 
841  /** Return an ordering of the elements that approximately follows the Hilbert
842  curve. The method performs a spatial (Hilbert) sort on the centers of all
843  elements and returns the resulting sequence, which can then be passed to
844  ReorderElements. This is a cheap alternative to GetGeckoElementOrdering.*/
845  void GetHilbertElementOrdering(Array<int> &ordering);
846 
847  /** Rebuilds the mesh with a different order of elements. For each element i,
848  the array ordering[i] contains its desired new index. Note that the method
849  reorders vertices, edges and faces along with the elements. */
850  void ReorderElements(const Array<int> &ordering, bool reorder_vertices = true);
851 
852  /// Deprecated: see @a MakeCartesian3D.
853  MFEM_DEPRECATED
854  Mesh(int nx, int ny, int nz, Element::Type type, bool generate_edges = false,
855  double sx = 1.0, double sy = 1.0, double sz = 1.0,
856  bool sfc_ordering = true)
857  {
858  Make3D(nx, ny, nz, type, sx, sy, sz, sfc_ordering);
859  Finalize(true); // refine = true
860  }
861 
862  /// Deprecated: see @a MakeCartesian2D.
863  MFEM_DEPRECATED
864  Mesh(int nx, int ny, Element::Type type, bool generate_edges = false,
865  double sx = 1.0, double sy = 1.0, bool sfc_ordering = true)
866  {
867  Make2D(nx, ny, type, sx, sy, generate_edges, sfc_ordering);
868  Finalize(true); // refine = true
869  }
870 
871  /// Deprecated: see @a MakeCartesian1D.
872  MFEM_DEPRECATED
873  explicit Mesh(int n, double sx = 1.0)
874  {
875  Make1D(n, sx);
876  // Finalize(); // reminder: not needed
877  }
878 
879  /** Creates mesh by reading a file in MFEM, Netgen, or VTK format. If
880  generate_edges = 0 (default) edges are not generated, if 1 edges are
881  generated. See also @a Mesh::LoadFromFile. */
882  explicit Mesh(const char *filename, int generate_edges = 0, int refine = 1,
883  bool fix_orientation = true);
884 
885  /** Creates mesh by reading data stream in MFEM, Netgen, or VTK format. If
886  generate_edges = 0 (default) edges are not generated, if 1 edges are
887  generated. */
888  explicit Mesh(std::istream &input, int generate_edges = 0, int refine = 1,
889  bool fix_orientation = true);
890 
891  /// Create a disjoint mesh from the given mesh array
892  Mesh(Mesh *mesh_array[], int num_pieces);
893 
894  /// Deprecated: see @a MakeRefined.
895  MFEM_DEPRECATED
896  Mesh(Mesh *orig_mesh, int ref_factor, int ref_type);
897 
898  /** This is similar to the mesh constructor with the same arguments, but here
899  the current mesh is destroyed and another one created based on the data
900  stream again given in MFEM, Netgen, or VTK format. If generate_edges = 0
901  (default) edges are not generated, if 1 edges are generated. */
902  /// \see mfem::ifgzstream() for on-the-fly decompression of compressed ascii
903  /// inputs.
904  virtual void Load(std::istream &input, int generate_edges = 0,
905  int refine = 1, bool fix_orientation = true)
906  {
907  Loader(input, generate_edges);
908  Finalize(refine, fix_orientation);
909  }
910 
911  /// Clear the contents of the Mesh.
912  void Clear() { Destroy(); SetEmpty(); }
913 
914  /** @brief Get the mesh generator/type.
915 
916  The purpose of this is to be able to quickly tell what type of elements
917  one has in the mesh. Examination of this bitmask along with knowledge
918  of the mesh dimension can be used to identify which element types are
919  present.
920 
921  @return A bitmask:
922  - bit 0 - simplices are present in the mesh (triangles, tets),
923  - bit 1 - tensor product elements are present in the mesh (quads, hexes),
924  - bit 2 - the mesh has wedge elements.
925  - bit 3 - the mesh has pyramid elements.
926 
927  In parallel, the result takes into account elements on all processors.
928  */
929  inline int MeshGenerator() { return meshgen; }
930 
931  /** @brief Returns number of vertices. Vertices are only at the corners of
932  elements, where you would expect them in the lowest-order mesh. */
933  inline int GetNV() const { return NumOfVertices; }
934 
935  /// Returns number of elements.
936  inline int GetNE() const { return NumOfElements; }
937 
938  /// Returns number of boundary elements.
939  inline int GetNBE() const { return NumOfBdrElements; }
940 
941  /// Return the number of edges.
942  inline int GetNEdges() const { return NumOfEdges; }
943 
944  /// Return the number of faces in a 3D mesh.
945  inline int GetNFaces() const { return NumOfFaces; }
946 
947  /// Return the number of faces (3D), edges (2D) or vertices (1D).
948  int GetNumFaces() const;
949 
950  /** @brief Return the number of faces (3D), edges (2D) or vertices (1D)
951  including ghost faces. */
952  int GetNumFacesWithGhost() const;
953 
954  /** @brief Returns the number of faces according to the requested type, does
955  not count master nonconforming faces.
956 
957  If type==Boundary returns only the number of true boundary faces
958  contrary to GetNBE() that returns all "boundary" elements which may
959  include actual interior faces.
960  Similarly, if type==Interior, only the true interior faces are counted
961  excluding all master nonconforming faces. */
962  virtual int GetNFbyType(FaceType type) const;
963 
964  /// Utility function: sum integers from all processors (Allreduce).
965  virtual long long ReduceInt(int value) const { return value; }
966 
967  /// Return the total (global) number of elements.
968  long long GetGlobalNE() const { return ReduceInt(NumOfElements); }
969 
970  /** Return vertex to vertex table. The connections stored in the table
971  are from smaller to bigger vertex index, i.e. if i<j and (i, j) is
972  in the table, then (j, i) is not stored. */
973  void GetVertexToVertexTable(DSTable &) const;
974 
975  /** @brief Return the mesh geometric factors corresponding to the given
976  integration rule.
977 
978  The IntegrationRule used with GetGeometricFactors needs to remain valid
979  until the internally stored GeometricFactors objects are destroyed (by
980  calling Mesh::DeleteGeometricFactors(), Mesh::NodesUpdated(), or the Mesh
981  destructor).
982 
983  If the device MemoryType parameter @a d_mt is specified, then the
984  returned object will use that type unless it was previously allocated
985  with a different type.
986 
987  The returned pointer points to an internal object that may be invalidated
988  by mesh operations such as refinement, vertex/node movement, etc. Since
989  not all such modifications can be tracked by the Mesh class (e.g. when
990  using the pointer returned by GetNodes() to change the nodes) one needs
991  to account for such changes by calling the method NodesUpdated() which,
992  in particular, will call DeleteGeometricFactors(). */
994  const IntegrationRule& ir,
995  const int flags,
997 
998  /** @brief Return the mesh geometric factors for the faces corresponding
999  to the given integration rule.
1000 
1001  The IntegrationRule used with GetFaceGeometricFactors needs to remain
1002  valid until the internally stored FaceGeometricFactors objects are
1003  destroyed (by either calling Mesh::DeleteGeometricFactors(),
1004  Mesh::NodesUpdated(), or the Mesh destructor).
1005 
1006  If the device MemoryType parameter @a d_mt is specified, then the
1007  returned object will use that type unless it was previously allocated
1008  with a different type.
1009 
1010  The returned pointer points to an internal object that may be invalidated
1011  by mesh operations such as refinement, vertex/node movement, etc. Since
1012  not all such modifications can be tracked by the Mesh class (e.g. when
1013  using the pointer returned by GetNodes() to change the nodes) one needs
1014  to account for such changes by calling the method NodesUpdated() which,
1015  in particular, will call DeleteGeometricFactors(). */
1017  const IntegrationRule& ir,
1018  const int flags,
1019  FaceType type,
1021 
1022  /// Destroy all GeometricFactors stored by the Mesh.
1023  /** This method can be used to force recomputation of the GeometricFactors,
1024  for example, after the mesh nodes are modified externally.
1025 
1026  @note In general, the preferred method for resetting the GeometricFactors
1027  should be to call NodesUpdated(). */
1028  void DeleteGeometricFactors();
1029 
1030  /** @brief This function should be called after the mesh node coordinates
1031  have been updated externally, e.g. by modifying the internal nodal
1032  GridFunction returned by GetNodes(). */
1033  /** It deletes internal quantities derived from the node coordinates,
1034  such as the (Face)GeometricFactors.
1035 
1036  @note Unlike the similarly named protected method UpdateNodes() this
1037  method does not modify the nodes. */
1039 
1040  /// Equals 1 + num_holes - num_loops
1041  inline int EulerNumber() const
1043  /// Equals 1 - num_holes
1044  inline int EulerNumber2D() const
1045  { return NumOfVertices - NumOfEdges + NumOfElements; }
1046 
1047  int Dimension() const { return Dim; }
1048  int SpaceDimension() const { return spaceDim; }
1049 
1050  /// @brief Return pointer to vertex i's coordinates.
1051  /// @warning For high-order meshes (when #Nodes != NULL) vertices may not be
1052  /// updated and should not be used!
1053  const double *GetVertex(int i) const { return vertices[i](); }
1054 
1055  /// @brief Return pointer to vertex i's coordinates.
1056  /// @warning For high-order meshes (when Nodes != NULL) vertices may not
1057  /// being updated and should not be used!
1058  double *GetVertex(int i) { return vertices[i](); }
1059 
1060  void GetElementData(int geom, Array<int> &elem_vtx, Array<int> &attr) const
1061  { GetElementData(elements, geom, elem_vtx, attr); }
1062 
1063  /// Checks if the mesh has boundary elements
1064  virtual bool HasBoundaryElements() const { return (NumOfBdrElements > 0); }
1065 
1066  void GetBdrElementData(int geom, Array<int> &bdr_elem_vtx,
1067  Array<int> &bdr_attr) const
1068  { GetElementData(boundary, geom, bdr_elem_vtx, bdr_attr); }
1069 
1070  /** @brief Set the internal Vertex array to point to the given @a vertices
1071  array without assuming ownership of the pointer. */
1072  /** If @a zerocopy is `true`, the vertices must be given as an array of 3
1073  doubles per vertex. If @a zerocopy is `false` then the current Vertex
1074  data is first copied to the @a vertices array. */
1075  void ChangeVertexDataOwnership(double *vertices, int len_vertices,
1076  bool zerocopy = false);
1077 
1078  const Element* const *GetElementsArray() const
1079  { return elements.GetData(); }
1080 
1081  const Element *GetElement(int i) const { return elements[i]; }
1082 
1083  Element *GetElement(int i) { return elements[i]; }
1084 
1085  const Element *GetBdrElement(int i) const { return boundary[i]; }
1086 
1087  Element *GetBdrElement(int i) { return boundary[i]; }
1088 
1089  const Element *GetFace(int i) const { return faces[i]; }
1090 
1091  /// Return the Geometry::Type associated with face @a i.
1092  Geometry::Type GetFaceGeometry(int i) const;
1093 
1095  {
1096  return elements[i]->GetGeometryType();
1097  }
1098 
1100  {
1101  return boundary[i]->GetGeometryType();
1102  }
1103 
1104  /// Deprecated in favor of Mesh::GetFaceGeometry
1105  MFEM_DEPRECATED Geometry::Type GetFaceBaseGeometry(int i) const
1106  { return GetFaceGeometry(i); }
1107 
1109  { return GetElementGeometry(i); }
1110 
1112  { return GetBdrElementGeometry(i); }
1113 
1114  /** @brief Return true iff the given @a geom is encountered in the mesh.
1115  Geometries of dimensions lower than Dimension() are counted as well. */
1116  bool HasGeometry(Geometry::Type geom) const
1117  { return mesh_geoms & (1 << geom); }
1118 
1119  /** @brief Return the number of geometries of the given dimension present in
1120  the mesh. */
1121  /** For a parallel mesh only the local geometries are counted. */
1122  int GetNumGeometries(int dim) const;
1123 
1124  /// Return all element geometries of the given dimension present in the mesh.
1125  /** For a parallel mesh only the local geometries are returned.
1126 
1127  The returned geometries are sorted. */
1128  void GetGeometries(int dim, Array<Geometry::Type> &el_geoms) const;
1129 
1130  /// List of mesh geometries stored as Array<Geometry::Type>.
1131  class GeometryList : public Array<Geometry::Type>
1132  {
1133  protected:
1135  public:
1136  /// Construct a GeometryList of all element geometries in @a mesh.
1138  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
1139  { mesh.GetGeometries(mesh.Dimension(), *this); }
1140  /** @brief Construct a GeometryList of all geometries of dimension @a dim
1141  in @a mesh. */
1142  GeometryList(Mesh &mesh, int dim)
1143  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
1144  { mesh.GetGeometries(dim, *this); }
1145  };
1146 
1147  /// Returns the indices of the vertices of element i.
1148  void GetElementVertices(int i, Array<int> &v) const
1149  { elements[i]->GetVertices(v); }
1150 
1151  /// Returns the indices of the vertices of boundary element i.
1152  void GetBdrElementVertices(int i, Array<int> &v) const
1153  { boundary[i]->GetVertices(v); }
1154 
1155  /// Return the indices and the orientations of all edges of element i.
1156  void GetElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
1157 
1158  /// Return the indices and the orientations of all edges of bdr element i.
1159  void GetBdrElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
1160 
1161  /** Return the indices and the orientations of all edges of face i.
1162  Works for both 2D (face=edge) and 3D faces. */
1163  void GetFaceEdges(int i, Array<int> &edges, Array<int> &o) const;
1164 
1165  /// Returns the indices of the vertices of face i.
1166  void GetFaceVertices(int i, Array<int> &vert) const
1167  {
1168  if (Dim == 1)
1169  {
1170  vert.SetSize(1); vert[0] = i;
1171  }
1172  else
1173  {
1174  faces[i]->GetVertices(vert);
1175  }
1176  }
1177 
1178  /// Returns the indices of the vertices of edge i.
1179  void GetEdgeVertices(int i, Array<int> &vert) const;
1180 
1181  /// Returns the face-to-edge Table (3D)
1182  Table *GetFaceEdgeTable() const;
1183 
1184  /// Returns the edge-to-vertex Table (3D)
1185  Table *GetEdgeVertexTable() const;
1186 
1187  /// Return the indices and the orientations of all faces of element i.
1188  void GetElementFaces(int i, Array<int> &faces, Array<int> &ori) const;
1189 
1190  /** @brief Returns the sorted, unique indices of elements sharing a face with
1191  element @a elem, including @a elem. */
1192  Array<int> FindFaceNeighbors(const int elem) const;
1193 
1194  /// Return the index and the orientation of the face of bdr element i. (3D)
1195  void GetBdrElementFace(int i, int *f, int *o) const;
1196 
1197  /** Return the vertex index of boundary element i. (1D)
1198  Return the edge index of boundary element i. (2D)
1199  Return the face index of boundary element i. (3D) */
1200  int GetBdrElementEdgeIndex(int i) const;
1201 
1202  /** @brief For the given boundary element, bdr_el, return its adjacent
1203  element and its info, i.e. 64*local_bdr_index+bdr_orientation.
1204 
1205  The returned bdr_orientation is that of the boundary element relative to
1206  the respective face element.
1207 
1208  @sa GetBdrElementAdjacentElement2() */
1209  void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const;
1210 
1211  /** @brief For the given boundary element, bdr_el, return its adjacent
1212  element and its info, i.e. 64*local_bdr_index+inverse_bdr_orientation.
1213 
1214  The returned inverse_bdr_orientation is the inverse of the orientation of
1215  the boundary element relative to the respective face element. In other
1216  words this is the orientation of the face element relative to the
1217  boundary element.
1218 
1219  @sa GetBdrElementAdjacentElement() */
1220  void GetBdrElementAdjacentElement2(int bdr_el, int &el, int &info) const;
1221 
1222  /// Returns the type of element i.
1223  Element::Type GetElementType(int i) const;
1224 
1225  /// Returns the type of boundary element i.
1226  Element::Type GetBdrElementType(int i) const;
1227 
1228  /* Return point matrix of element i of dimension Dim X #v, where for every
1229  vertex we give its coordinates in space of dimension Dim. */
1230  void GetPointMatrix(int i, DenseMatrix &pointmat) const;
1231 
1232  /* Return point matrix of boundary element i of dimension Dim X #v, where for
1233  every vertex we give its coordinates in space of dimension Dim. */
1234  void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const;
1235 
1237 
1238  /// Builds the transformation defining the i-th element in @a ElTr.
1239  /// @a ElTr must be allocated in advance and will be owned by the caller.
1241 
1242  /// Returns a pointer to the transformation defining the i-th element.
1243  /// Note that the pointer is owned by the class and is shared, i.e., calling
1244  /// this function resets pointers obtained from previous calls.
1246 
1247  /// Builds the transformation defining the i-th element in @a ElTr
1248  /// assuming position of the vertices/nodes are given by @a nodes.
1249  /// @a ElTr must be allocated in advance and will be owned by the caller.
1250  void GetElementTransformation(int i, const Vector &nodes,
1252 
1253  /// Returns a pointer to the transformation defining the i-th boundary
1254  /// element. Note that the pointer is owned by the class and is shared, i.e.,
1255  /// calling this function resets pointers obtained from previous calls.
1257 
1258  /// Builds the transformation defining the i-th boundary element in @a ElTr.
1259  /// @a ElTr must be allocated in advance and will be owned by the caller.
1261 
1262  /// Builds the transformation defining the i-th face element in @a FTr.
1263  /// @a FTr must be allocated in advance and will be owned by the caller.
1265 
1266  /** @brief A helper method that constructs a transformation from the
1267  reference space of a face to the reference space of an element. */
1268  /** The local index of the face as a face in the element and its orientation
1269  are given by the input parameter @a info, as @a info = 64*loc_face_idx +
1270  loc_face_orientation. */
1271  void GetLocalFaceTransformation(int face_type, int elem_type,
1273  int info);
1274 
1275  /// Returns a pointer to the transformation defining the given face element.
1276  /// Note that the pointer is owned by the class and is shared, i.e., calling
1277  /// this function resets pointers obtained from previous calls.
1279 
1280  /// Builds the transformation defining the i-th edge element in @a EdTr.
1281  /// @a EdTr must be allocated in advance and will be owned by the caller.
1283 
1284  /// Returns a pointer to the transformation defining the given edge element.
1285  /// Note that the pointer is owned by the class and is shared, i.e., calling
1286  /// this function resets pointers obtained from previous calls.
1288 
1289  /// Returns (a pointer to an object containing) the following data:
1290  ///
1291  /// 1) Elem1No - the index of the first element that contains this face this
1292  /// is the element that has the same outward unit normal vector as the
1293  /// face;
1294  ///
1295  /// 2) Elem2No - the index of the second element that contains this face this
1296  /// element has outward unit normal vector as the face multiplied with -1;
1297  ///
1298  /// 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first
1299  /// and the second element respectively;
1300  ///
1301  /// 4) Face - pointer to the ElementTransformation of the face;
1302  ///
1303  /// 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face
1304  /// coordinate system to the element coordinate system (both in their
1305  /// reference elements). Used to transform IntegrationPoints from face to
1306  /// element. More formally, let:
1307  /// TL1, TL2 be the transformations represented by Loc1, Loc2,
1308  /// TE1, TE2 - the transformations represented by Elem1, Elem2,
1309  /// TF - the transformation represented by Face, then
1310  /// TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
1311  ///
1312  /// 6) FaceGeom - the base geometry for the face.
1313  ///
1314  /// The mask specifies which fields in the structure to return:
1315  /// mask & 1 - Elem1, mask & 2 - Elem2
1316  /// mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.
1317  /// These mask values are defined in the ConfigMasks enum type as part of the
1318  /// FaceElementTransformations class in fem/eltrans.hpp.
1319  ///
1320  /// Note that the pointer is owned by the class and is shared, i.e., calling
1321  /// this function resets pointers obtained from previous calls.
1323  int FaceNo,
1324  int mask = 31);
1325 
1326  /// See GetFaceElementTransformations().
1328  {
1329  if (faces_info[FaceNo].Elem2No < 0) { return NULL; }
1330  return GetFaceElementTransformations (FaceNo);
1331  }
1332 
1333  /// Builds the transformation defining the given boundary face.
1334  /// The returned pointer is owned by the caller.
1336 
1337  /// Return the local face index for the given boundary face.
1338  int GetBdrFace(int BdrElemNo) const;
1339 
1340  /// Return true if the given face is interior. @sa FaceIsTrueInterior().
1341  bool FaceIsInterior(int FaceNo) const
1342  {
1343  return (faces_info[FaceNo].Elem2No >= 0);
1344  }
1345 
1346  /** This enumerated type describes the three main face topologies:
1347  - Boundary, for faces on the boundary of the computational domain,
1348  - Conforming, for conforming faces interior to the computational domain,
1349  - Nonconforming, for nonconforming faces interior to the computational
1350  domain. */
1351  enum class FaceTopology { Boundary,
1352  Conforming,
1353  Nonconforming,
1354  NA
1355  };
1356 
1357  /** This enumerated type describes the location of the two elements sharing a
1358  face, Local meaning that the element is local to the MPI rank, FaceNbr
1359  meaning that the element is distributed on a different MPI rank, this
1360  typically means that methods with FaceNbr should be used to access the
1361  relevant information, e.g., ParFiniteElementSpace::GetFaceNbrElementVDofs.
1362  */
1363  enum class ElementLocation { Local, FaceNbr, NA };
1364 
1365  /** This enumerated type describes the topological relation of an element to
1366  a face:
1367  - Coincident meaning that the element's face is topologically equal to
1368  the mesh face.
1369  - Superset meaning that the element's face is topologically coarser than
1370  the mesh face, i.e., the element's face contains the mesh face.
1371  - Subset meaning that the element's face is topologically finer than the
1372  mesh face, i.e., the element's face is contained in the mesh face.
1373  Superset and Subset are only relevant for nonconforming faces.
1374  Master nonconforming faces have a conforming element on one side, and a
1375  fine element on the other side. Slave nonconforming faces have a
1376  conforming element on one side, and a coarse element on the other side.
1377  */
1379 
1380  /** This enumerated type describes the corresponding FaceInfo internal
1381  representation (encoded cases), c.f. FaceInfo's documentation:
1382  Classification of a local (non-ghost) face based on its FaceInfo:
1383  - Elem2No >= 0 --> local interior face; can be either:
1384  - NCFace == -1 --> LocalConforming,
1385  - NCFace >= 0 --> LocalSlaveNonconforming,
1386  - Elem2No < 0 --> local "boundary" face; can be one of:
1387  - NCFace == -1 --> conforming face; can be either:
1388  - Elem2Inf < 0 --> Boundary,
1389  - Elem2Inf >= 0 --> SharedConforming,
1390  - NCFace >= 0 --> nonconforming face; can be one of:
1391  - Elem2Inf < 0 --> MasterNonconforming (shared or not shared),
1392  - Elem2Inf >= 0 --> SharedSlaveNonconforming.
1393  Classification of a ghost (non-local) face based on its FaceInfo:
1394  - Elem1No == -1 --> GhostMaster (includes other unused ghost faces),
1395  - Elem1No >= 0 --> GhostSlave.
1396  */
1397  enum class FaceInfoTag { Boundary,
1403  GhostSlave,
1404  GhostMaster
1405  };
1406 
1407  /** @brief This structure is used as a human readable output format that
1408  decipheres the information contained in Mesh::FaceInfo when using the
1409  Mesh::GetFaceInformation() method.
1410 
1411  The element indices in this structure don't need further processing,
1412  contrary to the ones obtained through Mesh::GetFacesElements and can
1413  directly be used, e.g., Elem1 and Elem2 indices.
1414  Likewise the orientations for Elem1 and Elem2 already take into account
1415  special cases and can be used as is.
1416  */
1418  {
1420 
1421  struct
1422  {
1425  int index;
1428  } element[2];
1429 
1431  int ncface;
1433 
1434  /** @brief Return true if the face is a local interior face which is NOT
1435  a master nonconforming face. */
1436  bool IsLocal() const
1437  {
1438  return element[1].location == Mesh::ElementLocation::Local;
1439  }
1440 
1441  /** @brief Return true if the face is a shared interior face which is NOT
1442  a master nonconforming face. */
1443  bool IsShared() const
1444  {
1445  return element[1].location == Mesh::ElementLocation::FaceNbr;
1446  }
1447 
1448  /** @brief return true if the face is an interior face to the computation
1449  domain, either a local or shared interior face (not a boundary face)
1450  which is NOT a master nonconforming face.
1451  */
1452  bool IsInterior() const
1453  {
1454  return topology == FaceTopology::Conforming ||
1456  }
1457 
1458  /** @brief Return true if the face is a boundary face. */
1459  bool IsBoundary() const
1460  {
1461  return topology == FaceTopology::Boundary;
1462  }
1463 
1464  /// @brief Return true if the face is of the same type as @a type.
1465  bool IsOfFaceType(FaceType type) const
1466  {
1467  switch (type)
1468  {
1469  case FaceType::Interior:
1470  return IsInterior();
1471  case FaceType::Boundary:
1472  return IsBoundary();
1473  default:
1474  return false;
1475  }
1476  }
1477 
1478  /// @brief Return true if the face is a conforming face.
1479  bool IsConforming() const
1480  {
1482  }
1483 
1484  /// @brief Return true if the face is a nonconforming fine face.
1485  bool IsNonconformingFine() const
1486  {
1488  (element[0].conformity == ElementConformity::Superset ||
1489  element[1].conformity == ElementConformity::Superset);
1490  }
1491 
1492  /// @brief Return true if the face is a nonconforming coarse face.
1493  /** Note that ghost nonconforming master faces cannot be clearly
1494  identified as such with the currently available information, so this
1495  method will return false for such faces. */
1497  {
1499  element[1].conformity == ElementConformity::Subset;
1500  }
1501 
1502  /// @brief cast operator from FaceInformation to FaceInfo.
1503  operator Mesh::FaceInfo() const;
1504  };
1505 
1506  /** This method aims to provide face information in a deciphered format, i.e.
1507  Mesh::FaceInformation, compared to the raw encoded information returned
1508  by Mesh::GetFaceElements() and Mesh::GetFaceInfos(). */
1509  FaceInformation GetFaceInformation(int f) const;
1510 
1511  void GetFaceElements (int Face, int *Elem1, int *Elem2) const;
1512  void GetFaceInfos (int Face, int *Inf1, int *Inf2) const;
1513  void GetFaceInfos (int Face, int *Inf1, int *Inf2, int *NCFace) const;
1514 
1515  /// Deprecated in favor of Mesh::GetFaceGeometry
1516  MFEM_DEPRECATED Geometry::Type GetFaceGeometryType(int Face) const
1517  { return GetFaceGeometry(Face); }
1518 
1519  Element::Type GetFaceElementType(int Face) const;
1520 
1521  Array<int> GetFaceToBdrElMap() const;
1522 
1523  /// Check (and optionally attempt to fix) the orientation of the elements
1524  /** @param[in] fix_it If `true`, attempt to fix the orientations of some
1525  elements: triangles, quads, and tets.
1526  @return The number of elements with wrong orientation.
1527 
1528  @note For meshes with nodes (e.g. high-order or periodic meshes), fixing
1529  the element orientations may require additional permutation of the nodal
1530  GridFunction of the mesh which is not performed by this method. Instead,
1531  the method Finalize() should be used with the parameter
1532  @a fix_orientation set to `true`.
1533 
1534  @note This method performs a simple check if an element is inverted, e.g.
1535  for most elements types, it checks if the Jacobian of the mapping from
1536  the reference element is non-negative at the center of the element. */
1537  int CheckElementOrientation(bool fix_it = true);
1538 
1539  /// Check the orientation of the boundary elements
1540  /** @return The number of boundary elements with wrong orientation. */
1541  int CheckBdrElementOrientation(bool fix_it = true);
1542 
1543  /// Return the attribute of element i.
1544  int GetAttribute(int i) const { return elements[i]->GetAttribute(); }
1545 
1546  /// Set the attribute of element i.
1547  void SetAttribute(int i, int attr) { elements[i]->SetAttribute(attr); }
1548 
1549  /// Return the attribute of boundary element i.
1550  int GetBdrAttribute(int i) const { return boundary[i]->GetAttribute(); }
1551 
1552  /// Set the attribute of boundary element i.
1553  void SetBdrAttribute(int i, int attr) { boundary[i]->SetAttribute(attr); }
1554 
1555  const Table &ElementToElementTable();
1556 
1557  const Table &ElementToFaceTable() const;
1558 
1559  const Table &ElementToEdgeTable() const;
1560 
1561  /// The returned Table must be destroyed by the caller
1563 
1564  /** Return the "face"-element Table. Here "face" refers to face (3D),
1565  edge (2D), or vertex (1D).
1566  The returned Table must be destroyed by the caller. */
1567  Table *GetFaceToElementTable() const;
1568 
1569  /** This method modifies a tetrahedral mesh so that Nedelec spaces of order
1570  greater than 1 can be defined on the mesh. Specifically, we
1571  1) rotate all tets in the mesh so that the vertices {v0, v1, v2, v3}
1572  satisfy: v0 < v1 < min(v2, v3).
1573  2) rotate all boundary triangles so that the vertices {v0, v1, v2}
1574  satisfy: v0 < min(v1, v2).
1575 
1576  @note Refinement does not work after a call to this method! */
1577  MFEM_DEPRECATED virtual void ReorientTetMesh();
1578 
1579  int *CartesianPartitioning(int nxyz[]);
1580  int *GeneratePartitioning(int nparts, int part_method = 1);
1581  void CheckPartitioning(int *partitioning_);
1582 
1583  void CheckDisplacements(const Vector &displacements, double &tmax);
1584 
1585  // Vertices are only at the corners of elements, where you would expect them
1586  // in the lowest-order mesh.
1587  void MoveVertices(const Vector &displacements);
1588  void GetVertices(Vector &vert_coord) const;
1589  void SetVertices(const Vector &vert_coord);
1590 
1591  // Nodes are only active for higher order meshes, and share locations with
1592  // the vertices, plus all the higher- order control points within the element
1593  // and along the edges and on the faces.
1594  void GetNode(int i, double *coord) const;
1595  void SetNode(int i, const double *coord);
1596 
1597  // Node operations for curved mesh.
1598  // They call the corresponding '...Vertices' method if the
1599  // mesh is not curved (i.e. Nodes == NULL).
1600  void MoveNodes(const Vector &displacements);
1601  void GetNodes(Vector &node_coord) const;
1602  /// Updates the vertex/node locations. Invokes NodesUpdated().
1603  void SetNodes(const Vector &node_coord);
1604 
1605  /// Return a pointer to the internal node GridFunction (may be NULL).
1606  /** If the mesh is straight-sided (low-order), it may not have a GridFunction
1607  for the nodes, in which case this function returns NULL. To ensure that
1608  the nodal GridFunction exists, call EnsureNodes().
1609  @sa SetCurvature(). */
1610  GridFunction *GetNodes() { return Nodes; }
1611  const GridFunction *GetNodes() const { return Nodes; }
1612  /// Return the mesh nodes ownership flag.
1613  bool OwnsNodes() const { return own_nodes; }
1614  /// Set the mesh nodes ownership flag.
1615  void SetNodesOwner(bool nodes_owner) { own_nodes = nodes_owner; }
1616  /// Replace the internal node GridFunction with the given GridFunction.
1617  /** Invokes NodesUpdated(). */
1618  void NewNodes(GridFunction &nodes, bool make_owner = false);
1619  /** @brief Swap the internal node GridFunction pointer and ownership flag
1620  members with the given ones. */
1621  /** Invokes NodesUpdated(). */
1622  void SwapNodes(GridFunction *&nodes, int &own_nodes_);
1623 
1624  /// Return the mesh nodes/vertices projected on the given GridFunction.
1625  void GetNodes(GridFunction &nodes) const;
1626  /** Replace the internal node GridFunction with a new GridFunction defined
1627  on the given FiniteElementSpace. The new node coordinates are projected
1628  (derived) from the current nodes/vertices. */
1629  virtual void SetNodalFESpace(FiniteElementSpace *nfes);
1630  /** Replace the internal node GridFunction with the given GridFunction. The
1631  given GridFunction is updated with node coordinates projected (derived)
1632  from the current nodes/vertices. */
1633  void SetNodalGridFunction(GridFunction *nodes, bool make_owner = false);
1634  /** Return the FiniteElementSpace on which the current mesh nodes are
1635  defined or NULL if the mesh does not have nodes. */
1636  const FiniteElementSpace *GetNodalFESpace() const;
1637  /** @brief Make sure that the mesh has valid nodes, i.e. its geometry is
1638  described by a vector finite element grid function (even if it is a
1639  low-order mesh with straight edges).
1640 
1641  @sa GetNodes(). */
1642  void EnsureNodes();
1643 
1644  /// Set the curvature of the mesh nodes using the given polynomial degree.
1645  /** Creates a nodal GridFunction if one doesn't already exist.
1646 
1647  @param[in] order Polynomial degree of the nodal FE space.
1648  @param[in] discont Whether to use a discontinuous or continuous
1649  finite element space (continuous is default).
1650  @param[in] space_dim The space dimension (optional).
1651  @param[in] ordering The Ordering of the finite element space
1652  (Ordering::byVDIM is the default). */
1653  virtual void SetCurvature(int order, bool discont = false, int space_dim = -1,
1654  int ordering = 1);
1655 
1656  /// Refine all mesh elements.
1657  /** @param[in] ref_algo %Refinement algorithm. Currently used only for pure
1658  tetrahedral meshes. If set to zero (default), a tet mesh will be refined
1659  using algorithm A, that produces elements with better quality compared to
1660  algorithm B used when the parameter is non-zero.
1661 
1662  For tetrahedral meshes, after using algorithm A, the mesh cannot be
1663  refined locally using methods like GeneralRefinement() unless it is
1664  re-finalized using Finalize() with the parameter @a refine set to true.
1665  Note that calling Finalize() in this way will generally invalidate any
1666  FiniteElementSpace%s and GridFunction%s defined on the mesh. */
1667  void UniformRefinement(int ref_algo = 0);
1668 
1669  /** Refine selected mesh elements. Refinement type can be specified for each
1670  element. The function can do conforming refinement of triangles and
1671  tetrahedra and nonconforming refinement (i.e., with hanging-nodes) of
1672  triangles, quadrilaterals and hexahedra. If 'nonconforming' = -1,
1673  suitable refinement method is selected automatically (namely, conforming
1674  refinement for triangles). Use nonconforming = 0/1 to force the method.
1675  For nonconforming refinements, nc_limit optionally specifies the maximum
1676  level of hanging nodes (unlimited by default). */
1677  void GeneralRefinement(const Array<Refinement> &refinements,
1678  int nonconforming = -1, int nc_limit = 0);
1679 
1680  /** Simplified version of GeneralRefinement taking a simple list of elements
1681  to refine, without refinement types. */
1682  void GeneralRefinement(const Array<int> &el_to_refine,
1683  int nonconforming = -1, int nc_limit = 0);
1684 
1685  /// Refine each element with given probability. Uses GeneralRefinement.
1686  void RandomRefinement(double prob, bool aniso = false,
1687  int nonconforming = -1, int nc_limit = 0);
1688 
1689  /// Refine elements sharing the specified vertex. Uses GeneralRefinement.
1690  void RefineAtVertex(const Vertex& vert,
1691  double eps = 0.0, int nonconforming = -1);
1692 
1693  /** Refine element i if elem_error[i] > threshold, for all i.
1694  Returns true if at least one element was refined, false otherwise. */
1695  bool RefineByError(const Array<double> &elem_error, double threshold,
1696  int nonconforming = -1, int nc_limit = 0);
1697 
1698  /** Refine element i if elem_error(i) > threshold, for all i.
1699  Returns true if at least one element was refined, false otherwise. */
1700  bool RefineByError(const Vector &elem_error, double threshold,
1701  int nonconforming = -1, int nc_limit = 0);
1702 
1703  /** Derefine the mesh based on an error measure associated with each
1704  element. A derefinement is performed if the sum of errors of its fine
1705  elements is smaller than 'threshold'. If 'nc_limit' > 0, derefinements
1706  that would increase the maximum level of hanging nodes of the mesh are
1707  skipped. Returns true if the mesh changed, false otherwise. */
1708  bool DerefineByError(Array<double> &elem_error, double threshold,
1709  int nc_limit = 0, int op = 1);
1710 
1711  /// Same as DerefineByError for an error vector.
1712  bool DerefineByError(const Vector &elem_error, double threshold,
1713  int nc_limit = 0, int op = 1);
1714 
1715  ///@{ @name NURBS mesh refinement methods
1716  void KnotInsert(Array<KnotVector *> &kv);
1717  void KnotInsert(Array<Vector *> &kv);
1718  /* For each knot vector:
1719  new_degree = max(old_degree, min(old_degree + rel_degree, degree)). */
1720  void DegreeElevate(int rel_degree, int degree = 16);
1721  ///@}
1722 
1723  /** Make sure that a quad/hex mesh is considered to be nonconforming (i.e.,
1724  has an associated NCMesh object). Simplex meshes can be both conforming
1725  (default) or nonconforming. */
1726  void EnsureNCMesh(bool simplices_nonconforming = false);
1727 
1728  bool Conforming() const { return ncmesh == NULL; }
1729  bool Nonconforming() const { return ncmesh != NULL; }
1730 
1731  /** Return fine element transformations following a mesh refinement.
1732  Space uses this to construct a global interpolation matrix. */
1734 
1735  /// Return type of last modification of the mesh.
1737 
1738  /** Return update counter. The counter starts at zero and is incremented
1739  each time refinement, derefinement, or rebalancing method is called.
1740  It is used for checking proper sequence of Space:: and GridFunction::
1741  Update() calls. */
1742  long GetSequence() const { return sequence; }
1743 
1744  /// Print the mesh to the given stream using Netgen/Truegrid format.
1745  virtual void PrintXG(std::ostream &os = mfem::out) const;
1746 
1747  /// Print the mesh to the given stream using the default MFEM mesh format.
1748  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1749  virtual void Print(std::ostream &os = mfem::out) const { Printer(os); }
1750 
1751  /// Save the mesh to a file using Mesh::Print. The given @a precision will be
1752  /// used for ASCII output.
1753  virtual void Save(const char *fname, int precision=16) const;
1754 
1755  /// Print the mesh to the given stream using the adios2 bp format
1756 #ifdef MFEM_USE_ADIOS2
1757  virtual void Print(adios2stream &os) const;
1758 #endif
1759  /// Print the mesh in VTK format (linear and quadratic meshes only).
1760  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1761  void PrintVTK(std::ostream &os);
1762  /** Print the mesh in VTK format. The parameter ref > 0 specifies an element
1763  subdivision number (useful for high order fields and curved meshes).
1764  If the optional field_data is set, we also add a FIELD section in the
1765  beginning of the file with additional dataset information. */
1766  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1767  void PrintVTK(std::ostream &os, int ref, int field_data=0);
1768  /** Print the mesh in VTU format. The parameter ref > 0 specifies an element
1769  subdivision number (useful for high order fields and curved meshes).
1770  If @a bdr_elements is true, then output (only) the boundary elements,
1771  otherwise output only the non-boundary elements. */
1772  void PrintVTU(std::ostream &os,
1773  int ref=1,
1774  VTKFormat format=VTKFormat::ASCII,
1775  bool high_order_output=false,
1776  int compression_level=0,
1777  bool bdr_elements=false);
1778  /** Print the mesh in VTU format with file name fname. */
1779  virtual void PrintVTU(std::string fname,
1780  VTKFormat format=VTKFormat::ASCII,
1781  bool high_order_output=false,
1782  int compression_level=0,
1783  bool bdr=false);
1784  /** Print the boundary elements of the mesh in VTU format, and output the
1785  boundary attributes as a data array (useful for boundary conditions). */
1786  void PrintBdrVTU(std::string fname,
1787  VTKFormat format=VTKFormat::ASCII,
1788  bool high_order_output=false,
1789  int compression_level=0);
1790 
1791  void GetElementColoring(Array<int> &colors, int el0 = 0);
1792 
1793  /** @brief Prints the mesh with boundary elements given by the boundary of
1794  the subdomains, so that the boundary of subdomain i has boundary
1795  attribute i+1. */
1796  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1797  void PrintWithPartitioning (int *partitioning,
1798  std::ostream &os, int elem_attr = 0) const;
1799 
1800  void PrintElementsWithPartitioning (int *partitioning,
1801  std::ostream &out,
1802  int interior_faces = 0);
1803 
1804  /// Print set of disjoint surfaces:
1805  /*!
1806  * If Aface_face(i,j) != 0, print face j as a boundary
1807  * element with attribute i+1.
1808  */
1809  void PrintSurfaces(const Table &Aface_face, std::ostream &out) const;
1810 
1811  void ScaleSubdomains (double sf);
1812  void ScaleElements (double sf);
1813 
1814  void Transform(void (*f)(const Vector&, Vector&));
1815  void Transform(VectorCoefficient &deformation);
1816 
1817  /// Remove unused vertices and rebuild mesh connectivity.
1818  void RemoveUnusedVertices();
1819 
1820  /** Remove boundary elements that lie in the interior of the mesh, i.e. that
1821  have two adjacent faces in 3D, or edges in 2D. */
1822  void RemoveInternalBoundaries();
1823 
1824  /** @brief Get the size of the i-th element relative to the perfect
1825  reference element. */
1826  double GetElementSize(int i, int type = 0);
1827 
1828  double GetElementSize(int i, const Vector &dir);
1829 
1830  double GetElementVolume(int i);
1831 
1832  void GetElementCenter(int i, Vector &center);
1833 
1834  /// Returns the minimum and maximum corners of the mesh bounding box.
1835  /** For high-order meshes, the geometry is first refined @a ref times. */
1836  void GetBoundingBox(Vector &min, Vector &max, int ref = 2);
1837 
1838  void GetCharacteristics(double &h_min, double &h_max,
1839  double &kappa_min, double &kappa_max,
1840  Vector *Vh = NULL, Vector *Vk = NULL);
1841 
1842  /// Auxiliary method used by PrintCharacteristics().
1843  /** It is also used in the `mesh-explorer` miniapp. */
1844  static void PrintElementsByGeometry(int dim,
1845  const Array<int> &num_elems_by_geom,
1846  std::ostream &out);
1847 
1848  /** @brief Compute and print mesh characteristics such as number of vertices,
1849  number of elements, number of boundary elements, minimal and maximal
1850  element sizes, minimal and maximal element aspect ratios, etc. */
1851  /** If @a Vh or @a Vk are not NULL, return the element sizes and aspect
1852  ratios for all elements in the given Vector%s. */
1853  void PrintCharacteristics(Vector *Vh = NULL, Vector *Vk = NULL,
1854  std::ostream &os = mfem::out);
1855 
1856  /** @brief In serial, this method calls PrintCharacteristics(). In parallel,
1857  additional information about the parallel decomposition is also printed.
1858  */
1859  virtual void PrintInfo(std::ostream &os = mfem::out)
1860  {
1861  PrintCharacteristics(NULL, NULL, os);
1862  }
1863 
1864  /** @brief Find the ids of the elements that contain the given points, and
1865  their corresponding reference coordinates.
1866 
1867  The DenseMatrix @a point_mat describes the given points - one point for
1868  each column; it should have SpaceDimension() rows.
1869 
1870  The InverseElementTransformation object, @a inv_trans, is used to attempt
1871  the element transformation inversion. If NULL pointer is given, the
1872  method will use a default constructed InverseElementTransformation. Note
1873  that the algorithms in the base class InverseElementTransformation can be
1874  completely overwritten by deriving custom classes that override the
1875  Transform() method.
1876 
1877  If no element is found for the i-th point, elem_ids[i] is set to -1.
1878 
1879  In the ParMesh implementation, the @a point_mat is expected to be the
1880  same on all ranks. If the i-th point is found by multiple ranks, only one
1881  of them will mark that point as found, i.e. set its elem_ids[i] to a
1882  non-negative number; the other ranks will set their elem_ids[i] to -2 to
1883  indicate that the point was found but assigned to another rank.
1884 
1885  @returns The total number of points that were found.
1886 
1887  @note This method is not 100 percent reliable, i.e. it is not guaranteed
1888  to find a point, even if it lies inside a mesh element. */
1889  virtual int FindPoints(DenseMatrix& point_mat, Array<int>& elem_ids,
1890  Array<IntegrationPoint>& ips, bool warn = true,
1891  InverseElementTransformation *inv_trans = NULL);
1892 
1893  /// Swaps internal data with another mesh. By default, non-geometry members
1894  /// like 'ncmesh' and 'NURBSExt' are only swapped when 'non_geometry' is set.
1895  void Swap(Mesh& other, bool non_geometry);
1896 
1897  /// Destroys Mesh.
1898  virtual ~Mesh() { DestroyPointers(); }
1899 
1900 #ifdef MFEM_DEBUG
1901  /// Output an NCMesh-compatible debug dump.
1902  void DebugDump(std::ostream &out) const;
1903 #endif
1904 };
1905 
1906 /** Overload operator<< for std::ostream and Mesh; valid also for the derived
1907  class ParMesh */
1908 std::ostream &operator<<(std::ostream &out, const Mesh &mesh);
1909 
1910 
1911 /** @brief Structure for storing mesh geometric factors: coordinates, Jacobians,
1912  and determinants of the Jacobians. */
1913 /** Typically objects of this type are constructed and owned by objects of class
1914  Mesh. See Mesh::GetGeometricFactors(). */
1916 {
1917 
1918 private:
1919  void Compute(const GridFunction &nodes,
1921 
1922 public:
1923  const Mesh *mesh;
1926 
1928  {
1929  COORDINATES = 1 << 0,
1930  JACOBIANS = 1 << 1,
1931  DETERMINANTS = 1 << 2,
1932  };
1933 
1934  GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags,
1936 
1937  GeometricFactors(const GridFunction &nodes, const IntegrationRule &ir,
1938  int flags,
1940 
1941  /// Mapped (physical) coordinates of all quadrature points.
1942  /** This array uses a column-major layout with dimensions (NQ x SDIM x NE)
1943  where
1944  - NQ = number of quadrature points per element,
1945  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1946  - NE = number of elements in the mesh. */
1948 
1949  /// Jacobians of the element transformations at all quadrature points.
1950  /** This array uses a column-major layout with dimensions (NQ x SDIM x DIM x
1951  NE) where
1952  - NQ = number of quadrature points per element,
1953  - SDIM = space dimension of the mesh = mesh.SpaceDimension(),
1954  - DIM = dimension of the mesh = mesh.Dimension(), and
1955  - NE = number of elements in the mesh. */
1957 
1958  /// Determinants of the Jacobians at all quadrature points.
1959  /** This array uses a column-major layout with dimensions (NQ x NE) where
1960  - NQ = number of quadrature points per element, and
1961  - NE = number of elements in the mesh. */
1963 };
1964 
1965 /** @brief Structure for storing face geometric factors: coordinates, Jacobians,
1966  determinants of the Jacobians, and normal vectors. */
1967 /** Typically objects of this type are constructed and owned by objects of class
1968  Mesh. See Mesh::GetFaceGeometricFactors(). */
1970 {
1971 public:
1972  const Mesh *mesh;
1976 
1978  {
1979  COORDINATES = 1 << 0,
1980  JACOBIANS = 1 << 1,
1981  DETERMINANTS = 1 << 2,
1982  NORMALS = 1 << 3,
1983  };
1984 
1985  FaceGeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags,
1987 
1988  /// Mapped (physical) coordinates of all quadrature points.
1989  /** This array uses a column-major layout with dimensions (NQ x SDIM x NF)
1990  where
1991  - NQ = number of quadrature points per face,
1992  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1993  - NF = number of faces in the mesh. */
1995 
1996  /// Jacobians of the element transformations at all quadrature points.
1997  /** This array uses a column-major layout with dimensions (NQ x SDIM x DIM x
1998  NF) where
1999  - NQ = number of quadrature points per face,
2000  - SDIM = space dimension of the mesh = mesh.SpaceDimension(),
2001  - DIM = dimension of the mesh = mesh.Dimension(), and
2002  - NF = number of faces in the mesh. */
2004 
2005  /// Determinants of the Jacobians at all quadrature points.
2006  /** This array uses a column-major layout with dimensions (NQ x NF) where
2007  - NQ = number of quadrature points per face, and
2008  - NF = number of faces in the mesh. */
2010 
2011  /// Normals at all quadrature points.
2012  /** This array uses a column-major layout with dimensions (NQ x DIM x NF) where
2013  - NQ = number of quadrature points per face,
2014  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
2015  - NF = number of faces in the mesh. */
2017 };
2018 
2019 /// Class used to extrude the nodes of a mesh
2021 {
2022 private:
2023  int n, layer;
2024  double p[2], s;
2025  Vector tip;
2026 public:
2027  NodeExtrudeCoefficient(const int dim, const int n_, const double s_);
2028  void SetLayer(const int l) { layer = l; }
2030  virtual void Eval(Vector &V, ElementTransformation &T,
2031  const IntegrationPoint &ip);
2033 };
2034 
2035 
2036 /// Extrude a 1D mesh
2037 Mesh *Extrude1D(Mesh *mesh, const int ny, const double sy,
2038  const bool closed = false);
2039 
2040 /// Extrude a 2D mesh
2041 Mesh *Extrude2D(Mesh *mesh, const int nz, const double sz);
2042 
2043 // shift cyclically 3 integers left-to-right
2044 inline void ShiftRight(int &a, int &b, int &c)
2045 {
2046  int t = a;
2047  a = c; c = b; b = t;
2048 }
2049 
2050 /// @brief Print function for Mesh::FaceInformation.
2051 std::ostream& operator<<(std::ostream& os, const Mesh::FaceInformation& info);
2052 
2053 }
2054 
2055 #endif
static Mesh MakePeriodic(const Mesh &orig_mesh, const std::vector< int > &v2v)
Create a periodic mesh by identifying vertices of orig_mesh.
Definition: mesh.cpp:4904
Abstract class for all finite elements.
Definition: fe_base.hpp:232
const IntegrationRule * IntRule
Definition: mesh.hpp:1973
MFEM_DEPRECATED Mesh(int n, double sx=1.0)
Deprecated: see MakeCartesian1D.
Definition: mesh.hpp:873
void Loader(std::istream &input, int generate_edges=0, std::string parse_tag="")
Definition: mesh.cpp:3979
double GetGeckoElementOrdering(Array< int > &ordering, int iterations=4, int window=4, int period=2, int seed=0, bool verbose=false, double time_limit=0)
Definition: mesh.cpp:2042
Table * GetEdgeVertexTable() const
Returns the edge-to-vertex Table (3D)
Definition: mesh.cpp:6066
Geometry::Constants< Geometry::PYRAMID > pyr_t
Definition: mesh.hpp:268
Array< int > GetFaceToBdrElMap() const
Definition: mesh.cpp:1437
static const int vtk_quadratic_hex[27]
Definition: mesh.hpp:251
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:7030
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5351
const DenseMatrix * PointMatrix
Definition: mesh.hpp:208
void ScaleElements(double sf)
Definition: mesh.cpp:11621
const GridFunction * GetNodes() const
Definition: mesh.hpp:1611
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
Definition: mesh.cpp:5945
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo)
void FreeElement(Element *E)
Definition: mesh.cpp:11937
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1686
int CheckElementOrientation(bool fix_it=true)
Check (and optionally attempt to fix) the orientation of the elements.
Definition: mesh.cpp:5433
static bool remove_unused_vertices
Definition: mesh.hpp:286
void SetVertices(const Vector &vert_coord)
Definition: mesh.cpp:7889
static const int vtk_quadratic_tet[10]
Definition: mesh.hpp:248
void Make2D(int nx, int ny, Element::Type type, double sx, double sy, bool generate_edges, bool sfc_ordering)
Definition: mesh.cpp:3351
Vector J
Jacobians of the element transformations at all quadrature points.
Definition: mesh.hpp:2003
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition: mesh.hpp:1947
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
Definition: mesh.cpp:12224
static Mesh MakeSimplicial(const Mesh &orig_mesh)
Definition: mesh.cpp:4575
Base class for vector Coefficients that optionally depend on time and space.
void UniformRefinement3D_base(Array< int > *f2qf=NULL, DSTable *v_to_v_p=NULL, bool update_nodes=true)
Definition: mesh.cpp:8206
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
Definition: mesh.cpp:840
int EulerNumber2D() const
Equals 1 - num_holes.
Definition: mesh.hpp:1044
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1108
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:6235
static void PrintElementsByGeometry(int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
Auxiliary method used by PrintCharacteristics().
Definition: mesh.cpp:228
void AddHexAsTets(const int *vi, int attr=1)
Definition: mesh.cpp:1765
static const int NumGeom
Definition: geom.hpp:42
void AddHexAsWedges(const int *vi, int attr=1)
Definition: mesh.cpp:1784
virtual void PrintInfo(std::ostream &os=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
Definition: mesh.hpp:1859
const Table & ElementToFaceTable() const
Definition: mesh.cpp:6472
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
Array< Element * > boundary
Definition: mesh.hpp:91
Geometry::Constants< Geometry::SQUARE > quad_t
Definition: mesh.hpp:264
int Dimension() const
Definition: mesh.hpp:1047
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:7074
CoarseFineTransformations CoarseFineTr
Definition: mesh.hpp:240
int own_nodes
Definition: mesh.hpp:246
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:5389
void MoveVertices(const Vector &displacements)
Definition: mesh.cpp:7869
bool HasGeometry(Geometry::Type geom) const
Return true iff the given geom is encountered in the mesh. Geometries of dimensions lower than Dimens...
Definition: mesh.hpp:1116
IsoparametricTransformation Transformation
Definition: mesh.hpp:234
void GetLocalPtToSegTransformation(IsoparametricTransformation &, int)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:641
void GetCharacteristics(double &h_min, double &h_max, double &kappa_min, double &kappa_max, Vector *Vh=NULL, Vector *Vk=NULL)
Definition: mesh.cpp:193
static FiniteElement * GetTransformationFEforElementType(Element::Type)
Definition: mesh.cpp:327
FaceGeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, FaceType type, MemoryType d_mt=MemoryType::DEFAULT)
Definition: mesh.cpp:12157
static Mesh MakeCartesian3D(int nx, int ny, int nz, Element::Type type, double sx=1.0, double sy=1.0, double sz=1.0, bool sfc_ordering=true)
Definition: mesh.cpp:3736
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:7976
int NumOfEdges
Definition: mesh.hpp:70
virtual void UniformRefinement2D()
Refine a mixed 2D mesh uniformly.
Definition: mesh.hpp:404
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones...
Definition: mesh.cpp:7998
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:6296
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:129
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Definition: mesh.cpp:1413
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
Definition: mesh.hpp:1341
void ReadNetgen2DMesh(std::istream &input, int &curved)
void ShiftRight(int &a, int &b, int &c)
Definition: mesh.hpp:2044
void GetVertexToVertexTable(DSTable &) const
Definition: mesh.cpp:6372
static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type)
Create a refined (by any factor) version of orig_mesh.
Definition: mesh.cpp:3746
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
void CreateVTKMesh(const Vector &points, const Array< int > &cell_data, const Array< int > &cell_offsets, const Array< int > &cell_types, const Array< int > &cell_attributes, int &curved, int &read_gf, bool &finalize_topo)
long GetSequence() const
Definition: mesh.hpp:1742
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1672
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition: mesh.hpp:2009
bool Nonconforming() const
Definition: mesh.hpp:1729
int AddBdrPoint(int v, int attr=1)
Definition: mesh.cpp:1892
void DeleteTables()
Definition: mesh.hpp:295
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
GridFunction * Nodes
Definition: mesh.hpp:245
int NumOfElements
Definition: mesh.hpp:69
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
Definition: mesh.hpp:1152
Mesh * Extrude1D(Mesh *mesh, const int ny, const double sy, const bool closed)
Extrude a 1D mesh.
Definition: mesh.cpp:12242
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:11691
bool IsSlaveFace(const FaceInfo &fi) const
Definition: mesh.cpp:1051
void FinalizeCheck()
Definition: mesh.cpp:1937
void ApplyLocalSlaveTransformation(FaceElementTransformations &FT, const FaceInfo &fi, bool is_ghost)
Definition: mesh.cpp:1056
void AverageVertices(const int *indexes, int n, int result)
Averages the vertices with given indexes and saves the result in vertices[result].
Definition: mesh.cpp:8011
void FinalizeWedgeMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a wedge Mesh.
Definition: mesh.cpp:2830
Element * GetElement(int i)
Definition: mesh.hpp:1083
void GenerateBoundaryElements()
Definition: mesh.cpp:1899
void GetElementData(const Array< Element *> &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.cpp:9403
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:942
const Element * GetElement(int i) const
Definition: mesh.hpp:1081
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &os=mfem::out)
Compute and print mesh characteristics such as number of vertices, number of elements, number of boundary elements, minimal and maximal element sizes, minimal and maximal element aspect ratios, etc.
Definition: mesh.cpp:242
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1915
NodeExtrudeCoefficient(const int dim, const int n_, const double s_)
Definition: mesh.cpp:12218
Data type for vertex.
Definition: vertex.hpp:22
void SetMeshGen()
Determine the mesh generator bitmask meshgen, see MeshGenerator().
Definition: mesh.cpp:3926
bool IsNonconformingFine() const
Return true if the face is a nonconforming fine face.
Definition: mesh.hpp:1485
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
Definition: mesh.cpp:5629
void ReadNetgen3DMesh(std::istream &input)
Data arrays will be written in ASCII format.
The inverse transformation of a given ElementTransformation.
Definition: eltrans.hpp:185
const Element *const * GetElementsArray() const
Definition: mesh.hpp:1078
void ReadInlineMesh(std::istream &input, bool generate_edges=false)
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
See GetFaceElementTransformations().
Definition: mesh.hpp:1327
ElementLocation
Definition: mesh.hpp:1363
void RemoveInternalBoundaries()
Definition: mesh.cpp:11849
ElementConformity
Definition: mesh.hpp:1378
Mesh * Extrude2D(Mesh *mesh, const int nz, const double sz)
Extrude a 2D mesh.
Definition: mesh.cpp:12402
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:961
void SetEmpty()
Definition: mesh.cpp:1473
Array< Element * > faces
Definition: mesh.hpp:92
std::vector< int > CreatePeriodicVertexMapping(const std::vector< Vector > &translations, double tol=1e-8) const
Creates a mapping v2v from the vertex indices of the mesh such that coincident vertices under the giv...
Definition: mesh.cpp:4938
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
Definition: mesh.cpp:2865
void GetPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:6301
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:939
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:359
Operation last_operation
Definition: mesh.hpp:289
void MarkTriMeshForRefinement()
Definition: mesh.cpp:2430
void ReadCubit(const char *filename, int &curved, int &read_gf)
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5921
void InitRefinementTransforms()
Definition: mesh.cpp:9968
IsoparametricTransformation FaceTransformation
Definition: mesh.hpp:236
void UniformRefinement2D_base(bool update_nodes=true)
Definition: mesh.cpp:8047
Array< FaceGeometricFactors * > face_geom_factors
Definition: mesh.hpp:280
Array< NCFaceInfo > nc_faces_info
Definition: mesh.hpp:218
Element * ReadElement(std::istream &)
Definition: mesh.cpp:3908
void GetLocalTriToPyrTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:746
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
Element * NewElement(int geom)
Definition: mesh.cpp:3854
FaceElementTransformations FaceElemTr
Definition: mesh.hpp:237
Table * el_to_face
Definition: mesh.hpp:221
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
Definition: mesh.cpp:574
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:1394
void MarkForRefinement()
Definition: mesh.cpp:2413
void SetVerticesFromNodes(const GridFunction *nodes)
Helper to set vertex coordinates given a high-order curvature function.
Definition: mesh.cpp:5375
Structure for storing face geometric factors: coordinates, Jacobians, determinants of the Jacobians...
Definition: mesh.hpp:1969
ElementTransformation * GetBdrElementTransformation(int i)
Definition: mesh.cpp:431
void DebugDump(std::ostream &out) const
Output an NCMesh-compatible debug dump.
Definition: mesh.cpp:12629
int AddBdrTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1849
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:933
void MakeRefined_(Mesh &orig_mesh, const Array< int > ref_factors, int ref_type)
Internal function used in Mesh::MakeRefined.
Definition: mesh.cpp:4346
Geometry::Constants< Geometry::SEGMENT > seg_t
Definition: mesh.hpp:262
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1544
void GetElementJacobian(int i, DenseMatrix &J)
Definition: mesh.cpp:60
bool IsLocal() const
Return true if the face is a local interior face which is NOT a master nonconforming face...
Definition: mesh.hpp:1436
static int GetTetOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
Definition: mesh.cpp:5677
void GetBdrElementAdjacentElement2(int bdr_el, int &el, int &info) const
For the given boundary element, bdr_el, return its adjacent element and its info, i...
Definition: mesh.cpp:6269
MFEM_DEPRECATED Mesh(int nx, int ny, int nz, Element::Type type, bool generate_edges=false, double sx=1.0, double sy=1.0, double sz=1.0, bool sfc_ordering=true)
Deprecated: see MakeCartesian3D.
Definition: mesh.hpp:854
int AddTri(const int *vi, int attr=1)
Definition: mesh.hpp:729
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1618
void GetBdrElementFace(int i, int *f, int *o) const
Return the index and the orientation of the face of bdr element i. (3D)
Definition: mesh.cpp:6212
Vector J
Jacobians of the element transformations at all quadrature points.
Definition: mesh.hpp:1956
FaceType
Definition: mesh.hpp:45
long long GetGlobalNE() const
Return the total (global) number of elements.
Definition: mesh.hpp:968
void MoveNodes(const Vector &displacements)
Definition: mesh.cpp:7937
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
Definition: mesh.cpp:2445
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
Definition: mesh.cpp:6161
NCFaceInfo(bool slave, int master, const DenseMatrix *pm)
Definition: mesh.hpp:213
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a triangular Mesh.
Definition: mesh.cpp:1951
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
void UpdateNURBS()
Definition: mesh.cpp:5137
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 GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:7878
void PrintVTK(std::ostream &os)
Definition: mesh.cpp:10365
void FinalizeTetMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a tetrahedral Mesh.
Definition: mesh.cpp:2789
void AddBdrQuadAsTriangles(const int *vi, int attr=1)
Definition: mesh.cpp:1877
void Make3D(int nx, int ny, int nz, Element::Type type, double sx, double sy, double sz, bool sfc_ordering)
Definition: mesh.cpp:3070
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
Definition: mesh.cpp:6522
virtual void SetAttributes()
Definition: mesh.cpp:1564
bool IsInterior() const
return true if the face is an interior face to the computation domain, either a local or shared inter...
Definition: mesh.hpp:1452
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:9535
GeometryList(Mesh &mesh, int dim)
Construct a GeometryList of all geometries of dimension dim in mesh.
Definition: mesh.hpp:1142
void GetNode(int i, double *coord) const
Definition: mesh.cpp:7898
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition: mesh.hpp:1962
void AddHexAsPyramids(const int *vi, int attr=1)
Definition: mesh.cpp:1802
int FindCoarseElement(int i)
Definition: mesh.cpp:9980
GeometryList(Mesh &mesh)
Construct a GeometryList of all element geometries in mesh.
Definition: mesh.hpp:1137
double b
Definition: lissajous.cpp:42
int AddTet(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1700
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Definition: mesh.cpp:497
virtual void PrintXG(std::ostream &os=mfem::out) const
Print the mesh to the given stream using Netgen/Truegrid format.
Definition: mesh.cpp:10053
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:656
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9878
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
Definition: mesh.cpp:6587
void CheckPartitioning(int *partitioning_)
Definition: mesh.cpp:7474
void KnotInsert(Array< KnotVector *> &kv)
Definition: mesh.cpp:5063
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
Definition: mesh.cpp:5202
void ReadGmshMesh(std::istream &input, int &curved, int &read_gf)
STable3D * GetElementToFaceTable(int ret_ftbl=0)
Definition: mesh.cpp:6839
bool Conforming() const
Definition: mesh.hpp:1728
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition: mesh.cpp:5356
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:6029
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:929
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:1166
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
void InitTables()
Definition: mesh.cpp:1466
void CheckDisplacements(const Vector &displacements, double &tmax)
Definition: mesh.cpp:7792
This structure stores the low level information necessary to interpret the configuration of elements ...
Definition: mesh.hpp:153
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1053
VTKFormat
Data array format for VTK and VTU files.
Definition: vtk.hpp:96
void GetElementData(int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.hpp:1060
int mesh_geoms
Definition: mesh.hpp:78
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1835
Geometry::Constants< Geometry::TETRAHEDRON > tet_t
Definition: mesh.hpp:265
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:1148
int AddElement(Element *elem)
The parameter elem should be allocated using the NewElement() method.
Definition: mesh.cpp:1821
void ReadLineMesh(std::istream &input)
virtual void Load(std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
Definition: mesh.hpp:904
void SetNodesOwner(bool nodes_owner)
Set the mesh nodes ownership flag.
Definition: mesh.hpp:1615
void SetLayer(const int l)
Definition: mesh.hpp:2028
Table * GetFaceToElementTable() const
Definition: mesh.cpp:6127
Array< int > FindFaceNeighbors(const int elem) const
Returns the sorted, unique indices of elements sharing a face with element elem, including elem...
Definition: mesh.cpp:6184
int GetElementToEdgeTable(Table &, Array< int > &)
Definition: mesh.cpp:6397
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:2209
int nbBoundaryFaces
Definition: mesh.hpp:75
void RandomRefinement(double prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
Refine each element with given probability. Uses GeneralRefinement.
Definition: mesh.cpp:9559
void GetBdrElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of bdr element i.
Definition: mesh.cpp:5967
void GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:771
void ResetLazyData()
Definition: mesh.cpp:1553
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:1400
static const int vtk_quadratic_pyramid[13]
Definition: mesh.hpp:249
bool OwnsNodes() const
Return the mesh nodes ownership flag.
Definition: mesh.hpp:1613
void GetElementCenter(int i, Vector &center)
Definition: mesh.cpp:68
MFEM_DEPRECATED Geometry::Type GetFaceBaseGeometry(int i) const
Deprecated in favor of Mesh::GetFaceGeometry.
Definition: mesh.hpp:1105
IsoparametricTransformation BdrTransformation
Definition: mesh.hpp:235
void GetLocalFaceTransformation(int face_type, int elem_type, IsoparametricTransformation &Transf, int info)
A helper method that constructs a transformation from the reference space of a face to the reference ...
Definition: mesh.cpp:896
IsoparametricTransformation Transformation2
Definition: mesh.hpp:234
void Init()
Definition: mesh.cpp:1448
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5298
Element * GetBdrElement(int i)
Definition: mesh.hpp:1087
Vector normal
Normals at all quadrature points.
Definition: mesh.hpp:2016
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
Definition: ncmesh.hpp:121
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:1095
bool IsOfFaceType(FaceType type) const
Return true if the face is of the same type as type.
Definition: mesh.hpp:1465
Geometry::Constants< Geometry::PRISM > pri_t
Definition: mesh.hpp:267
struct mfem::Mesh::FaceInformation::@13 element[2]
virtual MFEM_DEPRECATED void ReorientTetMesh()
Definition: mesh.cpp:6968
prob_type prob
Definition: ex25.cpp:153
int NumOfBdrElements
Definition: mesh.hpp:69
Table * el_to_edge
Definition: mesh.hpp:220
int GetBdrFace(int BdrElemNo) const
Return the local face index for the given boundary face.
Definition: mesh.cpp:1113
void ReadMFEMMesh(std::istream &input, int version, int &curved)
virtual void LocalRefinement(const Array< int > &marked_el, int type=3)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:8913
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
Definition: mesh.cpp:6038
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
virtual void MarkTetMeshForRefinement(DSTable &v_to_v)
Definition: mesh.cpp:2470
bool IsConforming() const
Return true if the face is a conforming face.
Definition: mesh.hpp:1479
Array< Triple< int, int, int > > tmp_vertex_parents
Definition: mesh.hpp:259
MFEM_DEPRECATED Mesh(int nx, int ny, Element::Type type, bool generate_edges=false, double sx=1.0, double sy=1.0, bool sfc_ordering=true)
Deprecated: see MakeCartesian2D.
Definition: mesh.hpp:864
List of mesh geometries stored as Array<Geometry::Type>.
Definition: mesh.hpp:1131
void Destroy()
Definition: mesh.cpp:1524
void PrintVTU(std::ostream &os, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr_elements=false)
Definition: mesh.cpp:10565
This structure is used as a human readable output format that decipheres the information contained in...
Definition: mesh.hpp:1417
int CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Definition: mesh.cpp:5810
Geometry::Type GetElementGeometry(int i) const
Definition: mesh.hpp:1094
GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, MemoryType d_mt=MemoryType::DEFAULT)
Definition: mesh.cpp:12077
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
Definition: mesh.cpp:11111
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1550
void GetLocalQuadToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:793
void AddPointFaceElement(int lf, int gf, int el)
Used in GenerateFaces()
Definition: mesh.cpp:6490
void Make1D(int n, double sx=1.0)
Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.
Definition: mesh.cpp:3529
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:9155
virtual ~Mesh()
Destroys Mesh.
Definition: mesh.hpp:1898
void FinalizeMesh(int refine=0, bool fix_orientation=true)
Finalize the construction of any type of Mesh.
Definition: mesh.cpp:2895
const Table & ElementToEdgeTable() const
Definition: mesh.cpp:6481
Table * el_to_el
Definition: mesh.hpp:222
int AddSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1658
bool FaceIsTrueInterior(int FaceNo) const
Definition: mesh.hpp:504
Geometry::Constants< Geometry::CUBE > hex_t
Definition: mesh.hpp:266
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
Geometry::Type geom_buf[Geometry::NumGeom]
Definition: mesh.hpp:1134
static Mesh LoadFromFile(const char *filename, int generate_edges=0, int refine=1, bool fix_orientation=true)
Definition: mesh.cpp:3708
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Definition: mesh.cpp:11742
static void PrintElement(const Element *, std::ostream &)
Definition: mesh.cpp:3920
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
Definition: mesh.cpp:5107
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
int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
Definition: mesh.cpp:1749
const IntegrationRule * IntRule
Definition: mesh.hpp:1924
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition: mesh.hpp:1994
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2902
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
Definition: mesh.hpp:1610
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:684
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
Definition: mesh.cpp:2493
int AddBdrQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1863
Array< Vertex > vertices
Definition: mesh.hpp:90
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
Definition: mesh.cpp:10309
static Mesh MakeCartesian2D(int nx, int ny, Element::Type type, bool generate_edges=false, double sx=1.0, double sy=1.0, bool sfc_ordering=true)
Definition: mesh.cpp:3726
void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Definition: mesh.cpp:1595
void SetBdrAttribute(int i, int attr)
Set the attribute of boundary element i.
Definition: mesh.hpp:1553
void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo, const std::string &xml_prefix="")
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
Definition: mesh.hpp:1066
void RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
Definition: mesh.cpp:9578
static Mesh MakeCartesian1D(int n, double sx=1.0)
Definition: mesh.cpp:3718
Array< Element * > elements
Definition: mesh.hpp:85
const Table & ElementToElementTable()
Definition: mesh.cpp:6436
int SpaceDimension() const
Definition: mesh.hpp:1048
int meshgen
Definition: mesh.hpp:77
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:9990
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9604
Table * bel_to_edge
Definition: mesh.hpp:224
void GetGeometries(int dim, Array< Geometry::Type > &el_geoms) const
Return all element geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5932
Geometry::Type GetBdrElementGeometry(int i) const
Definition: mesh.hpp:1099
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:936
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:696
void GetFaceEdges(int i, Array< int > &edges, Array< int > &o) const
Definition: mesh.cpp:5999
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:277
int nbInteriorFaces
Definition: mesh.hpp:75
void PrintBdrVTU(std::string fname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0)
Definition: mesh.cpp:10557
int EulerNumber() const
Equals 1 + num_holes - num_loops.
Definition: mesh.hpp:1041
Array< int > be_to_edge
Definition: mesh.hpp:223
Class for integration point with weight.
Definition: intrules.hpp:25
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:6291
void Swap(Mesh &other, bool non_geometry)
Definition: mesh.cpp:9346
Element * ReadElementWithoutAttr(std::istream &)
Definition: mesh.cpp:3878
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
static const int vtk_quadratic_wedge[18]
Definition: mesh.hpp:250
virtual long long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
Definition: mesh.hpp:965
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:676
virtual int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type, does not count master nonconforming face...
Definition: mesh.cpp:5405
Table * face_edge
Definition: mesh.hpp:231
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1980
const Element * GetFace(int i) const
Definition: mesh.hpp:1089
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2995
bool DerefineByError(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
Definition: mesh.cpp:9273
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:2261
void GetLocalQuadToPyrTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:817
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:5120
int dim
Definition: ex24.cpp:53
Table * edge_vertex
Definition: mesh.hpp:232
virtual bool HasBoundaryElements() const
Checks if the mesh has boundary elements.
Definition: mesh.hpp:1064
Element::Type GetFaceElementType(int Face) const
Definition: mesh.cpp:1432
long sequence
Definition: mesh.hpp:83
virtual ~NodeExtrudeCoefficient()
Definition: mesh.hpp:2032
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
Definition: hooke.cpp:45
const Mesh * mesh
Definition: mesh.hpp:1923
void MakeSimplicial_(const Mesh &orig_mesh, int *vglobal)
Definition: mesh.cpp:4582
void DestroyPointers()
Definition: mesh.cpp:1498
Array< FaceInfo > faces_info
Definition: mesh.hpp:217
bool IsBoundary() const
Return true if the face is a boundary face.
Definition: mesh.hpp:1459
STable3D * GetFacesTable()
Definition: mesh.cpp:6776
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
Definition: sparsemat.hpp:711
NCMesh * ncmesh
Optional nonconforming mesh extension.
Definition: mesh.hpp:278
virtual void UniformRefinement3D()
Refine a mixed 3D mesh uniformly.
Definition: mesh.hpp:414
int AddBdrElement(Element *elem)
Definition: mesh.cpp:1828
RefCoord t[3]
static void PrintElementWithoutAttr(const Element *, std::ostream &)
Definition: mesh.cpp:3896
Mesh(int Dim_, int NVert, int NElem, int NBdrElem=0, int spaceDim_=-1)
Init constructor: begin the construction of a Mesh object.
Definition: mesh.hpp:704
int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr=1)
Definition: mesh.cpp:1721
FaceInformation GetFaceInformation(int f) const
Definition: mesh.cpp:1131
virtual bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NC version of GeneralDerefinement.
Definition: mesh.cpp:9225
int Dim
Definition: mesh.hpp:66
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:912
ElementConformity conformity
Definition: mesh.hpp:1424
Geometry::Constants< Geometry::TRIANGLE > tri_t
Definition: mesh.hpp:263
double AggregateError(const Array< double > &elem_error, const int *fine, int nfine, int op)
Derefinement helper.
Definition: mesh.cpp:9205
double * GetVertex(int i)
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1058
virtual int FindPoints(DenseMatrix &point_mat, Array< int > &elem_ids, Array< IntegrationPoint > &ips, bool warn=true, InverseElementTransformation *inv_trans=NULL)
Find the ids of the elements that contain the given points, and their corresponding reference coordin...
Definition: mesh.cpp:11962
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:1547
Mesh & operator=(Mesh &&mesh)
Move assignment operstor.
Definition: mesh.cpp:3702
MFEM_DEPRECATED Geometry::Type GetFaceGeometryType(int Face) const
Deprecated in favor of Mesh::GetFaceGeometry.
Definition: mesh.hpp:1516
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:945
Vector data type.
Definition: vector.hpp:60
void ReadTrueGridMesh(std::istream &input)
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1749
const FaceGeometricFactors * GetFaceGeometricFactors(const IntegrationRule &ir, const int flags, FaceType type, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors for the faces corresponding to the given integration rule...
Definition: mesh.cpp:860
static void GetElementArrayEdgeTable(const Array< Element *> &elem_array, const DSTable &v_to_v, Table &el_to_edge)
Definition: mesh.cpp:6350
bool IsShared() const
Return true if the face is a shared interior face which is NOT a master nonconforming face...
Definition: mesh.hpp:1443
bool IsNonconformingCoarse() const
Return true if the face is a nonconforming coarse face.
Definition: mesh.hpp:1496
virtual void Save(const char *fname, int precision=16) const
Definition: mesh.cpp:10351
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
Definition: mesh.cpp:5581
void SetNode(int i, const double *coord)
Definition: mesh.cpp:7917
Table * face_to_elem
Definition: mesh.hpp:230
void DestroyTables()
Definition: mesh.cpp:1479
void PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
Definition: mesh.cpp:11484
void GetLocalTriToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:720
void GenerateFaces()
Definition: mesh.cpp:6614
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:10917
void ChangeVertexDataOwnership(double *vertices, int len_vertices, bool zerocopy=false)
Set the internal Vertex array to point to the given vertices array without assuming ownership of the ...
Definition: mesh.cpp:3788
void DeleteGeometricFactors()
Destroy all GeometricFactors stored by the Mesh.
Definition: mesh.cpp:882
int GetNumFacesWithGhost() const
Return the number of faces (3D), edges (2D) or vertices (1D) including ghost faces.
Definition: mesh.cpp:5400
int spaceDim
Definition: mesh.hpp:67
IsoparametricTransformation EdgeTransformation
Definition: mesh.hpp:236
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:70
void Printer(std::ostream &out=mfem::out, std::string section_delimiter="") const
Definition: mesh.cpp:10217
MemAlloc< Tetrahedron, 1024 > TetMemory
Definition: mesh.hpp:255
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
Definition: mesh.cpp:2559
void SetNodes(const Vector &node_coord)
Updates the vertex/node locations. Invokes NodesUpdated().
Definition: mesh.cpp:7961
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:1085
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:5345
Operation GetLastOperation() const
Return type of last modification of the mesh.
Definition: mesh.hpp:1736
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
int NumOfVertices
Definition: mesh.hpp:69
void ScaleSubdomains(double sf)
Definition: mesh.cpp:11551
Array< int > be_to_face
Definition: mesh.hpp:225
int AddPyramid(int v1, int v2, int v3, int v4, int v5, int attr=1)
Definition: mesh.cpp:1735
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
Definition: mesh.cpp:9302
void Bisection(int i, const DSTable &, int *, int *, int *)
Bisect a triangle: element with index i is bisected.
Definition: mesh.cpp:9633
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9467
ElementLocation location
Definition: mesh.hpp:1423
Class for parallel meshes.
Definition: pmesh.hpp:32
Geometry::Type GetBdrElementBaseGeometry(int i) const
Definition: mesh.hpp:1111
Abstract data type element.
Definition: element.hpp:28
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
Definition: mesh.cpp:6559
Array< GeometricFactors * > geom_factors
Optional geometric factors.
Definition: mesh.hpp:279
void BdrBisection(int i, const HashTable< Hashed2 > &)
Bisect a boundary triangle: boundary element with index i is bisected.
Definition: mesh.cpp:9841
void GenerateNCFaceInfo()
Definition: mesh.cpp:6715
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
Definition: mesh.cpp:6335
double GetElementVolume(int i)
Definition: mesh.cpp:112
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:273
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
Definition: mesh.cpp:6092
void EnsureNodes()
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element gr...
Definition: mesh.cpp:5304
void PrintWithPartitioning(int *partitioning, std::ostream &os, int elem_attr=0) const
Prints the mesh with boundary elements given by the boundary of the subdomains, so that the boundary ...
Definition: mesh.cpp:10992
const DenseMatrix * point_matrix
Definition: mesh.hpp:1432
void UpdateNodes()
Update the nodes of a curved mesh after the topological part of a Mesh::Operation, such as refinement, has been performed.
Definition: mesh.cpp:8032
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:6319
void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const
For the given boundary element, bdr_el, return its adjacent element and its info, i...
Definition: mesh.cpp:6247
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:365
void NodesUpdated()
This function should be called after the mesh node coordinates have been updated externally, e.g. by modifying the internal nodal GridFunction returned by GetNodes().
Definition: mesh.hpp:1038
Class used to extrude the nodes of a mesh.
Definition: mesh.hpp:2020
void AddVertexParents(int i, int p1, int p2)
Mark vertex i as nonconforming, with parent vertices p1 and p2.
Definition: mesh.cpp:1642
int NumOfFaces
Definition: mesh.hpp:70