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