MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_MESH
13 #define MFEM_MESH
14 
15 #include "../config/config.hpp"
16 #include "../general/stable3d.hpp"
17 #include "../general/globals.hpp"
18 #include "triangle.hpp"
19 #include "tetrahedron.hpp"
20 #include "vertex.hpp"
21 #include "ncmesh.hpp"
22 #include "../fem/eltrans.hpp"
23 #include "../fem/coefficient.hpp"
24 #include "../general/gzstream.hpp"
25 #include <iostream>
26 
27 namespace mfem
28 {
29 
30 // Data type mesh
31 
32 class KnotVector;
33 class NURBSExtension;
34 class FiniteElementSpace;
35 class GridFunction;
36 struct Refinement;
37 
38 #ifdef MFEM_USE_MPI
39 class ParMesh;
40 class ParNCMesh;
41 #endif
42 
43 
44 class Mesh
45 {
46 #ifdef MFEM_USE_MPI
47  friend class ParMesh;
48  friend class ParNCMesh;
49 #endif
50  friend class NURBSExtension;
51 
52 protected:
53  int Dim;
54  int spaceDim;
55 
58 
59  int BaseGeom, BaseBdrGeom; // element base geometries, -1 if not all the same
60 
61  int meshgen; // see MeshGenerator()
62 
63  // Counter for Mesh transformations: refinement, derefinement, rebalancing.
64  // Used for checking during Update operations on objects depending on the
65  // Mesh, such as FiniteElementSpace, GridFunction, etc.
66  long sequence;
67 
69  // Vertices are only at the corners of elements, where you would expect them
70  // in the lowest-order mesh. In some cases, e.g. in a Mesh that defines the
71  // patch topology for a NURBS mesh (see LoadPatchTopo()) the vertices may be
72  // empty while NumOfVertices is positive.
76 
77  struct FaceInfo
78  {
79  // Inf = 64 * LocalFaceIndex + FaceOrientation
81  int NCFace; /* -1 if this is a regular conforming/boundary face;
82  index into 'nc_faces_info' if >= 0. */
83  };
84  // NOTE: in NC meshes, master faces have Elem2No == -1. Slave faces on the
85  // other hand have Elem2No and Elem2Inf set to the master face's element and
86  // its local face number.
87  //
88  // A local face is one generated from a local element and has index i in
89  // faces_info such that i < GetNumFaces(). Also, Elem1No always refers to the
90  // element (slave or master, in the non-conforming case) that generated the
91  // face.
92  // Classification of a local (non-ghost) face based on its FaceInfo:
93  // - Elem2No >= 0 --> local internal face; can be either:
94  // - NCFace == -1 --> conforming face, or
95  // - NCFace >= 0 --> non-conforming slave face.
96  // - Elem2No < 0 --> local "boundary" face; can be one of:
97  // - NCFace == -1 --> conforming face; can be either:
98  // - Elem2Inf < 0 --> true boundary face (no element on side 2)
99  // - Elem2Inf >= 0 --> shared face where element 2 is a face-neighbor
100  // element with index -1-Elem2No. This state is initialized by
101  // ParMesh::ExchangeFaceNbrData().
102  // - NCFace >= 0 --> non-conforming master face. Elem2No is -1 or, in the
103  // case of a shared face, -1-Elem2No is the index of one of the adjacent
104  // (the last one?) slave ghost elements. Elem2Inf is -1.
105  //
106  // A ghost face is a non-conforming face that is generated by a non-local,
107  // i.e. ghost, element. A ghost face has index i in faces_info such that
108  // i >= GetNumFaces().
109  // Classification of a ghost (non-local) face based on its FaceInfo:
110  // - Elem1No == -1 --> master ghost face? These ghost faces also have:
111  // Elem2No == -1, Elem1Inf == Elem2Inf == -1, and NCFace == -1.
112  // - Elem1No >= 0 --> slave ghost face; Elem1No is the index of the local
113  // master side element, i.e. side 1 IS NOT the side that generated the
114  // face. Elem2No is < 0 and -1-Elem2No is the index of the ghost
115  // face-neighbor element that generated this slave ghost face. In this
116  // case, Elem2Inf >= 0.
117 
118  struct NCFaceInfo
119  {
120  bool Slave; // true if this is a slave face, false if master face
121  int MasterFace; // if Slave, this is the index of the master face
122  const DenseMatrix* PointMatrix; // if Slave, position within master face
123  // (NOTE: PointMatrix points to a matrix owned by NCMesh.)
124 
125  NCFaceInfo(bool slave, int master, const DenseMatrix* pm)
126  : Slave(slave), MasterFace(master), PointMatrix(pm) {}
127  };
128 
131 
136  Table *bel_to_edge; // for 3D
138  mutable Table *face_edge;
139  mutable Table *edge_vertex;
140 
144 
145  // refinement embeddings for forward compatibility with NCMesh
147 
148  // Nodes are only active for higher order meshes, and share locations with
149  // the vertices, plus all the higher- order control points within the
150  // element and along the edges and on the faces.
153 
154  static const int vtk_quadratic_tet[10];
155  static const int vtk_quadratic_hex[27];
156 
157 #ifdef MFEM_USE_MEMALLOC
158  friend class Tetrahedron;
160 #endif
161 
162 public:
168 
170 
171  /// A list of all unique element attributes used by the Mesh.
173  /// A list of all unique boundary attributes used by the Mesh.
175 
176  NURBSExtension *NURBSext; ///< Optional NURBS mesh extension.
177  NCMesh *ncmesh; ///< Optional non-conforming mesh extension.
178 
179  // Global parameter that can be used to control the removal of unused
180  // vertices performed when reading a mesh in MFEM format. The default value
181  // (true) is set in mesh_readers.cpp.
183 
184 protected:
186 
187  void Init();
188  void InitTables();
189  void SetEmpty(); // Init all data members with empty values
190  void DestroyTables();
192  void DestroyPointers(); // Delete data specifically allocated by class Mesh.
193  void Destroy(); // Delete all owned data.
194  void DeleteLazyTables();
195 
196  Element *ReadElementWithoutAttr(std::istream &);
197  static void PrintElementWithoutAttr(const Element *, std::ostream &);
198 
199  Element *ReadElement(std::istream &);
200  static void PrintElement(const Element *, std::ostream &);
201 
202  // Readers for different mesh formats, used in the Load() method.
203  // The implementations of these methods are in mesh_readers.cpp.
204  void ReadMFEMMesh(std::istream &input, bool mfem_v11, int &curved);
205  void ReadLineMesh(std::istream &input);
206  void ReadNetgen2DMesh(std::istream &input, int &curved);
207  void ReadNetgen3DMesh(std::istream &input);
208  void ReadTrueGridMesh(std::istream &input);
209  void ReadVTKMesh(std::istream &input, int &curved, int &read_gf);
210  void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf);
211  void ReadInlineMesh(std::istream &input, int generate_edges = 0);
212  void ReadGmshMesh(std::istream &input);
213  /* Note NetCDF (optional library) is used for reading cubit files */
214 #ifdef MFEM_USE_NETCDF
215  void ReadCubit(const char *filename, int &curved, int &read_gf);
216 #endif
217 
218  /// Determine the mesh generator bitmask #meshgen, see MeshGenerator().
219  void SetMeshGen();
220 
221  /// Return the length of the segment from node i to node j.
222  double GetLength(int i, int j) const;
223 
224  /** Compute the Jacobian of the transformation from the perfect
225  reference element at the center of the element. */
226  void GetElementJacobian(int i, DenseMatrix &J);
227 
228  void MarkForRefinement();
230  void GetEdgeOrdering(DSTable &v_to_v, Array<int> &order);
231  virtual void MarkTetMeshForRefinement(DSTable &v_to_v);
232 
233  void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert);
234  void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert);
235 
237  STable3D *GetElementToFaceTable(int ret_ftbl = 0);
238 
239  /** Red refinement. Element with index i is refined. The default
240  red refinement for now is Uniform. */
241  void RedRefinement(int i, const DSTable &v_to_v,
242  int *edge1, int *edge2, int *middle)
243  { UniformRefinement(i, v_to_v, edge1, edge2, middle); }
244 
245  /** Green refinement. Element with index i is refined. The default
246  refinement for now is Bisection. */
247  void GreenRefinement(int i, const DSTable &v_to_v,
248  int *edge1, int *edge2, int *middle)
249  { Bisection(i, v_to_v, edge1, edge2, middle); }
250 
251  /** Bisection. Element with index i is bisected. */
252  void Bisection(int i, const DSTable &, int *, int *, int *);
253 
254  /** Bisection. Boundary element with index i is bisected. */
255  void Bisection(int i, const DSTable &, int *);
256 
257  /** Uniform Refinement. Element with index i is refined uniformly. */
258  void UniformRefinement(int i, const DSTable &, int *, int *, int *);
259 
260  /** Averages the vertices with given indexes and saves the result in
261  vertices[result]. */
262  void AverageVertices (int * indexes, int n, int result);
263 
265  int FindCoarseElement(int i);
266 
267  /// Update the nodes of a curved mesh after refinement
268  void UpdateNodes();
269 
270  /// Refine quadrilateral mesh.
271  virtual void QuadUniformRefinement();
272 
273  /// Refine hexahedral mesh.
274  virtual void HexUniformRefinement();
275 
276  /// Refine NURBS mesh.
277  virtual void NURBSUniformRefinement();
278 
279  /// This function is not public anymore. Use GeneralRefinement instead.
280  virtual void LocalRefinement(const Array<int> &marked_el, int type = 3);
281 
282  /// This function is not public anymore. Use GeneralRefinement instead.
283  virtual void NonconformingRefinement(const Array<Refinement> &refinements,
284  int nc_limit = 0);
285 
286  /// NC version of GeneralDerefinement.
287  virtual bool NonconformingDerefinement(Array<double> &elem_error,
288  double threshold, int nc_limit = 0,
289  int op = 1);
290 
291  /// Derefine elements once a list of derefinements is known.
292  void DerefineMesh(const Array<int> &derefinements);
293 
294  /// Read NURBS patch/macro-element mesh
295  void LoadPatchTopo(std::istream &input, Array<int> &edge_to_knot);
296 
297  void UpdateNURBS();
298 
299  void PrintTopo(std::ostream &out, const Array<int> &e_to_k) const;
300 
301  /// Used in GetFaceElementTransformations (...)
304  int i);
306  int i);
307  /// Used in GetFaceElementTransformations (...)
309  int i);
310  /// Used in GetFaceElementTransformations (...)
312  int i);
313  /** Used in GetFaceElementTransformations to account for the fact that a
314  slave face occupies only a portion of its master face. */
316  const FaceInfo &fi);
317  bool IsSlaveFace(const FaceInfo &fi) const;
318 
319  /// Returns the orientation of "test" relative to "base"
320  static int GetTriOrientation (const int * base, const int * test);
321  /// Returns the orientation of "test" relative to "base"
322  static int GetQuadOrientation (const int * base, const int * test);
323 
324  static void GetElementArrayEdgeTable(const Array<Element*> &elem_array,
325  const DSTable &v_to_v,
326  Table &el_to_edge);
327 
328  /** Return vertex to vertex table. The connections stored in the table
329  are from smaller to bigger vertex index, i.e. if i<j and (i, j) is
330  in the table, then (j, i) is not stored. */
331  void GetVertexToVertexTable(DSTable &) const;
332 
333  /** Return element to edge table and the indices for the boundary edges.
334  The entries in the table are ordered according to the order of the
335  nodes in the elements. For example, if T is the element to edge table
336  T(i, 0) gives the index of edge in element i that connects vertex 0
337  to vertex 1, etc. Returns the number of the edges. */
339 
340  /// Used in GenerateFaces()
341  void AddPointFaceElement(int lf, int gf, int el);
342 
343  void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1);
344 
345  void AddTriangleFaceElement (int lf, int gf, int el,
346  int v0, int v1, int v2);
347 
348  void AddQuadFaceElement (int lf, int gf, int el,
349  int v0, int v1, int v2, int v3);
350  /** For a serial Mesh, return true if the face is interior. For a parallel
351  ParMesh return true if the face is interior or shared. In parallel, this
352  method only works if the face neighbor data is exchanged. */
353  bool FaceIsTrueInterior(int FaceNo) const
354  {
355  return FaceIsInterior(FaceNo) || (faces_info[FaceNo].Elem2Inf >= 0);
356  }
357 
358  // shift cyclically 3 integers left-to-right
359  inline static void ShiftL2R(int &, int &, int &);
360  // shift cyclically 3 integers so that the smallest is first
361  inline static void Rotate3(int &, int &, int &);
362 
363  void FreeElement(Element *E);
364 
365  void GenerateFaces();
366  void GenerateNCFaceInfo();
367 
368  /// Begin construction of a mesh
369  void InitMesh(int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem);
370 
371  void InitBaseGeom();
372 
373  // Used in the methods FinalizeXXXMesh() and FinalizeTopology()
374  void FinalizeCheck();
375 
376  void Loader(std::istream &input, int generate_edges = 0,
377  std::string parse_tag = "");
378 
379  // If NURBS mesh, write NURBS format. If NCMesh, write mfem v1.1 format.
380  // If section_delimiter is empty, write mfem v1.0 format. Otherwise, write
381  // mfem v1.2 format with the given section_delimiter at the end.
382  void Printer(std::ostream &out = mfem::out,
383  std::string section_delimiter = "") const;
384 
385  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
386  nx*ny*nz hexahedrals if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons
387  if type=TETRAHEDRON. If generate_edges = 0 (default) edges are not
388  generated, if 1 edges are generated. */
389  void Make3D(int nx, int ny, int nz, Element::Type type, int generate_edges,
390  double sx, double sy, double sz);
391 
392  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
393  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
394  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
395  if 1 edges are generated. */
396  void Make2D(int nx, int ny, Element::Type type, int generate_edges,
397  double sx, double sy);
398 
399  /// Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.
400  void Make1D(int n, double sx = 1.0);
401 
402  /// Initialize vertices/elements/boundary/tables from a nonconforming mesh.
403  void InitFromNCMesh(const NCMesh &ncmesh);
404 
405  /// Create from a nonconforming mesh.
406  explicit Mesh(const NCMesh &ncmesh);
407 
408  /// Swaps internal data with another mesh. By default, non-geometry members
409  /// like 'ncmesh' and 'NURBSExt' are only swapped when 'non_geometry' is set.
410  void Swap(Mesh& other, bool non_geometry = false);
411 
412  // used in GetElementData() and GetBdrElementData()
413  void GetElementData(const Array<Element*> &elem_array, int geom,
414  Array<int> &elem_vtx, Array<int> &attr) const;
415 
416 public:
417 
418  Mesh() { SetEmpty(); }
419 
420  /** Copy constructor. Performs a deep copy of (almost) all data, so that the
421  source mesh can be modified (e.g. deleted, refined) without affecting the
422  new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the
423  nodes, if present. */
424  explicit Mesh(const Mesh &mesh, bool copy_nodes = true);
425 
426  /// Construct a Mesh from the given primary data.
427  /** The array @a vertices is used as external data, i.e. the Mesh does not
428  copy the data and will not delete the pointer.
429 
430  The data from the other arrays is copied into the internal Mesh data
431  structures.
432 
433  This method calls the method FinalizeTopology(). The method Finalize()
434  may be called after this constructor and after optionally setting the
435  Mesh nodes. */
436  Mesh(double *vertices, int num_vertices,
437  int *element_indices, Geometry::Type element_type,
438  int *element_attributes, int num_elements,
439  int *boundary_indices, Geometry::Type boundary_type,
440  int *boundary_attributes, int num_boundary_elements,
441  int dimension, int space_dimension= -1);
442 
443  /** @anchor mfem_Mesh_init_ctor
444  @brief _Init_ constructor: begin the construction of a Mesh object. */
445  Mesh(int _Dim, int NVert, int NElem, int NBdrElem = 0, int _spaceDim = -1)
446  {
447  if (_spaceDim == -1)
448  {
449  _spaceDim = _Dim;
450  }
451  InitMesh(_Dim, _spaceDim, NVert, NElem, NBdrElem);
452  }
453 
454  /** @name Methods for Mesh construction.
455 
456  These methods are intended to be used with the @ref mfem_Mesh_init_ctor
457  "init constructor". */
458  ///@{
459 
460  Element *NewElement(int geom);
461 
462  void AddVertex(const double *);
463  void AddTri(const int *vi, int attr = 1);
464  void AddTriangle(const int *vi, int attr = 1);
465  void AddQuad(const int *vi, int attr = 1);
466  void AddTet(const int *vi, int attr = 1);
467  void AddHex(const int *vi, int attr = 1);
468  void AddHexAsTets(const int *vi, int attr = 1);
469  // 'elem' should be allocated using the NewElement method
470  void AddElement(Element *elem) { elements[NumOfElements++] = elem; }
471  void AddBdrElement(Element *elem) { boundary[NumOfBdrElements++] = elem; }
472  void AddBdrSegment(const int *vi, int attr = 1);
473  void AddBdrTriangle(const int *vi, int attr = 1);
474  void AddBdrQuad(const int *vi, int attr = 1);
475  void AddBdrQuadAsTriangles(const int *vi, int attr = 1);
476 
478  /// Finalize the construction of a triangular Mesh.
479  void FinalizeTriMesh(int generate_edges = 0, int refine = 0,
480  bool fix_orientation = true);
481  /// Finalize the construction of a quadrilateral Mesh.
482  void FinalizeQuadMesh(int generate_edges = 0, int refine = 0,
483  bool fix_orientation = true);
484  /// Finalize the construction of a tetrahedral Mesh.
485  void FinalizeTetMesh(int generate_edges = 0, int refine = 0,
486  bool fix_orientation = true);
487  /// Finalize the construction of a hexahedral Mesh.
488  void FinalizeHexMesh(int generate_edges = 0, int refine = 0,
489  bool fix_orientation = true);
490 
491  ///@}
492 
493  /** @brief Finalize the construction of the secondary topology (connectivity)
494  data of a Mesh. */
495  /** This method does not require any actual coordinate data (either vertex
496  coordinates for linear meshes or node coordinates for meshes with nodes)
497  to be available. However, the data generated by this method is generally
498  required by the FiniteElementSpace class.
499 
500  After calling this method, setting the Mesh vertices or nodes, it may be
501  appropriate to call the method Finalize(). */
502  void FinalizeTopology();
503 
504  /// Finalize the construction of a general Mesh.
505  /** This method will:
506  - check and optionally fix the orientation of regular elements
507  - check and fix the orientation of boundary elements
508  - assume that #vertices are defined, if #Nodes == NULL
509  - assume that #Nodes are defined, if #Nodes != NULL.
510  @param[in] refine If true, prepare the Mesh for conforming refinement of
511  triangular or tetrahedral meshes.
512  @param[in] fix_orientation
513  If true, fix the orientation of inverted mesh elements
514  by permuting their vertices.
515 
516  Before calling this method, call FinalizeTopology() and ensure that the
517  Mesh vertices or nodes are set. */
518  void Finalize(bool refine = false, bool fix_orientation = false);
519 
520  void SetAttributes();
521 
522 #ifdef MFEM_USE_GECKO
523  /** This is our integration with the Gecko library. This will call the
524  Gecko library to find an element ordering that will increase memory
525  coherency by putting elements that are in physical proximity closer in
526  memory. */
527  void GetGeckoElementReordering(Array<int> &ordering);
528 #endif
529 
530  /** Rebuilds the mesh with a different order of elements. The ordering
531  vector maps the old element number to the new element number. This also
532  reorders the vertices and nodes edges and faces along with the elements. */
533  void ReorderElements(const Array<int> &ordering, bool reorder_vertices = true);
534 
535  /** Creates mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into
536  nx*ny*nz hexahedrals if type=HEXAHEDRON or into 6*nx*ny*nz tetrahedrons
537  if type=TETRAHEDRON. If generate_edges = 0 (default) edges are not
538  generated, if 1 edges are generated. */
539  Mesh(int nx, int ny, int nz, Element::Type type, int generate_edges = 0,
540  double sx = 1.0, double sy = 1.0, double sz = 1.0)
541  {
542  Make3D(nx, ny, nz, type, generate_edges, sx, sy, sz);
543  }
544 
545  /** Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny
546  quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if
547  type = TRIANGLE. If generate_edges = 0 (default) edges are not generated,
548  if 1 edges are generated. */
549  Mesh(int nx, int ny, Element::Type type, int generate_edges = 0,
550  double sx = 1.0, double sy = 1.0)
551  {
552  Make2D(nx, ny, type, generate_edges, sx, sy);
553  }
554 
555  /** Creates 1D mesh , divided into n equal intervals. */
556  explicit Mesh(int n, double sx = 1.0)
557  {
558  Make1D(n, sx);
559  }
560 
561  /** Creates mesh by reading a file in MFEM, netgen, or VTK format. If
562  generate_edges = 0 (default) edges are not generated, if 1 edges are
563  generated. */
564  explicit Mesh(const char *filename, int generate_edges = 0, int refine = 1,
565  bool fix_orientation = true);
566 
567  /** Creates mesh by reading data stream in MFEM, netgen, or VTK format. If
568  generate_edges = 0 (default) edges are not generated, if 1 edges are
569  generated. */
570  explicit Mesh(std::istream &input, int generate_edges = 0, int refine = 1,
571  bool fix_orientation = true);
572 
573  /// Create a disjoint mesh from the given mesh array
574  Mesh(Mesh *mesh_array[], int num_pieces);
575 
576  /// Create a uniformly refined (by any factor) version of @a orig_mesh.
577  /** @param[in] orig_mesh The starting coarse mesh.
578  @param[in] ref_factor The refinement factor, an integer > 1.
579  @param[in] ref_type Specify the positions of the new vertices. The
580  options are BasisType::ClosedUniform or
581  BasisType::GaussLobatto.
582 
583  The refinement data which can be accessed with GetRefinementTransforms()
584  is set to reflect the performed refinements.
585 
586  @note The constructed Mesh is linear, i.e. it does not have nodes. */
587  Mesh(Mesh *orig_mesh, int ref_factor, int ref_type);
588 
589  /** This is similar to the mesh constructor with the same arguments, but here
590  the current mesh is destroyed and another one created based on the data
591  stream again given in MFEM, netgen, or VTK format. If generate_edges = 0
592  (default) edges are not generated, if 1 edges are generated. */
593  /// \see mfem::igzstream() for on-the-fly decompression of compressed ascii
594  /// inputs.
595  virtual void Load(std::istream &input, int generate_edges = 0,
596  int refine = 1, bool fix_orientation = true)
597  {
598  Loader(input, generate_edges);
599  Finalize(refine, fix_orientation);
600  }
601 
602  /// Clear the contents of the Mesh.
603  void Clear() { Destroy(); SetEmpty(); }
604 
605  /** @brief Get the mesh generator/type.
606 
607  @return A bitmask:
608  - bit 0 - simplices are present in the mesh (triangles, tets),
609  - bit 1 - tensor product elements are present in the mesh (quads, hexes).
610  */
611  inline int MeshGenerator() { return meshgen; }
612 
613  /** @brief Returns number of vertices. Vertices are only at the corners of
614  elements, where you would expect them in the lowest-order mesh. */
615  inline int GetNV() const { return NumOfVertices; }
616 
617  /// Returns number of elements.
618  inline int GetNE() const { return NumOfElements; }
619 
620  /// Returns number of boundary elements.
621  inline int GetNBE() const { return NumOfBdrElements; }
622 
623  /// Return the number of edges.
624  inline int GetNEdges() const { return NumOfEdges; }
625 
626  /// Return the number of faces in a 3D mesh.
627  inline int GetNFaces() const { return NumOfFaces; }
628 
629  /// Return the number of faces (3D), edges (2D) or vertices (1D).
630  int GetNumFaces() const;
631 
632  /// Utility function: sum integers from all processors (Allreduce).
633  virtual long ReduceInt(int value) const { return value; }
634 
635  /// Return the total (global) number of elements.
636  long GetGlobalNE() const { return ReduceInt(NumOfElements); }
637 
638  /// Equals 1 + num_holes - num_loops
639  inline int EulerNumber() const
641  /// Equals 1 - num_holes
642  inline int EulerNumber2D() const
643  { return NumOfVertices - NumOfEdges + NumOfElements; }
644 
645  int Dimension() const { return Dim; }
646  int SpaceDimension() const { return spaceDim; }
647 
648  /// @brief Return pointer to vertex i's coordinates.
649  /// @warning For high-order meshes (when #Nodes != NULL) vertices may not be
650  /// updated and should not be used!
651  const double *GetVertex(int i) const { return vertices[i](); }
652 
653  /// @brief Return pointer to vertex i's coordinates.
654  /// @warning For high-order meshes (when Nodes != NULL) vertices may not
655  /// being updated and should not be used!
656  double *GetVertex(int i) { return vertices[i](); }
657 
658  void GetElementData(int geom, Array<int> &elem_vtx, Array<int> &attr) const
659  { GetElementData(elements, geom, elem_vtx, attr); }
660 
661  void GetBdrElementData(int geom, Array<int> &bdr_elem_vtx,
662  Array<int> &bdr_attr) const
663  { GetElementData(boundary, geom, bdr_elem_vtx, bdr_attr); }
664 
665  /** @brief Set the internal Vertex array to point to the given @a vertices
666  array without assuming ownership of the pointer. */
667  /** If @a zerocopy is `true`, the vertices must be given as an array of 3
668  doubles per vertex. If @a zerocopy is `false` then the current Vertex
669  data is first copied to the @a vertices array. */
670  void ChangeVertexDataOwnership(double *vertices, int len_vertices,
671  bool zerocopy = false);
672 
673  const Element* const *GetElementsArray() const
674  { return elements.GetData(); }
675 
676  const Element *GetElement(int i) const { return elements[i]; }
677 
678  Element *GetElement(int i) { return elements[i]; }
679 
680  const Element *GetBdrElement(int i) const { return boundary[i]; }
681 
682  Element *GetBdrElement(int i) { return boundary[i]; }
683 
684  const Element *GetFace(int i) const { return faces[i]; }
685 
686  int GetFaceBaseGeometry(int i) const;
687 
688  int GetElementBaseGeometry(int i = 0) const
689  { return i < GetNE() ? elements[i]->GetGeometryType() : BaseGeom; }
690 
691  int GetBdrElementBaseGeometry(int i = 0) const
692  { return i < GetNBE() ? boundary[i]->GetGeometryType() : BaseBdrGeom; }
693 
694  /// Returns the indices of the vertices of element i.
695  void GetElementVertices(int i, Array<int> &v) const
696  { elements[i]->GetVertices(v); }
697 
698  /// Returns the indices of the vertices of boundary element i.
699  void GetBdrElementVertices(int i, Array<int> &v) const
700  { boundary[i]->GetVertices(v); }
701 
702  /// Return the indices and the orientations of all edges of element i.
703  void GetElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
704 
705  /// Return the indices and the orientations of all edges of bdr element i.
706  void GetBdrElementEdges(int i, Array<int> &edges, Array<int> &cor) const;
707 
708  /** Return the indices and the orientations of all edges of face i.
709  Works for both 2D (face=edge) and 3D faces. */
710  void GetFaceEdges(int i, Array<int> &, Array<int> &) const;
711 
712  /// Returns the indices of the vertices of face i.
713  void GetFaceVertices(int i, Array<int> &vert) const
714  {
715  if (Dim == 1)
716  {
717  vert.SetSize(1); vert[0] = i;
718  }
719  else
720  {
721  faces[i]->GetVertices(vert);
722  }
723  }
724 
725  /// Returns the indices of the vertices of edge i.
726  void GetEdgeVertices(int i, Array<int> &vert) const;
727 
728  /// Returns the face-to-edge Table (3D)
729  Table *GetFaceEdgeTable() const;
730 
731  /// Returns the edge-to-vertex Table (3D)
732  Table *GetEdgeVertexTable() const;
733 
734  /// Return the indices and the orientations of all faces of element i.
735  void GetElementFaces(int i, Array<int> &, Array<int> &) const;
736 
737  /// Return the index and the orientation of the face of bdr element i. (3D)
738  void GetBdrElementFace(int i, int *, int *) const;
739 
740  /** Return the vertex index of boundary element i. (1D)
741  Return the edge index of boundary element i. (2D)
742  Return the face index of boundary element i. (3D) */
743  int GetBdrElementEdgeIndex(int i) const;
744 
745  /** @brief For the given boundary element, bdr_el, return its adjacent
746  element and its info, i.e. 64*local_bdr_index+bdr_orientation. */
747  void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const;
748 
749  /// Returns the type of element i.
750  int GetElementType(int i) const;
751 
752  /// Returns the type of boundary element i.
753  int GetBdrElementType(int i) const;
754 
755  /* Return point matrix of element i of dimension Dim X #v, where for every
756  vertex we give its coordinates in space of dimension Dim. */
757  void GetPointMatrix(int i, DenseMatrix &pointmat) const;
758 
759  /* Return point matrix of boundary element i of dimension Dim X #v, where for
760  every vertex we give its coordinates in space of dimension Dim. */
761  void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const;
762 
764 
765  /** Builds the transformation defining the i-th element in the user-defined
766  variable. */
768 
769  /// Returns the transformation defining the i-th element
771 
772  /** Return the transformation defining the i-th element assuming
773  the position of the vertices/nodes are given by 'nodes'. */
774  void GetElementTransformation(int i, const Vector &nodes,
776 
777  /// Returns the transformation defining the i-th boundary element
780 
781  /** @brief Returns the transformation defining the given face element in a
782  user-defined variable. */
784 
785  /** @brief A helper method that constructs a transformation from the
786  reference space of a face to the reference space of an element. */
787  /** The local index of the face as a face in the element and its orientation
788  are given by the input parameter @a info, as @a info = 64*loc_face_idx +
789  loc_face_orientation. */
790  void GetLocalFaceTransformation(int face_type, int elem_type,
792  int info);
793 
794  /// Returns the transformation defining the given face element
796 
797  /** Returns the transformation defining the given edge element.
798  The transformation is stored in a user-defined variable. */
800 
801  /// Returns the transformation defining the given face element
803 
804  /// Returns (a pointer to a structure containing) the following data:
805  ///
806  /// 1) Elem1No - the index of the first element that contains this face this
807  /// is the element that has the same outward unit normal vector as the
808  /// face;
809  ///
810  /// 2) Elem2No - the index of the second element that contains this face this
811  /// element has outward unit normal vector as the face multiplied with -1;
812  ///
813  /// 3) Elem1, Elem2 - pointers to the ElementTransformation's of the first
814  /// and the second element respectively;
815  ///
816  /// 4) Face - pointer to the ElementTransformation of the face;
817  ///
818  /// 5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face
819  /// coordinate system to the element coordinate system (both in their
820  /// reference elements). Used to transform IntegrationPoints from face to
821  /// element. More formally, let:
822  /// TL1, TL2 be the transformations represented by Loc1, Loc2,
823  /// TE1, TE2 - the transformations represented by Elem1, Elem2,
824  /// TF - the transformation represented by Face, then
825  /// TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
826  ///
827  /// 6) FaceGeom - the base geometry for the face.
828  ///
829  /// The mask specifies which fields in the structure to return:
830  /// mask & 1 - Elem1, mask & 2 - Elem2
831  /// mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face.
833  int mask = 31);
834 
836  {
837  if (faces_info[FaceNo].Elem2No < 0) { return NULL; }
838  return GetFaceElementTransformations (FaceNo);
839  }
840 
842 
843  /// Return true if the given face is interior. @sa FaceIsTrueInterior().
844  bool FaceIsInterior(int FaceNo) const
845  {
846  return (faces_info[FaceNo].Elem2No >= 0);
847  }
848  void GetFaceElements (int Face, int *Elem1, int *Elem2);
849  void GetFaceInfos (int Face, int *Inf1, int *Inf2);
850 
851  int GetFaceGeometryType(int Face) const;
852  int GetFaceElementType(int Face) const;
853 
854  /// Check the orientation of the elements
855  /** @return The number of elements with wrong orientation. */
856  int CheckElementOrientation(bool fix_it = true);
857  /// Check the orientation of the boundary elements
858  /** @return The number of boundary elements with wrong orientation. */
859  int CheckBdrElementOrientation(bool fix_it = true);
860 
861  /// Return the attribute of element i.
862  int GetAttribute(int i) const { return elements[i]->GetAttribute(); }
863 
864  /// Set the attribute of element i.
865  void SetAttribute(int i, int attr) { elements[i]->SetAttribute(attr); }
866 
867  /// Return the attribute of boundary element i.
868  int GetBdrAttribute(int i) const { return boundary[i]->GetAttribute(); }
869 
870  const Table &ElementToElementTable();
871 
872  const Table &ElementToFaceTable() const;
873 
874  const Table &ElementToEdgeTable() const;
875 
876  /// The returned Table must be destroyed by the caller
878 
879  /** Return the "face"-element Table. Here "face" refers to face (3D),
880  edge (2D), or vertex (1D).
881  The returned Table must be destroyed by the caller. */
882  Table *GetFaceToElementTable() const;
883 
884  /** This method modifies a tetrahedral mesh so that Nedelec spaces of order
885  greater than 1 can be defined on the mesh. Specifically, we
886  1) rotate all tets in the mesh so that the vertices {v0, v1, v2, v3}
887  satisfy: v0 < v1 < min(v2, v3).
888  2) rotate all boundary triangles so that the vertices {v0, v1, v2}
889  satisfy: v0 < min(v1, v2).
890 
891  @note Refinement does not work after a call to this method! */
892  virtual void ReorientTetMesh();
893 
894  int *CartesianPartitioning(int nxyz[]);
895  int *GeneratePartitioning(int nparts, int part_method = 1);
896  void CheckPartitioning(int *partitioning);
897 
898  void CheckDisplacements(const Vector &displacements, double &tmax);
899 
900  // Vertices are only at the corners of elements, where you would expect them
901  // in the lowest-order mesh.
902  void MoveVertices(const Vector &displacements);
903  void GetVertices(Vector &vert_coord) const;
904  void SetVertices(const Vector &vert_coord);
905 
906  // Nodes are only active for higher order meshes, and share locations with
907  // the vertices, plus all the higher- order control points within the element
908  // and along the edges and on the faces.
909  void GetNode(int i, double *coord);
910  void SetNode(int i, const double *coord);
911 
912  // Node operations for curved mesh.
913  // They call the corresponding '...Vertices' method if the
914  // mesh is not curved (i.e. Nodes == NULL).
915  void MoveNodes(const Vector &displacements);
916  void GetNodes(Vector &node_coord) const;
917  void SetNodes(const Vector &node_coord);
918 
919  /// Return a pointer to the internal node GridFunction (may be NULL).
920  GridFunction *GetNodes() { return Nodes; }
921  const GridFunction *GetNodes() const { return Nodes; }
922  /// Return the mesh nodes ownership flag.
923  bool OwnsNodes() const { return own_nodes; }
924  /// Set the mesh nodes ownership flag.
925  void SetNodesOwner(bool nodes_owner) { own_nodes = nodes_owner; }
926  /// Replace the internal node GridFunction with the given GridFunction.
927  void NewNodes(GridFunction &nodes, bool make_owner = false);
928  /** Swap the internal node GridFunction pointer and ownership flag members
929  with the given ones. */
930  void SwapNodes(GridFunction *&nodes, int &own_nodes_);
931 
932  /// Return the mesh nodes/vertices projected on the given GridFunction.
933  void GetNodes(GridFunction &nodes) const;
934  /** Replace the internal node GridFunction with a new GridFunction defined
935  on the given FiniteElementSpace. The new node coordinates are projected
936  (derived) from the current nodes/vertices. */
938  /** Replace the internal node GridFunction with the given GridFunction. The
939  given GridFunction is updated with node coordinates projected (derived)
940  from the current nodes/vertices. */
941  void SetNodalGridFunction(GridFunction *nodes, bool make_owner = false);
942  /** Return the FiniteElementSpace on which the current mesh nodes are
943  defined or NULL if the mesh does not have nodes. */
944  const FiniteElementSpace *GetNodalFESpace() const;
945 
946  /** Set the curvature of the mesh nodes using the given polynomial degree,
947  'order', and optionally: discontinuous or continuous FE space, 'discont',
948  new space dimension, 'space_dim' (if != -1), and 'ordering'. */
949  void SetCurvature(int order, bool discont = false, int space_dim = -1,
950  int ordering = 1);
951 
952  /** Refine all mesh elements. */
953  void UniformRefinement();
954 
955  /** Refine selected mesh elements. Refinement type can be specified for each
956  element. The function can do conforming refinement of triangles and
957  tetrahedra and non-conforming refinement (i.e., with hanging-nodes) of
958  triangles, quadrilaterals and hexahedrons. If 'nonconforming' = -1,
959  suitable refinement method is selected automatically (namely, conforming
960  refinement for triangles). Use nonconforming = 0/1 to force the method.
961  For nonconforming refinements, nc_limit optionally specifies the maximum
962  level of hanging nodes (unlimited by default). */
963  void GeneralRefinement(const Array<Refinement> &refinements,
964  int nonconforming = -1, int nc_limit = 0);
965 
966  /** Simplified version of GeneralRefinement taking a simple list of elements
967  to refine, without refinement types. */
968  void GeneralRefinement(const Array<int> &el_to_refine,
969  int nonconforming = -1, int nc_limit = 0);
970 
971  /// Refine each element with given probability. Uses GeneralRefinement.
972  void RandomRefinement(double prob, bool aniso = false,
973  int nonconforming = -1, int nc_limit = 0);
974 
975  /// Refine elements sharing the specified vertex. Uses GeneralRefinement.
976  void RefineAtVertex(const Vertex& vert,
977  double eps = 0.0, int nonconforming = -1);
978 
979  /** Refine element i if elem_error[i] > threshold, for all i.
980  Returns true if at least one element was refined, false otherwise. */
981  bool RefineByError(const Array<double> &elem_error, double threshold,
982  int nonconforming = -1, int nc_limit = 0);
983 
984  /** Refine element i if elem_error(i) > threshold, for all i.
985  Returns true if at least one element was refined, false otherwise. */
986  bool RefineByError(const Vector &elem_error, double threshold,
987  int nonconforming = -1, int nc_limit = 0);
988 
989  /** Derefine the mesh based on an error measure associated with each
990  element. A derefinement is performed if the sum of errors of its fine
991  elements is smaller than 'threshold'. If 'nc_limit' > 0, derefinements
992  that would increase the maximum level of hanging nodes of the mesh are
993  skipped. Returns true if the mesh changed, false otherwise. */
994  bool DerefineByError(Array<double> &elem_error, double threshold,
995  int nc_limit = 0, int op = 1);
996 
997  /// Same as DerefineByError for an error vector.
998  bool DerefineByError(const Vector &elem_error, double threshold,
999  int nc_limit = 0, int op = 1);
1000 
1001  ///@{ @name NURBS mesh refinement methods
1002  void KnotInsert(Array<KnotVector *> &kv);
1003  /* For each knot vector:
1004  new_degree = max(old_degree, min(old_degree + rel_degree, degree)). */
1005  void DegreeElevate(int rel_degree, int degree = 16);
1006  ///@}
1007 
1008  /** Make sure that a quad/hex mesh is considered to be non-conforming (i.e.,
1009  has an associated NCMesh object). Triangles meshes can be both conforming
1010  (default) or non-conforming. */
1011  void EnsureNCMesh(bool triangles_nonconforming = false);
1012 
1013  bool Conforming() const { return ncmesh == NULL; }
1014  bool Nonconforming() const { return ncmesh != NULL; }
1015 
1016  /** Return fine element transformations following a mesh refinement.
1017  Space uses this to construct a global interpolation matrix. */
1019 
1020  /// Return type of last modification of the mesh.
1022 
1023  /** Return update counter. The counter starts at zero and is incremented
1024  each time refinement, derefinement, or rebalancing method is called.
1025  It is used for checking proper sequence of Space:: and GridFunction::
1026  Update() calls. */
1027  long GetSequence() const { return sequence; }
1028 
1029  /// Print the mesh to the given stream using Netgen/Truegrid format.
1030  virtual void PrintXG(std::ostream &out = mfem::out) const;
1031 
1032  /// Print the mesh to the given stream using the default MFEM mesh format.
1033  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1034  virtual void Print(std::ostream &out = mfem::out) const { Printer(out); }
1035 
1036  /// Print the mesh in VTK format (linear and quadratic meshes only).
1037  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1038  void PrintVTK(std::ostream &out);
1039 
1040  /** Print the mesh in VTK format. The parameter ref > 0 specifies an element
1041  subdivision number (useful for high order fields and curved meshes).
1042  If the optional field_data is set, we also add a FIELD section in the
1043  beginning of the file with additional dataset information. */
1044  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1045  void PrintVTK(std::ostream &out, int ref, int field_data=0);
1046 
1047  void GetElementColoring(Array<int> &colors, int el0 = 0);
1048 
1049  /** Prints the mesh with bdr elements given by the boundary of
1050  the subdomains, so that the boundary of subdomain i has bdr
1051  attribute i+1. */
1052  /// \see mfem::ogzstream() for on-the-fly compression of ascii outputs
1053  void PrintWithPartitioning (int *partitioning,
1054  std::ostream &out, int elem_attr = 0) const;
1055 
1056  void PrintElementsWithPartitioning (int *partitioning,
1057  std::ostream &out,
1058  int interior_faces = 0);
1059 
1060  /// Print set of disjoint surfaces:
1061  /*!
1062  * If Aface_face(i,j) != 0, print face j as a boundary
1063  * element with attribute i+1.
1064  */
1065  void PrintSurfaces(const Table &Aface_face, std::ostream &out) const;
1066 
1067  void ScaleSubdomains (double sf);
1068  void ScaleElements (double sf);
1069 
1070  void Transform(void (*f)(const Vector&, Vector&));
1071  void Transform(VectorCoefficient &deformation);
1072 
1073  /// Remove unused vertices and rebuild mesh connectivity.
1074  void RemoveUnusedVertices();
1075 
1076  /** Remove boundary elements that lie in the interior of the mesh, i.e. that
1077  have two adjacent faces in 3D, or edges in 2D. */
1078  void RemoveInternalBoundaries();
1079 
1080  /** Get the size of the i-th element relative to the perfect
1081  reference element. */
1082  double GetElementSize(int i, int type = 0);
1083 
1084  double GetElementSize(int i, const Vector &dir);
1085 
1086  double GetElementVolume(int i);
1087 
1088  /// Returns the minimum and maximum corners of the mesh bounding box. For
1089  /// high-order meshes, the geometry is refined first "ref" times.
1090  void GetBoundingBox(Vector &min, Vector &max, int ref = 2);
1091 
1092  void GetCharacteristics(double &h_min, double &h_max,
1093  double &kappa_min, double &kappa_max,
1094  Vector *Vh = NULL, Vector *Vk = NULL);
1095 
1096  void PrintCharacteristics(Vector *Vh = NULL, Vector *Vk = NULL,
1097  std::ostream &out = mfem::out);
1098 
1099  virtual void PrintInfo(std::ostream &out = mfem::out)
1100  {
1101  PrintCharacteristics(NULL, NULL, out);
1102  }
1103 
1104  void MesquiteSmooth(const int mesquite_option = 0);
1105 
1106  /** @brief Find the ids of the elements that contain the given points, and
1107  their corresponding reference coordinates.
1108 
1109  The DenseMatrix @a point_mat describes the given points - one point for
1110  each column; it should have SpaceDimension() rows.
1111 
1112  The InverseElementTransformation object, @a inv_trans, is used to attempt
1113  the element transformation inversion. If NULL pointer is given, the
1114  method will use a default constructed InverseElementTransformation. Note
1115  that the algorithms in the base class InverseElementTransformation can be
1116  completely overwritten by deriving custom classes that override the
1117  Transform() method.
1118 
1119  If no element is found for the i-th point, elem_ids[i] is set to -1.
1120 
1121  In the ParMesh implementation, the @a point_mat is expected to be the
1122  same on all ranks. If the i-th point is found by multiple ranks, only one
1123  of them will mark that point as found, i.e. set its elem_ids[i] to a
1124  non-negative number; the other ranks will set their elem_ids[i] to -2 to
1125  indicate that the point was found but assigned to another rank.
1126 
1127  @returns The total number of points that were found.
1128 
1129  @note This method is not 100 percent reliable, i.e. it is not guaranteed
1130  to find a point, even if it lies inside a mesh element. */
1131  virtual int FindPoints(DenseMatrix& point_mat, Array<int>& elem_ids,
1132  Array<IntegrationPoint>& ips, bool warn = true,
1133  InverseElementTransformation *inv_trans = NULL);
1134 
1135  /// Destroys Mesh.
1136  virtual ~Mesh() { DestroyPointers(); }
1137 };
1138 
1139 /** Overload operator<< for std::ostream and Mesh; valid also for the derived
1140  class ParMesh */
1141 std::ostream &operator<<(std::ostream &out, const Mesh &mesh);
1142 
1143 
1144 /// Class used to extrude the nodes of a mesh
1146 {
1147 private:
1148  int n, layer;
1149  double p[2], s;
1150  Vector tip;
1151 public:
1152  NodeExtrudeCoefficient(const int dim, const int _n, const double _s);
1153  void SetLayer(const int l) { layer = l; }
1155  virtual void Eval(Vector &V, ElementTransformation &T,
1156  const IntegrationPoint &ip);
1158 };
1159 
1160 
1161 /// Extrude a 1D mesh
1162 Mesh *Extrude1D(Mesh *mesh, const int ny, const double sy,
1163  const bool closed = false);
1164 
1165 
1166 /// Input file stream that remembers the input file name (useful for example
1167 /// when reading NetCDF meshes) and supports optional gzstream decompression.
1169 {
1170 public:
1171  const char *filename;
1172  named_ifgzstream(const char *mesh_name) :
1173  mfem::ifgzstream(mesh_name), filename(mesh_name) {}
1174 };
1175 
1176 
1177 // inline functions
1178 inline void Mesh::ShiftL2R(int &a, int &b, int &c)
1179 {
1180  int t = a;
1181  a = c; c = b; b = t;
1182 }
1183 
1184 inline void Mesh::Rotate3(int &a, int &b, int &c)
1185 {
1186  if (a < b)
1187  {
1188  if (a > c)
1189  {
1190  ShiftL2R(a, b, c);
1191  }
1192  }
1193  else
1194  {
1195  if (b < c)
1196  {
1197  ShiftL2R(c, b, a);
1198  }
1199  else
1200  {
1201  ShiftL2R(a, b, c);
1202  }
1203  }
1204 }
1205 
1206 }
1207 
1208 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:140
void AddHex(const int *vi, int attr=1)
Definition: mesh.cpp:1010
void Loader(std::istream &input, int generate_edges=0, std::string parse_tag="")
Definition: mesh.cpp:2618
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:8569
void GetFaceEdges(int i, Array< int > &, Array< int > &) const
Definition: mesh.cpp:3752
void PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
Definition: mesh.cpp:8100
void GetPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:4039
static const int vtk_quadratic_hex[27]
Definition: mesh.hpp:155
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1034
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:4601
const DenseMatrix * PointMatrix
Definition: mesh.hpp:122
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:868
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:651
void ScaleElements(double sf)
Definition: mesh.cpp:8236
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:3819
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:3720
void FreeElement(Element *E)
Definition: mesh.cpp:8550
void DerefineMesh(const Array< int > &derefinements)
Derefine elements once a list of derefinements is known.
Definition: mesh.cpp:6075
void GetFaceInfos(int Face, int *Inf1, int *Inf2)
Definition: mesh.cpp:773
int CheckElementOrientation(bool fix_it=true)
Check the orientation of the elements.
Definition: mesh.cpp:3389
static bool remove_unused_vertices
Definition: mesh.hpp:182
void SetVertices(const Vector &vert_coord)
Definition: mesh.cpp:5362
static const int vtk_quadratic_tet[10]
Definition: mesh.hpp:154
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:3782
bool FaceIsTrueInterior(int FaceNo) const
Definition: mesh.hpp:353
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Definition: mesh.cpp:8674
long GetGlobalNE() const
Return the total (global) number of elements.
Definition: mesh.hpp:636
void AddHexAsTets(const int *vi, int attr=1)
Definition: mesh.cpp:1015
void GetFaceElements(int Face, int *Elem1, int *Elem2)
Definition: mesh.cpp:767
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Array< Element * > boundary
Definition: mesh.hpp:74
Geometry::Constants< Geometry::SQUARE > quad_t
Definition: mesh.hpp:165
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:4640
CoarseFineTransformations CoarseFineTr
Definition: mesh.hpp:146
int own_nodes
Definition: mesh.hpp:152
bool Conforming() const
Definition: mesh.hpp:1013
void MoveVertices(const Vector &displacements)
Definition: mesh.cpp:5342
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:621
void DeleteLazyTables()
Definition: mesh.cpp:893
IsoparametricTransformation Transformation
Definition: mesh.hpp:141
void GetLocalPtToSegTransformation(IsoparametricTransformation &, int)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:523
void GetCharacteristics(double &h_min, double &h_max, double &kappa_min, double &kappa_max, Vector *Vh=NULL, Vector *Vk=NULL)
Definition: mesh.cpp:172
void PrintWithPartitioning(int *partitioning, std::ostream &out, int elem_attr=0) const
Definition: mesh.cpp:7610
int GetFaceGeometryType(int Face) const
Definition: mesh.cpp:779
const Geometry::Type geom
Definition: ex1.cpp:40
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:5446
int NumOfEdges
Definition: mesh.hpp:57
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:5460
void GetElementFaces(int i, Array< int > &, Array< int > &) const
Return the indices and the orientations of all faces of element i.
Definition: mesh.cpp:3914
void AddBdrElement(Element *elem)
Definition: mesh.hpp:471
int BaseGeom
Definition: mesh.hpp:59
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Definition: mesh.cpp:108
void ReadNetgen2DMesh(std::istream &input, int &curved)
int EulerNumber2D() const
Equals 1 - num_holes.
Definition: mesh.hpp:642
void AddTri(const int *vi, int attr=1)
Definition: mesh.cpp:982
void DeleteTables()
Definition: mesh.hpp:191
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:713
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
GridFunction * Nodes
Definition: mesh.hpp:151
int NumOfElements
Definition: mesh.hpp:56
Mesh * Extrude1D(Mesh *mesh, const int ny, const double sy, const bool closed)
Extrude a 1D mesh.
Definition: mesh.cpp:8692
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:8306
void FinalizeCheck()
Definition: mesh.cpp:1102
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:3997
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
const Element * GetFace(int i) const
Definition: mesh.hpp:684
Element * GetElement(int i)
Definition: mesh.hpp:678
void GenerateBoundaryElements()
Definition: mesh.cpp:1064
Data type for vertex.
Definition: vertex.hpp:21
void SetMeshGen()
Determine the mesh generator bitmask meshgen, see MeshGenerator().
Definition: mesh.cpp:2598
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3548
void GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:5351
void ReadNetgen3DMesh(std::istream &input)
int BaseBdrGeom
Definition: mesh.hpp:59
void AddBdrSegment(const int *vi, int attr=1)
Definition: mesh.cpp:1034
The inverse transformation of a given ElementTransformation.
Definition: eltrans.hpp:108
int GetFaceElementType(int Face) const
Definition: mesh.cpp:784
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
Definition: mesh.hpp:835
void RemoveInternalBoundaries()
Definition: mesh.cpp:8462
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:657
void SetEmpty()
Definition: mesh.cpp:812
Array< Element * > faces
Definition: mesh.hpp:75
virtual long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
Definition: mesh.hpp:633
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
Definition: mesh.cpp:1757
void GetVertexToVertexTable(DSTable &) const
Definition: mesh.cpp:4108
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:241
Operation last_operation
Definition: mesh.hpp:185
void MarkTriMeshForRefinement()
Definition: mesh.cpp:1389
void ReadCubit(const char *filename, int &curved, int &read_gf)
void InitRefinementTransforms()
Definition: mesh.cpp:6831
IsoparametricTransformation FaceTransformation
Definition: mesh.hpp:142
Array< NCFaceInfo > nc_faces_info
Definition: mesh.hpp:130
Element * ReadElement(std::istream &)
Definition: mesh.cpp:2580
void KnotInsert(Array< KnotVector * > &kv)
Definition: mesh.cpp:3085
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:61
Element * NewElement(int geom)
Definition: mesh.cpp:2530
FaceElementTransformations FaceElemTr
Definition: mesh.hpp:143
Table * el_to_face
Definition: mesh.hpp:133
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
Definition: mesh.cpp:460
const GridFunction * GetNodes() const
Definition: mesh.hpp:921
void MarkForRefinement()
Definition: mesh.cpp:1372
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:356
void AddVertex(const double *)
Definition: mesh.cpp:971
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:445
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
Definition: mesh.cpp:7157
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:3374
Geometry::Constants< Geometry::SEGMENT > seg_t
Definition: mesh.hpp:163
void AddElement(Element *elem)
Definition: mesh.hpp:470
void GetElementJacobian(int i, DenseMatrix &J)
Definition: mesh.cpp:56
void MesquiteSmooth(const int mesquite_option=0)
Definition: mesquite.cpp:1144
long GetSequence() const
Definition: mesh.hpp:1027
bool Nonconforming() const
Definition: mesh.hpp:1014
Mesh(int nx, int ny, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0)
Definition: mesh.hpp:549
void MoveNodes(const Vector &displacements)
Definition: mesh.cpp:5410
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
Definition: mesh.cpp:1404
void GetElementData(int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.hpp:658
int dim
Definition: ex3.cpp:47
NCFaceInfo(bool slave, int master, const DenseMatrix *pm)
Definition: mesh.hpp:125
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a triangular Mesh.
Definition: mesh.cpp:1115
void UpdateNURBS()
Definition: mesh.cpp:3137
Symmetric 3D Table.
Definition: stable3d.hpp:28
void FinalizeTetMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a tetrahedral Mesh.
Definition: mesh.cpp:1713
void AddBdrQuadAsTriangles(const int *vi, int attr=1)
Definition: mesh.cpp:1049
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
Definition: mesh.cpp:4243
void SetAttributes()
Definition: mesh.cpp:900
int FindCoarseElement(int i)
Definition: mesh.cpp:6843
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Returns the transformation defining the given face element in a user-defined variable.
Definition: mesh.cpp:391
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:538
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
Definition: mesh.cpp:4294
static FiniteElement * GetTransformationFEforElementType(int)
Definition: mesh.cpp:263
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
Definition: mesh.cpp:3200
STable3D * GetElementToFaceTable(int ret_ftbl=0)
Definition: mesh.cpp:4466
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3355
static void Rotate3(int &, int &, int &)
Definition: mesh.hpp:1184
void Make2D(int nx, int ny, Element::Type type, int generate_edges, double sx, double sy)
Definition: mesh.cpp:2116
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:611
const char * filename
Definition: mesh.hpp:1171
virtual void PrintInfo(std::ostream &out=mfem::out)
Definition: mesh.hpp:1099
Type
Constants for the classes derived from Element.
Definition: element.hpp:37
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:4055
void InitTables()
Definition: mesh.cpp:806
void CheckDisplacements(const Vector &displacements, double &tmax)
Definition: mesh.cpp:5265
void ReadInlineMesh(std::istream &input, int generate_edges=0)
Geometry::Constants< Geometry::TETRAHEDRON > tet_t
Definition: mesh.hpp:166
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:595
void SetNodesOwner(bool nodes_owner)
Set the mesh nodes ownership flag.
Definition: mesh.hpp:925
void SetLayer(const int l)
Definition: mesh.hpp:1153
int GetElementToEdgeTable(Table &, Array< int > &)
Definition: mesh.cpp:4133
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:4009
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:6424
int GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:4029
const Element * GetElement(int i) const
Definition: mesh.hpp:676
Mesh(int n, double sx=1.0)
Definition: mesh.hpp:556
void GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:602
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:624
IsoparametricTransformation Transformation2
Definition: mesh.hpp:141
void Init()
Definition: mesh.cpp:789
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3338
Element * GetBdrElement(int i)
Definition: mesh.hpp:682
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Definition: ncmesh.hpp:83
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:741
int Dimension() const
Definition: mesh.hpp:645
virtual void ReorientTetMesh()
Definition: mesh.cpp:4539
int NumOfBdrElements
Definition: mesh.hpp:56
Table * el_to_edge
Definition: mesh.hpp:132
void InitBaseGeom()
Definition: mesh.cpp:948
static void ShiftL2R(int &, int &, int &)
Definition: mesh.hpp:1178
virtual void LocalRefinement(const Array< int > &marked_el, int type=3)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:5767
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:4034
void GetGeckoElementReordering(Array< int > &ordering)
Definition: mesh.cpp:1180
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
virtual void MarkTetMeshForRefinement(DSTable &v_to_v)
Definition: mesh.cpp:1429
void Destroy()
Definition: mesh.cpp:865
double GetElementSize(int i, int type=0)
Definition: mesh.cpp:64
int CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Definition: mesh.cpp:3596
int SpaceDimension() const
Definition: mesh.hpp:646
void GetBdrElementFace(int i, int *, int *) const
Return the index and the orientation of the face of bdr element i. (3D)
Definition: mesh.cpp:3952
void CheckPartitioning(int *partitioning)
Definition: mesh.cpp:4947
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
Definition: mesh.cpp:7728
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
Definition: mesh.hpp:844
void ReadMFEMMesh(std::istream &input, bool mfem_v11, int &curved)
void AddPointFaceElement(int lf, int gf, int el)
Used in GenerateFaces()
Definition: mesh.cpp:4226
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:2276
const Element *const * GetElementsArray() const
Definition: mesh.hpp:673
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
void FinalizeTopology()
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:1790
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:6027
void PrintVTK(std::ostream &out)
Definition: mesh.cpp:7199
virtual ~Mesh()
Destroys Mesh.
Definition: mesh.hpp:1136
Table * el_to_el
Definition: mesh.hpp:134
Geometry::Constants< Geometry::CUBE > hex_t
Definition: mesh.hpp:167
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Definition: mesh.cpp:8355
void ApplyLocalSlaveTransformation(IsoparametricTransformation &transf, const FaceInfo &fi)
Definition: mesh.cpp:727
static void PrintElement(const Element *, std::ostream &)
Definition: mesh.cpp:2592
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
Definition: mesh.cpp:3107
static void GetElementArrayEdgeTable(const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
Definition: mesh.cpp:4086
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
Definition: mesh.hpp:920
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:569
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
Definition: mesh.cpp:1452
bool IsSlaveFace(const FaceInfo &fi) const
Definition: mesh.cpp:722
void AddQuad(const int *vi, int attr=1)
Definition: mesh.cpp:992
void ReadGmshMesh(std::istream &input)
Array< Vertex > vertices
Definition: mesh.hpp:73
virtual void PrintXG(std::ostream &out=mfem::out) const
Print the mesh to the given stream using Netgen/Truegrid format.
Definition: mesh.cpp:6911
virtual void QuadUniformRefinement()
Refine quadrilateral mesh.
Definition: mesh.cpp:5500
virtual void HexUniformRefinement()
Refine hexahedral mesh.
Definition: mesh.cpp:5602
void RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
Definition: mesh.cpp:6443
void Swap(Mesh &other, bool non_geometry=false)
Definition: mesh.cpp:6233
Array< Element * > elements
Definition: mesh.hpp:68
const Table & ElementToElementTable()
Definition: mesh.cpp:4172
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
Definition: mesh.cpp:4071
int meshgen
Definition: mesh.hpp:61
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:6853
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:6469
Table * bel_to_edge
Definition: mesh.hpp:136
Operation GetLastOperation() const
Return type of last modification of the mesh.
Definition: mesh.hpp:1021
int GetElementBaseGeometry(int i=0) const
Definition: mesh.hpp:688
void AddBdrTriangle(const int *vi, int attr=1)
Definition: mesh.cpp:1039
Table * GetFaceToElementTable() const
Definition: mesh.cpp:3880
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:578
int EulerNumber() const
Equals 1 + num_holes - num_loops.
Definition: mesh.hpp:639
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:176
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:615
void AddBdrQuad(const int *vi, int attr=1)
Definition: mesh.cpp:1044
Array< int > be_to_edge
Definition: mesh.hpp:135
const Table & ElementToFaceTable() const
Definition: mesh.cpp:4208
Class for integration point with weight.
Definition: intrules.hpp:25
void UniformRefinement()
Definition: mesh.cpp:6303
Element * ReadElementWithoutAttr(std::istream &)
Definition: mesh.cpp:2550
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:3698
void AddTriangle(const int *vi, int attr=1)
Definition: mesh.cpp:987
bool OwnsNodes() const
Return the mesh nodes ownership flag.
Definition: mesh.hpp:923
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:279
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:558
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:3350
Table * face_edge
Definition: mesh.hpp:138
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1147
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
Definition: mesh.hpp:661
void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:1868
bool DerefineByError(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
Definition: mesh.cpp:6146
NodeExtrudeCoefficient(const int dim, const int _n, const double _s)
Definition: mesh.cpp:8668
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:1225
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:3120
Table * edge_vertex
Definition: mesh.hpp:139
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf)
long sequence
Definition: mesh.hpp:66
virtual ~NodeExtrudeCoefficient()
Definition: mesh.hpp:1157
void DestroyPointers()
Definition: mesh.cpp:839
Array< FaceInfo > faces_info
Definition: mesh.hpp:129
STable3D * GetFacesTable()
Definition: mesh.cpp:4431
void Make3D(int nx, int ny, int nz, Element::Type type, int generate_edges, double sx, double sy, double sz)
Definition: mesh.cpp:1928
void EnsureNCMesh(bool triangles_nonconforming=false)
Definition: mesh.cpp:6407
NCMesh * ncmesh
Optional non-conforming mesh extension.
Definition: mesh.hpp:177
static void PrintElementWithoutAttr(const Element *, std::ostream &)
Definition: mesh.cpp:2568
virtual bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NC version of GeneralDerefinement.
Definition: mesh.cpp:6101
int Dim
Definition: mesh.hpp:53
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:603
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:624
Geometry::Constants< Geometry::TRIANGLE > tri_t
Definition: mesh.hpp:164
double * GetVertex(int i)
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:656
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:627
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:8575
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:865
Vector data type.
Definition: vector.hpp:48
void ReadTrueGridMesh(std::istream &input)
int GetFaceBaseGeometry(int i) const
Definition: mesh.cpp:3975
void AverageVertices(int *indexes, int n, int result)
Definition: mesh.cpp:5470
Mesh(int nx, int ny, int nz, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0, double sz=1.0)
Definition: mesh.hpp:539
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3500
void SetNode(int i, const double *coord)
Definition: mesh.cpp:5390
void GetNode(int i, double *coord)
Definition: mesh.cpp:5371
void DestroyTables()
Definition: mesh.cpp:824
void GenerateFaces()
Definition: mesh.cpp:4315
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:7535
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:2464
int spaceDim
Definition: mesh.hpp:54
IsoparametricTransformation EdgeTransformation
Definition: mesh.hpp:142
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:52
void AddTet(const int *vi, int attr=1)
Definition: mesh.cpp:997
MemAlloc< Tetrahedron, 1024 > TetMemory
Definition: mesh.hpp:159
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
Definition: mesh.cpp:3791
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
Definition: mesh.cpp:1494
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:5434
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:3344
named_ifgzstream(const char *mesh_name)
Definition: mesh.hpp:1172
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
int NumOfVertices
Definition: mesh.hpp:56
void ScaleSubdomains(double sf)
Definition: mesh.cpp:8166
Array< int > be_to_face
Definition: mesh.hpp:137
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:64
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
Definition: mesh.cpp:6175
void Bisection(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6498
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:6342
Class for parallel meshes.
Definition: pmesh.hpp:32
Abstract data type element.
Definition: element.hpp:27
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:862
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
Definition: mesh.cpp:4273
const Table & ElementToEdgeTable() const
Definition: mesh.cpp:4217
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
Definition: mesh.hpp:699
void GenerateNCFaceInfo()
Definition: mesh.cpp:4384
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:695
int GetBdrElementBaseGeometry(int i=0) const
Definition: mesh.hpp:691
double GetElementVolume(int i)
Definition: mesh.cpp:91
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:172
void InitMesh(int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Definition: mesh.cpp:931
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
Definition: mesh.cpp:3845
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:680
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.cpp:6275
void UpdateNodes()
Update the nodes of a curved mesh after refinement.
Definition: mesh.cpp:5491
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=mfem::out)
Definition: mesh.cpp:206
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:247
void Printer(std::ostream &out=mfem::out, std::string section_delimiter="") const
Definition: mesh.cpp:7075
Class used to extrude the nodes of a mesh.
Definition: mesh.hpp:1145
int NumOfFaces
Definition: mesh.hpp:57