MFEM  v4.3.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-2021, 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  struct FaceInfo
95  {
96  // Inf = 64 * LocalFaceIndex + FaceOrientation
98  int NCFace; /* -1 if this is a regular conforming/boundary face;
99  index into 'nc_faces_info' if >= 0. */
100  };
101  // NOTE: in NC meshes, master faces have Elem2No == -1. Slave faces on the
102  // other hand have Elem2No and Elem2Inf set to the master face's element and
103  // its local face number.
104  //
105  // A local face is one generated from a local element and has index i in
106  // faces_info such that i < GetNumFaces(). Also, Elem1No always refers to the
107  // element (slave or master, in the non-conforming case) that generated the
108  // face.
109  // Classification of a local (non-ghost) face based on its FaceInfo:
110  // - Elem2No >= 0 --> local internal face; can be either:
111  // - NCFace == -1 --> conforming face, or
112  // - NCFace >= 0 --> non-conforming slave face.
113  // - Elem2No < 0 --> local "boundary" face; can be one of:
114  // - NCFace == -1 --> conforming face; can be either:
115  // - Elem2Inf < 0 --> true boundary face (no element on side 2)
116  // - Elem2Inf >= 0 --> shared face where element 2 is a face-neighbor
117  // element with index -1-Elem2No. This state is initialized by
118  // ParMesh::ExchangeFaceNbrData().
119  // - NCFace >= 0 --> non-conforming master face. Elem2No is -1 or, in the
120  // case of a shared face, -1-Elem2No is the index of one of the adjacent
121  // (the last one?) slave ghost elements. Elem2Inf is -1.
122  //
123  // A ghost face is a non-conforming face that is generated by a non-local,
124  // i.e. ghost, element. A ghost face has index i in faces_info such that
125  // i >= GetNumFaces().
126  // Classification of a ghost (non-local) face based on its FaceInfo:
127  // - Elem1No == -1 --> master ghost face? These ghost faces also have:
128  // Elem2No == -1, Elem1Inf == Elem2Inf == -1, and NCFace == -1.
129  // - Elem1No >= 0 --> slave ghost face; Elem1No is the index of the local
130  // master side element, i.e. side 1 IS NOT the side that generated the
131  // face. Elem2No is < 0 and -1-Elem2No is the index of the ghost
132  // face-neighbor element that generated this slave ghost face. In this
133  // case, Elem2Inf >= 0.
134  // Relevant methods: GenerateFaces(), GenerateNCFaceInfo(),
135  // ParNCMesh::GetFaceNeighbors(),
136  // ParMesh::ExchangeFaceNbrData()
137 
138  struct NCFaceInfo
139  {
140  bool Slave; // true if this is a slave face, false if master face
141  int MasterFace; // if Slave, this is the index of the master face
142  // If not Slave, 'MasterFace' is the local face index of this master face
143  // as a face in the unique adjacent element.
144  const DenseMatrix* PointMatrix; // if Slave, position within master face
145  // (NOTE: PointMatrix points to a matrix owned by NCMesh.)
146 
147  NCFaceInfo() = default;
148 
149  NCFaceInfo(bool slave, int master, const DenseMatrix* pm)
150  : Slave(slave), MasterFace(master), PointMatrix(pm) {}
151  };
152 
155 
160  Table *bel_to_edge; // for 3D
162  mutable Table *face_edge;
163  mutable Table *edge_vertex;
164 
169 
170  // refinement embeddings for forward compatibility with NCMesh
172 
173  // Nodes are only active for higher order meshes, and share locations with
174  // the vertices, plus all the higher- order control points within the
175  // element and along the edges and on the faces.
178 
179  static const int vtk_quadratic_tet[10];
180  static const int vtk_quadratic_wedge[18];
181  static const int vtk_quadratic_hex[27];
182 
183 #ifdef MFEM_USE_MEMALLOC
184  friend class Tetrahedron;
186 #endif
187 
188  // used during NC mesh initialization only
190 
191 public:
198 
200 
201  /// A list of all unique element attributes used by the Mesh.
203  /// A list of all unique boundary attributes used by the Mesh.
205 
206  NURBSExtension *NURBSext; ///< Optional NURBS mesh extension.
207  NCMesh *ncmesh; ///< Optional non-conforming mesh extension.
208  Array<GeometricFactors*> geom_factors; ///< Optional geometric factors.
210  face_geom_factors; ///< Optional face geometric factors.
211 
212  // Global parameter that can be used to control the removal of unused
213  // vertices performed when reading a mesh in MFEM format. The default value
214  // (true) is set in mesh_readers.cpp.
216 
217 protected:
219 
220  void Init();
221  void InitTables();
222  void SetEmpty(); // Init all data members with empty values
223  void DestroyTables();
225  void DestroyPointers(); // Delete data specifically allocated by class Mesh.
226  void Destroy(); // Delete all owned data.
227  void ResetLazyData();
228 
229  Element *ReadElementWithoutAttr(std::istream &);
230  static void PrintElementWithoutAttr(const Element *, std::ostream &);
231 
232  Element *ReadElement(std::istream &);
233  static void PrintElement(const Element *, std::ostream &);
234 
235  // Readers for different mesh formats, used in the Load() method.
236  // The implementations of these methods are in mesh_readers.cpp.
237  void ReadMFEMMesh(std::istream &input, int version, int &curved);
238  void ReadLineMesh(std::istream &input);
239  void ReadNetgen2DMesh(std::istream &input, int &curved);
240  void ReadNetgen3DMesh(std::istream &input);
241  void ReadTrueGridMesh(std::istream &input);
242  void CreateVTKMesh(const Vector &points, const Array<int> &cell_data,
243  const Array<int> &cell_offsets,
244  const Array<int> &cell_types,
245  const Array<int> &cell_attributes,
246  int &curved, int &read_gf, bool &finalize_topo);
247  void ReadVTKMesh(std::istream &input, int &curved, int &read_gf,
248  bool &finalize_topo);
249  void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf,
250  bool &finalize_topo, const std::string &xml_prefix="");
251  void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf);
252  void ReadInlineMesh(std::istream &input, bool generate_edges = false);
253  void ReadGmshMesh(std::istream &input, int &curved, int &read_gf);
254  /* Note NetCDF (optional library) is used for reading cubit files */
255 #ifdef MFEM_USE_NETCDF
256  void ReadCubit(const char *filename, int &curved, int &read_gf);
257 #endif
258 
259  /// Determine the mesh generator bitmask #meshgen, see MeshGenerator().
260  /** Also, initializes #mesh_geoms. */
261  void SetMeshGen();
262 
263  /// Return the length of the segment from node i to node j.
264  double GetLength(int i, int j) const;
265 
266  /** Compute the Jacobian of the transformation from the perfect
267  reference element at the center of the element. */
268  void GetElementJacobian(int i, DenseMatrix &J);
269 
270  void MarkForRefinement();
272  void GetEdgeOrdering(DSTable &v_to_v, Array<int> &order);
273  virtual void MarkTetMeshForRefinement(DSTable &v_to_v);
274 
275  // Methods used to prepare and apply permutation of the mesh nodes assuming
276  // that the mesh elements may be rotated (e.g. to mark triangle or tet edges
277  // for refinement) between the two calls - PrepareNodeReorder() and
278  // DoNodeReorder(). The latter method assumes that the 'faces' have not been
279  // updated after the element rotations.
280  void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert);
281  void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert);
282 
284  STable3D *GetElementToFaceTable(int ret_ftbl = 0);
285 
286  /** Red refinement. Element with index i is refined. The default
287  red refinement for now is Uniform. */
288  void RedRefinement(int i, const DSTable &v_to_v,
289  int *edge1, int *edge2, int *middle)
290  { UniformRefinement(i, v_to_v, edge1, edge2, middle); }
291 
292  /** Green refinement. Element with index i is refined. The default
293  refinement for now is Bisection. */
294  void GreenRefinement(int i, const DSTable &v_to_v,
295  int *edge1, int *edge2, int *middle)
296  { Bisection(i, v_to_v, edge1, edge2, middle); }
297 
298  /// Bisect a triangle: element with index @a i is bisected.
299  void Bisection(int i, const DSTable &, int *, int *, int *);
300 
301  /// Bisect a tetrahedron: element with index @a i is bisected.
302  void Bisection(int i, HashTable<Hashed2> &);
303 
304  /// Bisect a boundary triangle: boundary element with index @a i is bisected.
305  void BdrBisection(int i, const HashTable<Hashed2> &);
306 
307  /** Uniform Refinement. Element with index i is refined uniformly. */
308  void UniformRefinement(int i, const DSTable &, int *, int *, int *);
309 
310  /** @brief Averages the vertices with given @a indexes and saves the result
311  in #vertices[result]. */
312  void AverageVertices(const int *indexes, int n, int result);
313 
315  int FindCoarseElement(int i);
316 
317  /// Update the nodes of a curved mesh after refinement
318  void UpdateNodes();
319 
320  /// Helper to set vertex coordinates given a high-order curvature function.
321  void SetVerticesFromNodes(const GridFunction *nodes);
322 
323  void UniformRefinement2D_base(bool update_nodes = true);
324 
325  /// Refine a mixed 2D mesh uniformly.
327 
328  /* If @a f2qf is not NULL, adds all quadrilateral faces to @a f2qf which
329  represents a "face-to-quad-face" index map. When all faces are quads, the
330  array @a f2qf is kept empty since it is not needed. */
331  void UniformRefinement3D_base(Array<int> *f2qf = NULL,
332  DSTable *v_to_v_p = NULL,
333  bool update_nodes = true);
334 
335  /// Refine a mixed 3D mesh uniformly.
337 
338  /// Refine NURBS mesh.
339  virtual void NURBSUniformRefinement();
340 
341  /// This function is not public anymore. Use GeneralRefinement instead.
342  virtual void LocalRefinement(const Array<int> &marked_el, int type = 3);
343 
344  /// This function is not public anymore. Use GeneralRefinement instead.
345  virtual void NonconformingRefinement(const Array<Refinement> &refinements,
346  int nc_limit = 0);
347 
348  /// NC version of GeneralDerefinement.
349  virtual bool NonconformingDerefinement(Array<double> &elem_error,
350  double threshold, int nc_limit = 0,
351  int op = 1);
352  /// Derefinement helper.
353  double AggregateError(const Array<double> &elem_error,
354  const int *fine, int nfine, int op);
355 
356  /// Read NURBS patch/macro-element mesh
357  void LoadPatchTopo(std::istream &input, Array<int> &edge_to_knot);
358 
359  void UpdateNURBS();
360 
361  void PrintTopo(std::ostream &out, const Array<int> &e_to_k) const;
362 
363  /// Used in GetFaceElementTransformations (...)
366  int i);
368  int i);
369  /// Used in GetFaceElementTransformations (...)
371  int i);
372  /// Used in GetFaceElementTransformations (...)
374  int i);
375  /// Used in GetFaceElementTransformations (...)
377  int i);
378  /// Used in GetFaceElementTransformations (...)
380  int i);
381 
382  /** Used in GetFaceElementTransformations to account for the fact that a
383  slave face occupies only a portion of its master face. */
385  const FaceInfo &fi, bool is_ghost);
386 
387  bool IsSlaveFace(const FaceInfo &fi) const;
388 
389  /// Returns the orientation of "test" relative to "base"
390  static int GetTriOrientation (const int * base, const int * test);
391  /// Returns the orientation of "test" relative to "base"
392  static int GetQuadOrientation (const int * base, const int * test);
393  /// Returns the orientation of "test" relative to "base"
394  static int GetTetOrientation (const int * base, const int * test);
395 
396  static void GetElementArrayEdgeTable(const Array<Element*> &elem_array,
397  const DSTable &v_to_v,
398  Table &el_to_edge);
399 
400  /** Return vertex to vertex table. The connections stored in the table
401  are from smaller to bigger vertex index, i.e. if i<j and (i, j) is
402  in the table, then (j, i) is not stored. */
403  void GetVertexToVertexTable(DSTable &) const;
404 
405  /** Return element to edge table and the indices for the boundary edges.
406  The entries in the table are ordered according to the order of the
407  nodes in the elements. For example, if T is the element to edge table
408  T(i, 0) gives the index of edge in element i that connects vertex 0
409  to vertex 1, etc. Returns the number of the edges. */
411 
412  /// Used in GenerateFaces()
413  void AddPointFaceElement(int lf, int gf, int el);
414 
415  void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1);
416 
417  void AddTriangleFaceElement (int lf, int gf, int el,
418  int v0, int v1, int v2);
419 
420  void AddQuadFaceElement (int lf, int gf, int el,
421  int v0, int v1, int v2, int v3);
422  /** For a serial Mesh, return true if the face is interior. For a parallel
423  ParMesh return true if the face is interior or shared. In parallel, this
424  method only works if the face neighbor data is exchanged. */
425  bool FaceIsTrueInterior(int FaceNo) const
426  {
427  return FaceIsInterior(FaceNo) || (faces_info[FaceNo].Elem2Inf >= 0);
428  }
429 
430  void FreeElement(Element *E);
431 
432  void GenerateFaces();
433  void GenerateNCFaceInfo();
434 
435  /// Begin construction of a mesh
436  void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem);
437 
438  // Used in the methods FinalizeXXXMesh() and FinalizeTopology()
439  void FinalizeCheck();
440 
441  void Loader(std::istream &input, int generate_edges = 0,
442  std::string parse_tag = "");
443 
444  // If NURBS mesh, write NURBS format. If NCMesh, write mfem v1.1 format.
445  // If section_delimiter is empty, write mfem v1.0 format. Otherwise, write
446  // mfem v1.2 format with the given section_delimiter at the end.
447  void Printer(std::ostream &out = mfem::out,
448  std::string section_delimiter = "") const;
449 
450  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
451  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
452  type=TETRAHEDRON. The parameter @a sfc_ordering controls how the elements
453  (when type=HEXAHEDRON) are ordered: true - use space-filling curve
454  ordering, or false - use lexicographic ordering. */
455  void Make3D(int nx, int ny, int nz, Element::Type type,
456  double sx, double sy, double sz, bool sfc_ordering);
457 
458  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
459  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
460  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
461  if 1 edges are generated. The parameter @a sfc_ordering controls how the
462  elements (when type=QUADRILATERAL) are ordered: true - use space-filling
463  curve ordering, or false - use lexicographic ordering. */
464  void Make2D(int nx, int ny, Element::Type type, double sx, double sy,
465  bool generate_edges, bool sfc_ordering);
466 
467  /// Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.
468  void Make1D(int n, double sx = 1.0);
469 
470  /// Internal function used in Mesh::MakeRefined
471  void MakeRefined_(Mesh &orig_mesh, const Array<int> ref_factors,
472  int ref_type);
473 
474  /// Initialize vertices/elements/boundary/tables from a nonconforming mesh.
475  void InitFromNCMesh(const NCMesh &ncmesh);
476 
477  /// Create from a nonconforming mesh.
478  explicit Mesh(const NCMesh &ncmesh);
479 
480  // used in GetElementData() and GetBdrElementData()
481  void GetElementData(const Array<Element*> &elem_array, int geom,
482  Array<int> &elem_vtx, Array<int> &attr) const;
483 
484  double GetElementSize(ElementTransformation *T, int type = 0);
485 
486  // Internal helper used in MakeSimplicial (and ParMesh::MakeSimplicial).
487  void MakeSimplicial_(const Mesh &orig_mesh, int *vglobal);
488 
489 public:
490 
491  Mesh() { SetEmpty(); }
492 
493  /** Copy constructor. Performs a deep copy of (almost) all data, so that the
494  source mesh can be modified (e.g. deleted, refined) without affecting the
495  new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the
496  nodes, if present. */
497  explicit Mesh(const Mesh &mesh, bool copy_nodes = true);
498 
499  /// Move constructor, useful for using a Mesh as a function return value.
500  Mesh(Mesh &&mesh);
501 
502  /// Move assignment operstor.
503  Mesh& operator=(Mesh &&mesh);
504 
505  /// Explicitly delete the copy assignment operator.
506  Mesh& operator=(Mesh &mesh) = delete;
507 
508  /** @name Named mesh constructors.
509 
510  Each of these constructors uses the move constructor, and can be used as
511  the right-hand side of an assignment when creating new meshes. */
512  ///@{
513 
514  /** Creates mesh by reading a file in MFEM, Netgen, or VTK format. If
515  generate_edges = 0 (default) edges are not generated, if 1 edges are
516  generated. */
517  static Mesh LoadFromFile(const char *filename,
518  int generate_edges = 0, int refine = 1,
519  bool fix_orientation = true);
520 
521  /** Creates 1D mesh , divided into n equal intervals. */
522  static Mesh MakeCartesian1D(int n, double sx = 1.0);
523 
524  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
525  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
526  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
527  if 1 edges are generated. If scf_ordering = true (default), elements are
528  ordered along a space-filling curve, instead of row by row. */
529  static Mesh MakeCartesian2D(
530  int nx, int ny, Element::Type type, bool generate_edges = false,
531  double sx = 1.0, double sy = 1.0, bool sfc_ordering = true);
532 
533  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
534  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
535  type=TETRAHEDRON. If sfc_ordering = true (default), elements are ordered
536  along a space-filling curve, instead of row by row and layer by layer. */
537  static Mesh MakeCartesian3D(
538  int nx, int ny, int nz, Element::Type type,
539  double sx = 1.0, double sy = 1.0, double sz = 1.0,
540  bool sfc_ordering = true);
541 
542  /// Create a refined (by any factor) version of @a orig_mesh.
543  /** @param[in] orig_mesh The starting coarse mesh.
544  @param[in] ref_factor The refinement factor, an integer > 1.
545  @param[in] ref_type Specify the positions of the new vertices. The
546  options are BasisType::ClosedUniform or
547  BasisType::GaussLobatto.
548 
549  The refinement data which can be accessed with GetRefinementTransforms()
550  is set to reflect the performed refinements.
551 
552  @note The constructed Mesh is straight-sided. */
553  static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type);
554 
555  /// Create a refined mesh, where each element of the original mesh may be
556  /// refined by a different factor.
557  /** @param[in] orig_mesh The starting coarse mesh.
558  @param[in] ref_factors An array of integers whose size is the number of
559  elements of @a orig_mesh. The @a ith element of
560  @a orig_mesh is refined by refinement factor
561  @a ref_factors[i].
562  @param[in] ref_type Specify the positions of the new vertices. The
563  options are BasisType::ClosedUniform or
564  BasisType::GaussLobatto.
565 
566  The refinement data which can be accessed with GetRefinementTransforms()
567  is set to reflect the performed refinements.
568 
569  @note The constructed Mesh is straight-sided. */
570  /// refined @a ref_factors[i] times in each dimension.
571  static Mesh MakeRefined(Mesh &orig_mesh, const Array<int> &ref_factors,
572  int ref_type);
573 
574  /** Create a mesh by splitting each element of @a orig_mesh into simplices.
575  Quadrilaterals are split into two triangles, prisms are split into
576  3 tetrahedra, and hexahedra are split into either 5 or 6 tetrahedra
577  depending on the configuration.
578  @warning The curvature of the original mesh is not carried over to the
579  new mesh. Periodic meshes are not supported. */
580  static Mesh MakeSimplicial(const Mesh &orig_mesh);
581 
582  /// Create a periodic mesh by identifying vertices of @a orig_mesh.
583  /** Each vertex @a i will be mapped to vertex @a v2v[i], such that all
584  vertices that are coincident under the periodic mapping get mapped to
585  the same index. The mapping @a v2v can be generated from translation
586  vectors using Mesh::CreatePeriodicVertexMapping.
587  @note MFEM requires that each edge of the resulting mesh be uniquely
588  identifiable by a pair of distinct vertices. As a consequence, periodic
589  boundaries must be connected by at least three edges. */
590  static Mesh MakePeriodic(const Mesh &orig_mesh, const std::vector<int> &v2v);
591 
592  ///@}
593 
594  /// @brief Creates a mapping @a v2v from the vertex indices of the mesh such
595  /// that coincident vertices under the given @a translations are identified.
596  /** Each Vector in @a translations should be of size @a sdim (the spatial
597  dimension of the mesh). Two vertices are considered coincident if the
598  translated coordinates of one vertex are within the given tolerance (@a
599  tol, relative to the mesh diameter) of the coordinates of the other
600  vertex.
601  @warning This algorithm does not scale well with the number of boundary
602  vertices in the mesh, and may run slowly on very large meshes. */
603  std::vector<int> CreatePeriodicVertexMapping(
604  const std::vector<Vector> &translations, double tol = 1e-8) const;
605 
606  /// Construct a Mesh from the given primary data.
607  /** The array @a vertices is used as external data, i.e. the Mesh does not
608  copy the data and will not delete the pointer.
609 
610  The data from the other arrays is copied into the internal Mesh data
611  structures.
612 
613  This method calls the method FinalizeTopology(). The method Finalize()
614  may be called after this constructor and after optionally setting the
615  Mesh nodes. */
616  Mesh(double *vertices, int num_vertices,
617  int *element_indices, Geometry::Type element_type,
618  int *element_attributes, int num_elements,
619  int *boundary_indices, Geometry::Type boundary_type,
620  int *boundary_attributes, int num_boundary_elements,
621  int dimension, int space_dimension = -1);
622 
623  /** @anchor mfem_Mesh_init_ctor
624  @brief _Init_ constructor: begin the construction of a Mesh object. */
625  Mesh(int Dim_, int NVert, int NElem, int NBdrElem = 0, int spaceDim_ = -1)
626  {
627  if (spaceDim_ == -1) { spaceDim_ = Dim_; }
628  InitMesh(Dim_, spaceDim_, NVert, NElem, NBdrElem);
629  }
630 
631  /** @name Methods for Mesh construction.
632 
633  These methods are intended to be used with the @ref mfem_Mesh_init_ctor
634  "init constructor". */
635  ///@{
636 
637  Element *NewElement(int geom);
638 
639  int AddVertex(double x, double y = 0.0, double z = 0.0);
640  int AddVertex(const double *coords);
641  /// Mark vertex @a i as non-conforming, with parent vertices @a p1 and @a p2.
642  void AddVertexParents(int i, int p1, int p2);
643 
644  int AddSegment(int v1, int v2, int attr = 1);
645  int AddSegment(const int *vi, int attr = 1);
646 
647  int AddTriangle(int v1, int v2, int v3, int attr = 1);
648  int AddTriangle(const int *vi, int attr = 1);
649  int AddTri(const int *vi, int attr = 1) { return AddTriangle(vi, attr); }
650 
651  int AddQuad(int v1, int v2, int v3, int v4, int attr = 1);
652  int AddQuad(const int *vi, int attr = 1);
653 
654  int AddTet(int v1, int v2, int v3, int v4, int attr = 1);
655  int AddTet(const int *vi, int attr = 1);
656 
657  int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr = 1);
658  int AddWedge(const int *vi, int attr = 1);
659 
660  int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8,
661  int attr = 1);
662  int AddHex(const int *vi, int attr = 1);
663  void AddHexAsTets(const int *vi, int attr = 1);
664  void AddHexAsWedges(const int *vi, int attr = 1);
665 
666  /// The parameter @a elem should be allocated using the NewElement() method
667  int AddElement(Element *elem);
668  int AddBdrElement(Element *elem);
669 
670  int AddBdrSegment(int v1, int v2, int attr = 1);
671  int AddBdrSegment(const int *vi, int attr = 1);
672 
673  int AddBdrTriangle(int v1, int v2, int v3, int attr = 1);
674  int AddBdrTriangle(const int *vi, int attr = 1);
675 
676  int AddBdrQuad(int v1, int v2, int v3, int v4, int attr = 1);
677  int AddBdrQuad(const int *vi, int attr = 1);
678  void AddBdrQuadAsTriangles(const int *vi, int attr = 1);
679 
680  int AddBdrPoint(int v, int attr = 1);
681 
683  /// Finalize the construction of a triangular Mesh.
684  void FinalizeTriMesh(int generate_edges = 0, int refine = 0,
685  bool fix_orientation = true);
686  /// Finalize the construction of a quadrilateral Mesh.
687  void FinalizeQuadMesh(int generate_edges = 0, int refine = 0,
688  bool fix_orientation = true);
689  /// Finalize the construction of a tetrahedral Mesh.
690  void FinalizeTetMesh(int generate_edges = 0, int refine = 0,
691  bool fix_orientation = true);
692  /// Finalize the construction of a wedge Mesh.
693  void FinalizeWedgeMesh(int generate_edges = 0, int refine = 0,
694  bool fix_orientation = true);
695  /// Finalize the construction of a hexahedral Mesh.
696  void FinalizeHexMesh(int generate_edges = 0, int refine = 0,
697  bool fix_orientation = true);
698  /// Finalize the construction of any type of Mesh.
699  /** This method calls FinalizeTopology() and Finalize(). */
700  void FinalizeMesh(int refine = 0, bool fix_orientation = true);
701 
702  ///@}
703 
704  /** @brief Finalize the construction of the secondary topology (connectivity)
705  data of a Mesh. */
706  /** This method does not require any actual coordinate data (either vertex
707  coordinates for linear meshes or node coordinates for meshes with nodes)
708  to be available. However, the data generated by this method is generally
709  required by the FiniteElementSpace class.
710 
711  After calling this method, setting the Mesh vertices or nodes, it may be
712  appropriate to call the method Finalize(). */
713  void FinalizeTopology(bool generate_bdr = true);
714 
715  /// Finalize the construction of a general Mesh.
716  /** This method will:
717  - check and optionally fix the orientation of regular elements
718  - check and fix the orientation of boundary elements
719  - assume that #vertices are defined, if #Nodes == NULL
720  - assume that #Nodes are defined, if #Nodes != NULL.
721  @param[in] refine If true, prepare the Mesh for conforming refinement of
722  triangular or tetrahedral meshes.
723  @param[in] fix_orientation
724  If true, fix the orientation of inverted mesh elements
725  by permuting their vertices.
726 
727  Before calling this method, call FinalizeTopology() and ensure that the
728  Mesh vertices or nodes are set. */
729  virtual void Finalize(bool refine = false, bool fix_orientation = false);
730 
731  virtual void SetAttributes();
732 
733  /** This is our integration with the Gecko library. The method finds an
734  element ordering that will increase memory coherency by putting elements
735  that are in physical proximity closer in memory. It can also be used to
736  obtain a space-filling curve ordering for ParNCMesh partitioning.
737  @param[out] ordering Output element ordering.
738  @param iterations Total number of V cycles. The ordering may improve with
739  more iterations. The best iteration is returned at the end.
740  @param window Initial window size. This determines the number of
741  permutations tested at each multigrid level and strongly influences the
742  quality of the result, but the cost of increasing 'window' is exponential.
743  @param period The window size is incremented every 'period' iterations.
744  @param seed Seed for initial random ordering (0 = skip random reorder).
745  @param verbose Print the progress of the optimization to mfem::out.
746  @param time_limit Optional time limit for the optimization, in seconds.
747  When reached, ordering from the best iteration so far is returned
748  (0 = no limit).
749  @return The final edge product cost of the ordering. The function may be
750  called in an external loop with different seeds, and the best ordering can
751  then be retained. */
752  double GetGeckoElementOrdering(Array<int> &ordering,
753  int iterations = 4, int window = 4,
754  int period = 2, int seed = 0,
755  bool verbose = false, double time_limit = 0);
756 
757  /** Return an ordering of the elements that approximately follows the Hilbert
758  curve. The method performs a spatial (Hilbert) sort on the centers of all
759  elements and returns the resulting sequence, which can then be passed to
760  ReorderElements. This is a cheap alternative to GetGeckoElementOrdering.*/
761  void GetHilbertElementOrdering(Array<int> &ordering);
762 
763  /** Rebuilds the mesh with a different order of elements. For each element i,
764  the array ordering[i] contains its desired new index. Note that the method
765  reorders vertices, edges and faces along with the elements. */
766  void ReorderElements(const Array<int> &ordering, bool reorder_vertices = true);
767 
768  /// Deprecated: see @a MakeCartesian3D.
769  MFEM_DEPRECATED
770  Mesh(int nx, int ny, int nz, Element::Type type, bool generate_edges = false,
771  double sx = 1.0, double sy = 1.0, double sz = 1.0,
772  bool sfc_ordering = true)
773  {
774  Make3D(nx, ny, nz, type, sx, sy, sz, sfc_ordering);
775  Finalize(true); // refine = true
776  }
777 
778  /// Deprecated: see @a MakeCartesian2D.
779  MFEM_DEPRECATED
780  Mesh(int nx, int ny, Element::Type type, bool generate_edges = false,
781  double sx = 1.0, double sy = 1.0, bool sfc_ordering = true)
782  {
783  Make2D(nx, ny, type, sx, sy, generate_edges, sfc_ordering);
784  Finalize(true); // refine = true
785  }
786 
787  /// Deprecated: see @a MakeCartesian1D.
788  MFEM_DEPRECATED
789  explicit Mesh(int n, double sx = 1.0)
790  {
791  Make1D(n, sx);
792  // Finalize(); // reminder: not needed
793  }
794 
795  /** Creates mesh by reading a file in MFEM, Netgen, or VTK format. If
796  generate_edges = 0 (default) edges are not generated, if 1 edges are
797  generated. See also @a Mesh::LoadFromFile. */
798  explicit Mesh(const char *filename, int generate_edges = 0, int refine = 1,
799  bool fix_orientation = true);
800 
801  /** Creates mesh by reading data stream in MFEM, Netgen, or VTK format. If
802  generate_edges = 0 (default) edges are not generated, if 1 edges are
803  generated. */
804  explicit Mesh(std::istream &input, int generate_edges = 0, int refine = 1,
805  bool fix_orientation = true);
806 
807  /// Create a disjoint mesh from the given mesh array
808  Mesh(Mesh *mesh_array[], int num_pieces);
809 
810  /// Deprecated: see @a MakeRefined.
811  MFEM_DEPRECATED
812  Mesh(Mesh *orig_mesh, int ref_factor, int ref_type);
813 
814  /** This is similar to the mesh constructor with the same arguments, but here
815  the current mesh is destroyed and another one created based on the data
816  stream again given in MFEM, Netgen, or VTK format. If generate_edges = 0
817  (default) edges are not generated, if 1 edges are generated. */
818  /// \see mfem::ifgzstream() for on-the-fly decompression of compressed ascii
819  /// inputs.
820  virtual void Load(std::istream &input, int generate_edges = 0,
821  int refine = 1, bool fix_orientation = true)
822  {
823  Loader(input, generate_edges);
824  Finalize(refine, fix_orientation);
825  }
826 
827  /// Clear the contents of the Mesh.
828  void Clear() { Destroy(); SetEmpty(); }
829 
830  /** @brief Get the mesh generator/type.
831 
832  @return A bitmask:
833  - bit 0 - simplices are present in the mesh (triangles, tets),
834  - bit 1 - tensor product elements are present in the mesh (quads, hexes),
835  - bit 2 - the mesh has wedge elements.
836 
837  In parallel, the result takes into account elements on all processors.
838  */
839  inline int MeshGenerator() { return meshgen; }
840 
841  /** @brief Returns number of vertices. Vertices are only at the corners of
842  elements, where you would expect them in the lowest-order mesh. */
843  inline int GetNV() const { return NumOfVertices; }
844 
845  /// Returns number of elements.
846  inline int GetNE() const { return NumOfElements; }
847 
848  /// Returns number of boundary elements.
849  inline int GetNBE() const { return NumOfBdrElements; }
850 
851  /// Return the number of edges.
852  inline int GetNEdges() const { return NumOfEdges; }
853 
854  /// Return the number of faces in a 3D mesh.
855  inline int GetNFaces() const { return NumOfFaces; }
856 
857  /// Return the number of faces (3D), edges (2D) or vertices (1D).
858  int GetNumFaces() const;
859 
860  /// Returns the number of faces according to the requested type.
861  /** If type==Boundary returns only the "true" number of boundary faces
862  contrary to GetNBE() that returns "fake" boundary faces associated to
863  visualization for GLVis.
864  Similarly, if type==Interior, the "fake" boundary faces associated to
865  visualization are counted as interior faces. */
866  int GetNFbyType(FaceType type) const;
867 
868  /// Utility function: sum integers from all processors (Allreduce).
869  virtual long ReduceInt(int value) const { return value; }
870 
871  /// Return the total (global) number of elements.
872  long GetGlobalNE() const { return ReduceInt(NumOfElements); }
873 
874  /** @brief Return the mesh geometric factors corresponding to the given
875  integration rule.
876 
877  The IntegrationRule used with GetGeometricFactors needs to remain valid
878  until the internally stored GeometricFactors objects are destroyed (by
879  either calling Mesh::DeleteGeometricFactors or the Mesh destructor). If
880  the device MemoryType parameter @a d_mt is specified, then the returned
881  object will use that type unless it was previously allocated with a
882  different type. */
884  const IntegrationRule& ir,
885  const int flags,
887 
888  /** @brief Return the mesh geometric factors for the faces corresponding
889  to the given integration rule.
890 
891  The IntegrationRule used with GetFaceGeometricFactors needs to remain
892  valid until the internally stored FaceGeometricFactors objects are
893  destroyed (by either calling Mesh::DeleteGeometricFactors or the Mesh
894  destructor). */
896  const int flags,
897  FaceType type);
898 
899  /// Destroy all GeometricFactors stored by the Mesh.
900  /** This method can be used to force recomputation of the GeometricFactors,
901  for example, after the mesh nodes are modified externally. */
902  void DeleteGeometricFactors();
903 
904  /// Equals 1 + num_holes - num_loops
905  inline int EulerNumber() const
907  /// Equals 1 - num_holes
908  inline int EulerNumber2D() const
909  { return NumOfVertices - NumOfEdges + NumOfElements; }
910 
911  int Dimension() const { return Dim; }
912  int SpaceDimension() const { return spaceDim; }
913 
914  /// @brief Return pointer to vertex i's coordinates.
915  /// @warning For high-order meshes (when #Nodes != NULL) vertices may not be
916  /// updated and should not be used!
917  const double *GetVertex(int i) const { return vertices[i](); }
918 
919  /// @brief Return pointer to vertex i's coordinates.
920  /// @warning For high-order meshes (when Nodes != NULL) vertices may not
921  /// being updated and should not be used!
922  double *GetVertex(int i) { return vertices[i](); }
923 
924  void GetElementData(int geom, Array<int> &elem_vtx, Array<int> &attr) const
925  { GetElementData(elements, geom, elem_vtx, attr); }
926 
927  void GetBdrElementData(int geom, Array<int> &bdr_elem_vtx,
928  Array<int> &bdr_attr) const
929  { GetElementData(boundary, geom, bdr_elem_vtx, bdr_attr); }
930 
931  /** @brief Set the internal Vertex array to point to the given @a vertices
932  array without assuming ownership of the pointer. */
933  /** If @a zerocopy is `true`, the vertices must be given as an array of 3
934  doubles per vertex. If @a zerocopy is `false` then the current Vertex
935  data is first copied to the @a vertices array. */
936  void ChangeVertexDataOwnership(double *vertices, int len_vertices,
937  bool zerocopy = false);
938 
939  const Element* const *GetElementsArray() const
940  { return elements.GetData(); }
941 
942  const Element *GetElement(int i) const { return elements[i]; }
943 
944  Element *GetElement(int i) { return elements[i]; }
945 
946  const Element *GetBdrElement(int i) const { return boundary[i]; }
947 
948  Element *GetBdrElement(int i) { return boundary[i]; }
949 
950  const Element *GetFace(int i) const { return faces[i]; }
951 
953  {
954  return faces[i]->GetGeometryType();
955  }
956 
958  {
959  return elements[i]->GetGeometryType();
960  }
961 
963  {
964  return boundary[i]->GetGeometryType();
965  }
966 
967  // deprecated: "base geometry" no longer means anything
969  { return GetFaceGeometry(i); }
970 
972  { return GetElementGeometry(i); }
973 
975  { return GetBdrElementGeometry(i); }
976 
977  /** @brief Return true iff the given @a geom is encountered in the mesh.
978  Geometries of dimensions lower than Dimension() are counted as well. */
980  { return mesh_geoms & (1 << geom); }
981 
982  /** @brief Return the number of geometries of the given dimension present in
983  the mesh. */
984  /** For a parallel mesh only the local geometries are counted. */
985  int GetNumGeometries(int dim) const;
986 
987  /// Return all element geometries of the given dimension present in the mesh.
988  /** For a parallel mesh only the local geometries are returned.
989 
990  The returned geometries are sorted. */
991  void GetGeometries(int dim, Array<Geometry::Type> &el_geoms) const;
992 
993  /// List of mesh geometries stored as Array<Geometry::Type>.
994  class GeometryList : public Array<Geometry::Type>
995  {
996  protected:
998  public:
999  /// Construct a GeometryList of all element geometries in @a mesh.
1001  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
1002  { mesh.GetGeometries(mesh.Dimension(), *this); }
1003  /** @brief Construct a GeometryList of all geometries of dimension @a dim
1004  in @a mesh. */
1005  GeometryList(Mesh &mesh, int dim)
1006  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
1007  { mesh.GetGeometries(dim, *this); }
1008  };
1009 
1010  /// Returns the indices of the vertices of element i.
1011  void GetElementVertices(int i, Array<int> &v) const
1012  { elements[i]->GetVertices(v); }
1013 
1014  /// Returns the indices of the vertices of boundary element i.
1015  void GetBdrElementVertices(int i, Array<int> &v) const
1016  { boundary[i]->GetVertices(v); }
1017 
1018  /// Return the indices and the orientations of all edges of element i.
1019  void GetElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
1020 
1021  /// Return the indices and the orientations of all edges of bdr element i.
1022  void GetBdrElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
1023 
1024  /** Return the indices and the orientations of all edges of face i.
1025  Works for both 2D (face=edge) and 3D faces. */
1026  void GetFaceEdges(int i, Array<int> &edges, Array<int> &o) const;
1027 
1028  /// Returns the indices of the vertices of face i.
1029  void GetFaceVertices(int i, Array<int> &vert) const
1030  {
1031  if (Dim == 1)
1032  {
1033  vert.SetSize(1); vert[0] = i;
1034  }
1035  else
1036  {
1037  faces[i]->GetVertices(vert);
1038  }
1039  }
1040 
1041  /// Returns the indices of the vertices of edge i.
1042  void GetEdgeVertices(int i, Array<int> &vert) const;
1043 
1044  /// Returns the face-to-edge Table (3D)
1045  Table *GetFaceEdgeTable() const;
1046 
1047  /// Returns the edge-to-vertex Table (3D)
1048  Table *GetEdgeVertexTable() const;
1049 
1050  /// Return the indices and the orientations of all faces of element i.
1051  void GetElementFaces(int i, Array<int> &faces, Array<int> &ori) const;
1052 
1053  /// Return the index and the orientation of the face of bdr element i. (3D)
1054  void GetBdrElementFace(int i, int *f, int *o) const;
1055 
1056  /** Return the vertex index of boundary element i. (1D)
1057  Return the edge index of boundary element i. (2D)
1058  Return the face index of boundary element i. (3D) */
1059  int GetBdrElementEdgeIndex(int i) const;
1060 
1061  /** @brief For the given boundary element, bdr_el, return its adjacent
1062  element and its info, i.e. 64*local_bdr_index+bdr_orientation. */
1063  void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const;
1064 
1065  /// Returns the type of element i.
1066  Element::Type GetElementType(int i) const;
1067 
1068  /// Returns the type of boundary element i.
1069  Element::Type GetBdrElementType(int i) const;
1070 
1071  /* Return point matrix of element i of dimension Dim X #v, where for every
1072  vertex we give its coordinates in space of dimension Dim. */
1073  void GetPointMatrix(int i, DenseMatrix &pointmat) const;
1074 
1075  /* Return point matrix of boundary element i of dimension Dim X #v, where for
1076  every vertex we give its coordinates in space of dimension Dim. */
1077  void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const;
1078 
1080 
1081  /** Builds the transformation defining the i-th element in the user-defined
1082  variable. */
1084 
1085  /// Returns the transformation defining the i-th element
1087 
1088  /** Return the transformation defining the i-th element assuming
1089  the position of the vertices/nodes are given by 'nodes'. */
1090  void GetElementTransformation(int i, const Vector &nodes,
1092 
1093  /// Returns the transformation defining the i-th boundary element
1096 
1097  /** @brief Returns the transformation defining the given face element in a
1098  user-defined variable. */
1100 
1101  /** @brief A helper method that constructs a transformation from the
1102  reference space of a face to the reference space of an element. */
1103  /** The local index of the face as a face in the element and its orientation
1104  are given by the input parameter @a info, as @a info = 64*loc_face_idx +
1105  loc_face_orientation. */
1106  void GetLocalFaceTransformation(int face_type, int elem_type,
1108  int info);
1109 
1110  /// Returns the transformation defining the given face element
1112 
1113  /** Returns the transformation defining the given edge element.
1114  The transformation is stored in a user-defined variable. */
1116 
1117  /// Returns the transformation defining the given face element
1119 
1120  /// Returns (a pointer to an object containing) the following data:
1121  ///
1122  /// 1) Elem1No - the index of the first element that contains this face this
1123  /// is the element that has the same outward unit normal vector as the
1124  /// face;
1125  ///
1126  /// 2) Elem2No - the index of the second element that contains this face this
1127  /// element has outward unit normal vector as the face multiplied with -1;
1128  ///
1129  /// 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first
1130  /// and the second element respectively;
1131  ///
1132  /// 4) Face - pointer to the ElementTransformation of the face;
1133  ///
1134  /// 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face
1135  /// coordinate system to the element coordinate system (both in their
1136  /// reference elements). Used to transform IntegrationPoints from face to
1137  /// element. More formally, let:
1138  /// TL1, TL2 be the transformations represented by Loc1, Loc2,
1139  /// TE1, TE2 - the transformations represented by Elem1, Elem2,
1140  /// TF - the transformation represented by Face, then
1141  /// TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
1142  ///
1143  /// 6) FaceGeom - the base geometry for the face.
1144  ///
1145  /// The mask specifies which fields in the structure to return:
1146  /// mask & 1 - Elem1, mask & 2 - Elem2
1147  /// mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.
1148  /// These mask values are defined in the ConfigMasks enum type as part of the
1149  /// FaceElementTransformations class in fem/eltrans.hpp.
1151  int mask = 31);
1152 
1154  {
1155  if (faces_info[FaceNo].Elem2No < 0) { return NULL; }
1156  return GetFaceElementTransformations (FaceNo);
1157  }
1158 
1160 
1161  /// Return the local face index for the given boundary face.
1162  int GetBdrFace(int BdrElemNo) const;
1163 
1164  /// Return true if the given face is interior. @sa FaceIsTrueInterior().
1165  bool FaceIsInterior(int FaceNo) const
1166  {
1167  return (faces_info[FaceNo].Elem2No >= 0);
1168  }
1169  void GetFaceElements (int Face, int *Elem1, int *Elem2) const;
1170  void GetFaceInfos (int Face, int *Inf1, int *Inf2) const;
1171  void GetFaceInfos (int Face, int *Inf1, int *Inf2, int *NCFace) const;
1172 
1173  Geometry::Type GetFaceGeometryType(int Face) const;
1174  Element::Type GetFaceElementType(int Face) const;
1175 
1176  /// Check (and optionally attempt to fix) the orientation of the elements
1177  /** @param[in] fix_it If `true`, attempt to fix the orientations of some
1178  elements: triangles, quads, and tets.
1179  @return The number of elements with wrong orientation.
1180 
1181  @note For meshes with nodes (e.g. high-order or periodic meshes), fixing
1182  the element orientations may require additional permutation of the nodal
1183  GridFunction of the mesh which is not performed by this method. Instead,
1184  the method Finalize() should be used with the parameter
1185  @a fix_orientation set to `true`.
1186 
1187  @note This method performs a simple check if an element is inverted, e.g.
1188  for most elements types, it checks if the Jacobian of the mapping from
1189  the reference element is non-negative at the center of the element. */
1190  int CheckElementOrientation(bool fix_it = true);
1191 
1192  /// Check the orientation of the boundary elements
1193  /** @return The number of boundary elements with wrong orientation. */
1194  int CheckBdrElementOrientation(bool fix_it = true);
1195 
1196  /// Return the attribute of element i.
1197  int GetAttribute(int i) const { return elements[i]->GetAttribute(); }
1198 
1199  /// Set the attribute of element i.
1200  void SetAttribute(int i, int attr) { elements[i]->SetAttribute(attr); }
1201 
1202  /// Return the attribute of boundary element i.
1203  int GetBdrAttribute(int i) const { return boundary[i]->GetAttribute(); }
1204 
1205  /// Set the attribute of boundary element i.
1206  void SetBdrAttribute(int i, int attr) { boundary[i]->SetAttribute(attr); }
1207 
1208  const Table &ElementToElementTable();
1209 
1210  const Table &ElementToFaceTable() const;
1211 
1212  const Table &ElementToEdgeTable() const;
1213 
1214  /// The returned Table must be destroyed by the caller
1216 
1217  /** Return the "face"-element Table. Here "face" refers to face (3D),
1218  edge (2D), or vertex (1D).
1219  The returned Table must be destroyed by the caller. */
1220  Table *GetFaceToElementTable() const;
1221 
1222  /** This method modifies a tetrahedral mesh so that Nedelec spaces of order
1223  greater than 1 can be defined on the mesh. Specifically, we
1224  1) rotate all tets in the mesh so that the vertices {v0, v1, v2, v3}
1225  satisfy: v0 < v1 < min(v2, v3).
1226  2) rotate all boundary triangles so that the vertices {v0, v1, v2}
1227  satisfy: v0 < min(v1, v2).
1228 
1229  @note Refinement does not work after a call to this method! */
1230  virtual void ReorientTetMesh();
1231 
1232  int *CartesianPartitioning(int nxyz[]);
1233  int *GeneratePartitioning(int nparts, int part_method = 1);
1234  void CheckPartitioning(int *partitioning_);
1235 
1236  void CheckDisplacements(const Vector &displacements, double &tmax);
1237 
1238  // Vertices are only at the corners of elements, where you would expect them
1239  // in the lowest-order mesh.
1240  void MoveVertices(const Vector &displacements);
1241  void GetVertices(Vector &vert_coord) const;
1242  void SetVertices(const Vector &vert_coord);
1243 
1244  // Nodes are only active for higher order meshes, and share locations with
1245  // the vertices, plus all the higher- order control points within the element
1246  // and along the edges and on the faces.
1247  void GetNode(int i, double *coord) const;
1248  void SetNode(int i, const double *coord);
1249 
1250  // Node operations for curved mesh.
1251  // They call the corresponding '...Vertices' method if the
1252  // mesh is not curved (i.e. Nodes == NULL).
1253  void MoveNodes(const Vector &displacements);
1254  void GetNodes(Vector &node_coord) const;
1255  void SetNodes(const Vector &node_coord);
1256 
1257  /// Return a pointer to the internal node GridFunction (may be NULL).
1258  GridFunction *GetNodes() { return Nodes; }
1259  const GridFunction *GetNodes() const { return Nodes; }
1260  /// Return the mesh nodes ownership flag.
1261  bool OwnsNodes() const { return own_nodes; }
1262  /// Set the mesh nodes ownership flag.
1263  void SetNodesOwner(bool nodes_owner) { own_nodes = nodes_owner; }
1264  /// Replace the internal node GridFunction with the given GridFunction.
1265  void NewNodes(GridFunction &nodes, bool make_owner = false);
1266  /** Swap the internal node GridFunction pointer and ownership flag members
1267  with the given ones. */
1268  void SwapNodes(GridFunction *&nodes, int &own_nodes_);
1269 
1270  /// Return the mesh nodes/vertices projected on the given GridFunction.
1271  void GetNodes(GridFunction &nodes) const;
1272  /** Replace the internal node GridFunction with a new GridFunction defined
1273  on the given FiniteElementSpace. The new node coordinates are projected
1274  (derived) from the current nodes/vertices. */
1275  void SetNodalFESpace(FiniteElementSpace *nfes);
1276  /** Replace the internal node GridFunction with the given GridFunction. The
1277  given GridFunction is updated with node coordinates projected (derived)
1278  from the current nodes/vertices. */
1279  void SetNodalGridFunction(GridFunction *nodes, bool make_owner = false);
1280  /** Return the FiniteElementSpace on which the current mesh nodes are
1281  defined or NULL if the mesh does not have nodes. */
1282  const FiniteElementSpace *GetNodalFESpace() const;
1283  /** Make sure that the mesh has valid nodes, i.e. its geometry is described
1284  by a vector finite element grid function (even if it is a low-order mesh
1285  with straight edges). */
1286  void EnsureNodes();
1287 
1288  /** Set the curvature of the mesh nodes using the given polynomial degree,
1289  'order', and optionally: discontinuous or continuous FE space, 'discont',
1290  new space dimension, 'space_dim' (if != -1), and 'ordering'. */
1291  virtual void SetCurvature(int order, bool discont = false, int space_dim = -1,
1292  int ordering = 1);
1293 
1294  /// Refine all mesh elements.
1295  /** @param[in] ref_algo %Refinement algorithm. Currently used only for pure
1296  tetrahedral meshes. If set to zero (default), a tet mesh will be refined
1297  using algorithm A, that produces elements with better quality compared to
1298  algorithm B used when the parameter is non-zero.
1299 
1300  For tetrahedral meshes, after using algorithm A, the mesh cannot be
1301  refined locally using methods like GeneralRefinement() unless it is
1302  re-finalized using Finalize() with the parameter @a refine set to true.
1303  Note that calling Finalize() in this way will generally invalidate any
1304  FiniteElementSpace%s and GridFunction%s defined on the mesh. */
1305  void UniformRefinement(int ref_algo = 0);
1306 
1307  /** Refine selected mesh elements. Refinement type can be specified for each
1308  element. The function can do conforming refinement of triangles and
1309  tetrahedra and non-conforming refinement (i.e., with hanging-nodes) of
1310  triangles, quadrilaterals and hexahedra. If 'nonconforming' = -1,
1311  suitable refinement method is selected automatically (namely, conforming
1312  refinement for triangles). Use nonconforming = 0/1 to force the method.
1313  For nonconforming refinements, nc_limit optionally specifies the maximum
1314  level of hanging nodes (unlimited by default). */
1315  void GeneralRefinement(const Array<Refinement> &refinements,
1316  int nonconforming = -1, int nc_limit = 0);
1317 
1318  /** Simplified version of GeneralRefinement taking a simple list of elements
1319  to refine, without refinement types. */
1320  void GeneralRefinement(const Array<int> &el_to_refine,
1321  int nonconforming = -1, int nc_limit = 0);
1322 
1323  /// Refine each element with given probability. Uses GeneralRefinement.
1324  void RandomRefinement(double prob, bool aniso = false,
1325  int nonconforming = -1, int nc_limit = 0);
1326 
1327  /// Refine elements sharing the specified vertex. Uses GeneralRefinement.
1328  void RefineAtVertex(const Vertex& vert,
1329  double eps = 0.0, int nonconforming = -1);
1330 
1331  /** Refine element i if elem_error[i] > threshold, for all i.
1332  Returns true if at least one element was refined, false otherwise. */
1333  bool RefineByError(const Array<double> &elem_error, double threshold,
1334  int nonconforming = -1, int nc_limit = 0);
1335 
1336  /** Refine element i if elem_error(i) > threshold, for all i.
1337  Returns true if at least one element was refined, false otherwise. */
1338  bool RefineByError(const Vector &elem_error, double threshold,
1339  int nonconforming = -1, int nc_limit = 0);
1340 
1341  /** Derefine the mesh based on an error measure associated with each
1342  element. A derefinement is performed if the sum of errors of its fine
1343  elements is smaller than 'threshold'. If 'nc_limit' > 0, derefinements
1344  that would increase the maximum level of hanging nodes of the mesh are
1345  skipped. Returns true if the mesh changed, false otherwise. */
1346  bool DerefineByError(Array<double> &elem_error, double threshold,
1347  int nc_limit = 0, int op = 1);
1348 
1349  /// Same as DerefineByError for an error vector.
1350  bool DerefineByError(const Vector &elem_error, double threshold,
1351  int nc_limit = 0, int op = 1);
1352 
1353  ///@{ @name NURBS mesh refinement methods
1354  void KnotInsert(Array<KnotVector *> &kv);
1355  void KnotInsert(Array<Vector *> &kv);
1356  /* For each knot vector:
1357  new_degree = max(old_degree, min(old_degree + rel_degree, degree)). */
1358  void DegreeElevate(int rel_degree, int degree = 16);
1359  ///@}
1360 
1361  /** Make sure that a quad/hex mesh is considered to be non-conforming (i.e.,
1362  has an associated NCMesh object). Simplex meshes can be both conforming
1363  (default) or non-conforming. */
1364  void EnsureNCMesh(bool simplices_nonconforming = false);
1365 
1366  bool Conforming() const { return ncmesh == NULL; }
1367  bool Nonconforming() const { return ncmesh != NULL; }
1368 
1369  /** Return fine element transformations following a mesh refinement.
1370  Space uses this to construct a global interpolation matrix. */
1372 
1373  /// Return type of last modification of the mesh.
1375 
1376  /** Return update counter. The counter starts at zero and is incremented
1377  each time refinement, derefinement, or rebalancing method is called.
1378  It is used for checking proper sequence of Space:: and GridFunction::
1379  Update() calls. */
1380  long GetSequence() const { return sequence; }
1381 
1382  /// Print the mesh to the given stream using Netgen/Truegrid format.
1383  virtual void PrintXG(std::ostream &out = mfem::out) const;
1384 
1385  /// Print the mesh to the given stream using the default MFEM mesh format.
1386  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1387  virtual void Print(std::ostream &out = mfem::out) const { Printer(out); }
1388 
1389  /// Save the mesh to a file using Mesh::Print. The given @a precision will be
1390  /// used for ASCII output.
1391  virtual void Save(const char *fname, int precision=16) const;
1392 
1393  /// Print the mesh to the given stream using the adios2 bp format
1394 #ifdef MFEM_USE_ADIOS2
1395  virtual void Print(adios2stream &out) const;
1396 #endif
1397  /// Print the mesh in VTK format (linear and quadratic meshes only).
1398  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1399  void PrintVTK(std::ostream &out);
1400  /** Print the mesh in VTK format. The parameter ref > 0 specifies an element
1401  subdivision number (useful for high order fields and curved meshes).
1402  If the optional field_data is set, we also add a FIELD section in the
1403  beginning of the file with additional dataset information. */
1404  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1405  void PrintVTK(std::ostream &out, int ref, int field_data=0);
1406  /** Print the mesh in VTU format. The parameter ref > 0 specifies an element
1407  subdivision number (useful for high order fields and curved meshes).
1408  If @a bdr_elements is true, then output (only) the boundary elements,
1409  otherwise output only the non-boundary elements. */
1410  void PrintVTU(std::ostream &out,
1411  int ref=1,
1412  VTKFormat format=VTKFormat::ASCII,
1413  bool high_order_output=false,
1414  int compression_level=0,
1415  bool bdr_elements=false);
1416  /** Print the mesh in VTU format with file name fname. */
1417  virtual void PrintVTU(std::string fname,
1418  VTKFormat format=VTKFormat::ASCII,
1419  bool high_order_output=false,
1420  int compression_level=0,
1421  bool bdr=false);
1422  /** Print the boundary elements of the mesh in VTU format, and output the
1423  boundary attributes as a data array (useful for boundary conditions). */
1424  void PrintBdrVTU(std::string fname,
1425  VTKFormat format=VTKFormat::ASCII,
1426  bool high_order_output=false,
1427  int compression_level=0);
1428 
1429  void GetElementColoring(Array<int> &colors, int el0 = 0);
1430 
1431  /** @brief Prints the mesh with boundary elements given by the boundary of
1432  the subdomains, so that the boundary of subdomain i has boundary
1433  attribute i+1. */
1434  /// \see mfem::ofgzstream() for on-the-fly compression of ascii outputs
1435  void PrintWithPartitioning (int *partitioning,
1436  std::ostream &out, int elem_attr = 0) const;
1437 
1438  void PrintElementsWithPartitioning (int *partitioning,
1439  std::ostream &out,
1440  int interior_faces = 0);
1441 
1442  /// Print set of disjoint surfaces:
1443  /*!
1444  * If Aface_face(i,j) != 0, print face j as a boundary
1445  * element with attribute i+1.
1446  */
1447  void PrintSurfaces(const Table &Aface_face, std::ostream &out) const;
1448 
1449  void ScaleSubdomains (double sf);
1450  void ScaleElements (double sf);
1451 
1452  void Transform(void (*f)(const Vector&, Vector&));
1453  void Transform(VectorCoefficient &deformation);
1454 
1455  /// Remove unused vertices and rebuild mesh connectivity.
1456  void RemoveUnusedVertices();
1457 
1458  /** Remove boundary elements that lie in the interior of the mesh, i.e. that
1459  have two adjacent faces in 3D, or edges in 2D. */
1460  void RemoveInternalBoundaries();
1461 
1462  /** @brief Get the size of the i-th element relative to the perfect
1463  reference element. */
1464  double GetElementSize(int i, int type = 0);
1465 
1466  double GetElementSize(int i, const Vector &dir);
1467 
1468  double GetElementVolume(int i);
1469 
1470  void GetElementCenter(int i, Vector &center);
1471 
1472  /// Returns the minimum and maximum corners of the mesh bounding box.
1473  /** For high-order meshes, the geometry is first refined @a ref times. */
1474  void GetBoundingBox(Vector &min, Vector &max, int ref = 2);
1475 
1476  void GetCharacteristics(double &h_min, double &h_max,
1477  double &kappa_min, double &kappa_max,
1478  Vector *Vh = NULL, Vector *Vk = NULL);
1479 
1480  /// Auxiliary method used by PrintCharacteristics().
1481  /** It is also used in the `mesh-explorer` miniapp. */
1482  static void PrintElementsByGeometry(int dim,
1483  const Array<int> &num_elems_by_geom,
1484  std::ostream &out);
1485 
1486  /** @brief Compute and print mesh characteristics such as number of vertices,
1487  number of elements, number of boundary elements, minimal and maximal
1488  element sizes, minimal and maximal element aspect ratios, etc. */
1489  /** If @a Vh or @a Vk are not NULL, return the element sizes and aspect
1490  ratios for all elements in the given Vector%s. */
1491  void PrintCharacteristics(Vector *Vh = NULL, Vector *Vk = NULL,
1492  std::ostream &out = mfem::out);
1493 
1494  /** @brief In serial, this method calls PrintCharacteristics(). In parallel,
1495  additional information about the parallel decomposition is also printed.
1496  */
1497  virtual void PrintInfo(std::ostream &out = mfem::out)
1498  {
1499  PrintCharacteristics(NULL, NULL, out);
1500  }
1501 
1502  void MesquiteSmooth(const int mesquite_option = 0);
1503 
1504  /** @brief Find the ids of the elements that contain the given points, and
1505  their corresponding reference coordinates.
1506 
1507  The DenseMatrix @a point_mat describes the given points - one point for
1508  each column; it should have SpaceDimension() rows.
1509 
1510  The InverseElementTransformation object, @a inv_trans, is used to attempt
1511  the element transformation inversion. If NULL pointer is given, the
1512  method will use a default constructed InverseElementTransformation. Note
1513  that the algorithms in the base class InverseElementTransformation can be
1514  completely overwritten by deriving custom classes that override the
1515  Transform() method.
1516 
1517  If no element is found for the i-th point, elem_ids[i] is set to -1.
1518 
1519  In the ParMesh implementation, the @a point_mat is expected to be the
1520  same on all ranks. If the i-th point is found by multiple ranks, only one
1521  of them will mark that point as found, i.e. set its elem_ids[i] to a
1522  non-negative number; the other ranks will set their elem_ids[i] to -2 to
1523  indicate that the point was found but assigned to another rank.
1524 
1525  @returns The total number of points that were found.
1526 
1527  @note This method is not 100 percent reliable, i.e. it is not guaranteed
1528  to find a point, even if it lies inside a mesh element. */
1529  virtual int FindPoints(DenseMatrix& point_mat, Array<int>& elem_ids,
1530  Array<IntegrationPoint>& ips, bool warn = true,
1531  InverseElementTransformation *inv_trans = NULL);
1532 
1533  /// Swaps internal data with another mesh. By default, non-geometry members
1534  /// like 'ncmesh' and 'NURBSExt' are only swapped when 'non_geometry' is set.
1535  void Swap(Mesh& other, bool non_geometry);
1536 
1537  /// Destroys Mesh.
1538  virtual ~Mesh() { DestroyPointers(); }
1539 
1540 #ifdef MFEM_DEBUG
1541  /// Output an NCMesh-compatible debug dump.
1542  void DebugDump(std::ostream &out) const;
1543 #endif
1544 };
1545 
1546 /** Overload operator<< for std::ostream and Mesh; valid also for the derived
1547  class ParMesh */
1548 std::ostream &operator<<(std::ostream &out, const Mesh &mesh);
1549 
1550 
1551 /** @brief Structure for storing mesh geometric factors: coordinates, Jacobians,
1552  and determinants of the Jacobians. */
1553 /** Typically objects of this type are constructed and owned by objects of class
1554  Mesh. See Mesh::GetGeometricFactors(). */
1556 {
1557 
1558 private:
1559  void Compute(const GridFunction &nodes,
1561 
1562 public:
1563  const Mesh *mesh;
1566 
1568  {
1569  COORDINATES = 1 << 0,
1570  JACOBIANS = 1 << 1,
1571  DETERMINANTS = 1 << 2,
1572  };
1573 
1574  GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags,
1576 
1577  GeometricFactors(const GridFunction &nodes, const IntegrationRule &ir,
1578  int flags,
1580 
1581  /// Mapped (physical) coordinates of all quadrature points.
1582  /** This array uses a column-major layout with dimensions (NQ x SDIM x NE)
1583  where
1584  - NQ = number of quadrature points per element,
1585  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1586  - NE = number of elements in the mesh. */
1588 
1589  /// Jacobians of the element transformations at all quadrature points.
1590  /** This array uses a column-major layout with dimensions (NQ x SDIM x DIM x
1591  NE) where
1592  - NQ = number of quadrature points per element,
1593  - SDIM = space dimension of the mesh = mesh.SpaceDimension(),
1594  - DIM = dimension of the mesh = mesh.Dimension(), and
1595  - NE = number of elements in the mesh. */
1597 
1598  /// Determinants of the Jacobians at all quadrature points.
1599  /** This array uses a column-major layout with dimensions (NQ x NE) where
1600  - NQ = number of quadrature points per element, and
1601  - NE = number of elements in the mesh. */
1603 };
1604 
1605 /** @brief Structure for storing face geometric factors: coordinates, Jacobians,
1606  determinants of the Jacobians, and normal vectors. */
1607 /** Typically objects of this type are constructed and owned by objects of class
1608  Mesh. See Mesh::GetFaceGeometricFactors(). */
1610 {
1611 public:
1612  const Mesh *mesh;
1616 
1618  {
1619  COORDINATES = 1 << 0,
1620  JACOBIANS = 1 << 1,
1621  DETERMINANTS = 1 << 2,
1622  NORMALS = 1 << 3,
1623  };
1624 
1625  FaceGeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags,
1626  FaceType type);
1627 
1628  /// Mapped (physical) coordinates of all quadrature points.
1629  /** This array uses a column-major layout with dimensions (NQ x SDIM x NF)
1630  where
1631  - NQ = number of quadrature points per face,
1632  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1633  - NF = number of faces in the mesh. */
1635 
1636  /// Jacobians of the element transformations at all quadrature points.
1637  /** This array uses a column-major layout with dimensions (NQ x SDIM x DIM x
1638  NF) where
1639  - NQ = number of quadrature points per face,
1640  - SDIM = space dimension of the mesh = mesh.SpaceDimension(),
1641  - DIM = dimension of the mesh = mesh.Dimension(), and
1642  - NF = number of faces in the mesh. */
1644 
1645  /// Determinants of the Jacobians at all quadrature points.
1646  /** This array uses a column-major layout with dimensions (NQ x NF) where
1647  - NQ = number of quadrature points per face, and
1648  - NF = number of faces in the mesh. */
1650 
1651  /// Normals at all quadrature points.
1652  /** This array uses a column-major layout with dimensions (NQ x DIM x NF) where
1653  - NQ = number of quadrature points per face,
1654  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1655  - NF = number of faces in the mesh. */
1657 };
1658 
1659 /// Class used to extrude the nodes of a mesh
1661 {
1662 private:
1663  int n, layer;
1664  double p[2], s;
1665  Vector tip;
1666 public:
1667  NodeExtrudeCoefficient(const int dim, const int n_, const double s_);
1668  void SetLayer(const int l) { layer = l; }
1670  virtual void Eval(Vector &V, ElementTransformation &T,
1671  const IntegrationPoint &ip);
1673 };
1674 
1675 
1676 /// Extrude a 1D mesh
1677 Mesh *Extrude1D(Mesh *mesh, const int ny, const double sy,
1678  const bool closed = false);
1679 
1680 /// Extrude a 2D mesh
1681 Mesh *Extrude2D(Mesh *mesh, const int nz, const double sz);
1682 
1683 // shift cyclically 3 integers left-to-right
1684 inline void ShiftRight(int &a, int &b, int &c)
1685 {
1686  int t = a;
1687  a = c; c = b; b = t;
1688 }
1689 
1690 }
1691 
1692 #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:4449
Abstract class for all finite elements.
Definition: fe.hpp:243
const IntegrationRule * IntRule
Definition: mesh.hpp:1613
MFEM_DEPRECATED Mesh(int n, double sx=1.0)
Deprecated: see MakeCartesian1D.
Definition: mesh.hpp:789
void Loader(std::istream &input, int generate_edges=0, std::string parse_tag="")
Definition: mesh.cpp:3530
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:1647
void PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
Definition: mesh.cpp:10801
virtual void Save(const char *fname, int precision=16) const
Definition: mesh.cpp:9615
void GetPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:5758
static const int vtk_quadratic_hex[27]
Definition: mesh.hpp:181
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1387
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:6438
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const DenseMatrix * PointMatrix
Definition: mesh.hpp:144
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1203
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:917
void ScaleElements(double sf)
Definition: mesh.cpp:10938
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:5573
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:5474
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo)
void FreeElement(Element *E)
Definition: mesh.cpp:11252
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1324
int CheckElementOrientation(bool fix_it=true)
Check (and optionally attempt to fix) the orientation of the elements.
Definition: mesh.cpp:4953
static bool remove_unused_vertices
Definition: mesh.hpp:215
void SetVertices(const Vector &vert_coord)
Definition: mesh.cpp:7283
static const int vtk_quadratic_tet[10]
Definition: mesh.hpp:179
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:5536
void Make2D(int nx, int ny, Element::Type type, double sx, double sy, bool generate_edges, bool sfc_ordering)
Definition: mesh.cpp:2913
bool FaceIsTrueInterior(int FaceNo) const
Definition: mesh.hpp:425
Vector J
Jacobians of the element transformations at all quadrature points.
Definition: mesh.hpp:1643
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition: mesh.hpp:1587
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:11527
static Mesh MakeSimplicial(const Mesh &orig_mesh)
Definition: mesh.cpp:4120
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:7588
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:784
long GetGlobalNE() const
Return the total (global) number of elements.
Definition: mesh.hpp:872
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:1389
static const int NumGeom
Definition: geom.hpp:41
void AddHexAsWedges(const int *vi, int attr=1)
Definition: mesh.cpp:1408
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:194
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:6477
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:5691
CoarseFineTransformations CoarseFineTr
Definition: mesh.hpp:171
int own_nodes
Definition: mesh.hpp:177
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5428
bool Conforming() const
Definition: mesh.hpp:1366
void MoveVertices(const Vector &displacements)
Definition: mesh.cpp:7263
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:849
IsoparametricTransformation Transformation
Definition: mesh.hpp:165
void GetLocalPtToSegTransformation(IsoparametricTransformation &, int)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:633
void GetCharacteristics(double &h_min, double &h_max, double &kappa_min, double &kappa_max, Vector *Vh=NULL, Vector *Vk=NULL)
Definition: mesh.cpp:193
void PrintWithPartitioning(int *partitioning, std::ostream &out, 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:10310
const Geometry::Type geom
Definition: ex1.cpp:40
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:3297
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:7367
int NumOfEdges
Definition: mesh.hpp:70
virtual void UniformRefinement2D()
Refine a mixed 2D mesh uniformly.
Definition: mesh.hpp:326
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:7386
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:5748
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:1684
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:5753
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:3307
int EulerNumber2D() const
Equals 1 - num_holes.
Definition: mesh.hpp:908
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
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:1310
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition: mesh.hpp:1649
int AddBdrPoint(int v, int attr=1)
Definition: mesh.cpp:1497
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:5439
void DeleteTables()
Definition: mesh.hpp:224
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:1029
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:804
GridFunction * Nodes
Definition: mesh.hpp:176
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:11545
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:11008
void FinalizeCheck()
Definition: mesh.cpp:1542
void ApplyLocalSlaveTransformation(FaceElementTransformations &FT, const FaceInfo &fi, bool is_ghost)
Definition: mesh.cpp:981
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:5714
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:7396
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
const Element * GetFace(int i) const
Definition: mesh.hpp:950
void FinalizeWedgeMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a wedge Mesh.
Definition: mesh.cpp:2435
Element * GetElement(int i)
Definition: mesh.hpp:944
void GenerateBoundaryElements()
Definition: mesh.cpp:1504
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:968
void PrintVTU(std::ostream &out, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr_elements=false)
Definition: mesh.cpp:9884
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1555
NodeExtrudeCoefficient(const int dim, const int n_, const double s_)
Definition: mesh.cpp:11521
Data type for vertex.
Definition: vertex.hpp:22
void SetMeshGen()
Determine the mesh generator bitmask meshgen, see MeshGenerator().
Definition: mesh.cpp:3486
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:5136
void GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:7272
void ReadNetgen3DMesh(std::istream &input)
The inverse transformation of a given ElementTransformation.
Definition: eltrans.hpp:172
void ReadInlineMesh(std::istream &input, bool generate_edges=false)
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
Definition: mesh.hpp:1153
void RemoveInternalBoundaries()
Definition: mesh.cpp:11164
void DebugDump(std::ostream &out) const
Output an NCMesh-compatible debug dump.
Definition: mesh.cpp:11932
Mesh * Extrude2D(Mesh *mesh, const int nz, const double sz)
Extrude a 2D mesh.
Definition: mesh.cpp:11705
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:886
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:971
void SetEmpty()
Definition: mesh.cpp:1122
Array< Element * > faces
Definition: mesh.hpp:92
virtual long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
Definition: mesh.hpp:869
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
Definition: mesh.cpp:2470
Geometry::Type GetBdrElementGeometry(int i) const
Definition: mesh.hpp:962
void GetVertexToVertexTable(DSTable &) const
Definition: mesh.cpp:5829
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:288
Operation last_operation
Definition: mesh.hpp:218
void MarkTriMeshForRefinement()
Definition: mesh.cpp:2035
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:979
void InitRefinementTransforms()
Definition: mesh.cpp:9233
IsoparametricTransformation FaceTransformation
Definition: mesh.hpp:167
void UniformRefinement2D_base(bool update_nodes=true)
Definition: mesh.cpp:7429
Array< FaceGeometricFactors * > face_geom_factors
Optional face geometric factors.
Definition: mesh.hpp:210
Array< NCFaceInfo > nc_faces_info
Definition: mesh.hpp:154
Element * ReadElement(std::istream &)
Definition: mesh.cpp:3468
void KnotInsert(Array< KnotVector * > &kv)
Definition: mesh.cpp:4608
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:4483
Element * NewElement(int geom)
Definition: mesh.cpp:3415
FaceElementTransformations FaceElemTr
Definition: mesh.hpp:168
Table * el_to_face
Definition: mesh.hpp:157
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
Definition: mesh.cpp:567
const GridFunction * GetNodes() const
Definition: mesh.hpp:1259
void MarkForRefinement()
Definition: mesh.cpp:2018
Element::Type GetFaceElementType(int Face) const
Definition: mesh.cpp:1093
void SetVerticesFromNodes(const GridFunction *nodes)
Helper to set vertex coordinates given a high-order curvature function.
Definition: mesh.cpp:4901
Structure for storing face geometric factors: coordinates, Jacobians, determinants of the Jacobians...
Definition: mesh.hpp:1609
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:428
int AddBdrTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1454
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
Definition: mesh.cpp:9573
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:4915
void MakeRefined_(Mesh &orig_mesh, const Array< int > ref_factors, int ref_type)
Internal function used in Mesh::MakeRefined.
Definition: mesh.cpp:3897
Geometry::Constants< Geometry::SEGMENT > seg_t
Definition: mesh.hpp:192
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:5184
long GetSequence() const
Definition: mesh.hpp:1380
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:770
int AddTri(const int *vi, int attr=1)
Definition: mesh.hpp:649
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1263
bool Nonconforming() const
Definition: mesh.hpp:1367
Vector J
Jacobians of the element transformations at all quadrature points.
Definition: mesh.hpp:1596
FaceType
Definition: mesh.hpp:45
void MoveNodes(const Vector &displacements)
Definition: mesh.cpp:7331
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
Definition: mesh.cpp:2050
void GetElementData(int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.hpp:924
NCFaceInfo(bool slave, int master, const DenseMatrix *pm)
Definition: mesh.hpp:149
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a triangular Mesh.
Definition: mesh.cpp:1556
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:4682
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 FinalizeTetMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a tetrahedral Mesh.
Definition: mesh.cpp:2394
void AddBdrQuadAsTriangles(const int *vi, int attr=1)
Definition: mesh.cpp:1482
void Make3D(int nx, int ny, int nz, Element::Type type, double sx, double sy, double sz, bool sfc_ordering)
Definition: mesh.cpp:2675
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
Definition: mesh.cpp:5979
virtual void SetAttributes()
Definition: mesh.cpp:1209
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:8800
GeometryList(Mesh &mesh, int dim)
Construct a GeometryList of all geometries of dimension dim in mesh.
Definition: mesh.hpp:1005
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition: mesh.hpp:1602
int FindCoarseElement(int i)
Definition: mesh.cpp:9245
GeometryList(Mesh &mesh)
Construct a GeometryList of all element geometries in mesh.
Definition: mesh.hpp:1000
double b
Definition: lissajous.cpp:42
int AddTet(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1338
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Returns the transformation defining the given face element in a user-defined variable.
Definition: mesh.cpp:492
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:648
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
Definition: mesh.cpp:6044
void CheckPartitioning(int *partitioning_)
Definition: mesh.cpp:6868
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
Definition: mesh.cpp:4747
void ReadGmshMesh(std::istream &input, int &curved, int &read_gf)
STable3D * GetElementToFaceTable(int ret_ftbl=0)
Definition: mesh.cpp:6263
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4882
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:839
virtual void PrintInfo(std::ostream &out=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
Definition: mesh.hpp:1497
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:5776
void InitTables()
Definition: mesh.cpp:1116
Geometry::Type GetFaceGeometry(int i) const
Definition: mesh.hpp:952
void CheckDisplacements(const Vector &displacements, double &tmax)
Definition: mesh.cpp:7186
VTKFormat
Definition: vtk.hpp:61
int mesh_geoms
Definition: mesh.hpp:78
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1440
Geometry::Constants< Geometry::TETRAHEDRON > tet_t
Definition: mesh.hpp:195
int AddElement(Element *elem)
The parameter elem should be allocated using the NewElement() method.
Definition: mesh.cpp:1426
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:820
void SetNodesOwner(bool nodes_owner)
Set the mesh nodes ownership flag.
Definition: mesh.hpp:1263
void SetLayer(const int l)
Definition: mesh.hpp:1668
int GetElementToEdgeTable(Table &, Array< int > &)
Definition: mesh.cpp:5854
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:5726
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:1814
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:8824
const Element * GetElement(int i) const
Definition: mesh.hpp:942
void GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:738
void ResetLazyData()
Definition: mesh.cpp:1199
void GetElementCenter(int i, Vector &center)
Definition: mesh.cpp:68
IsoparametricTransformation BdrTransformation
Definition: mesh.hpp:166
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:839
IsoparametricTransformation Transformation2
Definition: mesh.hpp:165
void Init()
Definition: mesh.cpp:1098
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:4843
Element * GetBdrElement(int i)
Definition: mesh.hpp:948
Vector normal
Normals at all quadrature points.
Definition: mesh.hpp:1656
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
Definition: ncmesh.hpp:109
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:1020
Geometry::Constants< Geometry::PRISM > pri_t
Definition: mesh.hpp:197
int Dimension() const
Definition: mesh.hpp:911
virtual void ReorientTetMesh()
Definition: mesh.cpp:6376
prob_type prob
Definition: ex25.cpp:153
int NumOfBdrElements
Definition: mesh.hpp:69
Table * el_to_edge
Definition: mesh.hpp:156
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:8180
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
virtual void MarkTetMeshForRefinement(DSTable &v_to_v)
Definition: mesh.cpp:2075
Array< Triple< int, int, int > > tmp_vertex_parents
Definition: mesh.hpp:189
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:780
List of mesh geometries stored as Array&lt;Geometry::Type&gt;.
Definition: mesh.hpp:994
void Destroy()
Definition: mesh.cpp:1170
int CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Definition: mesh.cpp:5317
int SpaceDimension() const
Definition: mesh.hpp:912
GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, MemoryType d_mt=MemoryType::DEFAULT)
Definition: mesh.cpp:11392
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
Definition: mesh.cpp:10429
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
Definition: mesh.hpp:1165
void GetLocalQuadToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:760
void AddPointFaceElement(int lf, int gf, int el)
Used in GenerateFaces()
Definition: mesh.cpp:5947
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:1055
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:3091
const Element *const * GetElementsArray() const
Definition: mesh.hpp:939
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:8422
void PrintVTK(std::ostream &out)
Definition: mesh.cpp:9629
virtual ~Mesh()
Destroys Mesh.
Definition: mesh.hpp:1538
void FinalizeMesh(int refine=0, bool fix_orientation=true)
Finalize the construction of any type of Mesh.
Definition: mesh.cpp:2500
Table * el_to_el
Definition: mesh.hpp:158
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:1061
int AddSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1296
Geometry::Constants< Geometry::CUBE > hex_t
Definition: mesh.hpp:196
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Geometry::Type geom_buf[Geometry::NumGeom]
Definition: mesh.hpp:997
static Mesh LoadFromFile(const char *filename, int generate_edges=0, int refine=1, bool fix_orientation=true)
Definition: mesh.cpp:3269
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Definition: mesh.cpp:11057
static void PrintElement(const Element *, std::ostream &)
Definition: mesh.cpp:3480
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
Definition: mesh.cpp:4652
int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
Definition: mesh.cpp:1373
const IntegrationRule * IntRule
Definition: mesh.hpp:1564
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition: mesh.hpp:1634
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2507
static void GetElementArrayEdgeTable(const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
Definition: mesh.cpp:5807
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
Definition: mesh.hpp:1258
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:674
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
Definition: mesh.cpp:2098
int AddBdrQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1468
bool IsSlaveFace(const FaceInfo &fi) const
Definition: mesh.cpp:976
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:3287
void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Definition: mesh.cpp:1240
void SetBdrAttribute(int i, int attr)
Set the attribute of boundary element i.
Definition: mesh.hpp:1206
void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo, const std::string &xml_prefix="")
virtual void PrintXG(std::ostream &out=mfem::out) const
Print the mesh to the given stream using Netgen/Truegrid format.
Definition: mesh.cpp:9318
Geometry::Type GetElementGeometry(int i) const
Definition: mesh.hpp:957
void RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
Definition: mesh.cpp:8843
static Mesh MakeCartesian1D(int n, double sx=1.0)
Definition: mesh.cpp:3279
Array< Element * > elements
Definition: mesh.hpp:85
const Table & ElementToElementTable()
Definition: mesh.cpp:5893
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
Definition: mesh.cpp:5792
int meshgen
Definition: mesh.hpp:77
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:9255
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:8869
Table * bel_to_edge
Definition: mesh.hpp:160
Operation GetLastOperation() const
Return type of last modification of the mesh.
Definition: mesh.hpp:1374
Table * GetFaceToElementTable() const
Definition: mesh.cpp:5634
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:688
double a
Definition: lissajous.cpp:41
int EulerNumber() const
Equals 1 + num_holes - num_loops.
Definition: mesh.hpp:905
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:206
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:843
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:9821
void GetNode(int i, double *coord) const
Definition: mesh.cpp:7292
Array< int > be_to_edge
Definition: mesh.hpp:159
const Table & ElementToFaceTable() const
Definition: mesh.cpp:5929
Class for integration point with weight.
Definition: intrules.hpp:25
void Swap(Mesh &other, bool non_geometry)
Definition: mesh.cpp:8612
Element * ReadElementWithoutAttr(std::istream &)
Definition: mesh.cpp:3438
FaceGeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, FaceType type)
Definition: mesh.cpp:11472
A standard isoparametric element transformation.
Definition: eltrans.hpp:348
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:5452
static const int vtk_quadratic_wedge[18]
Definition: mesh.hpp:180
bool OwnsNodes() const
Return the mesh nodes ownership flag.
Definition: mesh.hpp:1261
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:347
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:668
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:4877
Table * face_edge
Definition: mesh.hpp:162
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1585
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
Definition: mesh.hpp:927
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2600
bool DerefineByError(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
Definition: mesh.cpp:8539
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:1866
void GetFaceEdges(int i, Array< int > &edges, Array< int > &o) const
Definition: mesh.cpp:5506
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:4665
int dim
Definition: ex24.cpp:53
Table * edge_vertex
Definition: mesh.hpp:163
long sequence
Definition: mesh.hpp:83
virtual ~NodeExtrudeCoefficient()
Definition: mesh.hpp:1672
const Mesh * mesh
Definition: mesh.hpp:1563
void MakeSimplicial_(const Mesh &orig_mesh, int *vglobal)
Definition: mesh.cpp:4127
void DestroyPointers()
Definition: mesh.cpp:1144
Array< FaceInfo > faces_info
Definition: mesh.hpp:153
STable3D * GetFacesTable()
Definition: mesh.cpp:6214
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
Definition: sparsemat.hpp:655
NCMesh * ncmesh
Optional non-conforming mesh extension.
Definition: mesh.hpp:207
Geometry::Type GetBdrElementBaseGeometry(int i) const
Definition: mesh.hpp:974
virtual void UniformRefinement3D()
Refine a mixed 3D mesh uniformly.
Definition: mesh.hpp:336
int AddBdrElement(Element *elem)
Definition: mesh.cpp:1433
RefCoord t[3]
static void PrintElementWithoutAttr(const Element *, std::ostream &)
Definition: mesh.cpp:3456
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:625
int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr=1)
Definition: mesh.cpp:1359
virtual bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NC version of GeneralDerefinement.
Definition: mesh.cpp:8491
int Dim
Definition: mesh.hpp:66
Geometry::Type GetFaceGeometryType(int Face) const
Definition: mesh.cpp:1074
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:828
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:852
Geometry::Constants< Geometry::TRIANGLE > tri_t
Definition: mesh.hpp:193
double AggregateError(const Array< double > &elem_error, const int *fine, int nfine, int op)
Derefinement helper.
Definition: mesh.cpp:8472
double * GetVertex(int i)
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:922
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:855
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:11277
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:1200
Mesh & operator=(Mesh &&mesh)
Move assignment operstor.
Definition: mesh.cpp:3263
Vector data type.
Definition: vector.hpp:60
void ReadTrueGridMesh(std::istream &input)
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: mesh.cpp:4941
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:5088
void SetNode(int i, const double *coord)
Definition: mesh.cpp:7311
int GetBdrFace(int BdrElemNo) const
Return the local face index for the given boundary face.
Definition: mesh.cpp:1037
void DestroyTables()
Definition: mesh.cpp:1128
void GetLocalTriToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:712
void GenerateFaces()
Definition: mesh.cpp:6071
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:10235
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:3349
void DeleteGeometricFactors()
Destroy all GeometricFactors stored by the Mesh.
Definition: mesh.cpp:825
int spaceDim
Definition: mesh.hpp:67
IsoparametricTransformation EdgeTransformation
Definition: mesh.hpp:167
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:60
MemAlloc< Tetrahedron, 1024 > TetMemory
Definition: mesh.hpp:185
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
Definition: mesh.cpp:5545
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
Definition: mesh.cpp:2164
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:7355
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:4871
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
int NumOfVertices
Definition: mesh.hpp:69
void ScaleSubdomains(double sf)
Definition: mesh.cpp:10868
Array< int > be_to_face
Definition: mesh.hpp:161
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:8568
void Bisection(int i, const DSTable &, int *, int *, int *)
Bisect a triangle: element with index i is bisected.
Definition: mesh.cpp:8898
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:8732
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:1197
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
Definition: mesh.cpp:6016
const Table & ElementToEdgeTable() const
Definition: mesh.cpp:5938
Array< GeometricFactors * > geom_factors
Optional geometric factors.
Definition: mesh.hpp:208
void BdrBisection(int i, const HashTable< Hashed2 > &)
Bisect a boundary triangle: boundary element with index i is bisected.
Definition: mesh.cpp:9106
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
Definition: mesh.hpp:1015
void GenerateNCFaceInfo()
Definition: mesh.cpp:6156
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:1011
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:202
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
Definition: mesh.cpp:5599
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:946
void EnsureNodes()
Definition: mesh.cpp:4849
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.cpp:8668
void UpdateNodes()
Update the nodes of a curved mesh after refinement.
Definition: mesh.cpp:7417
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=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
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:5668
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:294
void Printer(std::ostream &out=mfem::out, std::string section_delimiter="") const
Definition: mesh.cpp:9482
Class used to extrude the nodes of a mesh.
Definition: mesh.hpp:1660
void AddVertexParents(int i, int p1, int p2)
Mark vertex i as non-conforming, with parent vertices p1 and p2.
Definition: mesh.cpp:1280
int NumOfFaces
Definition: mesh.hpp:70