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