MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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  mutable Table *face_edge;
227  mutable Table *edge_vertex;
228 
233 
234  // refinement embeddings for forward compatibility with NCMesh
236 
237  // Nodes are only active for higher order meshes, and share locations with
238  // the vertices, plus all the higher- order control points within the
239  // element and along the edges and on the faces.
242 
243  static const int vtk_quadratic_tet[10];
244  static const int vtk_quadratic_pyramid[13];
245  static const int vtk_quadratic_wedge[18];
246  static const int vtk_quadratic_hex[27];
247 
248 #ifdef MFEM_USE_MEMALLOC
249  friend class Tetrahedron;
251 #endif
252 
253  // used during NC mesh initialization only
255 
256 public:
264 
266 
267  /// A list of all unique element attributes used by the Mesh.
269  /// A list of all unique boundary attributes used by the Mesh.
271 
272  NURBSExtension *NURBSext; ///< Optional NURBS mesh extension.
273  NCMesh *ncmesh; ///< Optional nonconforming mesh extension.
274  Array<GeometricFactors*> geom_factors; ///< Optional geometric factors.
276  face_geom_factors; ///< Optional face geometric factors.
277 
278  // Global parameter that can be used to control the removal of unused
279  // vertices performed when reading a mesh in MFEM format. The default value
280  // (true) is set in mesh_readers.cpp.
282 
283 protected:
285 
286  void Init();
287  void InitTables();
288  void SetEmpty(); // Init all data members with empty values
289  void DestroyTables();
291  void DestroyPointers(); // Delete data specifically allocated by class Mesh.
292  void Destroy(); // Delete all owned data.
293  void ResetLazyData();
294 
295  Element *ReadElementWithoutAttr(std::istream &);
296  static void PrintElementWithoutAttr(const Element *, std::ostream &);
297 
298  Element *ReadElement(std::istream &);
299  static void PrintElement(const Element *, std::ostream &);
300 
301  // Readers for different mesh formats, used in the Load() method.
302  // The implementations of these methods are in mesh_readers.cpp.
303  void ReadMFEMMesh(std::istream &input, int version, int &curved);
304  void ReadLineMesh(std::istream &input);
305  void ReadNetgen2DMesh(std::istream &input, int &curved);
306  void ReadNetgen3DMesh(std::istream &input);
307  void ReadTrueGridMesh(std::istream &input);
308  void CreateVTKMesh(const Vector &points, const Array<int> &cell_data,
309  const Array<int> &cell_offsets,
310  const Array<int> &cell_types,
311  const Array<int> &cell_attributes,
312  int &curved, int &read_gf, bool &finalize_topo);
313  void ReadVTKMesh(std::istream &input, int &curved, int &read_gf,
314  bool &finalize_topo);
315  void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf,
316  bool &finalize_topo, const std::string &xml_prefix="");
317  void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf);
318  void ReadInlineMesh(std::istream &input, bool generate_edges = false);
319  void ReadGmshMesh(std::istream &input, int &curved, int &read_gf);
320  /* Note NetCDF (optional library) is used for reading cubit files */
321 #ifdef MFEM_USE_NETCDF
322  void ReadCubit(const char *filename, int &curved, int &read_gf);
323 #endif
324 
325  /// Determine the mesh generator bitmask #meshgen, see MeshGenerator().
326  /** Also, initializes #mesh_geoms. */
327  void SetMeshGen();
328 
329  /// Return the length of the segment from node i to node j.
330  double GetLength(int i, int j) const;
331 
332  /** Compute the Jacobian of the transformation from the perfect
333  reference element at the center of the element. */
334  void GetElementJacobian(int i, DenseMatrix &J);
335 
336  void MarkForRefinement();
338  void GetEdgeOrdering(DSTable &v_to_v, Array<int> &order);
339  virtual void MarkTetMeshForRefinement(DSTable &v_to_v);
340 
341  // Methods used to prepare and apply permutation of the mesh nodes assuming
342  // that the mesh elements may be rotated (e.g. to mark triangle or tet edges
343  // for refinement) between the two calls - PrepareNodeReorder() and
344  // DoNodeReorder(). The latter method assumes that the 'faces' have not been
345  // updated after the element rotations.
346  void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert);
347  void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert);
348 
350  STable3D *GetElementToFaceTable(int ret_ftbl = 0);
351 
352  /** Red refinement. Element with index i is refined. The default
353  red refinement for now is Uniform. */
354  void RedRefinement(int i, const DSTable &v_to_v,
355  int *edge1, int *edge2, int *middle)
356  { UniformRefinement(i, v_to_v, edge1, edge2, middle); }
357 
358  /** Green refinement. Element with index i is refined. The default
359  refinement for now is Bisection. */
360  void GreenRefinement(int i, const DSTable &v_to_v,
361  int *edge1, int *edge2, int *middle)
362  { Bisection(i, v_to_v, edge1, edge2, middle); }
363 
364  /// Bisect a triangle: element with index @a i is bisected.
365  void Bisection(int i, const DSTable &, int *, int *, int *);
366 
367  /// Bisect a tetrahedron: element with index @a i is bisected.
368  void Bisection(int i, HashTable<Hashed2> &);
369 
370  /// Bisect a boundary triangle: boundary element with index @a i is bisected.
371  void BdrBisection(int i, const HashTable<Hashed2> &);
372 
373  /** Uniform Refinement. Element with index i is refined uniformly. */
374  void UniformRefinement(int i, const DSTable &, int *, int *, int *);
375 
376  /** @brief Averages the vertices with given @a indexes and saves the result
377  in #vertices[result]. */
378  void AverageVertices(const int *indexes, int n, int result);
379 
381  int FindCoarseElement(int i);
382 
383  /// Update the nodes of a curved mesh after refinement
384  void UpdateNodes();
385 
386  /// Helper to set vertex coordinates given a high-order curvature function.
387  void SetVerticesFromNodes(const GridFunction *nodes);
388 
389  void UniformRefinement2D_base(bool update_nodes = true);
390 
391  /// Refine a mixed 2D mesh uniformly.
393 
394  /* If @a f2qf is not NULL, adds all quadrilateral faces to @a f2qf which
395  represents a "face-to-quad-face" index map. When all faces are quads, the
396  array @a f2qf is kept empty since it is not needed. */
397  void UniformRefinement3D_base(Array<int> *f2qf = NULL,
398  DSTable *v_to_v_p = NULL,
399  bool update_nodes = true);
400 
401  /// Refine a mixed 3D mesh uniformly.
403 
404  /// Refine NURBS mesh.
405  virtual void NURBSUniformRefinement();
406 
407  /// This function is not public anymore. Use GeneralRefinement instead.
408  virtual void LocalRefinement(const Array<int> &marked_el, int type = 3);
409 
410  /// This function is not public anymore. Use GeneralRefinement instead.
411  virtual void NonconformingRefinement(const Array<Refinement> &refinements,
412  int nc_limit = 0);
413 
414  /// NC version of GeneralDerefinement.
415  virtual bool NonconformingDerefinement(Array<double> &elem_error,
416  double threshold, int nc_limit = 0,
417  int op = 1);
418  /// Derefinement helper.
419  double AggregateError(const Array<double> &elem_error,
420  const int *fine, int nfine, int op);
421 
422  /// Read NURBS patch/macro-element mesh
423  void LoadPatchTopo(std::istream &input, Array<int> &edge_to_knot);
424 
425  void UpdateNURBS();
426 
427  void PrintTopo(std::ostream &out, const Array<int> &e_to_k) const;
428 
429  /// Used in GetFaceElementTransformations (...)
432  int i);
434  int i);
435  /// Used in GetFaceElementTransformations (...)
437  int i);
438  /// Used in GetFaceElementTransformations (...)
440  int i);
441  /// Used in GetFaceElementTransformations (...)
443  int i);
444  /// Used in GetFaceElementTransformations (...)
446  int i);
447  /// Used in GetFaceElementTransformations (...)
449  int i);
450  /// Used in GetFaceElementTransformations (...)
452  int i);
453 
454  /** Used in GetFaceElementTransformations to account for the fact that a
455  slave face occupies only a portion of its master face. */
457  const FaceInfo &fi, bool is_ghost);
458 
459  bool IsSlaveFace(const FaceInfo &fi) const;
460 
461  /// Returns the orientation of "test" relative to "base"
462  static int GetTriOrientation (const int * base, const int * test);
463  /// Returns the orientation of "test" relative to "base"
464  static int GetQuadOrientation (const int * base, const int * test);
465  /// Returns the orientation of "test" relative to "base"
466  static int GetTetOrientation (const int * base, const int * test);
467 
468  static void GetElementArrayEdgeTable(const Array<Element*> &elem_array,
469  const DSTable &v_to_v,
470  Table &el_to_edge);
471 
472  /** Return vertex to vertex table. The connections stored in the table
473  are from smaller to bigger vertex index, i.e. if i<j and (i, j) is
474  in the table, then (j, i) is not stored. */
475  void GetVertexToVertexTable(DSTable &) const;
476 
477  /** Return element to edge table and the indices for the boundary edges.
478  The entries in the table are ordered according to the order of the
479  nodes in the elements. For example, if T is the element to edge table
480  T(i, 0) gives the index of edge in element i that connects vertex 0
481  to vertex 1, etc. Returns the number of the edges. */
483 
484  /// Used in GenerateFaces()
485  void AddPointFaceElement(int lf, int gf, int el);
486 
487  void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1);
488 
489  void AddTriangleFaceElement (int lf, int gf, int el,
490  int v0, int v1, int v2);
491 
492  void AddQuadFaceElement (int lf, int gf, int el,
493  int v0, int v1, int v2, int v3);
494  /** For a serial Mesh, return true if the face is interior. For a parallel
495  ParMesh return true if the face is interior or shared. In parallel, this
496  method only works if the face neighbor data is exchanged. */
497  bool FaceIsTrueInterior(int FaceNo) const
498  {
499  return FaceIsInterior(FaceNo) || (faces_info[FaceNo].Elem2Inf >= 0);
500  }
501 
502  void FreeElement(Element *E);
503 
504  void GenerateFaces();
505  void GenerateNCFaceInfo();
506 
507  /// Begin construction of a mesh
508  void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem);
509 
510  // Used in the methods FinalizeXXXMesh() and FinalizeTopology()
511  void FinalizeCheck();
512 
513  void Loader(std::istream &input, int generate_edges = 0,
514  std::string parse_tag = "");
515 
516  // If NURBS mesh, write NURBS format. If NCMesh, write mfem v1.1 format.
517  // If section_delimiter is empty, write mfem v1.0 format. Otherwise, write
518  // mfem v1.2 format with the given section_delimiter at the end.
519  void Printer(std::ostream &out = mfem::out,
520  std::string section_delimiter = "") const;
521 
522  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
523  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
524  type=TETRAHEDRON. The parameter @a sfc_ordering controls how the elements
525  (when type=HEXAHEDRON) are ordered: true - use space-filling curve
526  ordering, or false - use lexicographic ordering. */
527  void Make3D(int nx, int ny, int nz, Element::Type type,
528  double sx, double sy, double sz, bool sfc_ordering);
529 
530  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
531  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
532  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
533  if 1 edges are generated. The parameter @a sfc_ordering controls how the
534  elements (when type=QUADRILATERAL) are ordered: true - use space-filling
535  curve ordering, or false - use lexicographic ordering. */
536  void Make2D(int nx, int ny, Element::Type type, double sx, double sy,
537  bool generate_edges, bool sfc_ordering);
538 
539  /// Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.
540  void Make1D(int n, double sx = 1.0);
541 
542  /// Internal function used in Mesh::MakeRefined
543  void MakeRefined_(Mesh &orig_mesh, const Array<int> ref_factors,
544  int ref_type);
545 
546  /// Initialize vertices/elements/boundary/tables from a nonconforming mesh.
547  void InitFromNCMesh(const NCMesh &ncmesh);
548 
549  /// Create from a nonconforming mesh.
550  explicit Mesh(const NCMesh &ncmesh);
551 
552  // used in GetElementData() and GetBdrElementData()
553  void GetElementData(const Array<Element*> &elem_array, int geom,
554  Array<int> &elem_vtx, Array<int> &attr) const;
555 
556  double GetElementSize(ElementTransformation *T, int type = 0);
557 
558  // Internal helper used in MakeSimplicial (and ParMesh::MakeSimplicial).
559  void MakeSimplicial_(const Mesh &orig_mesh, int *vglobal);
560 
561 public:
562 
563  Mesh() { SetEmpty(); }
564 
565  /** Copy constructor. Performs a deep copy of (almost) all data, so that the
566  source mesh can be modified (e.g. deleted, refined) without affecting the
567  new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the
568  nodes, if present. */
569  explicit Mesh(const Mesh &mesh, bool copy_nodes = true);
570 
571  /// Move constructor, useful for using a Mesh as a function return value.
572  Mesh(Mesh &&mesh);
573 
574  /// Move assignment operstor.
575  Mesh& operator=(Mesh &&mesh);
576 
577  /// Explicitly delete the copy assignment operator.
578  Mesh& operator=(const Mesh &mesh) = delete;
579 
580  /** @name Named mesh constructors.
581 
582  Each of these constructors uses the move constructor, and can be used as
583  the right-hand side of an assignment when creating new meshes. */
584  ///@{
585 
586  /** Creates mesh by reading a file in MFEM, Netgen, or VTK format. If
587  generate_edges = 0 (default) edges are not generated, if 1 edges are
588  generated. */
589  static Mesh LoadFromFile(const char *filename,
590  int generate_edges = 0, int refine = 1,
591  bool fix_orientation = true);
592 
593  /** Creates 1D mesh , divided into n equal intervals. */
594  static Mesh MakeCartesian1D(int n, double sx = 1.0);
595 
596  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
597  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
598  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
599  if 1 edges are generated. If scf_ordering = true (default), elements are
600  ordered along a space-filling curve, instead of row by row. */
601  static Mesh MakeCartesian2D(
602  int nx, int ny, Element::Type type, bool generate_edges = false,
603  double sx = 1.0, double sy = 1.0, bool sfc_ordering = true);
604 
605  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
606  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
607  type=TETRAHEDRON. If sfc_ordering = true (default), elements are ordered
608  along a space-filling curve, instead of row by row and layer by layer. */
609  static Mesh MakeCartesian3D(
610  int nx, int ny, int nz, Element::Type type,
611  double sx = 1.0, double sy = 1.0, double sz = 1.0,
612  bool sfc_ordering = true);
613 
614  /// Create a refined (by any factor) version of @a orig_mesh.
615  /** @param[in] orig_mesh The starting coarse mesh.
616  @param[in] ref_factor The refinement factor, an integer > 1.
617  @param[in] ref_type Specify the positions of the new vertices. The
618  options are BasisType::ClosedUniform or
619  BasisType::GaussLobatto.
620 
621  The refinement data which can be accessed with GetRefinementTransforms()
622  is set to reflect the performed refinements.
623 
624  @note The constructed Mesh is straight-sided. */
625  static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type);
626 
627  /// Create a refined mesh, where each element of the original mesh may be
628  /// refined by a different factor.
629  /** @param[in] orig_mesh The starting coarse mesh.
630  @param[in] ref_factors An array of integers whose size is the number of
631  elements of @a orig_mesh. The @a ith element of
632  @a orig_mesh is refined by refinement factor
633  @a ref_factors[i].
634  @param[in] ref_type Specify the positions of the new vertices. The
635  options are BasisType::ClosedUniform or
636  BasisType::GaussLobatto.
637 
638  The refinement data which can be accessed with GetRefinementTransforms()
639  is set to reflect the performed refinements.
640 
641  @note The constructed Mesh is straight-sided. */
642  /// refined @a ref_factors[i] times in each dimension.
643  static Mesh MakeRefined(Mesh &orig_mesh, const Array<int> &ref_factors,
644  int ref_type);
645 
646  /** Create a mesh by splitting each element of @a orig_mesh into simplices.
647  Quadrilaterals are split into two triangles, prisms are split into
648  3 tetrahedra, and hexahedra are split into either 5 or 6 tetrahedra
649  depending on the configuration.
650  @warning The curvature of the original mesh is not carried over to the
651  new mesh. Periodic meshes are not supported. */
652  static Mesh MakeSimplicial(const Mesh &orig_mesh);
653 
654  /// Create a periodic mesh by identifying vertices of @a orig_mesh.
655  /** Each vertex @a i will be mapped to vertex @a v2v[i], such that all
656  vertices that are coincident under the periodic mapping get mapped to
657  the same index. The mapping @a v2v can be generated from translation
658  vectors using Mesh::CreatePeriodicVertexMapping.
659  @note MFEM requires that each edge of the resulting mesh be uniquely
660  identifiable by a pair of distinct vertices. As a consequence, periodic
661  boundaries must be connected by at least three edges. */
662  static Mesh MakePeriodic(const Mesh &orig_mesh, const std::vector<int> &v2v);
663 
664  ///@}
665 
666  /// @brief Creates a mapping @a v2v from the vertex indices of the mesh such
667  /// that coincident vertices under the given @a translations are identified.
668  /** Each Vector in @a translations should be of size @a sdim (the spatial
669  dimension of the mesh). Two vertices are considered coincident if the
670  translated coordinates of one vertex are within the given tolerance (@a
671  tol, relative to the mesh diameter) of the coordinates of the other
672  vertex.
673  @warning This algorithm does not scale well with the number of boundary
674  vertices in the mesh, and may run slowly on very large meshes. */
675  std::vector<int> CreatePeriodicVertexMapping(
676  const std::vector<Vector> &translations, double tol = 1e-8) const;
677 
678  /// Construct a Mesh from the given primary data.
679  /** The array @a vertices is used as external data, i.e. the Mesh does not
680  copy the data and will not delete the pointer.
681 
682  The data from the other arrays is copied into the internal Mesh data
683  structures.
684 
685  This method calls the method FinalizeTopology(). The method Finalize()
686  may be called after this constructor and after optionally setting the
687  Mesh nodes. */
688  Mesh(double *vertices, int num_vertices,
689  int *element_indices, Geometry::Type element_type,
690  int *element_attributes, int num_elements,
691  int *boundary_indices, Geometry::Type boundary_type,
692  int *boundary_attributes, int num_boundary_elements,
693  int dimension, int space_dimension = -1);
694 
695  /** @anchor mfem_Mesh_init_ctor
696  @brief _Init_ constructor: begin the construction of a Mesh object. */
697  Mesh(int Dim_, int NVert, int NElem, int NBdrElem = 0, int spaceDim_ = -1)
698  {
699  if (spaceDim_ == -1) { spaceDim_ = Dim_; }
700  InitMesh(Dim_, spaceDim_, NVert, NElem, NBdrElem);
701  }
702 
703  /** @name Methods for Mesh construction.
704 
705  These methods are intended to be used with the @ref mfem_Mesh_init_ctor
706  "init constructor". */
707  ///@{
708 
709  Element *NewElement(int geom);
710 
711  int AddVertex(double x, double y = 0.0, double z = 0.0);
712  int AddVertex(const double *coords);
713  /// Mark vertex @a i as nonconforming, with parent vertices @a p1 and @a p2.
714  void AddVertexParents(int i, int p1, int p2);
715 
716  int AddSegment(int v1, int v2, int attr = 1);
717  int AddSegment(const int *vi, int attr = 1);
718 
719  int AddTriangle(int v1, int v2, int v3, int attr = 1);
720  int AddTriangle(const int *vi, int attr = 1);
721  int AddTri(const int *vi, int attr = 1) { return AddTriangle(vi, attr); }
722 
723  int AddQuad(int v1, int v2, int v3, int v4, int attr = 1);
724  int AddQuad(const int *vi, int attr = 1);
725 
726  int AddTet(int v1, int v2, int v3, int v4, int attr = 1);
727  int AddTet(const int *vi, int attr = 1);
728 
729  int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr = 1);
730  int AddWedge(const int *vi, int attr = 1);
731 
732  int AddPyramid(int v1, int v2, int v3, int v4, int v5, int attr = 1);
733  int AddPyramid(const int *vi, int attr = 1);
734 
735  int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8,
736  int attr = 1);
737  int AddHex(const int *vi, int attr = 1);
738  void AddHexAsTets(const int *vi, int attr = 1);
739  void AddHexAsWedges(const int *vi, int attr = 1);
740  void AddHexAsPyramids(const int *vi, int attr = 1);
741 
742  /// The parameter @a elem should be allocated using the NewElement() method
743  int AddElement(Element *elem);
744  int AddBdrElement(Element *elem);
745 
746  int AddBdrSegment(int v1, int v2, int attr = 1);
747  int AddBdrSegment(const int *vi, int attr = 1);
748 
749  int AddBdrTriangle(int v1, int v2, int v3, int attr = 1);
750  int AddBdrTriangle(const int *vi, int attr = 1);
751 
752  int AddBdrQuad(int v1, int v2, int v3, int v4, int attr = 1);
753  int AddBdrQuad(const int *vi, int attr = 1);
754  void AddBdrQuadAsTriangles(const int *vi, int attr = 1);
755 
756  int AddBdrPoint(int v, int attr = 1);
757 
759  /// Finalize the construction of a triangular Mesh.
760  void FinalizeTriMesh(int generate_edges = 0, int refine = 0,
761  bool fix_orientation = true);
762  /// Finalize the construction of a quadrilateral Mesh.
763  void FinalizeQuadMesh(int generate_edges = 0, int refine = 0,
764  bool fix_orientation = true);
765  /// Finalize the construction of a tetrahedral Mesh.
766  void FinalizeTetMesh(int generate_edges = 0, int refine = 0,
767  bool fix_orientation = true);
768  /// Finalize the construction of a wedge Mesh.
769  void FinalizeWedgeMesh(int generate_edges = 0, int refine = 0,
770  bool fix_orientation = true);
771  /// Finalize the construction of a hexahedral Mesh.
772  void FinalizeHexMesh(int generate_edges = 0, int refine = 0,
773  bool fix_orientation = true);
774  /// Finalize the construction of any type of Mesh.
775  /** This method calls FinalizeTopology() and Finalize(). */
776  void FinalizeMesh(int refine = 0, bool fix_orientation = true);
777 
778  ///@}
779 
780  /** @brief Finalize the construction of the secondary topology (connectivity)
781  data of a Mesh. */
782  /** This method does not require any actual coordinate data (either vertex
783  coordinates for linear meshes or node coordinates for meshes with nodes)
784  to be available. However, the data generated by this method is generally
785  required by the FiniteElementSpace class.
786 
787  After calling this method, setting the Mesh vertices or nodes, it may be
788  appropriate to call the method Finalize(). */
789  void FinalizeTopology(bool generate_bdr = true);
790 
791  /// Finalize the construction of a general Mesh.
792  /** This method will:
793  - check and optionally fix the orientation of regular elements
794  - check and fix the orientation of boundary elements
795  - assume that #vertices are defined, if #Nodes == NULL
796  - assume that #Nodes are defined, if #Nodes != NULL.
797  @param[in] refine If true, prepare the Mesh for conforming refinement of
798  triangular or tetrahedral meshes.
799  @param[in] fix_orientation
800  If true, fix the orientation of inverted mesh elements
801  by permuting their vertices.
802 
803  Before calling this method, call FinalizeTopology() and ensure that the
804  Mesh vertices or nodes are set. */
805  virtual void Finalize(bool refine = false, bool fix_orientation = false);
806 
807  virtual void SetAttributes();
808 
809  /** This is our integration with the Gecko library. The method finds an
810  element ordering that will increase memory coherency by putting elements
811  that are in physical proximity closer in memory. It can also be used to
812  obtain a space-filling curve ordering for ParNCMesh partitioning.
813  @param[out] ordering Output element ordering.
814  @param iterations Total number of V cycles. The ordering may improve with
815  more iterations. The best iteration is returned at the end.
816  @param window Initial window size. This determines the number of
817  permutations tested at each multigrid level and strongly influences the
818  quality of the result, but the cost of increasing 'window' is exponential.
819  @param period The window size is incremented every 'period' iterations.
820  @param seed Seed for initial random ordering (0 = skip random reorder).
821  @param verbose Print the progress of the optimization to mfem::out.
822  @param time_limit Optional time limit for the optimization, in seconds.
823  When reached, ordering from the best iteration so far is returned
824  (0 = no limit).
825  @return The final edge product cost of the ordering. The function may be
826  called in an external loop with different seeds, and the best ordering can
827  then be retained. */
828  double GetGeckoElementOrdering(Array<int> &ordering,
829  int iterations = 4, int window = 4,
830  int period = 2, int seed = 0,
831  bool verbose = false, double time_limit = 0);
832 
833  /** Return an ordering of the elements that approximately follows the Hilbert
834  curve. The method performs a spatial (Hilbert) sort on the centers of all
835  elements and returns the resulting sequence, which can then be passed to
836  ReorderElements. This is a cheap alternative to GetGeckoElementOrdering.*/
837  void GetHilbertElementOrdering(Array<int> &ordering);
838 
839  /** Rebuilds the mesh with a different order of elements. For each element i,
840  the array ordering[i] contains its desired new index. Note that the method
841  reorders vertices, edges and faces along with the elements. */
842  void ReorderElements(const Array<int> &ordering, bool reorder_vertices = true);
843 
844  /// Deprecated: see @a MakeCartesian3D.
845  MFEM_DEPRECATED
846  Mesh(int nx, int ny, int nz, Element::Type type, bool generate_edges = false,
847  double sx = 1.0, double sy = 1.0, double sz = 1.0,
848  bool sfc_ordering = true)
849  {
850  Make3D(nx, ny, nz, type, sx, sy, sz, sfc_ordering);
851  Finalize(true); // refine = true
852  }
853 
854  /// Deprecated: see @a MakeCartesian2D.
855  MFEM_DEPRECATED
856  Mesh(int nx, int ny, Element::Type type, bool generate_edges = false,
857  double sx = 1.0, double sy = 1.0, bool sfc_ordering = true)
858  {
859  Make2D(nx, ny, type, sx, sy, generate_edges, sfc_ordering);
860  Finalize(true); // refine = true
861  }
862 
863  /// Deprecated: see @a MakeCartesian1D.
864  MFEM_DEPRECATED
865  explicit Mesh(int n, double sx = 1.0)
866  {
867  Make1D(n, sx);
868  // Finalize(); // reminder: not needed
869  }
870 
871  /** Creates mesh by reading a file in MFEM, Netgen, or VTK format. If
872  generate_edges = 0 (default) edges are not generated, if 1 edges are
873  generated. See also @a Mesh::LoadFromFile. */
874  explicit Mesh(const char *filename, int generate_edges = 0, int refine = 1,
875  bool fix_orientation = true);
876 
877  /** Creates mesh by reading data stream in MFEM, Netgen, or VTK format. If
878  generate_edges = 0 (default) edges are not generated, if 1 edges are
879  generated. */
880  explicit Mesh(std::istream &input, int generate_edges = 0, int refine = 1,
881  bool fix_orientation = true);
882 
883  /// Create a disjoint mesh from the given mesh array
884  Mesh(Mesh *mesh_array[], int num_pieces);
885 
886  /// Deprecated: see @a MakeRefined.
887  MFEM_DEPRECATED
888  Mesh(Mesh *orig_mesh, int ref_factor, int ref_type);
889 
890  /** This is similar to the mesh constructor with the same arguments, but here
891  the current mesh is destroyed and another one created based on the data
892  stream again given in MFEM, Netgen, or VTK format. If generate_edges = 0
893  (default) edges are not generated, if 1 edges are generated. */
894  /// \see mfem::ifgzstream() for on-the-fly decompression of compressed ascii
895  /// inputs.
896  virtual void Load(std::istream &input, int generate_edges = 0,
897  int refine = 1, bool fix_orientation = true)
898  {
899  Loader(input, generate_edges);
900  Finalize(refine, fix_orientation);
901  }
902 
903  /// Clear the contents of the Mesh.
904  void Clear() { Destroy(); SetEmpty(); }
905 
906  /** @brief Get the mesh generator/type.
907 
908  The purpose of this is to be able to quickly tell what type of elements
909  one has in the mesh. Examination of this bitmask along with knowledge
910  of the mesh dimension can be used to identify which element types are
911  present.
912 
913  @return A bitmask:
914  - bit 0 - simplices are present in the mesh (triangles, tets),
915  - bit 1 - tensor product elements are present in the mesh (quads, hexes),
916  - bit 2 - the mesh has wedge elements.
917  - bit 3 - the mesh has pyramid elements.
918 
919  In parallel, the result takes into account elements on all processors.
920  */
921  inline int MeshGenerator() { return meshgen; }
922 
923  /** @brief Returns number of vertices. Vertices are only at the corners of
924  elements, where you would expect them in the lowest-order mesh. */
925  inline int GetNV() const { return NumOfVertices; }
926 
927  /// Returns number of elements.
928  inline int GetNE() const { return NumOfElements; }
929 
930  /// Returns number of boundary elements.
931  inline int GetNBE() const { return NumOfBdrElements; }
932 
933  /// Return the number of edges.
934  inline int GetNEdges() const { return NumOfEdges; }
935 
936  /// Return the number of faces in a 3D mesh.
937  inline int GetNFaces() const { return NumOfFaces; }
938 
939  /// Return the number of faces (3D), edges (2D) or vertices (1D).
940  int GetNumFaces() const;
941 
942  /** @brief Return the number of faces (3D), edges (2D) or vertices (1D)
943  including ghost faces. */
944  int GetNumFacesWithGhost() const;
945 
946  /** @brief Returns the number of faces according to the requested type, does
947  not count master nonconforming faces.
948 
949  If type==Boundary returns only the number of true boundary faces
950  contrary to GetNBE() that returns all "boundary" elements which may
951  include actual interior faces.
952  Similarly, if type==Interior, only the true interior faces are counted
953  excluding all master nonconforming faces. */
954  virtual int GetNFbyType(FaceType type) const;
955 
956  /// Utility function: sum integers from all processors (Allreduce).
957  virtual long ReduceInt(int value) const { return value; }
958 
959  /// Return the total (global) number of elements.
960  long GetGlobalNE() const { return ReduceInt(NumOfElements); }
961 
962  /** @brief Return the mesh geometric factors corresponding to the given
963  integration rule.
964 
965  The IntegrationRule used with GetGeometricFactors needs to remain valid
966  until the internally stored GeometricFactors objects are destroyed (by
967  either calling Mesh::DeleteGeometricFactors or the Mesh destructor). If
968  the device MemoryType parameter @a d_mt is specified, then the returned
969  object will use that type unless it was previously allocated with a
970  different type. */
972  const IntegrationRule& ir,
973  const int flags,
975 
976  /** @brief Return the mesh geometric factors for the faces corresponding
977  to the given integration rule.
978 
979  The IntegrationRule used with GetFaceGeometricFactors needs to remain
980  valid until the internally stored FaceGeometricFactors objects are
981  destroyed (by either calling Mesh::DeleteGeometricFactors or the Mesh
982  destructor). */
984  const int flags,
985  FaceType type);
986 
987  /// Destroy all GeometricFactors stored by the Mesh.
988  /** This method can be used to force recomputation of the GeometricFactors,
989  for example, after the mesh nodes are modified externally. */
990  void DeleteGeometricFactors();
991 
992  /// Equals 1 + num_holes - num_loops
993  inline int EulerNumber() const
995  /// Equals 1 - num_holes
996  inline int EulerNumber2D() const
997  { return NumOfVertices - NumOfEdges + NumOfElements; }
998 
999  int Dimension() const { return Dim; }
1000  int SpaceDimension() const { return spaceDim; }
1001 
1002  /// @brief Return pointer to vertex i's coordinates.
1003  /// @warning For high-order meshes (when #Nodes != NULL) vertices may not be
1004  /// updated and should not be used!
1005  const double *GetVertex(int i) const { return vertices[i](); }
1006 
1007  /// @brief Return pointer to vertex i's coordinates.
1008  /// @warning For high-order meshes (when Nodes != NULL) vertices may not
1009  /// being updated and should not be used!
1010  double *GetVertex(int i) { return vertices[i](); }
1011 
1012  void GetElementData(int geom, Array<int> &elem_vtx, Array<int> &attr) const
1013  { GetElementData(elements, geom, elem_vtx, attr); }
1014 
1015  /// Checks if the mesh has boundary elements
1016  virtual bool HasBoundaryElements() const { return (NumOfBdrElements > 0); }
1017 
1018  void GetBdrElementData(int geom, Array<int> &bdr_elem_vtx,
1019  Array<int> &bdr_attr) const
1020  { GetElementData(boundary, geom, bdr_elem_vtx, bdr_attr); }
1021 
1022  /** @brief Set the internal Vertex array to point to the given @a vertices
1023  array without assuming ownership of the pointer. */
1024  /** If @a zerocopy is `true`, the vertices must be given as an array of 3
1025  doubles per vertex. If @a zerocopy is `false` then the current Vertex
1026  data is first copied to the @a vertices array. */
1027  void ChangeVertexDataOwnership(double *vertices, int len_vertices,
1028  bool zerocopy = false);
1029 
1030  const Element* const *GetElementsArray() const
1031  { return elements.GetData(); }
1032 
1033  const Element *GetElement(int i) const { return elements[i]; }
1034 
1035  Element *GetElement(int i) { return elements[i]; }
1036 
1037  const Element *GetBdrElement(int i) const { return boundary[i]; }
1038 
1039  Element *GetBdrElement(int i) { return boundary[i]; }
1040 
1041  const Element *GetFace(int i) const { return faces[i]; }
1042 
1044  {
1045  return faces[i]->GetGeometryType();
1046  }
1047 
1049  {
1050  return elements[i]->GetGeometryType();
1051  }
1052 
1054  {
1055  return boundary[i]->GetGeometryType();
1056  }
1057 
1058  // deprecated: "base geometry" no longer means anything
1060  { return GetFaceGeometry(i); }
1061 
1063  { return GetElementGeometry(i); }
1064 
1066  { return GetBdrElementGeometry(i); }
1067 
1068  /** @brief Return true iff the given @a geom is encountered in the mesh.
1069  Geometries of dimensions lower than Dimension() are counted as well. */
1070  bool HasGeometry(Geometry::Type geom) const
1071  { return mesh_geoms & (1 << geom); }
1072 
1073  /** @brief Return the number of geometries of the given dimension present in
1074  the mesh. */
1075  /** For a parallel mesh only the local geometries are counted. */
1076  int GetNumGeometries(int dim) const;
1077 
1078  /// Return all element geometries of the given dimension present in the mesh.
1079  /** For a parallel mesh only the local geometries are returned.
1080 
1081  The returned geometries are sorted. */
1082  void GetGeometries(int dim, Array<Geometry::Type> &el_geoms) const;
1083 
1084  /// List of mesh geometries stored as Array<Geometry::Type>.
1085  class GeometryList : public Array<Geometry::Type>
1086  {
1087  protected:
1089  public:
1090  /// Construct a GeometryList of all element geometries in @a mesh.
1092  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
1093  { mesh.GetGeometries(mesh.Dimension(), *this); }
1094  /** @brief Construct a GeometryList of all geometries of dimension @a dim
1095  in @a mesh. */
1096  GeometryList(Mesh &mesh, int dim)
1097  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
1098  { mesh.GetGeometries(dim, *this); }
1099  };
1100 
1101  /// Returns the indices of the vertices of element i.
1102  void GetElementVertices(int i, Array<int> &v) const
1103  { elements[i]->GetVertices(v); }
1104 
1105  /// Returns the indices of the vertices of boundary element i.
1106  void GetBdrElementVertices(int i, Array<int> &v) const
1107  { boundary[i]->GetVertices(v); }
1108 
1109  /// Return the indices and the orientations of all edges of element i.
1110  void GetElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
1111 
1112  /// Return the indices and the orientations of all edges of bdr element i.
1113  void GetBdrElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
1114 
1115  /** Return the indices and the orientations of all edges of face i.
1116  Works for both 2D (face=edge) and 3D faces. */
1117  void GetFaceEdges(int i, Array<int> &edges, Array<int> &o) const;
1118 
1119  /// Returns the indices of the vertices of face i.
1120  void GetFaceVertices(int i, Array<int> &vert) const
1121  {
1122  if (Dim == 1)
1123  {
1124  vert.SetSize(1); vert[0] = i;
1125  }
1126  else
1127  {
1128  faces[i]->GetVertices(vert);
1129  }
1130  }
1131 
1132  /// Returns the indices of the vertices of edge i.
1133  void GetEdgeVertices(int i, Array<int> &vert) const;
1134 
1135  /// Returns the face-to-edge Table (3D)
1136  Table *GetFaceEdgeTable() const;
1137 
1138  /// Returns the edge-to-vertex Table (3D)
1139  Table *GetEdgeVertexTable() const;
1140 
1141  /// Return the indices and the orientations of all faces of element i.
1142  void GetElementFaces(int i, Array<int> &faces, Array<int> &ori) const;
1143 
1144  /// Return the index and the orientation of the face of bdr element i. (3D)
1145  void GetBdrElementFace(int i, int *f, int *o) const;
1146 
1147  /** Return the vertex index of boundary element i. (1D)
1148  Return the edge index of boundary element i. (2D)
1149  Return the face index of boundary element i. (3D) */
1150  int GetBdrElementEdgeIndex(int i) const;
1151 
1152  /** @brief For the given boundary element, bdr_el, return its adjacent
1153  element and its info, i.e. 64*local_bdr_index+bdr_orientation. */
1154  void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const;
1155 
1156  /// Returns the type of element i.
1157  Element::Type GetElementType(int i) const;
1158 
1159  /// Returns the type of boundary element i.
1160  Element::Type GetBdrElementType(int i) const;
1161 
1162  /* Return point matrix of element i of dimension Dim X #v, where for every
1163  vertex we give its coordinates in space of dimension Dim. */
1164  void GetPointMatrix(int i, DenseMatrix &pointmat) const;
1165 
1166  /* Return point matrix of boundary element i of dimension Dim X #v, where for
1167  every vertex we give its coordinates in space of dimension Dim. */
1168  void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const;
1169 
1171 
1172  /** Builds the transformation defining the i-th element in the user-defined
1173  variable. */
1175 
1176  /// Returns the transformation defining the i-th element
1178 
1179  /** Return the transformation defining the i-th element assuming
1180  the position of the vertices/nodes are given by 'nodes'. */
1181  void GetElementTransformation(int i, const Vector &nodes,
1183 
1184  /// Returns the transformation defining the i-th boundary element
1187 
1188  /** @brief Returns the transformation defining the given face element in a
1189  user-defined variable. */
1191 
1192  /** @brief A helper method that constructs a transformation from the
1193  reference space of a face to the reference space of an element. */
1194  /** The local index of the face as a face in the element and its orientation
1195  are given by the input parameter @a info, as @a info = 64*loc_face_idx +
1196  loc_face_orientation. */
1197  void GetLocalFaceTransformation(int face_type, int elem_type,
1199  int info);
1200 
1201  /// Returns the transformation defining the given face element
1203 
1204  /** Returns the transformation defining the given edge element.
1205  The transformation is stored in a user-defined variable. */
1207 
1208  /// Returns the transformation defining the given face element
1210 
1211  /// Returns (a pointer to an object containing) the following data:
1212  ///
1213  /// 1) Elem1No - the index of the first element that contains this face this
1214  /// is the element that has the same outward unit normal vector as the
1215  /// face;
1216  ///
1217  /// 2) Elem2No - the index of the second element that contains this face this
1218  /// element has outward unit normal vector as the face multiplied with -1;
1219  ///
1220  /// 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first
1221  /// and the second element respectively;
1222  ///
1223  /// 4) Face - pointer to the ElementTransformation of the face;
1224  ///
1225  /// 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face
1226  /// coordinate system to the element coordinate system (both in their
1227  /// reference elements). Used to transform IntegrationPoints from face to
1228  /// element. More formally, let:
1229  /// TL1, TL2 be the transformations represented by Loc1, Loc2,
1230  /// TE1, TE2 - the transformations represented by Elem1, Elem2,
1231  /// TF - the transformation represented by Face, then
1232  /// TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
1233  ///
1234  /// 6) FaceGeom - the base geometry for the face.
1235  ///
1236  /// The mask specifies which fields in the structure to return:
1237  /// mask & 1 - Elem1, mask & 2 - Elem2
1238  /// mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.
1239  /// These mask values are defined in the ConfigMasks enum type as part of the
1240  /// FaceElementTransformations class in fem/eltrans.hpp.
1242  int FaceNo,
1243  int mask = 31);
1244 
1246  {
1247  if (faces_info[FaceNo].Elem2No < 0) { return NULL; }
1248  return GetFaceElementTransformations (FaceNo);
1249  }
1250 
1252 
1253  /// Return the local face index for the given boundary face.
1254  int GetBdrFace(int BdrElemNo) const;
1255 
1256  /// Return true if the given face is interior. @sa FaceIsTrueInterior().
1257  bool FaceIsInterior(int FaceNo) const
1258  {
1259  return (faces_info[FaceNo].Elem2No >= 0);
1260  }
1261 
1262  /** This enumerated type describes the three main face topologies:
1263  - Boundary, for faces on the boundary of the computational domain,
1264  - Conforming, for conforming faces interior to the computational domain,
1265  - Nonconforming, for nonconforming faces interior to the computational
1266  domain. */
1267  enum class FaceTopology { Boundary,
1268  Conforming,
1269  Nonconforming,
1270  NA
1271  };
1272 
1273  /** This enumerated type describes the location of the two elements sharing a
1274  face, Local meaning that the element is local to the MPI rank, FaceNbr
1275  meaning that the element is distributed on a different MPI rank, this
1276  typically means that methods with FaceNbr should be used to access the
1277  relevant information, e.g., ParFiniteElementSpace::GetFaceNbrElementVDofs.
1278  */
1279  enum class ElementLocation { Local, FaceNbr, NA };
1280 
1281  /** This enumerated type describes the topological relation of an element to
1282  a face:
1283  - Coincident meaning that the element's face is topologically equal to
1284  the mesh face.
1285  - Superset meaning that the element's face is topologically coarser than
1286  the mesh face, i.e., the element's face contains the mesh face.
1287  - Subset meaning that the element's face is topologically finer than the
1288  mesh face, i.e., the element's face is contained in the mesh face.
1289  Superset and Subset are only relevant for nonconforming faces.
1290  Master nonconforming faces have a conforming element on one side, and a
1291  fine element on the other side. Slave nonconforming faces have a
1292  conforming element on one side, and a coarse element on the other side.
1293  */
1295 
1296  /** This enumerated type describes the corresponding FaceInfo internal
1297  representation (encoded cases), c.f. FaceInfo's documentation:
1298  Classification of a local (non-ghost) face based on its FaceInfo:
1299  - Elem2No >= 0 --> local interior face; can be either:
1300  - NCFace == -1 --> LocalConforming,
1301  - NCFace >= 0 --> LocalSlaveNonconforming,
1302  - Elem2No < 0 --> local "boundary" face; can be one of:
1303  - NCFace == -1 --> conforming face; can be either:
1304  - Elem2Inf < 0 --> Boundary,
1305  - Elem2Inf >= 0 --> SharedConforming,
1306  - NCFace >= 0 --> nonconforming face; can be one of:
1307  - Elem2Inf < 0 --> MasterNonconforming (shared or not shared),
1308  - Elem2Inf >= 0 --> SharedSlaveNonconforming.
1309  Classification of a ghost (non-local) face based on its FaceInfo:
1310  - Elem1No == -1 --> GhostMaster (includes other unused ghost faces),
1311  - Elem1No >= 0 --> GhostSlave.
1312  */
1313  enum class FaceInfoTag { Boundary,
1319  GhostSlave,
1320  GhostMaster
1321  };
1322 
1323  /** @brief This structure is used as a human readable output format that
1324  decipheres the information contained in Mesh::FaceInfo when using the
1325  Mesh::GetFaceInformation() method.
1326 
1327  The element indices in this structure don't need further processing,
1328  contrary to the ones obtained through Mesh::GetFacesElements and can
1329  directly be used, e.g., Elem1 and Elem2 indices.
1330  Likewise the orientations for Elem1 and Elem2 already take into account
1331  special cases and can be used as is.
1332  */
1334  {
1336 
1337  struct
1338  {
1341  int index;
1344  } element[2];
1345 
1347  int ncface;
1349 
1350  /** @brief Return true if the face is a local interior face which is NOT
1351  a master nonconforming face. */
1352  bool IsLocal() const
1353  {
1354  return element[1].location == Mesh::ElementLocation::Local;
1355  }
1356 
1357  /** @brief Return true if the face is a shared interior face which is NOT
1358  a master nonconforming face. */
1359  bool IsShared() const
1360  {
1361  return element[1].location == Mesh::ElementLocation::FaceNbr;
1362  }
1363 
1364  /** @brief return true if the face is an interior face to the computation
1365  domain, either a local or shared interior face (not a boundary face)
1366  which is NOT a master nonconforming face.
1367  */
1368  bool IsInterior() const
1369  {
1370  return topology == FaceTopology::Conforming ||
1372  }
1373 
1374  /** @brief Return true if the face is a boundary face. */
1375  bool IsBoundary() const
1376  {
1377  return topology == FaceTopology::Boundary;
1378  }
1379 
1380  /// @brief Return true if the face is of the same type as @a type.
1381  bool IsOfFaceType(FaceType type) const
1382  {
1383  switch (type)
1384  {
1385  case FaceType::Interior:
1386  return IsInterior();
1387  case FaceType::Boundary:
1388  return IsBoundary();
1389  default:
1390  return false;
1391  }
1392  }
1393 
1394  /// @brief Return true if the face is a conforming face.
1395  bool IsConforming() const
1396  {
1398  }
1399 
1400  /// @brief Return true if the face is a nonconforming fine face.
1401  bool IsNonconformingFine() const
1402  {
1404  (element[0].conformity == ElementConformity::Superset ||
1405  element[1].conformity == ElementConformity::Superset);
1406  }
1407 
1408  /// @brief Return true if the face is a nonconforming coarse face.
1409  /** Note that ghost nonconforming master faces cannot be clearly
1410  identified as such with the currently available information, so this
1411  method will return false for such faces. */
1413  {
1415  element[1].conformity == ElementConformity::Subset;
1416  }
1417 
1418  /// @brief cast operator from FaceInformation to FaceInfo.
1419  operator Mesh::FaceInfo() const;
1420  };
1421 
1422  /** This method aims to provide face information in a deciphered format, i.e.
1423  Mesh::FaceInformation, compared to the raw encoded information returned
1424  by Mesh::GetFaceElements() and Mesh::GetFaceInfos(). */
1425  FaceInformation GetFaceInformation(int f) const;
1426 
1427  void GetFaceElements (int Face, int *Elem1, int *Elem2) const;
1428  void GetFaceInfos (int Face, int *Inf1, int *Inf2) const;
1429  void GetFaceInfos (int Face, int *Inf1, int *Inf2, int *NCFace) const;
1430 
1431  Geometry::Type GetFaceGeometryType(int Face) const;
1432  Element::Type GetFaceElementType(int Face) const;
1433 
1434  /// Check (and optionally attempt to fix) the orientation of the elements
1435  /** @param[in] fix_it If `true`, attempt to fix the orientations of some
1436  elements: triangles, quads, and tets.
1437  @return The number of elements with wrong orientation.
1438 
1439  @note For meshes with nodes (e.g. high-order or periodic meshes), fixing
1440  the element orientations may require additional permutation of the nodal
1441  GridFunction of the mesh which is not performed by this method. Instead,
1442  the method Finalize() should be used with the parameter
1443  @a fix_orientation set to `true`.
1444 
1445  @note This method performs a simple check if an element is inverted, e.g.
1446  for most elements types, it checks if the Jacobian of the mapping from
1447  the reference element is non-negative at the center of the element. */
1448  int CheckElementOrientation(bool fix_it = true);
1449 
1450  /// Check the orientation of the boundary elements
1451  /** @return The number of boundary elements with wrong orientation. */
1452  int CheckBdrElementOrientation(bool fix_it = true);
1453 
1454  /// Return the attribute of element i.
1455  int GetAttribute(int i) const { return elements[i]->GetAttribute(); }
1456 
1457  /// Set the attribute of element i.
1458  void SetAttribute(int i, int attr) { elements[i]->SetAttribute(attr); }
1459 
1460  /// Return the attribute of boundary element i.
1461  int GetBdrAttribute(int i) const { return boundary[i]->GetAttribute(); }
1462 
1463  /// Set the attribute of boundary element i.
1464  void SetBdrAttribute(int i, int attr) { boundary[i]->SetAttribute(attr); }
1465 
1466  const Table &ElementToElementTable();
1467 
1468  const Table &ElementToFaceTable() const;
1469 
1470  const Table &ElementToEdgeTable() const;
1471 
1472  /// The returned Table must be destroyed by the caller
1474 
1475  /** Return the "face"-element Table. Here "face" refers to face (3D),
1476  edge (2D), or vertex (1D).
1477  The returned Table must be destroyed by the caller. */
1478  Table *GetFaceToElementTable() const;
1479 
1480  /** This method modifies a tetrahedral mesh so that Nedelec spaces of order
1481  greater than 1 can be defined on the mesh. Specifically, we
1482  1) rotate all tets in the mesh so that the vertices {v0, v1, v2, v3}
1483  satisfy: v0 < v1 < min(v2, v3).
1484  2) rotate all boundary triangles so that the vertices {v0, v1, v2}
1485  satisfy: v0 < min(v1, v2).
1486 
1487  @note Refinement does not work after a call to this method! */
1488  MFEM_DEPRECATED virtual void ReorientTetMesh();
1489 
1490  int *CartesianPartitioning(int nxyz[]);
1491  int *GeneratePartitioning(int nparts, int part_method = 1);
1492  void CheckPartitioning(int *partitioning_);
1493 
1494  void CheckDisplacements(const Vector &displacements, double &tmax);
1495 
1496  // Vertices are only at the corners of elements, where you would expect them
1497  // in the lowest-order mesh.
1498  void MoveVertices(const Vector &displacements);
1499  void GetVertices(Vector &vert_coord) const;
1500  void SetVertices(const Vector &vert_coord);
1501 
1502  // Nodes are only active for higher order meshes, and share locations with
1503  // the vertices, plus all the higher- order control points within the element
1504  // and along the edges and on the faces.
1505  void GetNode(int i, double *coord) const;
1506  void SetNode(int i, const double *coord);
1507 
1508  // Node operations for curved mesh.
1509  // They call the corresponding '...Vertices' method if the
1510  // mesh is not curved (i.e. Nodes == NULL).
1511  void MoveNodes(const Vector &displacements);
1512  void GetNodes(Vector &node_coord) const;
1513  void SetNodes(const Vector &node_coord);
1514 
1515  /// Return a pointer to the internal node GridFunction (may be NULL).
1516  GridFunction *GetNodes() { return Nodes; }
1517  const GridFunction *GetNodes() const { return Nodes; }
1518  /// Return the mesh nodes ownership flag.
1519  bool OwnsNodes() const { return own_nodes; }
1520  /// Set the mesh nodes ownership flag.
1521  void SetNodesOwner(bool nodes_owner) { own_nodes = nodes_owner; }
1522  /// Replace the internal node GridFunction with the given GridFunction.
1523  void NewNodes(GridFunction &nodes, bool make_owner = false);
1524  /** Swap the internal node GridFunction pointer and ownership flag members
1525  with the given ones. */
1526  void SwapNodes(GridFunction *&nodes, int &own_nodes_);
1527 
1528  /// Return the mesh nodes/vertices projected on the given GridFunction.
1529  void GetNodes(GridFunction &nodes) const;
1530  /** Replace the internal node GridFunction with a new GridFunction defined
1531  on the given FiniteElementSpace. The new node coordinates are projected
1532  (derived) from the current nodes/vertices. */
1533  void SetNodalFESpace(FiniteElementSpace *nfes);
1534  /** Replace the internal node GridFunction with the given GridFunction. The
1535  given GridFunction is updated with node coordinates projected (derived)
1536  from the current nodes/vertices. */
1537  void SetNodalGridFunction(GridFunction *nodes, bool make_owner = false);
1538  /** Return the FiniteElementSpace on which the current mesh nodes are
1539  defined or NULL if the mesh does not have nodes. */
1540  const FiniteElementSpace *GetNodalFESpace() const;
1541  /** Make sure that the mesh has valid nodes, i.e. its geometry is described
1542  by a vector finite element grid function (even if it is a low-order mesh
1543  with straight edges). */
1544  void EnsureNodes();
1545 
1546  /** Set the curvature of the mesh nodes using the given polynomial degree,
1547  'order', and optionally: discontinuous or continuous FE space, 'discont',
1548  new space dimension, 'space_dim' (if != -1), and 'ordering' (byVDim by
1549  default). */
1550  virtual void SetCurvature(int order, bool discont = false, int space_dim = -1,
1551  int ordering = 1);
1552 
1553  /// Refine all mesh elements.
1554  /** @param[in] ref_algo %Refinement algorithm. Currently used only for pure
1555  tetrahedral meshes. If set to zero (default), a tet mesh will be refined
1556  using algorithm A, that produces elements with better quality compared to
1557  algorithm B used when the parameter is non-zero.
1558 
1559  For tetrahedral meshes, after using algorithm A, the mesh cannot be
1560  refined locally using methods like GeneralRefinement() unless it is
1561  re-finalized using Finalize() with the parameter @a refine set to true.
1562  Note that calling Finalize() in this way will generally invalidate any
1563  FiniteElementSpace%s and GridFunction%s defined on the mesh. */
1564  void UniformRefinement(int ref_algo = 0);
1565 
1566  /** Refine selected mesh elements. Refinement type can be specified for each
1567  element. The function can do conforming refinement of triangles and
1568  tetrahedra and nonconforming refinement (i.e., with hanging-nodes) of
1569  triangles, quadrilaterals and hexahedra. If 'nonconforming' = -1,
1570  suitable refinement method is selected automatically (namely, conforming
1571  refinement for triangles). Use nonconforming = 0/1 to force the method.
1572  For nonconforming refinements, nc_limit optionally specifies the maximum
1573  level of hanging nodes (unlimited by default). */
1574  void GeneralRefinement(const Array<Refinement> &refinements,
1575  int nonconforming = -1, int nc_limit = 0);
1576 
1577  /** Simplified version of GeneralRefinement taking a simple list of elements
1578  to refine, without refinement types. */
1579  void GeneralRefinement(const Array<int> &el_to_refine,
1580  int nonconforming = -1, int nc_limit = 0);
1581 
1582  /// Refine each element with given probability. Uses GeneralRefinement.
1583  void RandomRefinement(double prob, bool aniso = false,
1584  int nonconforming = -1, int nc_limit = 0);
1585 
1586  /// Refine elements sharing the specified vertex. Uses GeneralRefinement.
1587  void RefineAtVertex(const Vertex& vert,
1588  double eps = 0.0, int nonconforming = -1);
1589 
1590  /** Refine element i if elem_error[i] > threshold, for all i.
1591  Returns true if at least one element was refined, false otherwise. */
1592  bool RefineByError(const Array<double> &elem_error, double threshold,
1593  int nonconforming = -1, int nc_limit = 0);
1594 
1595  /** Refine element i if elem_error(i) > threshold, for all i.
1596  Returns true if at least one element was refined, false otherwise. */
1597  bool RefineByError(const Vector &elem_error, double threshold,
1598  int nonconforming = -1, int nc_limit = 0);
1599 
1600  /** Derefine the mesh based on an error measure associated with each
1601  element. A derefinement is performed if the sum of errors of its fine
1602  elements is smaller than 'threshold'. If 'nc_limit' > 0, derefinements
1603  that would increase the maximum level of hanging nodes of the mesh are
1604  skipped. Returns true if the mesh changed, false otherwise. */
1605  bool DerefineByError(Array<double> &elem_error, double threshold,
1606  int nc_limit = 0, int op = 1);
1607 
1608  /// Same as DerefineByError for an error vector.
1609  bool DerefineByError(const Vector &elem_error, double threshold,
1610  int nc_limit = 0, int op = 1);
1611 
1612  ///@{ @name NURBS mesh refinement methods
1613  void KnotInsert(Array<KnotVector *> &kv);
1614  void KnotInsert(Array<Vector *> &kv);
1615  /* For each knot vector:
1616  new_degree = max(old_degree, min(old_degree + rel_degree, degree)). */
1617  void DegreeElevate(int rel_degree, int degree = 16);
1618  ///@}
1619 
1620  /** Make sure that a quad/hex mesh is considered to be nonconforming (i.e.,
1621  has an associated NCMesh object). Simplex meshes can be both conforming
1622  (default) or nonconforming. */
1623  void EnsureNCMesh(bool simplices_nonconforming = false);
1624 
1625  bool Conforming() const { return ncmesh == NULL; }
1626  bool Nonconforming() const { return ncmesh != NULL; }
1627 
1628  /** Return fine element transformations following a mesh refinement.
1629  Space uses this to construct a global interpolation matrix. */
1631 
1632  /// Return type of last modification of the mesh.
1634 
1635  /** Return update counter. The counter starts at zero and is incremented
1636  each time refinement, derefinement, or rebalancing method is called.
1637  It is used for checking proper sequence of Space:: and GridFunction::
1638  Update() calls. */
1639  long GetSequence() const { return sequence; }
1640 
1641  /// Print the mesh to the given stream using Netgen/Truegrid format.
1642  virtual void PrintXG(std::ostream &os = mfem::out) const;
1643 
1644  /// Print the mesh to the given stream using the default MFEM mesh format.
1645  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1646  virtual void Print(std::ostream &os = mfem::out) const { Printer(os); }
1647 
1648  /// Save the mesh to a file using Mesh::Print. The given @a precision will be
1649  /// used for ASCII output.
1650  virtual void Save(const char *fname, int precision=16) const;
1651 
1652  /// Print the mesh to the given stream using the adios2 bp format
1653 #ifdef MFEM_USE_ADIOS2
1654  virtual void Print(adios2stream &os) const;
1655 #endif
1656  /// Print the mesh in VTK format (linear and quadratic meshes only).
1657  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1658  void PrintVTK(std::ostream &os);
1659  /** Print the mesh in VTK format. The parameter ref > 0 specifies an element
1660  subdivision number (useful for high order fields and curved meshes).
1661  If the optional field_data is set, we also add a FIELD section in the
1662  beginning of the file with additional dataset information. */
1663  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1664  void PrintVTK(std::ostream &os, int ref, int field_data=0);
1665  /** Print the mesh in VTU format. The parameter ref > 0 specifies an element
1666  subdivision number (useful for high order fields and curved meshes).
1667  If @a bdr_elements is true, then output (only) the boundary elements,
1668  otherwise output only the non-boundary elements. */
1669  void PrintVTU(std::ostream &os,
1670  int ref=1,
1671  VTKFormat format=VTKFormat::ASCII,
1672  bool high_order_output=false,
1673  int compression_level=0,
1674  bool bdr_elements=false);
1675  /** Print the mesh in VTU format with file name fname. */
1676  virtual void PrintVTU(std::string fname,
1677  VTKFormat format=VTKFormat::ASCII,
1678  bool high_order_output=false,
1679  int compression_level=0,
1680  bool bdr=false);
1681  /** Print the boundary elements of the mesh in VTU format, and output the
1682  boundary attributes as a data array (useful for boundary conditions). */
1683  void PrintBdrVTU(std::string fname,
1684  VTKFormat format=VTKFormat::ASCII,
1685  bool high_order_output=false,
1686  int compression_level=0);
1687 
1688  void GetElementColoring(Array<int> &colors, int el0 = 0);
1689 
1690  /** @brief Prints the mesh with boundary elements given by the boundary of
1691  the subdomains, so that the boundary of subdomain i has boundary
1692  attribute i+1. */
1693  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1694  void PrintWithPartitioning (int *partitioning,
1695  std::ostream &os, int elem_attr = 0) const;
1696 
1697  void PrintElementsWithPartitioning (int *partitioning,
1698  std::ostream &out,
1699  int interior_faces = 0);
1700 
1701  /// Print set of disjoint surfaces:
1702  /*!
1703  * If Aface_face(i,j) != 0, print face j as a boundary
1704  * element with attribute i+1.
1705  */
1706  void PrintSurfaces(const Table &Aface_face, std::ostream &out) const;
1707 
1708  void ScaleSubdomains (double sf);
1709  void ScaleElements (double sf);
1710 
1711  void Transform(void (*f)(const Vector&, Vector&));
1712  void Transform(VectorCoefficient &deformation);
1713 
1714  /// Remove unused vertices and rebuild mesh connectivity.
1715  void RemoveUnusedVertices();
1716 
1717  /** Remove boundary elements that lie in the interior of the mesh, i.e. that
1718  have two adjacent faces in 3D, or edges in 2D. */
1719  void RemoveInternalBoundaries();
1720 
1721  /** @brief Get the size of the i-th element relative to the perfect
1722  reference element. */
1723  double GetElementSize(int i, int type = 0);
1724 
1725  double GetElementSize(int i, const Vector &dir);
1726 
1727  double GetElementVolume(int i);
1728 
1729  void GetElementCenter(int i, Vector &center);
1730 
1731  /// Returns the minimum and maximum corners of the mesh bounding box.
1732  /** For high-order meshes, the geometry is first refined @a ref times. */
1733  void GetBoundingBox(Vector &min, Vector &max, int ref = 2);
1734 
1735  void GetCharacteristics(double &h_min, double &h_max,
1736  double &kappa_min, double &kappa_max,
1737  Vector *Vh = NULL, Vector *Vk = NULL);
1738 
1739  /// Auxiliary method used by PrintCharacteristics().
1740  /** It is also used in the `mesh-explorer` miniapp. */
1741  static void PrintElementsByGeometry(int dim,
1742  const Array<int> &num_elems_by_geom,
1743  std::ostream &out);
1744 
1745  /** @brief Compute and print mesh characteristics such as number of vertices,
1746  number of elements, number of boundary elements, minimal and maximal
1747  element sizes, minimal and maximal element aspect ratios, etc. */
1748  /** If @a Vh or @a Vk are not NULL, return the element sizes and aspect
1749  ratios for all elements in the given Vector%s. */
1750  void PrintCharacteristics(Vector *Vh = NULL, Vector *Vk = NULL,
1751  std::ostream &os = mfem::out);
1752 
1753  /** @brief In serial, this method calls PrintCharacteristics(). In parallel,
1754  additional information about the parallel decomposition is also printed.
1755  */
1756  virtual void PrintInfo(std::ostream &os = mfem::out)
1757  {
1758  PrintCharacteristics(NULL, NULL, os);
1759  }
1760 
1761  void MesquiteSmooth(const int mesquite_option = 0);
1762 
1763  /** @brief Find the ids of the elements that contain the given points, and
1764  their corresponding reference coordinates.
1765 
1766  The DenseMatrix @a point_mat describes the given points - one point for
1767  each column; it should have SpaceDimension() rows.
1768 
1769  The InverseElementTransformation object, @a inv_trans, is used to attempt
1770  the element transformation inversion. If NULL pointer is given, the
1771  method will use a default constructed InverseElementTransformation. Note
1772  that the algorithms in the base class InverseElementTransformation can be
1773  completely overwritten by deriving custom classes that override the
1774  Transform() method.
1775 
1776  If no element is found for the i-th point, elem_ids[i] is set to -1.
1777 
1778  In the ParMesh implementation, the @a point_mat is expected to be the
1779  same on all ranks. If the i-th point is found by multiple ranks, only one
1780  of them will mark that point as found, i.e. set its elem_ids[i] to a
1781  non-negative number; the other ranks will set their elem_ids[i] to -2 to
1782  indicate that the point was found but assigned to another rank.
1783 
1784  @returns The total number of points that were found.
1785 
1786  @note This method is not 100 percent reliable, i.e. it is not guaranteed
1787  to find a point, even if it lies inside a mesh element. */
1788  virtual int FindPoints(DenseMatrix& point_mat, Array<int>& elem_ids,
1789  Array<IntegrationPoint>& ips, bool warn = true,
1790  InverseElementTransformation *inv_trans = NULL);
1791 
1792  /// Swaps internal data with another mesh. By default, non-geometry members
1793  /// like 'ncmesh' and 'NURBSExt' are only swapped when 'non_geometry' is set.
1794  void Swap(Mesh& other, bool non_geometry);
1795 
1796  /// Destroys Mesh.
1797  virtual ~Mesh() { DestroyPointers(); }
1798 
1799 #ifdef MFEM_DEBUG
1800  /// Output an NCMesh-compatible debug dump.
1801  void DebugDump(std::ostream &out) const;
1802 #endif
1803 };
1804 
1805 /** Overload operator<< for std::ostream and Mesh; valid also for the derived
1806  class ParMesh */
1807 std::ostream &operator<<(std::ostream &out, const Mesh &mesh);
1808 
1809 
1810 /** @brief Structure for storing mesh geometric factors: coordinates, Jacobians,
1811  and determinants of the Jacobians. */
1812 /** Typically objects of this type are constructed and owned by objects of class
1813  Mesh. See Mesh::GetGeometricFactors(). */
1815 {
1816 
1817 private:
1818  void Compute(const GridFunction &nodes,
1820 
1821 public:
1822  const Mesh *mesh;
1825 
1827  {
1828  COORDINATES = 1 << 0,
1829  JACOBIANS = 1 << 1,
1830  DETERMINANTS = 1 << 2,
1831  };
1832 
1833  GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags,
1835 
1836  GeometricFactors(const GridFunction &nodes, const IntegrationRule &ir,
1837  int flags,
1839 
1840  /// Mapped (physical) coordinates of all quadrature points.
1841  /** This array uses a column-major layout with dimensions (NQ x SDIM x NE)
1842  where
1843  - NQ = number of quadrature points per element,
1844  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1845  - NE = number of elements in the mesh. */
1847 
1848  /// Jacobians of the element transformations at all quadrature points.
1849  /** This array uses a column-major layout with dimensions (NQ x SDIM x DIM x
1850  NE) where
1851  - NQ = number of quadrature points per element,
1852  - SDIM = space dimension of the mesh = mesh.SpaceDimension(),
1853  - DIM = dimension of the mesh = mesh.Dimension(), and
1854  - NE = number of elements in the mesh. */
1856 
1857  /// Determinants of the Jacobians at all quadrature points.
1858  /** This array uses a column-major layout with dimensions (NQ x NE) where
1859  - NQ = number of quadrature points per element, and
1860  - NE = number of elements in the mesh. */
1862 };
1863 
1864 /** @brief Structure for storing face geometric factors: coordinates, Jacobians,
1865  determinants of the Jacobians, and normal vectors. */
1866 /** Typically objects of this type are constructed and owned by objects of class
1867  Mesh. See Mesh::GetFaceGeometricFactors(). */
1869 {
1870 public:
1871  const Mesh *mesh;
1875 
1877  {
1878  COORDINATES = 1 << 0,
1879  JACOBIANS = 1 << 1,
1880  DETERMINANTS = 1 << 2,
1881  NORMALS = 1 << 3,
1882  };
1883 
1884  FaceGeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags,
1885  FaceType type);
1886 
1887  /// Mapped (physical) coordinates of all quadrature points.
1888  /** This array uses a column-major layout with dimensions (NQ x SDIM x NF)
1889  where
1890  - NQ = number of quadrature points per face,
1891  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1892  - NF = number of faces in the mesh. */
1894 
1895  /// Jacobians of the element transformations at all quadrature points.
1896  /** This array uses a column-major layout with dimensions (NQ x SDIM x DIM x
1897  NF) where
1898  - NQ = number of quadrature points per face,
1899  - SDIM = space dimension of the mesh = mesh.SpaceDimension(),
1900  - DIM = dimension of the mesh = mesh.Dimension(), and
1901  - NF = number of faces in the mesh. */
1903 
1904  /// Determinants of the Jacobians at all quadrature points.
1905  /** This array uses a column-major layout with dimensions (NQ x NF) where
1906  - NQ = number of quadrature points per face, and
1907  - NF = number of faces in the mesh. */
1909 
1910  /// Normals at all quadrature points.
1911  /** This array uses a column-major layout with dimensions (NQ x DIM x NF) where
1912  - NQ = number of quadrature points per face,
1913  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1914  - NF = number of faces in the mesh. */
1916 };
1917 
1918 /// Class used to extrude the nodes of a mesh
1920 {
1921 private:
1922  int n, layer;
1923  double p[2], s;
1924  Vector tip;
1925 public:
1926  NodeExtrudeCoefficient(const int dim, const int n_, const double s_);
1927  void SetLayer(const int l) { layer = l; }
1929  virtual void Eval(Vector &V, ElementTransformation &T,
1930  const IntegrationPoint &ip);
1932 };
1933 
1934 
1935 /// Extrude a 1D mesh
1936 Mesh *Extrude1D(Mesh *mesh, const int ny, const double sy,
1937  const bool closed = false);
1938 
1939 /// Extrude a 2D mesh
1940 Mesh *Extrude2D(Mesh *mesh, const int nz, const double sz);
1941 
1942 // shift cyclically 3 integers left-to-right
1943 inline void ShiftRight(int &a, int &b, int &c)
1944 {
1945  int t = a;
1946  a = c; c = b; b = t;
1947 }
1948 
1949 /// @brief Print function for Mesh::FaceInformation.
1950 std::ostream& operator<<(std::ostream& os, const Mesh::FaceInformation& info);
1951 
1952 }
1953 
1954 #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:4879
Abstract class for all finite elements.
Definition: fe_base.hpp:235
const IntegrationRule * IntRule
Definition: mesh.hpp:1872
MFEM_DEPRECATED Mesh(int n, double sx=1.0)
Deprecated: see MakeCartesian1D.
Definition: mesh.hpp:865
void Loader(std::istream &input, int generate_edges=0, std::string parse_tag="")
Definition: mesh.cpp:3954
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:2018
Geometry::Constants< Geometry::PYRAMID > pyr_t
Definition: mesh.hpp:263
void PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
Definition: mesh.cpp:11344
virtual void Save(const char *fname, int precision=16) const
Definition: mesh.cpp:10210
void GetPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:6207
bool IsNonconformingFine() const
Return true if the face is a nonconforming fine face.
Definition: mesh.hpp:1401
static const int vtk_quadratic_hex[27]
Definition: mesh.hpp:246
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:6936
FaceInformation GetFaceInformation(int f) const
Definition: mesh.cpp:1130
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const DenseMatrix * PointMatrix
Definition: mesh.hpp:208
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1461
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1005
void ScaleElements(double sf)
Definition: mesh.cpp:11481
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Table * GetEdgeVertexTable() const
Returns the edge-to-vertex Table (3D)
Definition: mesh.cpp:6022
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:5923
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo)
void FreeElement(Element *E)
Definition: mesh.cpp:11795
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1662
int CheckElementOrientation(bool fix_it=true)
Check (and optionally attempt to fix) the orientation of the elements.
Definition: mesh.cpp:5389
static bool remove_unused_vertices
Definition: mesh.hpp:281
void SetVertices(const Vector &vert_coord)
Definition: mesh.cpp:7781
static const int vtk_quadratic_tet[10]
Definition: mesh.hpp:243
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:5985
void Make2D(int nx, int ny, Element::Type type, double sx, double sy, bool generate_edges, bool sfc_ordering)
Definition: mesh.cpp:3327
bool FaceIsTrueInterior(int FaceNo) const
Definition: mesh.hpp:497
Vector J
Jacobians of the element transformations at all quadrature points.
Definition: mesh.hpp:1902
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition: mesh.hpp:1846
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:12070
static Mesh MakeSimplicial(const Mesh &orig_mesh)
Definition: mesh.cpp:4550
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:8086
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
bool IsNonconformingCoarse() const
Return true if the face is a nonconforming coarse face.
Definition: mesh.hpp:1412
long GetGlobalNE() const
Return the total (global) number of elements.
Definition: mesh.hpp:960
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:1741
static const int NumGeom
Definition: geom.hpp:42
void AddHexAsWedges(const int *vi, int attr=1)
Definition: mesh.cpp:1760
virtual void PrintInfo(std::ostream &os=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
Definition: mesh.hpp:1756
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:259
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:6975
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:6140
CoarseFineTransformations CoarseFineTr
Definition: mesh.hpp:235
int own_nodes
Definition: mesh.hpp:241
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5877
virtual bool HasBoundaryElements() const
Checks if the mesh has boundary elements.
Definition: mesh.hpp:1016
bool Conforming() const
Definition: mesh.hpp:1625
void MoveVertices(const Vector &displacements)
Definition: mesh.cpp:7761
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:931
IsoparametricTransformation Transformation
Definition: mesh.hpp:229
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
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:3711
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:7865
int NumOfEdges
Definition: mesh.hpp:70
virtual void UniformRefinement2D()
Refine a mixed 2D mesh uniformly.
Definition: mesh.hpp:392
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:7884
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:6197
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:129
void ReadNetgen2DMesh(std::istream &input, int &curved)
void ShiftRight(int &a, int &b, int &c)
Definition: mesh.hpp:1943
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:6202
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:3721
int EulerNumber2D() const
Equals 1 - num_holes.
Definition: mesh.hpp:996
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)
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1648
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition: mesh.hpp:1908
int AddBdrPoint(int v, int attr=1)
Definition: mesh.cpp:1868
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:5888
void DeleteTables()
Definition: mesh.hpp:290
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:1120
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
const FaceGeometricFactors * GetFaceGeometricFactors(const IntegrationRule &ir, const int flags, FaceType type)
Return the mesh geometric factors for the faces corresponding to the given integration rule...
Definition: mesh.cpp:860
GridFunction * Nodes
Definition: mesh.hpp:240
int NumOfElements
Definition: mesh.hpp:69
Mesh * Extrude1D(Mesh *mesh, const int ny, const double sy, const bool closed)
Extrude a 1D mesh.
Definition: mesh.cpp:12088
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:11551
void FinalizeCheck()
Definition: mesh.cpp:1913
void ApplyLocalSlaveTransformation(FaceElementTransformations &FT, const FaceInfo &fi, bool is_ghost)
Definition: mesh.cpp:1055
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:6163
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:7894
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
const Element * GetFace(int i) const
Definition: mesh.hpp:1041
void FinalizeWedgeMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a wedge Mesh.
Definition: mesh.cpp:2806
Element * GetElement(int i)
Definition: mesh.hpp:1035
void GenerateBoundaryElements()
Definition: mesh.cpp:1875
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:1059
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:1814
NodeExtrudeCoefficient(const int dim, const int n_, const double s_)
Definition: mesh.cpp:12064
Data type for vertex.
Definition: vertex.hpp:22
void SetMeshGen()
Determine the mesh generator bitmask meshgen, see MeshGenerator().
Definition: mesh.cpp:3901
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:5585
void GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:7770
void ReadNetgen3DMesh(std::istream &input)
Data arrays will be written in ASCII format.
The inverse transformation of a given ElementTransformation.
Definition: eltrans.hpp:185
void ReadInlineMesh(std::istream &input, bool generate_edges=false)
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
Definition: mesh.hpp:1245
ElementLocation
Definition: mesh.hpp:1279
void RemoveInternalBoundaries()
Definition: mesh.cpp:11707
ElementConformity
Definition: mesh.hpp:1294
void DebugDump(std::ostream &out) const
Output an NCMesh-compatible debug dump.
Definition: mesh.cpp:12475
Mesh * Extrude2D(Mesh *mesh, const int nz, const double sz)
Extrude a 2D mesh.
Definition: mesh.cpp:12248
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:960
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1062
void SetEmpty()
Definition: mesh.cpp:1460
Array< Element * > faces
Definition: mesh.hpp:92
virtual long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
Definition: mesh.hpp:957
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
Definition: mesh.cpp:2841
Geometry::Type GetBdrElementGeometry(int i) const
Definition: mesh.hpp:1053
void GetVertexToVertexTable(DSTable &) const
Definition: mesh.cpp:6278
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:354
Operation last_operation
Definition: mesh.hpp:284
void MarkTriMeshForRefinement()
Definition: mesh.cpp:2406
void ReadCubit(const char *filename, int &curved, int &read_gf)
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:1070
void InitRefinementTransforms()
Definition: mesh.cpp:9827
IsoparametricTransformation FaceTransformation
Definition: mesh.hpp:231
void UniformRefinement2D_base(bool update_nodes=true)
Definition: mesh.cpp:7927
Array< FaceGeometricFactors * > face_geom_factors
Optional face geometric factors.
Definition: mesh.hpp:276
Array< NCFaceInfo > nc_faces_info
Definition: mesh.hpp:218
Element * ReadElement(std::istream &)
Definition: mesh.cpp:3883
void KnotInsert(Array< KnotVector * > &kv)
Definition: mesh.cpp:5038
void GetLocalTriToPyrTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:746
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
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:4913
Element * NewElement(int geom)
Definition: mesh.cpp:3829
FaceElementTransformations FaceElemTr
Definition: mesh.hpp:232
Table * el_to_face
Definition: mesh.hpp:221
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
Definition: mesh.cpp:574
const GridFunction * GetNodes() const
Definition: mesh.hpp:1517
void MarkForRefinement()
Definition: mesh.cpp:2389
Element::Type GetFaceElementType(int Face) const
Definition: mesh.cpp:1431
void SetVerticesFromNodes(const GridFunction *nodes)
Helper to set vertex coordinates given a high-order curvature function.
Definition: mesh.cpp:5331
Structure for storing face geometric factors: coordinates, Jacobians, determinants of the Jacobians...
Definition: mesh.hpp:1868
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:431
int AddBdrTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1825
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
Definition: mesh.cpp:10168
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:5345
void MakeRefined_(Mesh &orig_mesh, const Array< int > ref_factors, int ref_type)
Internal function used in Mesh::MakeRefined.
Definition: mesh.cpp:4321
Geometry::Constants< Geometry::SEGMENT > seg_t
Definition: mesh.hpp:257
bool IsShared() const
Return true if the face is a shared interior face which is NOT a master nonconforming face...
Definition: mesh.hpp:1359
void GetElementJacobian(int i, DenseMatrix &J)
Definition: mesh.cpp:60
void MesquiteSmooth(const int mesquite_option=0)
Definition: mesquite.cpp:1144
static int GetTetOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:5633
long GetSequence() const
Definition: mesh.hpp:1639
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:846
int AddTri(const int *vi, int attr=1)
Definition: mesh.hpp:721
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1601
bool Nonconforming() const
Definition: mesh.hpp:1626
Vector J
Jacobians of the element transformations at all quadrature points.
Definition: mesh.hpp:1855
FaceType
Definition: mesh.hpp:45
void MoveNodes(const Vector &displacements)
Definition: mesh.cpp:7829
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
Definition: mesh.cpp:2421
void GetElementData(int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.hpp:1012
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:10852
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:1368
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:1927
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:5112
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 PrintVTK(std::ostream &os)
Definition: mesh.cpp:10224
void FinalizeTetMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a tetrahedral Mesh.
Definition: mesh.cpp:2765
void AddBdrQuadAsTriangles(const int *vi, int attr=1)
Definition: mesh.cpp:1853
void Make3D(int nx, int ny, int nz, Element::Type type, double sx, double sy, double sz, bool sfc_ordering)
Definition: mesh.cpp:3046
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
Definition: mesh.cpp:6428
virtual void SetAttributes()
Definition: mesh.cpp:1547
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:9394
GeometryList(Mesh &mesh, int dim)
Construct a GeometryList of all geometries of dimension dim in mesh.
Definition: mesh.hpp:1096
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition: mesh.hpp:1861
void AddHexAsPyramids(const int *vi, int attr=1)
Definition: mesh.cpp:1778
int FindCoarseElement(int i)
Definition: mesh.cpp:9839
GeometryList(Mesh &mesh)
Construct a GeometryList of all element geometries in mesh.
Definition: mesh.hpp:1091
double b
Definition: lissajous.cpp:42
int AddTet(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1676
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Returns the transformation defining the given face element in a user-defined variable.
Definition: mesh.cpp:497
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:656
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
Definition: mesh.cpp:6493
void CheckPartitioning(int *partitioning_)
Definition: mesh.cpp:7366
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
Definition: mesh.cpp:5177
void ReadGmshMesh(std::istream &input, int &curved, int &read_gf)
STable3D * GetElementToFaceTable(int ret_ftbl=0)
Definition: mesh.cpp:6745
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:5312
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:921
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:6225
void InitTables()
Definition: mesh.cpp:1454
Geometry::Type GetFaceGeometry(int i) const
Definition: mesh.hpp:1043
void CheckDisplacements(const Vector &displacements, double &tmax)
Definition: mesh.cpp:7684
This structure stores the low level information necessary to interpret the configuration of elements ...
Definition: mesh.hpp:153
VTKFormat
Data array format for VTK and VTU files.
Definition: vtk.hpp:93
int mesh_geoms
Definition: mesh.hpp:78
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1811
Geometry::Constants< Geometry::TETRAHEDRON > tet_t
Definition: mesh.hpp:260
int AddElement(Element *elem)
The parameter elem should be allocated using the NewElement() method.
Definition: mesh.cpp:1797
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:896
void SetNodesOwner(bool nodes_owner)
Set the mesh nodes ownership flag.
Definition: mesh.hpp:1521
void SetLayer(const int l)
Definition: mesh.hpp:1927
int GetElementToEdgeTable(Table &, Array< int > &)
Definition: mesh.cpp:6303
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:6175
bool IsBoundary() const
Return true if the face is a boundary face.
Definition: mesh.hpp:1375
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:2185
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:9418
const Element * GetElement(int i) const
Definition: mesh.hpp:1033
void GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:771
void ResetLazyData()
Definition: mesh.cpp:1537
static const int vtk_quadratic_pyramid[13]
Definition: mesh.hpp:244
void GetElementCenter(int i, Vector &center)
Definition: mesh.cpp:68
IsoparametricTransformation BdrTransformation
Definition: mesh.hpp:230
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:895
IsoparametricTransformation Transformation2
Definition: mesh.hpp:229
void Init()
Definition: mesh.cpp:1436
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5273
Element * GetBdrElement(int i)
Definition: mesh.hpp:1039
Vector normal
Normals at all quadrature points.
Definition: mesh.hpp:1915
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:1094
Geometry::Constants< Geometry::PRISM > pri_t
Definition: mesh.hpp:262
int Dimension() const
Definition: mesh.hpp:999
virtual MFEM_DEPRECATED void ReorientTetMesh()
Definition: mesh.cpp:6874
prob_type prob
Definition: ex25.cpp:153
int NumOfBdrElements
Definition: mesh.hpp:69
Table * el_to_edge
Definition: mesh.hpp:220
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:8774
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
virtual void MarkTetMeshForRefinement(DSTable &v_to_v)
Definition: mesh.cpp:2446
Array< Triple< int, int, int > > tmp_vertex_parents
Definition: mesh.hpp:254
bool IsOfFaceType(FaceType type) const
Return true if the face is of the same type as type.
Definition: mesh.hpp:1381
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:856
List of mesh geometries stored as Array&lt;Geometry::Type&gt;.
Definition: mesh.hpp:1085
void Destroy()
Definition: mesh.cpp:1508
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:10424
This structure is used as a human readable output format that decipheres the information contained in...
Definition: mesh.hpp:1333
int CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Definition: mesh.cpp:5766
int SpaceDimension() const
Definition: mesh.hpp:1000
GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, MemoryType d_mt=MemoryType::DEFAULT)
Definition: mesh.cpp:11935
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
Definition: mesh.cpp:10971
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
Definition: mesh.hpp:1257
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:6396
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:1393
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:3505
const Element *const * GetElementsArray() const
Definition: mesh.hpp:1030
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:9016
virtual ~Mesh()
Destroys Mesh.
Definition: mesh.hpp:1797
void FinalizeMesh(int refine=0, bool fix_orientation=true)
Finalize the construction of any type of Mesh.
Definition: mesh.cpp:2871
Table * el_to_el
Definition: mesh.hpp:222
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:1399
int AddSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1634
struct mfem::Mesh::FaceInformation::@3 element[2]
Geometry::Constants< Geometry::CUBE > hex_t
Definition: mesh.hpp:261
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
Geometry::Type geom_buf[Geometry::NumGeom]
Definition: mesh.hpp:1088
static Mesh LoadFromFile(const char *filename, int generate_edges=0, int refine=1, bool fix_orientation=true)
Definition: mesh.cpp:3683
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Definition: mesh.cpp:11600
static void PrintElement(const Element *, std::ostream &)
Definition: mesh.cpp:3895
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
Definition: mesh.cpp:5082
int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
Definition: mesh.cpp:1725
const IntegrationRule * IntRule
Definition: mesh.hpp:1823
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition: mesh.hpp:1893
bool IsLocal() const
Return true if the face is a local interior face which is NOT a master nonconforming face...
Definition: mesh.hpp:1352
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2878
static void GetElementArrayEdgeTable(const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
Definition: mesh.cpp:6256
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
Definition: mesh.hpp:1516
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
Definition: mesh.cpp:2469
int AddBdrQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1839
bool IsSlaveFace(const FaceInfo &fi) const
Definition: mesh.cpp:1050
Array< Vertex > vertices
Definition: mesh.hpp:90
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:3701
void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Definition: mesh.cpp:1578
void SetBdrAttribute(int i, int attr)
Set the attribute of boundary element i.
Definition: mesh.hpp:1464
void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo, const std::string &xml_prefix="")
Geometry::Type GetElementGeometry(int i) const
Definition: mesh.hpp:1048
void RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
Definition: mesh.cpp:9437
static Mesh MakeCartesian1D(int n, double sx=1.0)
Definition: mesh.cpp:3693
Array< Element * > elements
Definition: mesh.hpp:85
const Table & ElementToElementTable()
Definition: mesh.cpp:6342
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
Definition: mesh.cpp:6241
int meshgen
Definition: mesh.hpp:77
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:9849
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9463
Table * bel_to_edge
Definition: mesh.hpp:224
Operation GetLastOperation() const
Return type of last modification of the mesh.
Definition: mesh.hpp:1633
Table * GetFaceToElementTable() const
Definition: mesh.cpp:6083
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:696
double a
Definition: lissajous.cpp:41
int EulerNumber() const
Equals 1 + num_holes - num_loops.
Definition: mesh.hpp:993
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:272
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:925
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:10416
void GetNode(int i, double *coord) const
Definition: mesh.cpp:7790
Array< int > be_to_edge
Definition: mesh.hpp:223
const Table & ElementToFaceTable() const
Definition: mesh.cpp:6378
Class for integration point with weight.
Definition: intrules.hpp:25
void Swap(Mesh &other, bool non_geometry)
Definition: mesh.cpp:9206
Element * ReadElementWithoutAttr(std::istream &)
Definition: mesh.cpp:3853
FaceGeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, FaceType type)
Definition: mesh.cpp:12015
int GetNumFacesWithGhost() const
Return the number of faces (3D), edges (2D) or vertices (1D) including ghost faces.
Definition: mesh.cpp:5356
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
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:5901
static const int vtk_quadratic_wedge[18]
Definition: mesh.hpp:245
bool OwnsNodes() const
Return the mesh nodes ownership flag.
Definition: mesh.hpp:1519
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:676
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5307
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1646
Table * face_edge
Definition: mesh.hpp:226
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1956
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
Definition: mesh.hpp:1018
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2971
bool DerefineByError(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
Definition: mesh.cpp:9133
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:2237
void GetFaceEdges(int i, Array< int > &edges, Array< int > &o) const
Definition: mesh.cpp:5955
void GetLocalQuadToPyrTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:817
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:5095
int dim
Definition: ex24.cpp:53
Table * edge_vertex
Definition: mesh.hpp:227
long sequence
Definition: mesh.hpp:83
virtual ~NodeExtrudeCoefficient()
Definition: mesh.hpp:1931
const Mesh * mesh
Definition: mesh.hpp:1822
void MakeSimplicial_(const Mesh &orig_mesh, int *vglobal)
Definition: mesh.cpp:4557
void DestroyPointers()
Definition: mesh.cpp:1482
Array< FaceInfo > faces_info
Definition: mesh.hpp:217
STable3D * GetFacesTable()
Definition: mesh.cpp:6682
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
Definition: sparsemat.hpp:680
NCMesh * ncmesh
Optional nonconforming mesh extension.
Definition: mesh.hpp:273
Geometry::Type GetBdrElementBaseGeometry(int i) const
Definition: mesh.hpp:1065
virtual void UniformRefinement3D()
Refine a mixed 3D mesh uniformly.
Definition: mesh.hpp:402
int AddBdrElement(Element *elem)
Definition: mesh.cpp:1804
RefCoord t[3]
static void PrintElementWithoutAttr(const Element *, std::ostream &)
Definition: mesh.cpp:3871
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:697
int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr=1)
Definition: mesh.cpp:1697
virtual bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NC version of GeneralDerefinement.
Definition: mesh.cpp:9085
int Dim
Definition: mesh.hpp:66
Geometry::Type GetFaceGeometryType(int Face) const
Definition: mesh.cpp:1412
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:904
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:934
ElementConformity conformity
Definition: mesh.hpp:1340
Geometry::Constants< Geometry::TRIANGLE > tri_t
Definition: mesh.hpp:258
double AggregateError(const Array< double > &elem_error, const int *fine, int nfine, int op)
Derefinement helper.
Definition: mesh.cpp:9066
double * GetVertex(int i)
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1010
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:937
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:11820
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:1458
Mesh & operator=(Mesh &&mesh)
Move assignment operstor.
Definition: mesh.cpp:3677
virtual void PrintXG(std::ostream &os=mfem::out) const
Print the mesh to the given stream using Netgen/Truegrid format.
Definition: mesh.cpp:9912
Vector data type.
Definition: vector.hpp:60
void ReadTrueGridMesh(std::istream &input)
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:5361
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:5537
void SetNode(int i, const double *coord)
Definition: mesh.cpp:7809
int GetBdrFace(int BdrElemNo) const
Return the local face index for the given boundary face.
Definition: mesh.cpp:1112
void DestroyTables()
Definition: mesh.cpp:1466
void GetLocalTriToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:720
void GenerateFaces()
Definition: mesh.cpp:6520
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:10777
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:3763
void DeleteGeometricFactors()
Destroy all GeometricFactors stored by the Mesh.
Definition: mesh.cpp:881
int spaceDim
Definition: mesh.hpp:67
IsoparametricTransformation EdgeTransformation
Definition: mesh.hpp:231
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:70
MemAlloc< Tetrahedron, 1024 > TetMemory
Definition: mesh.hpp:250
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
Definition: mesh.cpp:5994
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
Definition: mesh.cpp:2535
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:7853
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:5301
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
int NumOfVertices
Definition: mesh.hpp:69
void ScaleSubdomains(double sf)
Definition: mesh.cpp:11411
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:1711
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
Definition: mesh.cpp:9162
void Bisection(int i, const DSTable &, int *, int *, int *)
Bisect a triangle: element with index i is bisected.
Definition: mesh.cpp:9492
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9326
ElementLocation location
Definition: mesh.hpp:1339
Class for parallel meshes.
Definition: pmesh.hpp:32
Abstract data type element.
Definition: element.hpp:28
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1455
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
Definition: mesh.cpp:6465
const Table & ElementToEdgeTable() const
Definition: mesh.cpp:6387
Array< GeometricFactors * > geom_factors
Optional geometric factors.
Definition: mesh.hpp:274
void BdrBisection(int i, const HashTable< Hashed2 > &)
Bisect a boundary triangle: boundary element with index i is bisected.
Definition: mesh.cpp:9700
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
Definition: mesh.hpp:1106
void GenerateNCFaceInfo()
Definition: mesh.cpp:6621
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:1102
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:268
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
Definition: mesh.cpp:6048
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:1037
void EnsureNodes()
Definition: mesh.cpp:5279
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.cpp:9262
const DenseMatrix * point_matrix
Definition: mesh.hpp:1348
void UpdateNodes()
Update the nodes of a curved mesh after refinement.
Definition: mesh.cpp:7915
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:6117
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:360
bool IsConforming() const
Return true if the face is a conforming face.
Definition: mesh.hpp:1395
void Printer(std::ostream &out=mfem::out, std::string section_delimiter="") const
Definition: mesh.cpp:10076
Class used to extrude the nodes of a mesh.
Definition: mesh.hpp:1919
void AddVertexParents(int i, int p1, int p2)
Mark vertex i as nonconforming, with parent vertices p1 and p2.
Definition: mesh.cpp:1618
int NumOfFaces
Definition: mesh.hpp:70