MFEM  v4.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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
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 "ncmesh.hpp"
22 #include "../fem/eltrans.hpp"
23 #include "../fem/coefficient.hpp"
24 #include "../general/gzstream.hpp"
25 #include <iostream>
26 
27 namespace mfem
28 {
29 
30 // Data type mesh
31 
32 class GeometricFactors;
33 class KnotVector;
34 class NURBSExtension;
35 class FiniteElementSpace;
36 class GridFunction;
37 struct Refinement;
38 
39 #ifdef MFEM_USE_MPI
40 class ParMesh;
41 class ParNCMesh;
42 #endif
43 
44 
45 class Mesh
46 {
47 #ifdef MFEM_USE_MPI
48  friend class ParMesh;
49  friend class ParNCMesh;
50 #endif
51  friend class NURBSExtension;
52 
53 protected:
54  int Dim;
55  int spaceDim;
56 
59 
60  int meshgen; // see MeshGenerator()
61  int mesh_geoms; // sum of (1 << geom) for all geom of all dimensions
62 
63  // Counter for Mesh transformations: refinement, derefinement, rebalancing.
64  // Used for checking during Update operations on objects depending on the
65  // Mesh, such as FiniteElementSpace, GridFunction, etc.
66  long sequence;
67 
69  // Vertices are only at the corners of elements, where you would expect them
70  // in the lowest-order mesh. In some cases, e.g. in a Mesh that defines the
71  // patch topology for a NURBS mesh (see LoadPatchTopo()) the vertices may be
72  // empty while NumOfVertices is positive.
76 
77  struct FaceInfo
78  {
79  // Inf = 64 * LocalFaceIndex + FaceOrientation
81  int NCFace; /* -1 if this is a regular conforming/boundary face;
82  index into 'nc_faces_info' if >= 0. */
83  };
84  // NOTE: in NC meshes, master faces have Elem2No == -1. Slave faces on the
85  // other hand have Elem2No and Elem2Inf set to the master face's element and
86  // its local face number.
87  //
88  // A local face is one generated from a local element and has index i in
89  // faces_info such that i < GetNumFaces(). Also, Elem1No always refers to the
90  // element (slave or master, in the non-conforming case) that generated the
91  // face.
92  // Classification of a local (non-ghost) face based on its FaceInfo:
93  // - Elem2No >= 0 --> local internal face; can be either:
94  // - NCFace == -1 --> conforming face, or
95  // - NCFace >= 0 --> non-conforming slave face.
96  // - Elem2No < 0 --> local "boundary" face; can be one of:
97  // - NCFace == -1 --> conforming face; can be either:
98  // - Elem2Inf < 0 --> true boundary face (no element on side 2)
99  // - Elem2Inf >= 0 --> shared face where element 2 is a face-neighbor
100  // element with index -1-Elem2No. This state is initialized by
101  // ParMesh::ExchangeFaceNbrData().
102  // - NCFace >= 0 --> non-conforming master face. Elem2No is -1 or, in the
103  // case of a shared face, -1-Elem2No is the index of one of the adjacent
104  // (the last one?) slave ghost elements. Elem2Inf is -1.
105  //
106  // A ghost face is a non-conforming face that is generated by a non-local,
107  // i.e. ghost, element. A ghost face has index i in faces_info such that
108  // i >= GetNumFaces().
109  // Classification of a ghost (non-local) face based on its FaceInfo:
110  // - Elem1No == -1 --> master ghost face? These ghost faces also have:
111  // Elem2No == -1, Elem1Inf == Elem2Inf == -1, and NCFace == -1.
112  // - Elem1No >= 0 --> slave ghost face; Elem1No is the index of the local
113  // master side element, i.e. side 1 IS NOT the side that generated the
114  // face. Elem2No is < 0 and -1-Elem2No is the index of the ghost
115  // face-neighbor element that generated this slave ghost face. In this
116  // case, Elem2Inf >= 0.
117 
118  struct NCFaceInfo
119  {
120  bool Slave; // true if this is a slave face, false if master face
121  int MasterFace; // if Slave, this is the index of the master face
122  const DenseMatrix* PointMatrix; // if Slave, position within master face
123  // (NOTE: PointMatrix points to a matrix owned by NCMesh.)
124 
125  NCFaceInfo() = default;
126 
127  NCFaceInfo(bool slave, int master, const DenseMatrix* pm)
128  : Slave(slave), MasterFace(master), PointMatrix(pm) {}
129  };
130 
133 
138  Table *bel_to_edge; // for 3D
140  mutable Table *face_edge;
141  mutable Table *edge_vertex;
142 
147 
148  // refinement embeddings for forward compatibility with NCMesh
150 
151  // Nodes are only active for higher order meshes, and share locations with
152  // the vertices, plus all the higher- order control points within the
153  // element and along the edges and on the faces.
156 
157  static const int vtk_quadratic_tet[10];
158  static const int vtk_quadratic_wedge[18];
159  static const int vtk_quadratic_hex[27];
160 
161 #ifdef MFEM_USE_MEMALLOC
162  friend class Tetrahedron;
164 #endif
165 
166 public:
173 
175 
176  /// A list of all unique element attributes used by the Mesh.
178  /// A list of all unique boundary attributes used by the Mesh.
180 
181  NURBSExtension *NURBSext; ///< Optional NURBS mesh extension.
182  NCMesh *ncmesh; ///< Optional non-conforming mesh extension.
183  Array<GeometricFactors*> geom_factors; ///< Optional geometric factors.
184 
185  // Global parameter that can be used to control the removal of unused
186  // vertices performed when reading a mesh in MFEM format. The default value
187  // (true) is set in mesh_readers.cpp.
189 
190 protected:
192 
193  void Init();
194  void InitTables();
195  void SetEmpty(); // Init all data members with empty values
196  void DestroyTables();
198  void DestroyPointers(); // Delete data specifically allocated by class Mesh.
199  void Destroy(); // Delete all owned data.
200  void DeleteLazyTables();
201 
202  Element *ReadElementWithoutAttr(std::istream &);
203  static void PrintElementWithoutAttr(const Element *, std::ostream &);
204 
205  Element *ReadElement(std::istream &);
206  static void PrintElement(const Element *, std::ostream &);
207 
208  // Readers for different mesh formats, used in the Load() method.
209  // The implementations of these methods are in mesh_readers.cpp.
210  void ReadMFEMMesh(std::istream &input, bool mfem_v11, int &curved);
211  void ReadLineMesh(std::istream &input);
212  void ReadNetgen2DMesh(std::istream &input, int &curved);
213  void ReadNetgen3DMesh(std::istream &input);
214  void ReadTrueGridMesh(std::istream &input);
215  void ReadVTKMesh(std::istream &input, int &curved, int &read_gf,
216  bool &finalize_topo);
217  void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf);
218  void ReadInlineMesh(std::istream &input, bool generate_edges = false);
219  void ReadGmshMesh(std::istream &input);
220  /* Note NetCDF (optional library) is used for reading cubit files */
221 #ifdef MFEM_USE_NETCDF
222  void ReadCubit(const char *filename, int &curved, int &read_gf);
223 #endif
224 
225  /// Determine the mesh generator bitmask #meshgen, see MeshGenerator().
226  /** Also, initializes #mesh_geoms. */
227  void SetMeshGen();
228 
229  /// Return the length of the segment from node i to node j.
230  double GetLength(int i, int j) const;
231 
232  /** Compute the Jacobian of the transformation from the perfect
233  reference element at the center of the element. */
234  void GetElementJacobian(int i, DenseMatrix &J);
235 
236  void GetElementCenter(int i, Vector &c);
237 
238  void MarkForRefinement();
240  void GetEdgeOrdering(DSTable &v_to_v, Array<int> &order);
241  virtual void MarkTetMeshForRefinement(DSTable &v_to_v);
242 
243  // Methods used to prepare and apply permutation of the mesh nodes assuming
244  // that the mesh elements may be rotated (e.g. to mark triangle or tet edges
245  // for refinement) between the two calls - PrepareNodeReorder() and
246  // DoNodeReorder(). The latter method assumes that the 'faces' have not been
247  // updated after the element rotations.
248  void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert);
249  void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert);
250 
252  STable3D *GetElementToFaceTable(int ret_ftbl = 0);
253 
254  /** Red refinement. Element with index i is refined. The default
255  red refinement for now is Uniform. */
256  void RedRefinement(int i, const DSTable &v_to_v,
257  int *edge1, int *edge2, int *middle)
258  { UniformRefinement(i, v_to_v, edge1, edge2, middle); }
259 
260  /** Green refinement. Element with index i is refined. The default
261  refinement for now is Bisection. */
262  void GreenRefinement(int i, const DSTable &v_to_v,
263  int *edge1, int *edge2, int *middle)
264  { Bisection(i, v_to_v, edge1, edge2, middle); }
265 
266  /// Bisect a triangle: element with index @a i is bisected.
267  void Bisection(int i, const DSTable &, int *, int *, int *);
268 
269  /// Bisect a tetrahedron: element with index @a i is bisected.
270  void Bisection(int i, HashTable<Hashed2> &);
271 
272  /// Bisect a boundary triangle: boundary element with index @a i is bisected.
273  void BdrBisection(int i, const HashTable<Hashed2> &);
274 
275  /** Uniform Refinement. Element with index i is refined uniformly. */
276  void UniformRefinement(int i, const DSTable &, int *, int *, int *);
277 
278  /** @brief Averages the vertices with given @a indexes and saves the result
279  in #vertices[result]. */
280  void AverageVertices(const int *indexes, int n, int result);
281 
283  int FindCoarseElement(int i);
284 
285  /// Update the nodes of a curved mesh after refinement
286  void UpdateNodes();
287 
288  /// Refine a mixed 2D mesh uniformly.
289  virtual void UniformRefinement2D();
290 
291  /* If @a f2qf is not NULL, adds all quadrilateral faces to @a f2qf which
292  represents a "face-to-quad-face" index map. When all faces are quads, the
293  array @a f2qf is kept empty since it is not needed. */
294  void UniformRefinement3D_base(Array<int> *f2qf = NULL,
295  DSTable *v_to_v_p = NULL);
296 
297  /// Refine a mixed 3D mesh uniformly.
299 
300  /// Refine NURBS mesh.
301  virtual void NURBSUniformRefinement();
302 
303  /// This function is not public anymore. Use GeneralRefinement instead.
304  virtual void LocalRefinement(const Array<int> &marked_el, int type = 3);
305 
306  /// This function is not public anymore. Use GeneralRefinement instead.
307  virtual void NonconformingRefinement(const Array<Refinement> &refinements,
308  int nc_limit = 0);
309 
310  /// NC version of GeneralDerefinement.
311  virtual bool NonconformingDerefinement(Array<double> &elem_error,
312  double threshold, int nc_limit = 0,
313  int op = 1);
314  /// Derefinement helper.
315  double AggregateError(const Array<double> &elem_error,
316  const int *fine, int nfine, int op);
317 
318  /// Read NURBS patch/macro-element mesh
319  void LoadPatchTopo(std::istream &input, Array<int> &edge_to_knot);
320 
321  void UpdateNURBS();
322 
323  void PrintTopo(std::ostream &out, const Array<int> &e_to_k) const;
324 
325  /// Used in GetFaceElementTransformations (...)
328  int i);
330  int i);
331  /// Used in GetFaceElementTransformations (...)
333  int i);
334  /// Used in GetFaceElementTransformations (...)
336  int i);
337  /// Used in GetFaceElementTransformations (...)
339  int i);
340  /// Used in GetFaceElementTransformations (...)
342  int i);
343 
344  /** Used in GetFaceElementTransformations to account for the fact that a
345  slave face occupies only a portion of its master face. */
347  const FaceInfo &fi);
348  bool IsSlaveFace(const FaceInfo &fi) const;
349 
350  /// Returns the orientation of "test" relative to "base"
351  static int GetTriOrientation (const int * base, const int * test);
352  /// Returns the orientation of "test" relative to "base"
353  static int GetQuadOrientation (const int * base, const int * test);
354 
355  static void GetElementArrayEdgeTable(const Array<Element*> &elem_array,
356  const DSTable &v_to_v,
357  Table &el_to_edge);
358 
359  /** Return vertex to vertex table. The connections stored in the table
360  are from smaller to bigger vertex index, i.e. if i<j and (i, j) is
361  in the table, then (j, i) is not stored. */
362  void GetVertexToVertexTable(DSTable &) const;
363 
364  /** Return element to edge table and the indices for the boundary edges.
365  The entries in the table are ordered according to the order of the
366  nodes in the elements. For example, if T is the element to edge table
367  T(i, 0) gives the index of edge in element i that connects vertex 0
368  to vertex 1, etc. Returns the number of the edges. */
370 
371  /// Used in GenerateFaces()
372  void AddPointFaceElement(int lf, int gf, int el);
373 
374  void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1);
375 
376  void AddTriangleFaceElement (int lf, int gf, int el,
377  int v0, int v1, int v2);
378 
379  void AddQuadFaceElement (int lf, int gf, int el,
380  int v0, int v1, int v2, int v3);
381  /** For a serial Mesh, return true if the face is interior. For a parallel
382  ParMesh return true if the face is interior or shared. In parallel, this
383  method only works if the face neighbor data is exchanged. */
384  bool FaceIsTrueInterior(int FaceNo) const
385  {
386  return FaceIsInterior(FaceNo) || (faces_info[FaceNo].Elem2Inf >= 0);
387  }
388 
389  // shift cyclically 3 integers left-to-right
390  inline static void ShiftL2R(int &, int &, int &);
391  // shift cyclically 3 integers so that the smallest is first
392  inline static void Rotate3(int &, int &, int &);
393 
394  void FreeElement(Element *E);
395 
396  void GenerateFaces();
397  void GenerateNCFaceInfo();
398 
399  /// Begin construction of a mesh
400  void InitMesh(int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem);
401 
402  // Used in the methods FinalizeXXXMesh() and FinalizeTopology()
403  void FinalizeCheck();
404 
405  void Loader(std::istream &input, int generate_edges = 0,
406  std::string parse_tag = "");
407 
408  // If NURBS mesh, write NURBS format. If NCMesh, write mfem v1.1 format.
409  // If section_delimiter is empty, write mfem v1.0 format. Otherwise, write
410  // mfem v1.2 format with the given section_delimiter at the end.
411  void Printer(std::ostream &out = mfem::out,
412  std::string section_delimiter = "") const;
413 
414  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
415  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
416  type=TETRAHEDRON. The parameter @a sfc_ordering controls how the elements
417  (when type=HEXAHEDRON) are ordered: true - use space-filling curve
418  ordering, or false - use lexicographic ordering. */
419  void Make3D(int nx, int ny, int nz, Element::Type type,
420  double sx, double sy, double sz, bool sfc_ordering);
421 
422  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
423  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
424  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
425  if 1 edges are generated. The parameter @a sfc_ordering controls how the
426  elements (when type=QUADRILATERAL) are ordered: true - use space-filling
427  curve ordering, or false - use lexicographic ordering. */
428  void Make2D(int nx, int ny, Element::Type type, double sx, double sy,
429  bool generate_edges, bool sfc_ordering);
430 
431  /// Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.
432  void Make1D(int n, double sx = 1.0);
433 
434  /// Initialize vertices/elements/boundary/tables from a nonconforming mesh.
435  void InitFromNCMesh(const NCMesh &ncmesh);
436 
437  /// Create from a nonconforming mesh.
438  explicit Mesh(const NCMesh &ncmesh);
439 
440  /// Swaps internal data with another mesh. By default, non-geometry members
441  /// like 'ncmesh' and 'NURBSExt' are only swapped when 'non_geometry' is set.
442  void Swap(Mesh& other, bool non_geometry);
443 
444  // used in GetElementData() and GetBdrElementData()
445  void GetElementData(const Array<Element*> &elem_array, int geom,
446  Array<int> &elem_vtx, Array<int> &attr) const;
447 
448 public:
449 
450  Mesh() { SetEmpty(); }
451 
452  /** Copy constructor. Performs a deep copy of (almost) all data, so that the
453  source mesh can be modified (e.g. deleted, refined) without affecting the
454  new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the
455  nodes, if present. */
456  explicit Mesh(const Mesh &mesh, bool copy_nodes = true);
457 
458  /// Construct a Mesh from the given primary data.
459  /** The array @a vertices is used as external data, i.e. the Mesh does not
460  copy the data and will not delete the pointer.
461 
462  The data from the other arrays is copied into the internal Mesh data
463  structures.
464 
465  This method calls the method FinalizeTopology(). The method Finalize()
466  may be called after this constructor and after optionally setting the
467  Mesh nodes. */
468  Mesh(double *vertices, int num_vertices,
469  int *element_indices, Geometry::Type element_type,
470  int *element_attributes, int num_elements,
471  int *boundary_indices, Geometry::Type boundary_type,
472  int *boundary_attributes, int num_boundary_elements,
473  int dimension, int space_dimension= -1);
474 
475  /** @anchor mfem_Mesh_init_ctor
476  @brief _Init_ constructor: begin the construction of a Mesh object. */
477  Mesh(int _Dim, int NVert, int NElem, int NBdrElem = 0, int _spaceDim = -1)
478  {
479  if (_spaceDim == -1)
480  {
481  _spaceDim = _Dim;
482  }
483  InitMesh(_Dim, _spaceDim, NVert, NElem, NBdrElem);
484  }
485 
486  /** @name Methods for Mesh construction.
487 
488  These methods are intended to be used with the @ref mfem_Mesh_init_ctor
489  "init constructor". */
490  ///@{
491 
492  Element *NewElement(int geom);
493 
494  void AddVertex(const double *);
495  void AddTri(const int *vi, int attr = 1);
496  void AddTriangle(const int *vi, int attr = 1);
497  void AddQuad(const int *vi, int attr = 1);
498  void AddTet(const int *vi, int attr = 1);
499  void AddWedge(const int *vi, int attr = 1);
500  void AddHex(const int *vi, int attr = 1);
501  void AddHexAsTets(const int *vi, int attr = 1);
502  void AddHexAsWedges(const int *vi, int attr = 1);
503  /// The parameter @a elem should be allocated using the NewElement() method
504  void AddElement(Element *elem) { elements[NumOfElements++] = elem; }
505  void AddBdrElement(Element *elem) { boundary[NumOfBdrElements++] = elem; }
506  void AddBdrSegment(const int *vi, int attr = 1);
507  void AddBdrTriangle(const int *vi, int attr = 1);
508  void AddBdrQuad(const int *vi, int attr = 1);
509  void AddBdrQuadAsTriangles(const int *vi, int attr = 1);
510 
512  /// Finalize the construction of a triangular Mesh.
513  void FinalizeTriMesh(int generate_edges = 0, int refine = 0,
514  bool fix_orientation = true);
515  /// Finalize the construction of a quadrilateral Mesh.
516  void FinalizeQuadMesh(int generate_edges = 0, int refine = 0,
517  bool fix_orientation = true);
518  /// Finalize the construction of a tetrahedral Mesh.
519  void FinalizeTetMesh(int generate_edges = 0, int refine = 0,
520  bool fix_orientation = true);
521  /// Finalize the construction of a wedge Mesh.
522  void FinalizeWedgeMesh(int generate_edges = 0, int refine = 0,
523  bool fix_orientation = true);
524  /// Finalize the construction of a hexahedral Mesh.
525  void FinalizeHexMesh(int generate_edges = 0, int refine = 0,
526  bool fix_orientation = true);
527  /// Finalize the construction of any type of Mesh.
528  /** This method calls FinalizeTopology() and Finalize(). */
529  void FinalizeMesh(int refine = 0, bool fix_orientation = true);
530 
531  ///@}
532 
533  /** @brief Finalize the construction of the secondary topology (connectivity)
534  data of a Mesh. */
535  /** This method does not require any actual coordinate data (either vertex
536  coordinates for linear meshes or node coordinates for meshes with nodes)
537  to be available. However, the data generated by this method is generally
538  required by the FiniteElementSpace class.
539 
540  After calling this method, setting the Mesh vertices or nodes, it may be
541  appropriate to call the method Finalize(). */
542  void FinalizeTopology();
543 
544  /// Finalize the construction of a general Mesh.
545  /** This method will:
546  - check and optionally fix the orientation of regular elements
547  - check and fix the orientation of boundary elements
548  - assume that #vertices are defined, if #Nodes == NULL
549  - assume that #Nodes are defined, if #Nodes != NULL.
550  @param[in] refine If true, prepare the Mesh for conforming refinement of
551  triangular or tetrahedral meshes.
552  @param[in] fix_orientation
553  If true, fix the orientation of inverted mesh elements
554  by permuting their vertices.
555 
556  Before calling this method, call FinalizeTopology() and ensure that the
557  Mesh vertices or nodes are set. */
558  virtual void Finalize(bool refine = false, bool fix_orientation = false);
559 
560  void SetAttributes();
561 
562 #ifdef MFEM_USE_GECKO
563  /** This is our integration with the Gecko library. This will call the
564  Gecko library to find an element ordering that will increase memory
565  coherency by putting elements that are in physical proximity closer in
566  memory. It can also be used to get a space-filling curve ordering for
567  ParNCMesh partitioning.
568  @param[out] ordering Output element ordering.
569  @param[in] iterations Number of V cycles (default 1).
570  @param[in] window Initial window size (default 2).
571  @param[in] period Iterations between window increment (default 1).
572  @param[in] seed Random number seed (default 0). */
573  void GetGeckoElementReordering(Array<int> &ordering,
574  int iterations = 1, int window = 2,
575  int period = 1, int seed = 0);
576 #endif
577 
578  /** Rebuilds the mesh with a different order of elements. The ordering
579  vector maps the old element number to the new element number. This also
580  reorders the vertices and nodes edges and faces along with the elements. */
581  void ReorderElements(const Array<int> &ordering, bool reorder_vertices = true);
582 
583  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
584  nx*ny*nz hexahedra if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if
585  type=TETRAHEDRON. If sfc_ordering = true (default), elements are ordered
586  along a space-filling curve, instead of row by row and layer by layer.
587  The parameter @a generate_edges is ignored (for now, it is kept for
588  backward compatibility). */
589  Mesh(int nx, int ny, int nz, Element::Type type, bool generate_edges = false,
590  double sx = 1.0, double sy = 1.0, double sz = 1.0,
591  bool sfc_ordering = true)
592  {
593  Make3D(nx, ny, nz, type, sx, sy, sz, sfc_ordering);
594  Finalize(true); // refine = true
595  }
596 
597  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
598  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
599  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
600  if 1 edges are generated. If scf_ordering = true (default), elements are
601  ordered along a space-filling curve, instead of row by row. */
602  Mesh(int nx, int ny, Element::Type type, bool generate_edges = false,
603  double sx = 1.0, double sy = 1.0, bool sfc_ordering = true)
604  {
605  Make2D(nx, ny, type, sx, sy, generate_edges, sfc_ordering);
606  Finalize(true); // refine = true
607  }
608 
609  /** Creates 1D mesh , divided into n equal intervals. */
610  explicit Mesh(int n, double sx = 1.0)
611  {
612  Make1D(n, sx);
613  // Finalize(); // reminder: not needed
614  }
615 
616  /** Creates mesh by reading a file in MFEM, Netgen, or VTK format. If
617  generate_edges = 0 (default) edges are not generated, if 1 edges are
618  generated. */
619  explicit Mesh(const char *filename, int generate_edges = 0, int refine = 1,
620  bool fix_orientation = true);
621 
622  /** Creates mesh by reading data stream in MFEM, Netgen, or VTK format. If
623  generate_edges = 0 (default) edges are not generated, if 1 edges are
624  generated. */
625  explicit Mesh(std::istream &input, int generate_edges = 0, int refine = 1,
626  bool fix_orientation = true);
627 
628  /// Create a disjoint mesh from the given mesh array
629  Mesh(Mesh *mesh_array[], int num_pieces);
630 
631  /// Create a uniformly refined (by any factor) version of @a orig_mesh.
632  /** @param[in] orig_mesh The starting coarse mesh.
633  @param[in] ref_factor The refinement factor, an integer > 1.
634  @param[in] ref_type Specify the positions of the new vertices. The
635  options are BasisType::ClosedUniform or
636  BasisType::GaussLobatto.
637 
638  The refinement data which can be accessed with GetRefinementTransforms()
639  is set to reflect the performed refinements.
640 
641  @note The constructed Mesh is linear, i.e. it does not have nodes. */
642  Mesh(Mesh *orig_mesh, int ref_factor, int ref_type);
643 
644  /** This is similar to the mesh constructor with the same arguments, but here
645  the current mesh is destroyed and another one created based on the data
646  stream again given in MFEM, Netgen, or VTK format. If generate_edges = 0
647  (default) edges are not generated, if 1 edges are generated. */
648  /// \see mfem::igzstream() for on-the-fly decompression of compressed ascii
649  /// inputs.
650  virtual void Load(std::istream &input, int generate_edges = 0,
651  int refine = 1, bool fix_orientation = true)
652  {
653  Loader(input, generate_edges);
654  Finalize(refine, fix_orientation);
655  }
656 
657  /// Clear the contents of the Mesh.
658  void Clear() { Destroy(); SetEmpty(); }
659 
660  /** @brief Get the mesh generator/type.
661 
662  @return A bitmask:
663  - bit 0 - simplices are present in the mesh (triangles, tets),
664  - bit 1 - tensor product elements are present in the mesh (quads, hexes),
665  - bit 2 - the mesh has wedge elements.
666 
667  In parallel, the result takes into account elements on all processors.
668  */
669  inline int MeshGenerator() { return meshgen; }
670 
671  /** @brief Returns number of vertices. Vertices are only at the corners of
672  elements, where you would expect them in the lowest-order mesh. */
673  inline int GetNV() const { return NumOfVertices; }
674 
675  /// Returns number of elements.
676  inline int GetNE() const { return NumOfElements; }
677 
678  /// Returns number of boundary elements.
679  inline int GetNBE() const { return NumOfBdrElements; }
680 
681  /// Return the number of edges.
682  inline int GetNEdges() const { return NumOfEdges; }
683 
684  /// Return the number of faces in a 3D mesh.
685  inline int GetNFaces() const { return NumOfFaces; }
686 
687  /// Return the number of faces (3D), edges (2D) or vertices (1D).
688  int GetNumFaces() const;
689 
690  /// Utility function: sum integers from all processors (Allreduce).
691  virtual long ReduceInt(int value) const { return value; }
692 
693  /// Return the total (global) number of elements.
694  long GetGlobalNE() const { return ReduceInt(NumOfElements); }
695 
696  /** @brief Return the mesh geometric factors corresponding to the given
697  integration rule. */
699  const int flags);
700 
701  /// Destroy all GeometricFactors stored by the Mesh.
702  /** This method can be used to force recomputation of the GeometricFactors,
703  for example, after the mesh nodes are modified externally. */
704  void DeleteGeometricFactors();
705 
706  /// Equals 1 + num_holes - num_loops
707  inline int EulerNumber() const
709  /// Equals 1 - num_holes
710  inline int EulerNumber2D() const
711  { return NumOfVertices - NumOfEdges + NumOfElements; }
712 
713  int Dimension() const { return Dim; }
714  int SpaceDimension() const { return spaceDim; }
715 
716  /// @brief Return pointer to vertex i's coordinates.
717  /// @warning For high-order meshes (when #Nodes != NULL) vertices may not be
718  /// updated and should not be used!
719  const double *GetVertex(int i) const { return vertices[i](); }
720 
721  /// @brief Return pointer to vertex i's coordinates.
722  /// @warning For high-order meshes (when Nodes != NULL) vertices may not
723  /// being updated and should not be used!
724  double *GetVertex(int i) { return vertices[i](); }
725 
726  void GetElementData(int geom, Array<int> &elem_vtx, Array<int> &attr) const
727  { GetElementData(elements, geom, elem_vtx, attr); }
728 
729  void GetBdrElementData(int geom, Array<int> &bdr_elem_vtx,
730  Array<int> &bdr_attr) const
731  { GetElementData(boundary, geom, bdr_elem_vtx, bdr_attr); }
732 
733  /** @brief Set the internal Vertex array to point to the given @a vertices
734  array without assuming ownership of the pointer. */
735  /** If @a zerocopy is `true`, the vertices must be given as an array of 3
736  doubles per vertex. If @a zerocopy is `false` then the current Vertex
737  data is first copied to the @a vertices array. */
738  void ChangeVertexDataOwnership(double *vertices, int len_vertices,
739  bool zerocopy = false);
740 
741  const Element* const *GetElementsArray() const
742  { return elements.GetData(); }
743 
744  const Element *GetElement(int i) const { return elements[i]; }
745 
746  Element *GetElement(int i) { return elements[i]; }
747 
748  const Element *GetBdrElement(int i) const { return boundary[i]; }
749 
750  Element *GetBdrElement(int i) { return boundary[i]; }
751 
752  const Element *GetFace(int i) const { return faces[i]; }
753 
755  {
756  return faces[i]->GetGeometryType();
757  }
758 
760  {
761  return elements[i]->GetGeometryType();
762  }
763 
765  {
766  return boundary[i]->GetGeometryType();
767  }
768 
769  /** @brief Return true iff the given @a geom is encountered in the mesh.
770  Geometries of dimensions lower than Dimension() are counted as well. */
772  { return mesh_geoms & (1 << geom); }
773 
774  /** @brief Return the number of geometries of the given dimension present in
775  the mesh. */
776  /** For a parallel mesh only the local geometries are counted. */
777  int GetNumGeometries(int dim) const;
778 
779  /// Return all element geometries of the given dimension present in the mesh.
780  /** For a parallel mesh only the local geometries are returned.
781 
782  The returned geometries are sorted. */
783  void GetGeometries(int dim, Array<Geometry::Type> &el_geoms) const;
784 
785  /// List of mesh geometries stored as Array<Geometry::Type>.
786  class GeometryList : public Array<Geometry::Type>
787  {
788  protected:
790  public:
791  /// Construct a GeometryList of all element geometries in @a mesh.
793  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
794  { mesh.GetGeometries(mesh.Dimension(), *this); }
795  /** @brief Construct a GeometryList of all geometries of dimension @a dim
796  in @a mesh. */
797  GeometryList(Mesh &mesh, int dim)
798  : Array<Geometry::Type>(geom_buf, Geometry::NumGeom)
799  { mesh.GetGeometries(dim, *this); }
800  };
801 
802  /// Returns the indices of the vertices of element i.
803  void GetElementVertices(int i, Array<int> &v) const
804  { elements[i]->GetVertices(v); }
805 
806  /// Returns the indices of the vertices of boundary element i.
807  void GetBdrElementVertices(int i, Array<int> &v) const
808  { boundary[i]->GetVertices(v); }
809 
810  /// Return the indices and the orientations of all edges of element i.
811  void GetElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
812 
813  /// Return the indices and the orientations of all edges of bdr element i.
814  void GetBdrElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
815 
816  /** Return the indices and the orientations of all edges of face i.
817  Works for both 2D (face=edge) and 3D faces. */
818  void GetFaceEdges(int i, Array<int> &, Array<int> &) const;
819 
820  /// Returns the indices of the vertices of face i.
821  void GetFaceVertices(int i, Array<int> &vert) const
822  {
823  if (Dim == 1)
824  {
825  vert.SetSize(1); vert[0] = i;
826  }
827  else
828  {
829  faces[i]->GetVertices(vert);
830  }
831  }
832 
833  /// Returns the indices of the vertices of edge i.
834  void GetEdgeVertices(int i, Array<int> &vert) const;
835 
836  /// Returns the face-to-edge Table (3D)
837  Table *GetFaceEdgeTable() const;
838 
839  /// Returns the edge-to-vertex Table (3D)
840  Table *GetEdgeVertexTable() const;
841 
842  /// Return the indices and the orientations of all faces of element i.
843  void GetElementFaces(int i, Array<int> &, Array<int> &) const;
844 
845  /// Return the index and the orientation of the face of bdr element i. (3D)
846  void GetBdrElementFace(int i, int *, int *) const;
847 
848  /** Return the vertex index of boundary element i. (1D)
849  Return the edge index of boundary element i. (2D)
850  Return the face index of boundary element i. (3D) */
851  int GetBdrElementEdgeIndex(int i) const;
852 
853  /** @brief For the given boundary element, bdr_el, return its adjacent
854  element and its info, i.e. 64*local_bdr_index+bdr_orientation. */
855  void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const;
856 
857  /// Returns the type of element i.
858  Element::Type GetElementType(int i) const;
859 
860  /// Returns the type of boundary element i.
861  Element::Type GetBdrElementType(int i) const;
862 
863  /* Return point matrix of element i of dimension Dim X #v, where for every
864  vertex we give its coordinates in space of dimension Dim. */
865  void GetPointMatrix(int i, DenseMatrix &pointmat) const;
866 
867  /* Return point matrix of boundary element i of dimension Dim X #v, where for
868  every vertex we give its coordinates in space of dimension Dim. */
869  void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const;
870 
872 
873  /** Builds the transformation defining the i-th element in the user-defined
874  variable. */
876 
877  /// Returns the transformation defining the i-th element
879 
880  /** Return the transformation defining the i-th element assuming
881  the position of the vertices/nodes are given by 'nodes'. */
882  void GetElementTransformation(int i, const Vector &nodes,
884 
885  /// Returns the transformation defining the i-th boundary element
888 
889  /** @brief Returns the transformation defining the given face element in a
890  user-defined variable. */
892 
893  /** @brief A helper method that constructs a transformation from the
894  reference space of a face to the reference space of an element. */
895  /** The local index of the face as a face in the element and its orientation
896  are given by the input parameter @a info, as @a info = 64*loc_face_idx +
897  loc_face_orientation. */
898  void GetLocalFaceTransformation(int face_type, int elem_type,
900  int info);
901 
902  /// Returns the transformation defining the given face element
904 
905  /** Returns the transformation defining the given edge element.
906  The transformation is stored in a user-defined variable. */
908 
909  /// Returns the transformation defining the given face element
911 
912  /// Returns (a pointer to a structure containing) the following data:
913  ///
914  /// 1) Elem1No - the index of the first element that contains this face this
915  /// is the element that has the same outward unit normal vector as the
916  /// face;
917  ///
918  /// 2) Elem2No - the index of the second element that contains this face this
919  /// element has outward unit normal vector as the face multiplied with -1;
920  ///
921  /// 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first
922  /// and the second element respectively;
923  ///
924  /// 4) Face - pointer to the ElementTransformation of the face;
925  ///
926  /// 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face
927  /// coordinate system to the element coordinate system (both in their
928  /// reference elements). Used to transform IntegrationPoints from face to
929  /// element. More formally, let:
930  /// TL1, TL2 be the transformations represented by Loc1, Loc2,
931  /// TE1, TE2 - the transformations represented by Elem1, Elem2,
932  /// TF - the transformation represented by Face, then
933  /// TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
934  ///
935  /// 6) FaceGeom - the base geometry for the face.
936  ///
937  /// The mask specifies which fields in the structure to return:
938  /// mask & 1 - Elem1, mask & 2 - Elem2
939  /// mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.
941  int mask = 31);
942 
944  {
945  if (faces_info[FaceNo].Elem2No < 0) { return NULL; }
946  return GetFaceElementTransformations (FaceNo);
947  }
948 
950 
951  /// Return true if the given face is interior. @sa FaceIsTrueInterior().
952  bool FaceIsInterior(int FaceNo) const
953  {
954  return (faces_info[FaceNo].Elem2No >= 0);
955  }
956  void GetFaceElements (int Face, int *Elem1, int *Elem2) const;
957  void GetFaceInfos (int Face, int *Inf1, int *Inf2) const;
958 
959  Geometry::Type GetFaceGeometryType(int Face) const;
960  Element::Type GetFaceElementType(int Face) const;
961 
962  /// Check the orientation of the elements
963  /** @return The number of elements with wrong orientation. */
964  int CheckElementOrientation(bool fix_it = true);
965  /// Check the orientation of the boundary elements
966  /** @return The number of boundary elements with wrong orientation. */
967  int CheckBdrElementOrientation(bool fix_it = true);
968 
969  /// Return the attribute of element i.
970  int GetAttribute(int i) const { return elements[i]->GetAttribute(); }
971 
972  /// Set the attribute of element i.
973  void SetAttribute(int i, int attr) { elements[i]->SetAttribute(attr); }
974 
975  /// Return the attribute of boundary element i.
976  int GetBdrAttribute(int i) const { return boundary[i]->GetAttribute(); }
977 
978  const Table &ElementToElementTable();
979 
980  const Table &ElementToFaceTable() const;
981 
982  const Table &ElementToEdgeTable() const;
983 
984  /// The returned Table must be destroyed by the caller
986 
987  /** Return the "face"-element Table. Here "face" refers to face (3D),
988  edge (2D), or vertex (1D).
989  The returned Table must be destroyed by the caller. */
990  Table *GetFaceToElementTable() const;
991 
992  /** This method modifies a tetrahedral mesh so that Nedelec spaces of order
993  greater than 1 can be defined on the mesh. Specifically, we
994  1) rotate all tets in the mesh so that the vertices {v0, v1, v2, v3}
995  satisfy: v0 < v1 < min(v2, v3).
996  2) rotate all boundary triangles so that the vertices {v0, v1, v2}
997  satisfy: v0 < min(v1, v2).
998 
999  @note Refinement does not work after a call to this method! */
1000  virtual void ReorientTetMesh();
1001 
1002  int *CartesianPartitioning(int nxyz[]);
1003  int *GeneratePartitioning(int nparts, int part_method = 1);
1004  void CheckPartitioning(int *partitioning);
1005 
1006  void CheckDisplacements(const Vector &displacements, double &tmax);
1007 
1008  // Vertices are only at the corners of elements, where you would expect them
1009  // in the lowest-order mesh.
1010  void MoveVertices(const Vector &displacements);
1011  void GetVertices(Vector &vert_coord) const;
1012  void SetVertices(const Vector &vert_coord);
1013 
1014  // Nodes are only active for higher order meshes, and share locations with
1015  // the vertices, plus all the higher- order control points within the element
1016  // and along the edges and on the faces.
1017  void GetNode(int i, double *coord);
1018  void SetNode(int i, const double *coord);
1019 
1020  // Node operations for curved mesh.
1021  // They call the corresponding '...Vertices' method if the
1022  // mesh is not curved (i.e. Nodes == NULL).
1023  void MoveNodes(const Vector &displacements);
1024  void GetNodes(Vector &node_coord) const;
1025  void SetNodes(const Vector &node_coord);
1026 
1027  /// Return a pointer to the internal node GridFunction (may be NULL).
1028  GridFunction *GetNodes() { return Nodes; }
1029  const GridFunction *GetNodes() const { return Nodes; }
1030  /// Return the mesh nodes ownership flag.
1031  bool OwnsNodes() const { return own_nodes; }
1032  /// Set the mesh nodes ownership flag.
1033  void SetNodesOwner(bool nodes_owner) { own_nodes = nodes_owner; }
1034  /// Replace the internal node GridFunction with the given GridFunction.
1035  void NewNodes(GridFunction &nodes, bool make_owner = false);
1036  /** Swap the internal node GridFunction pointer and ownership flag members
1037  with the given ones. */
1038  void SwapNodes(GridFunction *&nodes, int &own_nodes_);
1039 
1040  /// Return the mesh nodes/vertices projected on the given GridFunction.
1041  void GetNodes(GridFunction &nodes) const;
1042  /** Replace the internal node GridFunction with a new GridFunction defined
1043  on the given FiniteElementSpace. The new node coordinates are projected
1044  (derived) from the current nodes/vertices. */
1045  void SetNodalFESpace(FiniteElementSpace *nfes);
1046  /** Replace the internal node GridFunction with the given GridFunction. The
1047  given GridFunction is updated with node coordinates projected (derived)
1048  from the current nodes/vertices. */
1049  void SetNodalGridFunction(GridFunction *nodes, bool make_owner = false);
1050  /** Return the FiniteElementSpace on which the current mesh nodes are
1051  defined or NULL if the mesh does not have nodes. */
1052  const FiniteElementSpace *GetNodalFESpace() const;
1053  /** Make sure that the mesh has valid nodes, i.e. its geometry is described
1054  by a vector finite element grid function (even if it is a low-order mesh
1055  with straight edges). */
1056  void EnsureNodes();
1057 
1058  /** Set the curvature of the mesh nodes using the given polynomial degree,
1059  'order', and optionally: discontinuous or continuous FE space, 'discont',
1060  new space dimension, 'space_dim' (if != -1), and 'ordering'. */
1061  void SetCurvature(int order, bool discont = false, int space_dim = -1,
1062  int ordering = 1);
1063 
1064  /// Refine all mesh elements.
1065  /** @param[in] ref_algo %Refinement algorithm. Currently used only for pure
1066  tetrahedral meshes. If set to zero (default), a tet mesh will be refined
1067  using algorithm A, that produces elements with better quality compared to
1068  algorithm B used when the parameter is non-zero.
1069 
1070  For tetrahedral meshes, after using algorithm A, the mesh cannot be
1071  refined locally using methods like GeneralRefinement() unless it is
1072  re-finalized using Finalize() with the parameter @a refine set to true.
1073  Note that calling Finalize() in this way will generally invalidate any
1074  FiniteElementSpace%s and GridFunction%s defined on the mesh. */
1075  void UniformRefinement(int ref_algo = 0);
1076 
1077  /** Refine selected mesh elements. Refinement type can be specified for each
1078  element. The function can do conforming refinement of triangles and
1079  tetrahedra and non-conforming refinement (i.e., with hanging-nodes) of
1080  triangles, quadrilaterals and hexahedra. If 'nonconforming' = -1,
1081  suitable refinement method is selected automatically (namely, conforming
1082  refinement for triangles). Use nonconforming = 0/1 to force the method.
1083  For nonconforming refinements, nc_limit optionally specifies the maximum
1084  level of hanging nodes (unlimited by default). */
1085  void GeneralRefinement(const Array<Refinement> &refinements,
1086  int nonconforming = -1, int nc_limit = 0);
1087 
1088  /** Simplified version of GeneralRefinement taking a simple list of elements
1089  to refine, without refinement types. */
1090  void GeneralRefinement(const Array<int> &el_to_refine,
1091  int nonconforming = -1, int nc_limit = 0);
1092 
1093  /// Refine each element with given probability. Uses GeneralRefinement.
1094  void RandomRefinement(double prob, bool aniso = false,
1095  int nonconforming = -1, int nc_limit = 0);
1096 
1097  /// Refine elements sharing the specified vertex. Uses GeneralRefinement.
1098  void RefineAtVertex(const Vertex& vert,
1099  double eps = 0.0, int nonconforming = -1);
1100 
1101  /** Refine element i if elem_error[i] > threshold, for all i.
1102  Returns true if at least one element was refined, false otherwise. */
1103  bool RefineByError(const Array<double> &elem_error, double threshold,
1104  int nonconforming = -1, int nc_limit = 0);
1105 
1106  /** Refine element i if elem_error(i) > threshold, for all i.
1107  Returns true if at least one element was refined, false otherwise. */
1108  bool RefineByError(const Vector &elem_error, double threshold,
1109  int nonconforming = -1, int nc_limit = 0);
1110 
1111  /** Derefine the mesh based on an error measure associated with each
1112  element. A derefinement is performed if the sum of errors of its fine
1113  elements is smaller than 'threshold'. If 'nc_limit' > 0, derefinements
1114  that would increase the maximum level of hanging nodes of the mesh are
1115  skipped. Returns true if the mesh changed, false otherwise. */
1116  bool DerefineByError(Array<double> &elem_error, double threshold,
1117  int nc_limit = 0, int op = 1);
1118 
1119  /// Same as DerefineByError for an error vector.
1120  bool DerefineByError(const Vector &elem_error, double threshold,
1121  int nc_limit = 0, int op = 1);
1122 
1123  ///@{ @name NURBS mesh refinement methods
1124  void KnotInsert(Array<KnotVector *> &kv);
1125  /* For each knot vector:
1126  new_degree = max(old_degree, min(old_degree + rel_degree, degree)). */
1127  void DegreeElevate(int rel_degree, int degree = 16);
1128  ///@}
1129 
1130  /** Make sure that a quad/hex mesh is considered to be non-conforming (i.e.,
1131  has an associated NCMesh object). Triangles meshes can be both conforming
1132  (default) or non-conforming. */
1133  void EnsureNCMesh(bool triangles_nonconforming = false);
1134 
1135  bool Conforming() const { return ncmesh == NULL; }
1136  bool Nonconforming() const { return ncmesh != NULL; }
1137 
1138  /** Return fine element transformations following a mesh refinement.
1139  Space uses this to construct a global interpolation matrix. */
1141 
1142  /// Return type of last modification of the mesh.
1144 
1145  /** Return update counter. The counter starts at zero and is incremented
1146  each time refinement, derefinement, or rebalancing method is called.
1147  It is used for checking proper sequence of Space:: and GridFunction::
1148  Update() calls. */
1149  long GetSequence() const { return sequence; }
1150 
1151  /// Print the mesh to the given stream using Netgen/Truegrid format.
1152  virtual void PrintXG(std::ostream &out = mfem::out) const;
1153 
1154  /// Print the mesh to the given stream using the default MFEM mesh format.
1155  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1156  virtual void Print(std::ostream &out = mfem::out) const { Printer(out); }
1157 
1158  /// Print the mesh in VTK format (linear and quadratic meshes only).
1159  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1160  void PrintVTK(std::ostream &out);
1161 
1162  /** Print the mesh in VTK format. The parameter ref > 0 specifies an element
1163  subdivision number (useful for high order fields and curved meshes).
1164  If the optional field_data is set, we also add a FIELD section in the
1165  beginning of the file with additional dataset information. */
1166  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1167  void PrintVTK(std::ostream &out, int ref, int field_data=0);
1168 
1169  void GetElementColoring(Array<int> &colors, int el0 = 0);
1170 
1171  /** @brief Prints the mesh with boundary elements given by the boundary of
1172  the subdomains, so that the boundary of subdomain i has boundary
1173  attribute i+1. */
1174  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1175  void PrintWithPartitioning (int *partitioning,
1176  std::ostream &out, int elem_attr = 0) const;
1177 
1178  void PrintElementsWithPartitioning (int *partitioning,
1179  std::ostream &out,
1180  int interior_faces = 0);
1181 
1182  /// Print set of disjoint surfaces:
1183  /*!
1184  * If Aface_face(i,j) != 0, print face j as a boundary
1185  * element with attribute i+1.
1186  */
1187  void PrintSurfaces(const Table &Aface_face, std::ostream &out) const;
1188 
1189  void ScaleSubdomains (double sf);
1190  void ScaleElements (double sf);
1191 
1192  void Transform(void (*f)(const Vector&, Vector&));
1193  void Transform(VectorCoefficient &deformation);
1194 
1195  /// Remove unused vertices and rebuild mesh connectivity.
1196  void RemoveUnusedVertices();
1197 
1198  /** Remove boundary elements that lie in the interior of the mesh, i.e. that
1199  have two adjacent faces in 3D, or edges in 2D. */
1200  void RemoveInternalBoundaries();
1201 
1202  /** @brief Get the size of the i-th element relative to the perfect
1203  reference element. */
1204  double GetElementSize(int i, int type = 0);
1205 
1206  double GetElementSize(int i, const Vector &dir);
1207 
1208  double GetElementVolume(int i);
1209 
1210  /// Returns the minimum and maximum corners of the mesh bounding box.
1211  /** For high-order meshes, the geometry is first refined @a ref times. */
1212  void GetBoundingBox(Vector &min, Vector &max, int ref = 2);
1213 
1214  void GetCharacteristics(double &h_min, double &h_max,
1215  double &kappa_min, double &kappa_max,
1216  Vector *Vh = NULL, Vector *Vk = NULL);
1217 
1218  /// Auxiliary method used by PrintCharacteristics().
1219  /** It is also used in the `mesh-explorer` miniapp. */
1220  static void PrintElementsByGeometry(int dim,
1221  const Array<int> &num_elems_by_geom,
1222  std::ostream &out);
1223 
1224  /** @brief Compute and print mesh characteristics such as number of vertices,
1225  number of elements, number of boundary elements, minimal and maximal
1226  element sizes, minimal and maximal element aspect ratios, etc. */
1227  /** If @a Vh or @a Vk are not NULL, return the element sizes and aspect
1228  ratios for all elements in the given Vector%s. */
1229  void PrintCharacteristics(Vector *Vh = NULL, Vector *Vk = NULL,
1230  std::ostream &out = mfem::out);
1231 
1232  /** @brief In serial, this method calls PrintCharacteristics(). In parallel,
1233  additional information about the parallel decomposition is also printed.
1234  */
1235  virtual void PrintInfo(std::ostream &out = mfem::out)
1236  {
1237  PrintCharacteristics(NULL, NULL, out);
1238  }
1239 
1240  void MesquiteSmooth(const int mesquite_option = 0);
1241 
1242  /** @brief Find the ids of the elements that contain the given points, and
1243  their corresponding reference coordinates.
1244 
1245  The DenseMatrix @a point_mat describes the given points - one point for
1246  each column; it should have SpaceDimension() rows.
1247 
1248  The InverseElementTransformation object, @a inv_trans, is used to attempt
1249  the element transformation inversion. If NULL pointer is given, the
1250  method will use a default constructed InverseElementTransformation. Note
1251  that the algorithms in the base class InverseElementTransformation can be
1252  completely overwritten by deriving custom classes that override the
1253  Transform() method.
1254 
1255  If no element is found for the i-th point, elem_ids[i] is set to -1.
1256 
1257  In the ParMesh implementation, the @a point_mat is expected to be the
1258  same on all ranks. If the i-th point is found by multiple ranks, only one
1259  of them will mark that point as found, i.e. set its elem_ids[i] to a
1260  non-negative number; the other ranks will set their elem_ids[i] to -2 to
1261  indicate that the point was found but assigned to another rank.
1262 
1263  @returns The total number of points that were found.
1264 
1265  @note This method is not 100 percent reliable, i.e. it is not guaranteed
1266  to find a point, even if it lies inside a mesh element. */
1267  virtual int FindPoints(DenseMatrix& point_mat, Array<int>& elem_ids,
1268  Array<IntegrationPoint>& ips, bool warn = true,
1269  InverseElementTransformation *inv_trans = NULL);
1270 
1271  /// Destroys Mesh.
1272  virtual ~Mesh() { DestroyPointers(); }
1273 };
1274 
1275 /** Overload operator<< for std::ostream and Mesh; valid also for the derived
1276  class ParMesh */
1277 std::ostream &operator<<(std::ostream &out, const Mesh &mesh);
1278 
1279 
1280 /** @brief Structure for storing mesh geometric factors: coordinates, Jacobians,
1281  and determinants of the Jacobians. */
1282 /** Typically objects of this type are constructed and owned by objects of class
1283  Mesh. See Mesh::GetGeometricFactors(). */
1285 {
1286 public:
1287  const Mesh *mesh;
1290 
1292  {
1293  COORDINATES = 1 << 0,
1294  JACOBIANS = 1 << 1,
1295  DETERMINANTS = 1 << 2,
1296  };
1297 
1298  GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags);
1299 
1300  /// Mapped (physical) coordinates of all quadrature points.
1301  /** This array uses a column-major layout with dimensions (NQ x SDIM x NE)
1302  where
1303  - NQ = number of quadrature points per element,
1304  - SDIM = space dimension of the mesh = mesh.SpaceDimension(), and
1305  - NE = number of elements in the mesh. */
1307 
1308  /// Jacobians of the element transformations at all quadrature points.
1309  /** This array uses a column-major layout with dimensions (NQ x SDIM x DIM x
1310  NE) where
1311  - NQ = number of quadrature points per element,
1312  - SDIM = space dimension of the mesh = mesh.SpaceDimension(),
1313  - DIM = dimension of the mesh = mesh.Dimension(), and
1314  - NE = number of elements in the mesh. */
1316 
1317  /// Determinants of the Jacobians at all quadrature points.
1318  /** This array uses a column-major layout with dimensions (NQ x NE) where
1319  - NQ = number of quadrature points per element, and
1320  - NE = number of elements in the mesh. */
1322 };
1323 
1324 
1325 /// Class used to extrude the nodes of a mesh
1327 {
1328 private:
1329  int n, layer;
1330  double p[2], s;
1331  Vector tip;
1332 public:
1333  NodeExtrudeCoefficient(const int dim, const int _n, const double _s);
1334  void SetLayer(const int l) { layer = l; }
1336  virtual void Eval(Vector &V, ElementTransformation &T,
1337  const IntegrationPoint &ip);
1339 };
1340 
1341 
1342 /// Extrude a 1D mesh
1343 Mesh *Extrude1D(Mesh *mesh, const int ny, const double sy,
1344  const bool closed = false);
1345 
1346 /// Extrude a 2D mesh
1347 Mesh *Extrude2D(Mesh *mesh, const int nz, const double sz);
1348 
1349 
1350 /// Input file stream that remembers the input file name (useful for example
1351 /// when reading NetCDF meshes) and supports optional gzstream decompression.
1353 {
1354 public:
1355  const char *filename;
1356  named_ifgzstream(const char *mesh_name) :
1357  mfem::ifgzstream(mesh_name), filename(mesh_name) {}
1358 };
1359 
1360 
1361 // inline functions
1362 inline void Mesh::ShiftL2R(int &a, int &b, int &c)
1363 {
1364  int t = a;
1365  a = c; c = b; b = t;
1366 }
1367 
1368 inline void Mesh::Rotate3(int &a, int &b, int &c)
1369 {
1370  if (a < b)
1371  {
1372  if (a > c)
1373  {
1374  ShiftL2R(a, b, c);
1375  }
1376  }
1377  else
1378  {
1379  if (b < c)
1380  {
1381  ShiftL2R(c, b, a);
1382  }
1383  else
1384  {
1385  ShiftL2R(a, b, c);
1386  }
1387  }
1388 }
1389 
1390 }
1391 
1392 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:229
void AddHex(const int *vi, int attr=1)
Definition: mesh.cpp:1161
void Loader(std::istream &input, int generate_edges=0, std::string parse_tag="")
Definition: mesh.cpp:2959
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:9443
void GetFaceEdges(int i, Array< int > &, Array< int > &) const
Definition: mesh.cpp:4096
void PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
Definition: mesh.cpp:8973
void GetPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:4361
static const int vtk_quadratic_hex[27]
Definition: mesh.hpp:159
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1156
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:4977
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:85
const DenseMatrix * PointMatrix
Definition: mesh.hpp:122
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:976
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:719
void ScaleElements(double sf)
Definition: mesh.cpp:9110
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Table * GetEdgeVertexTable() const
Returns the edge-to-vertex Table (3D)
Definition: mesh.cpp:4163
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:4064
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo)
void FreeElement(Element *E)
Definition: mesh.cpp:9424
int CheckElementOrientation(bool fix_it=true)
Check the orientation of the elements.
Definition: mesh.cpp:3678
void UniformRefinement3D_base(Array< int > *f2qf=NULL, DSTable *v_to_v_p=NULL)
Definition: mesh.cpp:6064
static bool remove_unused_vertices
Definition: mesh.hpp:188
void SetVertices(const Vector &vert_coord)
Definition: mesh.cpp:5777
static const int vtk_quadratic_tet[10]
Definition: mesh.hpp:157
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:4126
void Make2D(int nx, int ny, Element::Type type, double sx, double sy, bool generate_edges, bool sfc_ordering)
Definition: mesh.cpp:2410
bool FaceIsTrueInterior(int FaceNo) const
Definition: mesh.hpp:384
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition: mesh.hpp:1306
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:9615
long GetGlobalNE() const
Return the total (global) number of elements.
Definition: mesh.hpp:694
static void PrintElementsByGeometry(int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
Auxiliary method used by PrintCharacteristics().
Definition: mesh.cpp:217
void AddHexAsTets(const int *vi, int attr=1)
Definition: mesh.cpp:1166
static const int NumGeom
Definition: geom.hpp:41
void AddHexAsWedges(const int *vi, int attr=1)
Definition: mesh.cpp:1185
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:74
Geometry::Constants< Geometry::SQUARE > quad_t
Definition: mesh.hpp:169
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:5016
CoarseFineTransformations CoarseFineTr
Definition: mesh.hpp:149
int own_nodes
Definition: mesh.hpp:155
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:4018
bool Conforming() const
Definition: mesh.hpp:1135
void MoveVertices(const Vector &displacements)
Definition: mesh.cpp:5757
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:679
void DeleteLazyTables()
Definition: mesh.cpp:1061
IsoparametricTransformation Transformation
Definition: mesh.hpp:143
void GetLocalPtToSegTransformation(IsoparametricTransformation &, int)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:605
void GetCharacteristics(double &h_min, double &h_max, double &kappa_min, double &kappa_max, Vector *Vh=NULL, Vector *Vk=NULL)
Definition: mesh.cpp:182
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:8482
const Geometry::Type geom
Definition: ex1.cpp:40
static FiniteElement * GetTransformationFEforElementType(Element::Type)
Definition: mesh.cpp:316
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:5861
int NumOfEdges
Definition: mesh.hpp:58
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:5875
void GetElementFaces(int i, Array< int > &, Array< int > &) const
Return the indices and the orientations of all faces of element i.
Definition: mesh.cpp:4258
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:4351
void AddBdrElement(Element *elem)
Definition: mesh.hpp:505
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:118
void ReadNetgen2DMesh(std::istream &input, int &curved)
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:4356
int EulerNumber2D() const
Equals 1 - num_holes.
Definition: mesh.hpp:710
void AddTri(const int *vi, int attr=1)
Definition: mesh.cpp:1128
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:4029
void DeleteTables()
Definition: mesh.hpp:197
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:821
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
GridFunction * Nodes
Definition: mesh.hpp:154
int NumOfElements
Definition: mesh.hpp:57
Mesh * Extrude1D(Mesh *mesh, const int ny, const double sy, const bool closed)
Extrude a 1D mesh.
Definition: mesh.cpp:9633
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:9180
void FinalizeCheck()
Definition: mesh.cpp:1271
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:4319
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:5885
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
const Element * GetFace(int i) const
Definition: mesh.hpp:752
void FinalizeWedgeMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a wedge Mesh.
Definition: mesh.cpp:1952
Element * GetElement(int i)
Definition: mesh.hpp:746
void GenerateBoundaryElements()
Definition: mesh.cpp:1233
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)
Definition: mesh.hpp:589
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:754
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1284
Data type for vertex.
Definition: vertex.hpp:21
void SetMeshGen()
Determine the mesh generator bitmask meshgen, see MeshGenerator().
Definition: mesh.cpp:2915
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3859
void GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:5766
void ReadNetgen3DMesh(std::istream &input)
void AddBdrSegment(const int *vi, int attr=1)
Definition: mesh.cpp:1203
The inverse transformation of a given ElementTransformation.
Definition: eltrans.hpp:109
void ReadInlineMesh(std::istream &input, bool generate_edges=false)
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
Definition: mesh.hpp:943
void RemoveInternalBoundaries()
Definition: mesh.cpp:9336
Mesh * Extrude2D(Mesh *mesh, const int nz, const double sz)
Extrude a 2D mesh.
Definition: mesh.cpp:9793
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:831
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:759
void SetEmpty()
Definition: mesh.cpp:985
Array< Element * > faces
Definition: mesh.hpp:75
virtual long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
Definition: mesh.hpp:691
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
Definition: mesh.cpp:1987
void GetVertexToVertexTable(DSTable &) const
Definition: mesh.cpp:4430
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:256
Operation last_operation
Definition: mesh.hpp:191
void MarkTriMeshForRefinement()
Definition: mesh.cpp:1555
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:771
void InitRefinementTransforms()
Definition: mesh.cpp:7681
IsoparametricTransformation FaceTransformation
Definition: mesh.hpp:145
Array< NCFaceInfo > nc_faces_info
Definition: mesh.hpp:132
Element * ReadElement(std::istream &)
Definition: mesh.cpp:2897
void KnotInsert(Array< KnotVector * > &kv)
Definition: mesh.cpp:3410
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
Element * NewElement(int geom)
Definition: mesh.cpp:2844
FaceElementTransformations FaceElemTr
Definition: mesh.hpp:146
Table * el_to_face
Definition: mesh.hpp:135
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
Definition: mesh.cpp:542
const GridFunction * GetNodes() const
Definition: mesh.hpp:1029
void MarkForRefinement()
Definition: mesh.cpp:1538
Element::Type GetFaceElementType(int Face) const
Definition: mesh.cpp:958
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:414
void AddVertex(const double *)
Definition: mesh.cpp:1117
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:477
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
Definition: mesh.cpp:8018
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:3663
Geometry::Constants< Geometry::SEGMENT > seg_t
Definition: mesh.hpp:167
void AddElement(Element *elem)
The parameter elem should be allocated using the NewElement() method.
Definition: mesh.hpp:504
void GetElementJacobian(int i, DenseMatrix &J)
Definition: mesh.cpp:58
void MesquiteSmooth(const int mesquite_option=0)
Definition: mesquite.cpp:1144
long GetSequence() const
Definition: mesh.hpp:1149
bool Nonconforming() const
Definition: mesh.hpp:1136
Vector J
Jacobians of the element transformations at all quadrature points.
Definition: mesh.hpp:1315
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags)
Return the mesh geometric factors corresponding to the given integration rule.
Definition: mesh.cpp:756
void MoveNodes(const Vector &displacements)
Definition: mesh.cpp:5825
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
Definition: mesh.cpp:1570
void GetElementData(int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.hpp:726
int dim
Definition: ex3.cpp:48
NCFaceInfo(bool slave, int master, const DenseMatrix *pm)
Definition: mesh.hpp:127
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a triangular Mesh.
Definition: mesh.cpp:1285
void UpdateNURBS()
Definition: mesh.cpp:3462
Symmetric 3D Table.
Definition: stable3d.hpp:29
void FinalizeTetMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a tetrahedral Mesh.
Definition: mesh.cpp:1911
void AddBdrQuadAsTriangles(const int *vi, int attr=1)
Definition: mesh.cpp:1218
void Make3D(int nx, int ny, int nz, Element::Type type, double sx, double sy, double sz, bool sfc_ordering)
Definition: mesh.cpp:2172
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
Definition: mesh.cpp:4565
void SetAttributes()
Definition: mesh.cpp:1069
GeometryList(Mesh &mesh, int dim)
Construct a GeometryList of all geometries of dimension dim in mesh.
Definition: mesh.hpp:797
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition: mesh.hpp:1321
int FindCoarseElement(int i)
Definition: mesh.cpp:7698
GeometryList(Mesh &mesh)
Construct a GeometryList of all element geometries in mesh.
Definition: mesh.hpp:792
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Returns the transformation defining the given face element in a user-defined variable.
Definition: mesh.cpp:473
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:620
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
Definition: mesh.cpp:4622
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
Definition: mesh.cpp:3525
STable3D * GetElementToFaceTable(int ret_ftbl=0)
Definition: mesh.cpp:4826
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3644
static void Rotate3(int &, int &, int &)
Definition: mesh.hpp:1368
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:669
const char * filename
Definition: mesh.hpp:1355
virtual void PrintInfo(std::ostream &out=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
Definition: mesh.hpp:1235
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:4377
void InitTables()
Definition: mesh.cpp:979
void CheckDisplacements(const Vector &displacements, double &tmax)
Definition: mesh.cpp:5680
void GetElementCenter(int i, Vector &c)
Definition: mesh.cpp:66
int mesh_geoms
Definition: mesh.hpp:61
Geometry::Constants< Geometry::TETRAHEDRON > tet_t
Definition: mesh.hpp:170
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:650
void SetNodesOwner(bool nodes_owner)
Set the mesh nodes ownership flag.
Definition: mesh.hpp:1033
void SetLayer(const int l)
Definition: mesh.hpp:1334
int GetElementToEdgeTable(Table &, Array< int > &)
Definition: mesh.cpp:4455
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:4331
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:7272
const Element * GetElement(int i) const
Definition: mesh.hpp:744
Mesh(int n, double sx=1.0)
Definition: mesh.hpp:610
void GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:710
IsoparametricTransformation BdrTransformation
Definition: mesh.hpp:144
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:784
IsoparametricTransformation Transformation2
Definition: mesh.hpp:143
void Init()
Definition: mesh.cpp:963
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3621
Element * GetBdrElement(int i)
Definition: mesh.hpp:750
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Definition: ncmesh.hpp:100
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:915
Geometry::Constants< Geometry::PRISM > pri_t
Definition: mesh.hpp:172
int Dimension() const
Definition: mesh.hpp:713
virtual void ReorientTetMesh()
Definition: mesh.cpp:4915
int NumOfBdrElements
Definition: mesh.hpp:57
Table * el_to_edge
Definition: mesh.hpp:134
static void ShiftL2R(int &, int &, int &)
Definition: mesh.hpp:1362
virtual void LocalRefinement(const Array< int > &marked_el, int type=3)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:6637
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
virtual void MarkTetMeshForRefinement(DSTable &v_to_v)
Definition: mesh.cpp:1595
List of mesh geometries stored as Array&lt;Geometry::Type&gt;.
Definition: mesh.hpp:786
void Destroy()
Definition: mesh.cpp:1033
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition: mesh.cpp:74
int CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Definition: mesh.cpp:3907
int SpaceDimension() const
Definition: mesh.hpp:714
void GetBdrElementFace(int i, int *, int *) const
Return the index and the orientation of the face of bdr element i. (3D)
Definition: mesh.cpp:4296
void CheckPartitioning(int *partitioning)
Definition: mesh.cpp:5362
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
Definition: mesh.cpp:8601
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
Definition: mesh.hpp:952
void ReadMFEMMesh(std::istream &input, bool mfem_v11, int &curved)
void GetLocalQuadToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:732
void AddPointFaceElement(int lf, int gf, int el)
Used in GenerateFaces()
Definition: mesh.cpp:4548
GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags)
Definition: mesh.cpp:9564
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:941
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:2588
const Element *const * GetElementsArray() const
Definition: mesh.hpp:741
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
void FinalizeTopology()
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2024
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:6879
void PrintVTK(std::ostream &out)
Definition: mesh.cpp:8060
virtual ~Mesh()
Destroys Mesh.
Definition: mesh.hpp:1272
void FinalizeMesh(int refine=0, bool fix_orientation=true)
Finalize the construction of any type of Mesh.
Definition: mesh.cpp:2017
Mesh(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.hpp:602
Table * el_to_el
Definition: mesh.hpp:136
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:947
Geometry::Constants< Geometry::CUBE > hex_t
Definition: mesh.hpp:171
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
Geometry::Type geom_buf[Geometry::NumGeom]
Definition: mesh.hpp:789
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Definition: mesh.cpp:9229
void ApplyLocalSlaveTransformation(IsoparametricTransformation &transf, const FaceInfo &fi)
Definition: mesh.cpp:901
static void PrintElement(const Element *, std::ostream &)
Definition: mesh.cpp:2909
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
Definition: mesh.cpp:3432
const IntegrationRule * IntRule
Definition: mesh.hpp:1288
static void GetElementArrayEdgeTable(const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
Definition: mesh.cpp:4408
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
Definition: mesh.hpp:1028
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:618
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
Definition: mesh.cpp:1618
bool IsSlaveFace(const FaceInfo &fi) const
Definition: mesh.cpp:896
void AddQuad(const int *vi, int attr=1)
Definition: mesh.cpp:1138
void ReadGmshMesh(std::istream &input)
Array< Vertex > vertices
Definition: mesh.hpp:73
void AddWedge(const int *vi, int attr=1)
Definition: mesh.cpp:1156
virtual void PrintXG(std::ostream &out=mfem::out) const
Print the mesh to the given stream using Netgen/Truegrid format.
Definition: mesh.cpp:7771
void RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
Definition: mesh.cpp:7291
Array< Element * > elements
Definition: mesh.hpp:68
const Table & ElementToElementTable()
Definition: mesh.cpp:4494
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
Definition: mesh.cpp:4393
int meshgen
Definition: mesh.hpp:60
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:7708
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:7317
Table * bel_to_edge
Definition: mesh.hpp:138
Operation GetLastOperation() const
Return type of last modification of the mesh.
Definition: mesh.hpp:1143
void AddBdrTriangle(const int *vi, int attr=1)
Definition: mesh.cpp:1208
Table * GetFaceToElementTable() const
Definition: mesh.cpp:4224
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:660
int EulerNumber() const
Equals 1 + num_holes - num_loops.
Definition: mesh.hpp:707
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:181
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:673
void AddBdrQuad(const int *vi, int attr=1)
Definition: mesh.cpp:1213
Array< int > be_to_edge
Definition: mesh.hpp:137
const Table & ElementToFaceTable() const
Definition: mesh.cpp:4530
void GetGeckoElementReordering(Array< int > &ordering, int iterations=1, int window=2, int period=1, int seed=0)
Definition: mesh.cpp:1344
Class for integration point with weight.
Definition: intrules.hpp:25
void Swap(Mesh &other, bool non_geometry)
Definition: mesh.cpp:7072
Element * ReadElementWithoutAttr(std::istream &)
Definition: mesh.cpp:2867
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:4042
static const int vtk_quadratic_wedge[18]
Definition: mesh.hpp:158
void AddTriangle(const int *vi, int attr=1)
Definition: mesh.cpp:1133
bool OwnsNodes() const
Return the mesh nodes ownership flag.
Definition: mesh.hpp:1031
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:336
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:640
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:3639
Table * face_edge
Definition: mesh.hpp:140
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1314
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
Definition: mesh.hpp:729
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2098
bool DerefineByError(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
Definition: mesh.cpp:7000
NodeExtrudeCoefficient(const int dim, const int _n, const double _s)
Definition: mesh.cpp:9609
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:1386
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:3445
Table * edge_vertex
Definition: mesh.hpp:141
virtual void UniformRefinement2D()
Refine a mixed 2D mesh uniformly.
Definition: mesh.cpp:5915
long sequence
Definition: mesh.hpp:66
virtual ~NodeExtrudeCoefficient()
Definition: mesh.hpp:1338
const Mesh * mesh
Definition: mesh.hpp:1287
void DestroyPointers()
Definition: mesh.cpp:1007
Array< FaceInfo > faces_info
Definition: mesh.hpp:131
STable3D * GetFacesTable()
Definition: mesh.cpp:4777
void EnsureNCMesh(bool triangles_nonconforming=false)
Definition: mesh.cpp:7253
NCMesh * ncmesh
Optional non-conforming mesh extension.
Definition: mesh.hpp:182
Geometry::Type GetBdrElementBaseGeometry(int i) const
Definition: mesh.hpp:764
virtual void UniformRefinement3D()
Refine a mixed 3D mesh uniformly.
Definition: mesh.hpp:298
static void PrintElementWithoutAttr(const Element *, std::ostream &)
Definition: mesh.cpp:2885
virtual bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NC version of GeneralDerefinement.
Definition: mesh.cpp:6948
int Dim
Definition: mesh.hpp:54
Geometry::Type GetFaceGeometryType(int Face) const
Definition: mesh.cpp:953
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:658
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:682
Geometry::Constants< Geometry::TRIANGLE > tri_t
Definition: mesh.hpp:168
double AggregateError(const Array< double > &elem_error, const int *fine, int nfine, int op)
Derefinement helper.
Definition: mesh.cpp:6929
double * GetVertex(int i)
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:724
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:685
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:9449
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:973
Vector data type.
Definition: vector.hpp:48
void ReadTrueGridMesh(std::istream &input)
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3811
void SetNode(int i, const double *coord)
Definition: mesh.cpp:5805
void GetNode(int i, double *coord)
Definition: mesh.cpp:5786
void DestroyTables()
Definition: mesh.cpp:991
void GetLocalTriToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:684
void GenerateFaces()
Definition: mesh.cpp:4645
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:8407
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:2778
void DeleteGeometricFactors()
Destroy all GeometricFactors stored by the Mesh.
Definition: mesh.cpp:775
int spaceDim
Definition: mesh.hpp:55
IsoparametricTransformation EdgeTransformation
Definition: mesh.hpp:145
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:59
void AddTet(const int *vi, int attr=1)
Definition: mesh.cpp:1143
MemAlloc< Tetrahedron, 1024 > TetMemory
Definition: mesh.hpp:163
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
Definition: mesh.cpp:4135
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
Definition: mesh.cpp:1684
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:5849
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:3633
named_ifgzstream(const char *mesh_name)
Definition: mesh.hpp:1356
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
int NumOfVertices
Definition: mesh.hpp:57
void ScaleSubdomains(double sf)
Definition: mesh.cpp:9040
Array< int > be_to_face
Definition: mesh.hpp:139
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:64
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
Definition: mesh.cpp:7029
void Bisection(int i, const DSTable &, int *, int *, int *)
Bisect a triangle: element with index i is bisected.
Definition: mesh.cpp:7346
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:7185
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:970
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
Definition: mesh.cpp:4598
const Table & ElementToEdgeTable() const
Definition: mesh.cpp:4539
Array< GeometricFactors * > geom_factors
Optional geometric factors.
Definition: mesh.hpp:183
void BdrBisection(int i, const HashTable< Hashed2 > &)
Bisect a boundary triangle: boundary element with index i is bisected.
Definition: mesh.cpp:7554
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
Definition: mesh.hpp:807
void GenerateNCFaceInfo()
Definition: mesh.cpp:4730
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:803
double GetElementVolume(int i)
Definition: mesh.cpp:101
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:177
void InitMesh(int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Definition: mesh.cpp:1100
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
Definition: mesh.cpp:4189
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:748
void EnsureNodes()
Definition: mesh.cpp:3627
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.cpp:7117
void UpdateNodes()
Update the nodes of a curved mesh after refinement.
Definition: mesh.cpp:5906
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:231
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:262
void Printer(std::ostream &out=mfem::out, std::string section_delimiter="") const
Definition: mesh.cpp:7935
Class used to extrude the nodes of a mesh.
Definition: mesh.hpp:1326
int NumOfFaces
Definition: mesh.hpp:58