MFEM  v3.0
 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.googlecode.com.
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 "triangle.hpp"
18 #include "tetrahedron.hpp"
19 #include "vertex.hpp"
20 #include "ncmesh.hpp"
21 #include "../fem/eltrans.hpp"
22 #include "../fem/coefficient.hpp"
23 #include <iostream>
24 
25 namespace mfem
26 {
27 
28 // Data type mesh
29 
30 class KnotVector;
31 class NURBSExtension;
32 class FiniteElementSpace;
33 class GridFunction;
34 struct Refinement;
35 
36 #ifdef MFEM_USE_MPI
37 class ParMesh;
38 #endif
39 
40 class Mesh
41 {
42 #ifdef MFEM_USE_MPI
43  friend class ParMesh;
44 #endif
45  friend class NURBSExtension;
46 
47 protected:
48  int Dim;
49  int spaceDim;
50 
53 
55 
56  // 0 = Empty, 1 = Standard (NetGen), 2 = TrueGrid
57  int meshgen;
58 
63 
65  // Vertices are only at the corners of elements, where you would expect them
66  // in the lowest-order mesh.
70 
71  class FaceInfo { public: int Elem1No, Elem2No, Elem1Inf, Elem2Inf; };
73 
78  Table *bel_to_edge; // for 3D
80  mutable Table *face_edge;
81  mutable Table *edge_vertex;
82 
84  Array<int> fc_be_to_edge; // swapped with be_to_edge when switching state
85  Table *c_el_to_face, *f_el_to_face; // for 3D two-level state
86  Array<FaceInfo> fc_faces_info; // for 3D two-level state
87 
91 
92  // Nodes are only active for higher order meshes, and share locations with
93  // the vertecies, plus all the higher- order control points within the
94  // element and along the edges and on the faces.
96  int own_nodes;
97 
98  // Backup of the coarse mesh. Only used if WantTwoLevelState == 1 and
99  // nonconforming refinements are used.
101 
102  static const int tet_faces[4][3];
103  static const int hex_faces[6][4]; // same as Hexahedron::faces
104 
105  static const int tri_orientations[6][3];
106  static const int quad_orientations[8][4];
107 
108 #ifdef MFEM_USE_MEMALLOC
109  friend class Tetrahedron;
110 
113 #endif
114 
115  void Init();
116 
117  void InitTables();
118 
119  void DeleteTables();
120 
123  void DeleteCoarseTables();
124 
125  Element *ReadElementWithoutAttr(std::istream &);
126  static void PrintElementWithoutAttr(const Element *, std::ostream &);
127 
128  Element *ReadElement(std::istream &);
129  static void PrintElement(const Element *, std::ostream &);
130 
131  void SetMeshGen(); // set 'meshgen'
132 
134  double GetLength(int i, int j) const;
135 
138  void GetElementJacobian(int i, DenseMatrix &J);
139 
140  void MarkForRefinement();
142  void GetEdgeOrdering(DSTable &v_to_v, Array<int> &order);
144 
145  void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert);
146  void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert);
147 
149  STable3D *GetElementToFaceTable(int ret_ftbl = 0);
150 
153  void RedRefinement(int i, const DSTable &v_to_v,
154  int *edge1, int *edge2, int *middle)
155  { UniformRefinement(i, v_to_v, edge1, edge2, middle); }
156 
159  void GreenRefinement(int i, const DSTable &v_to_v,
160  int *edge1, int *edge2, int *middle)
161  { Bisection(i, v_to_v, edge1, edge2, middle); }
162 
164  void Bisection(int i, const DSTable &, int *, int *, int *);
165 
167  void Bisection(int i, const DSTable &, int *);
168 
170  void UniformRefinement(int i, const DSTable &, int *, int *, int *);
171 
174  void AverageVertices (int * indexes, int n, int result);
175 
177  void UpdateNodes();
178 
180  virtual void QuadUniformRefinement();
181 
183  virtual void HexUniformRefinement();
184 
186  virtual void NURBSUniformRefinement();
187 
189  virtual void LocalRefinement(const Array<int> &marked_el, int type = 3);
190 
192  void NonconformingRefinement(const Array<Refinement> &refinements,
193  int nc_limit = 0);
194 
196  void LoadPatchTopo(std::istream &input, Array<int> &edge_to_knot);
197 
198  void UpdateNURBS();
199 
200  void PrintTopo(std::ostream &out, const Array<int> &e_to_k) const;
201 
202  void BisectTriTrans (DenseMatrix &pointmat, Triangle *tri,
203  int child);
204 
205  void BisectTetTrans (DenseMatrix &pointmat, Tetrahedron *tet,
206  int child);
207 
208  int GetFineElemPath (int i, int j);
209 
211 
215  int i);
217  int i);
220  int i);
223  int i);
224 
226  static int GetTriOrientation (const int * base, const int * test);
228  static int GetQuadOrientation (const int * base, const int * test);
229 
230  static void GetElementArrayEdgeTable(const Array<Element*> &elem_array,
231  const DSTable &v_to_v,
232  Table &el_to_edge);
233 
237  void GetVertexToVertexTable(DSTable &) const;
238 
245 
247  void AddPointFaceElement(int lf, int gf, int el);
248 
249  void AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1);
250 
251  void AddTriangleFaceElement (int lf, int gf, int el,
252  int v0, int v1, int v2);
253 
254  void AddQuadFaceElement (int lf, int gf, int el,
255  int v0, int v1, int v2, int v3);
259  bool FaceIsTrueInterior(int FaceNo) const
260  {
261  return FaceIsInterior(FaceNo) || (faces_info[FaceNo].Elem2Inf >= 0);
262  }
263 
264  // shift cyclically 3 integers left-to-right
265  inline static void ShiftL2R(int &, int &, int &);
266  // shift cyclically 3 integers so that the smallest is first
267  inline static void Rotate3(int &, int &, int &);
268 
269  void FreeElement (Element *E);
270 
271  void GenerateFaces();
272 
274  void InitMesh(int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem);
275 
280  void Make3D(int nx, int ny, int nz, Element::Type type, int generate_edges,
281  double sx, double sy, double sz);
282 
287  void Make2D(int nx, int ny, Element::Type type, int generate_edges,
288  double sx, double sy);
289 
291  void Make1D(int n, double sx = 1.0);
292 
294  Mesh(NCMesh& ncmesh);
295 
298  void Swap(Mesh& other, bool non_geometry = false);
299 
300 public:
301 
303 
306 
309 
310  Mesh() { Init(); InitTables(); meshgen = 0; Dim = 0; }
311 
312  Mesh(int _Dim, int NVert, int NElem, int NBdrElem = 0, int _spaceDim= -1)
313  {
314  if (_spaceDim == -1)
315  _spaceDim = _Dim;
316  InitMesh(_Dim, _spaceDim, NVert, NElem, NBdrElem);
317  }
318 
319  Element *NewElement(int geom);
320 
321  void AddVertex(const double *);
322  void AddTri(const int *vi, int attr = 1);
323  void AddTriangle(const int *vi, int attr = 1);
324  void AddQuad(const int *vi, int attr = 1);
325  void AddTet(const int *vi, int attr = 1);
326  void AddHex(const int *vi, int attr = 1);
327  void AddHexAsTets(const int *vi, int attr = 1);
328  // 'elem' should be allocated using the NewElement method
329  void AddElement(Element *elem) { elements[NumOfElements++] = elem; }
330  void AddBdrSegment(const int *vi, int attr = 1);
331  void AddBdrTriangle(const int *vi, int attr = 1);
332  void AddBdrQuad(const int *vi, int attr = 1);
333  void AddBdrQuadAsTriangles(const int *vi, int attr = 1);
335  void FinalizeTriMesh(int generate_edges = 0, int refine = 0,
336  bool fix_orientation = true);
337  void FinalizeQuadMesh(int generate_edges = 0, int refine = 0,
338  bool fix_orientation = true);
339  void FinalizeTetMesh(int generate_edges = 0, int refine = 0,
340  bool fix_orientation = true);
341  void FinalizeHexMesh(int generate_edges = 0, int refine = 0,
342  bool fix_orientation = true);
343 
344  void SetAttributes();
345 
350  Mesh(int nx, int ny, int nz, Element::Type type, int generate_edges = 0,
351  double sx = 1.0, double sy = 1.0, double sz = 1.0)
352  {
353  Make3D(nx, ny, nz, type, generate_edges, sx, sy, sz);
354  }
355 
360  Mesh(int nx, int ny, Element::Type type, int generate_edges = 0,
361  double sx = 1.0, double sy = 1.0)
362  {
363  Make2D(nx, ny, type, generate_edges, sx, sy);
364  }
365 
367  explicit Mesh(int n, double sx = 1.0)
368  {
369  Make1D(n, sx);
370  }
371 
375  Mesh(std::istream &input, int generate_edges = 0, int refine = 1,
376  bool fix_orientation = true);
377 
379  Mesh(Mesh *mesh_array[], int num_pieces);
380 
381  /* This is similar to the above mesh constructor, but here the current
382  mesh is destroyed and another one created based on the data stream
383  again given in MFEM, netgen, or VTK format. If generate_edges = 0
384  (default) edges are not generated, if 1 edges are generated. */
385  void Load(std::istream &input, int generate_edges = 0, int refine = 1,
386  bool fix_orientation = true);
387 
389  inline int MeshGenerator() { return meshgen; }
390 
393  inline int GetNV() const { return NumOfVertices; }
394 
396  inline int GetNE() const { return NumOfElements; }
397 
399  inline int GetNBE() const { return NumOfBdrElements; }
400 
402  inline int GetNEdges() const { return NumOfEdges; }
403 
405  inline int GetNFaces() const { return NumOfFaces; }
406 
408  int GetNumFaces() const;
409 
411  inline int EulerNumber() const
414  inline int EulerNumber2D() const
415  { return NumOfVertices - NumOfEdges + NumOfElements; }
416 
417  int Dimension() const { return Dim; }
418  int SpaceDimension() const { return spaceDim; }
419 
421  const double *GetVertex(int i) const { return vertices[i](); }
422  double *GetVertex(int i) { return vertices[i](); }
423 
424  const Element *GetElement(int i) const { return elements[i]; }
425 
426  Element *GetElement(int i) { return elements[i]; }
427 
428  const Element *GetBdrElement(int i) const { return boundary[i]; }
429 
430  Element *GetBdrElement(int i) { return boundary[i]; }
431 
432  const Element *GetFace(int i) const { return faces[i]; }
433 
434  int GetFaceBaseGeometry(int i) const;
435 
436  int GetElementBaseGeometry(int i) const
437  { return elements[i]->GetGeometryType(); }
438 
439  int GetBdrElementBaseGeometry(int i) const
440  { return boundary[i]->GetGeometryType(); }
441 
443  void GetElementVertices(int i, Array<int> &dofs) const
444  { elements[i]->GetVertices(dofs); }
445 
447  void GetBdrElementVertices(int i, Array<int> &dofs) const
448  { boundary[i]->GetVertices(dofs); }
449 
451  void GetElementEdges(int i, Array<int> &, Array<int> &) const;
452 
454  void GetBdrElementEdges(int i, Array<int> &, Array<int> &) const;
455 
458  void GetFaceEdges(int i, Array<int> &, Array<int> &) const;
459 
461  void GetFaceVertices(int i, Array<int> &vert) const
462  {
463  if (Dim == 1)
464  {
465  vert.SetSize(1); vert[0] = i;
466  }
467  else
468  faces[i]->GetVertices(vert);
469  }
470 
472  void GetEdgeVertices(int i, Array<int> &vert) const;
473 
475  Table *GetFaceEdgeTable() const;
476 
478  Table *GetEdgeVertexTable() const;
479 
481  void GetElementFaces(int i, Array<int> &, Array<int> &) const;
482 
484  void GetBdrElementFace(int i, int *, int *) const;
485 
488  int GetBdrElementEdgeIndex(int i) const;
489 
491  int GetElementType(int i) const;
492 
494  int GetBdrElementType(int i) const;
495 
496  /* Return point matrix of element i of dimension Dim X #dofs, where for
497  every degree of freedom we give its coordinates in space of dimension
498  Dim. */
499  void GetPointMatrix(int i, DenseMatrix &pointmat) const;
500 
501  /* Return point matrix of boundary element i of dimension Dim X #dofs,
502  where for every degree of freedom we give its coordinates in space
503  of dimension Dim. */
504  void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const;
505 
507 
511 
514 
517  void GetElementTransformation(int i, const Vector &nodes,
519 
523 
527 
530 
534 
537 
561  int mask = 31);
562 
564  { if (faces_info[FaceNo].Elem2No < 0) return NULL;
565  return GetFaceElementTransformations (FaceNo); }
566 
568 
570  bool FaceIsInterior(int FaceNo) const
571  {
572  return (faces_info[FaceNo].Elem2No >= 0);
573  }
574  void GetFaceElements (int Face, int *Elem1, int *Elem2);
575  void GetFaceInfos (int Face, int *Inf1, int *Inf2);
576 
578  void CheckElementOrientation(bool fix_it = true);
580  void CheckBdrElementOrientation(bool fix_it = true);
581 
583  int GetAttribute(int i) const { return elements[i]->GetAttribute();}
584 
586  int GetBdrAttribute(int i) const { return boundary[i]->GetAttribute(); }
587 
588  const Table &ElementToElementTable();
589 
590  const Table &ElementToFaceTable() const;
591 
592  const Table &ElementToEdgeTable() const;
593 
596 
600  Table *GetFaceToElementTable() const;
601 
610  virtual void ReorientTetMesh();
611 
612  int *CartesianPartitioning(int nxyz[]);
613  int *GeneratePartitioning(int nparts, int part_method = 1);
614  void CheckPartitioning(int *partitioning);
615 
616  void CheckDisplacements(const Vector &displacements, double &tmax);
617 
618  // Vertices are only at the corners of elements, where you would expect them
619  // in the lowest-order mesh.
620  void MoveVertices(const Vector &displacements);
621  void GetVertices(Vector &vert_coord) const;
622  void SetVertices(const Vector &vert_coord);
623 
624  // Nodes are only active for higher order meshes, and share locations with
625  // the vertecies, plus all the higher- order control points within the
626  // element and along the edges and on the faces.
627  void GetNode(int i, double *coord);
628  void SetNode(int i, const double *coord);
629 
630  // Node operations for curved mesh.
631  // They call the corresponding '...Vertices' method if the
632  // mesh is not curved (i.e. Nodes == NULL).
633  void MoveNodes(const Vector &displacements);
634  void GetNodes(Vector &node_coord) const;
635  void SetNodes(const Vector &node_coord);
636 
638  GridFunction *GetNodes() { return Nodes; }
640  void NewNodes(GridFunction &nodes, bool make_owner = false);
643  void SwapNodes(GridFunction *&nodes, int &own_nodes_);
644 
646  void GetNodes(GridFunction &nodes) const;
654  void SetNodalGridFunction(GridFunction *nodes, bool make_owner = false);
658 
660  void UniformRefinement();
661 
670  void GeneralRefinement(Array<Refinement> &refinements,
671  int nonconforming = -1, int nc_limit = 0);
672 
675  void GeneralRefinement(Array<int> &el_to_refine,
676  int nonconforming = -1, int nc_limit = 0);
677 
678  // NURBS mesh refinement methods
680  void DegreeElevate(int t);
681 
684  void UseTwoLevelState (int use)
685  {
686  if (!use && State != Mesh::NORMAL)
688  WantTwoLevelState = use;
689  }
690 
692  void SetState (int s);
693 
694  int GetState() const { return State; }
695 
698  int GetNumFineElems (int i);
699 
701  int GetRefinementType (int i);
702 
705  int GetFineElem (int i, int j);
706 
710  ElementTransformation * GetFineElemTrans (int i, int j);
711 
713  virtual void PrintXG(std::ostream &out = std::cout) const;
714 
716  virtual void Print(std::ostream &out = std::cout) const;
717 
719  void PrintVTK(std::ostream &out);
720 
725  void PrintVTK(std::ostream &out, int ref, int field_data=0);
726 
727  void GetElementColoring(Array<int> &colors, int el0 = 0);
728 
732  void PrintWithPartitioning (int *partitioning,
733  std::ostream &out, int elem_attr = 0) const;
734 
735  void PrintElementsWithPartitioning (int *partitioning,
736  std::ostream &out,
737  int interior_faces = 0);
738 
740 
744  void PrintSurfaces(const Table &Aface_face, std::ostream &out) const;
745 
746  void ScaleSubdomains (double sf);
747  void ScaleElements (double sf);
748 
749  void Transform(void (*f)(const Vector&, Vector&));
750 
753  double GetElementSize(int i, int type = 0);
754 
755  double GetElementSize(int i, const Vector &dir);
756 
757  double GetElementVolume(int i);
758 
759  void PrintCharacteristics(Vector *Vh = NULL, Vector *Vk = NULL);
760 
761  void MesquiteSmooth(const int mesquite_option = 0);
762 
764  virtual ~Mesh();
765 };
766 
769 std::ostream &operator<<(std::ostream &out, const Mesh &mesh);
770 
771 
774 {
775 private:
776  int n, layer;
777  double p[2], s;
778  Vector tip;
779 public:
780  NodeExtrudeCoefficient(const int dim, const int _n, const double _s);
781  void SetLayer(const int l) { layer = l; }
783  virtual void Eval(Vector &V, ElementTransformation &T,
784  const IntegrationPoint &ip);
786 };
787 
788 
790 Mesh *Extrude1D(Mesh *mesh, const int ny, const double sy,
791  const bool closed = false);
792 
793 
794 // inline functions
795 inline void Mesh::ShiftL2R(int &a, int &b, int &c)
796 {
797  int t = a;
798  a = c; c = b; b = t;
799 }
800 
801 inline void Mesh::Rotate3(int &a, int &b, int &c)
802 {
803  if (a < b)
804  {
805  if (a > c)
806  ShiftL2R(a, b, c);
807  }
808  else
809  {
810  if (b < c)
811  ShiftL2R(c, b, a);
812  else
813  ShiftL2R(a, b, c);
814  }
815 }
816 
817 }
818 
819 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:42
void AddHex(const int *vi, int attr=1)
Definition: mesh.cpp:797
void SetState(int s)
Change the mesh state to NORMAL, TWO_LEVEL_COARSE, TWO_LEVEL_FINE.
Definition: mesh.cpp:6298
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:8243
void GetFaceEdges(int i, Array< int > &, Array< int > &) const
Definition: mesh.cpp:3393
void PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
Definition: mesh.cpp:8000
Table * f_el_to_face
Definition: mesh.hpp:85
void GetPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:3660
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:4209
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:586
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:421
void ScaleElements(double sf)
Definition: mesh.cpp:8121
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:26
Table * GetEdgeVertexTable() const
Returns the edge-to-vertex Table (3D)
Definition: mesh.cpp:3451
void FreeElement(Element *E)
Definition: mesh.cpp:8206
void GetFaceInfos(int Face, int *Inf1, int *Inf2)
Definition: mesh.cpp:621
void SetVertices(const Vector &vert_coord)
Definition: mesh.cpp:4902
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:3421
int f_NumOfEdges
Definition: mesh.hpp:62
bool FaceIsTrueInterior(int FaceNo) const
Definition: mesh.hpp:259
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Definition: mesh.cpp:8256
int GetBdrElementBaseGeometry(int i) const
Definition: mesh.hpp:439
int f_NumOfElements
Definition: mesh.hpp:60
Table * c_bel_to_edge
Definition: mesh.hpp:83
void AddHexAsTets(const int *vi, int attr=1)
Definition: mesh.cpp:802
void GetFaceElements(int Face, int *Elem1, int *Elem2)
Definition: mesh.cpp:615
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Array< Element * > boundary
Definition: mesh.hpp:68
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:4252
void BisectTetTrans(DenseMatrix &pointmat, Tetrahedron *tet, int child)
Definition: mesh.cpp:6722
int own_nodes
Definition: mesh.hpp:96
void MoveVertices(const Vector &displacements)
Definition: mesh.cpp:4886
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:399
IsoparametricTransformation Transformation
Definition: mesh.hpp:88
void GetLocalPtToSegTransformation(IsoparametricTransformation &, int)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:374
void PrintWithPartitioning(int *partitioning, std::ostream &out, int elem_attr=0) const
Definition: mesh.cpp:7582
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:4964
int NumOfEdges
Definition: mesh.hpp:52
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:4978
void GetElementFaces(int i, Array< int > &, Array< int > &) const
Return the indices and the orientations of all faces of element i.
Definition: mesh.cpp:3534
int GetElementBaseGeometry(int i) const
Definition: mesh.hpp:436
void DegreeElevate(int t)
Definition: mesh.cpp:2878
int EulerNumber2D() const
Equals 1 - num_holes.
Definition: mesh.hpp:414
void AddTri(const int *vi, int attr=1)
Definition: mesh.cpp:769
int WantTwoLevelState
Definition: mesh.hpp:54
void DeleteTables()
Definition: mesh.cpp:645
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:461
Data type dense matrix.
Definition: densemat.hpp:22
GridFunction * Nodes
Definition: mesh.hpp:95
int NumOfElements
Definition: mesh.hpp:51
Mesh * Extrude1D(Mesh *mesh, const int ny, const double sy, const bool closed)
Extrude a 1D mesh.
Definition: mesh.cpp:8274
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:8183
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:3625
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:396
const Element * GetFace(int i) const
Definition: mesh.hpp:432
Element * GetElement(int i)
Definition: mesh.hpp:426
void BisectTriTrans(DenseMatrix &pointmat, Triangle *tri, int child)
Definition: mesh.cpp:6698
void GenerateBoundaryElements()
Definition: mesh.cpp:845
void GetElementEdges(int i, Array< int > &, Array< int > &) const
Return the indices and the orientations of all edges of element i.
Definition: mesh.cpp:3345
void GetElementVertices(int i, Array< int > &dofs) const
Returns the indices of the dofs of element i.
Definition: mesh.hpp:443
void DeleteCoarseTables()
Definition: mesh.cpp:671
void SetMeshGen()
Definition: mesh.cpp:1733
ElementTransformation * GetFineElemTrans(int i, int j)
Definition: mesh.cpp:6802
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3230
void GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:4893
virtual ~Mesh()
Destroys mesh.
Definition: mesh.cpp:8221
void AddBdrSegment(const int *vi, int attr=1)
Definition: mesh.cpp:817
Table * c_el_to_face
Definition: mesh.hpp:85
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
Definition: mesh.hpp:563
void CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Definition: mesh.cpp:3266
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:499
Array< Element * > faces
Definition: mesh.hpp:69
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Definition: mesh.cpp:1292
void GetVertexToVertexTable(DSTable &) const
Definition: mesh.cpp:3723
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:153
void MarkTriMeshForRefinement()
Definition: mesh.cpp:950
IsoparametricTransformation FaceTransformation
Definition: mesh.hpp:89
Element * ReadElement(std::istream &)
Definition: mesh.cpp:1715
void KnotInsert(Array< KnotVector * > &kv)
Definition: mesh.cpp:2853
Array< FaceInfo > fc_faces_info
Definition: mesh.hpp:86
int GetFineElemPath(int i, int j)
Definition: mesh.cpp:6769
Element * NewElement(int geom)
Definition: mesh.cpp:1670
FaceElementTransformations FaceElemTr
Definition: mesh.hpp:90
Table * el_to_face
Definition: mesh.hpp:75
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
Definition: mesh.cpp:329
void MarkForRefinement()
Definition: mesh.cpp:939
static const int tet_faces[4][3]
Definition: mesh.hpp:102
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:218
int GetState() const
Definition: mesh.hpp:694
void AddVertex(const double *)
Definition: mesh.cpp:760
Mesh(int _Dim, int NVert, int NElem, int NBdrElem=0, int _spaceDim=-1)
Definition: mesh.hpp:312
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
Definition: mesh.cpp:7189
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:3083
void AddElement(Element *elem)
Definition: mesh.hpp:329
void Load(std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
Definition: mesh.cpp:1772
void GetElementJacobian(int i, DenseMatrix &J)
Definition: mesh.cpp:31
void MesquiteSmooth(const int mesquite_option=0)
Definition: mesquite.cpp:1051
Mesh(int nx, int ny, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0)
Definition: mesh.hpp:360
void MoveNodes(const Vector &displacements)
Definition: mesh.cpp:4940
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
Definition: mesh.cpp:963
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Definition: mesh.cpp:892
void UpdateNURBS()
Definition: mesh.cpp:2896
Symmetric 3D Table.
Definition: stable3d.hpp:28
static const int quad_orientations[8][4]
Definition: mesh.hpp:106
void FinalizeTetMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Definition: mesh.cpp:1254
void AddBdrQuadAsTriangles(const int *vi, int attr=1)
Definition: mesh.cpp:832
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
Definition: mesh.cpp:3884
void SetAttributes()
Definition: mesh.cpp:683
void MarkTetMeshForRefinement()
Definition: mesh.cpp:987
int State
Definition: mesh.hpp:54
Table * f_el_to_edge
Definition: mesh.hpp:83
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Definition: mesh.cpp:248
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:388
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
Definition: mesh.cpp:3930
int GetRefinementType(int i)
For a given coarse element i returns its refinement type.
Definition: mesh.cpp:6541
static FiniteElement * GetTransformationFEforElementType(int)
Definition: mesh.cpp:141
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
Definition: mesh.cpp:2951
STable3D * GetElementToFaceTable(int ret_ftbl=0)
Definition: mesh.cpp:4053
static void Rotate3(int &, int &, int &)
Definition: mesh.hpp:801
void Make2D(int nx, int ny, Element::Type type, int generate_edges, double sx, double sy)
Definition: mesh.cpp:1473
int MeshGenerator()
Truegrid or NetGen?
Definition: mesh.hpp:389
Type
Constants for the classes derived from Element.
Definition: element.hpp:37
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
Definition: mesh.cpp:3674
void InitTables()
Definition: mesh.cpp:639
void GetBdrElementEdges(int i, Array< int > &, Array< int > &) const
Return the indices and the orientations of all edges of bdr element i.
Definition: mesh.cpp:3364
void CheckDisplacements(const Vector &displacements, double &tmax)
Definition: mesh.cpp:4819
void CheckElementOrientation(bool fix_it=true)
Check the orientation of the elements.
Definition: mesh.cpp:3098
Table * f_bel_to_edge
Definition: mesh.hpp:83
void SetLayer(const int l)
Definition: mesh.hpp:781
int GetElementToEdgeTable(Table &, Array< int > &)
Definition: mesh.cpp:3748
int c_NumOfEdges
Definition: mesh.hpp:61
int GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:3632
Data type triangle element.
Definition: triangle.hpp:22
const Element * GetElement(int i) const
Definition: mesh.hpp:424
Mesh(int n, double sx=1.0)
Definition: mesh.hpp:367
void GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:478
IsoparametricTransformation Transformation2
Definition: mesh.hpp:88
void Init()
Definition: mesh.cpp:627
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3066
Element * GetBdrElement(int i)
Definition: mesh.hpp:430
int f_NumOfVertices
Definition: mesh.hpp:60
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Definition: ncmesh.hpp:65
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:598
int Dimension() const
Definition: mesh.hpp:417
virtual void ReorientTetMesh()
Definition: mesh.cpp:4118
int NumOfBdrElements
Definition: mesh.hpp:51
Table * el_to_edge
Definition: mesh.hpp:74
static void ShiftL2R(int &, int &, int &)
Definition: mesh.hpp:795
virtual void LocalRefinement(const Array< int > &marked_el, int type=3)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:5371
int c_NumOfBdrElements
Definition: mesh.hpp:59
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:3647
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
double GetElementSize(int i, int type=0)
Definition: mesh.cpp:39
int SpaceDimension() const
Definition: mesh.hpp:418
void GetBdrElementFace(int i, int *, int *) const
Return the index and the orientation of the face of bdr element i. (3D)
Definition: mesh.cpp:3560
int c_NumOfVertices
Definition: mesh.hpp:59
void CheckPartitioning(int *partitioning)
Definition: mesh.cpp:4533
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Definition: mesh.cpp:7136
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
Definition: mesh.cpp:7678
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
Definition: mesh.hpp:570
void AddPointFaceElement(int lf, int gf, int el)
Used in GenerateFaces()
Definition: mesh.cpp:3862
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:1621
int GetNumFineElems(int i)
Definition: mesh.cpp:6445
Array< int > bdr_attributes
Definition: mesh.hpp:305
void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
Definition: mesh.cpp:5738
void PrintVTK(std::ostream &out)
Print the mesh in VTK format (linear and quadratic meshes only).
Definition: mesh.cpp:7225
void UseTwoLevelState(int use)
Definition: mesh.hpp:684
Table * el_to_el
Definition: mesh.hpp:76
int GetBisectionHierarchy(Element *E)
Definition: mesh.cpp:6502
Abstract finite element space.
Definition: fespace.hpp:61
static void PrintElement(const Element *, std::ostream &)
Definition: mesh.cpp:1727
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
Definition: mesh.cpp:2868
static void GetElementArrayEdgeTable(const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
Definition: mesh.cpp:3701
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
Definition: mesh.hpp:638
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:293
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
Definition: mesh.cpp:1006
void AddQuad(const int *vi, int attr=1)
Definition: mesh.cpp:779
Array< Vertex > vertices
Definition: mesh.hpp:67
MemAlloc< BisectedElement, 1024 > BEMemory
Definition: mesh.hpp:112
virtual void QuadUniformRefinement()
Refine quadrilateral mesh.
Definition: mesh.cpp:5009
virtual void HexUniformRefinement()
Refine hexahedral mesh.
Definition: mesh.cpp:5156
void Swap(Mesh &other, bool non_geometry=false)
Definition: mesh.cpp:5832
Array< Element * > elements
Definition: mesh.hpp:64
const Table & ElementToElementTable()
Definition: mesh.cpp:3783
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
Definition: mesh.cpp:3688
int meshgen
Definition: mesh.hpp:57
Table * bel_to_edge
Definition: mesh.hpp:78
void AddBdrTriangle(const int *vi, int attr=1)
Definition: mesh.cpp:822
Array< int > fc_be_to_edge
Definition: mesh.hpp:84
Table * GetFaceToElementTable() const
Definition: mesh.cpp:3506
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:455
int EulerNumber() const
Equals 1 + num_holes - num_loops.
Definition: mesh.hpp:411
NURBSExtension * NURBSext
Definition: mesh.hpp:307
int GetNV() const
Definition: mesh.hpp:393
void AddBdrQuad(const int *vi, int attr=1)
Definition: mesh.cpp:827
Array< int > be_to_edge
Definition: mesh.hpp:77
const Table & ElementToFaceTable() const
Definition: mesh.cpp:3848
Table * c_el_to_edge
Definition: mesh.hpp:83
Class for integration point with weight.
Definition: intrules.hpp:25
void UniformRefinement()
Definition: mesh.cpp:5874
Element * ReadElementWithoutAttr(std::istream &)
Definition: mesh.cpp:1690
void AddTriangle(const int *vi, int attr=1)
Definition: mesh.cpp:774
int f_NumOfFaces
Definition: mesh.hpp:62
static const int hex_faces[6][4]
Definition: mesh.hpp:103
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:157
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
Definition: mesh.cpp:411
const FiniteElementSpace * GetNodalFESpace()
Definition: mesh.cpp:3078
int c_NumOfElements
Definition: mesh.hpp:59
Table * face_edge
Definition: mesh.hpp:80
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Definition: mesh.cpp:916
int c_NumOfFaces
Definition: mesh.hpp:61
NodeExtrudeCoefficient(const int dim, const int _n, const double _s)
Definition: mesh.cpp:8250
Table * edge_vertex
Definition: mesh.hpp:81
void GetBdrElementVertices(int i, Array< int > &dofs) const
Returns the indices of the dofs of boundary element i.
Definition: mesh.hpp:447
virtual ~NodeExtrudeCoefficient()
Definition: mesh.hpp:785
Array< FaceInfo > faces_info
Definition: mesh.hpp:72
STable3D * GetFacesTable()
Definition: mesh.cpp:4020
void Make3D(int nx, int ny, int nz, Element::Type type, int generate_edges, double sx, double sy, double sz)
Definition: mesh.cpp:1317
NCMesh * ncmesh
Definition: mesh.hpp:308
static void PrintElementWithoutAttr(const Element *, std::ostream &)
Definition: mesh.cpp:1705
int GetFineElem(int i, int j)
Definition: mesh.cpp:6601
int Dim
Definition: mesh.hpp:48
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:402
double * GetVertex(int i)
Definition: mesh.hpp:422
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:405
Vector data type.
Definition: vector.hpp:29
int GetFaceBaseGeometry(int i) const
Definition: mesh.cpp:3589
void AverageVertices(int *indexes, int n, int result)
Definition: mesh.cpp:4988
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:350
int f_NumOfBdrElements
Definition: mesh.hpp:60
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3200
void SetNode(int i, const double *coord)
Definition: mesh.cpp:4924
void GetNode(int i, double *coord)
Definition: mesh.cpp:4909
void GenerateFaces()
Definition: mesh.cpp:3954
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:7515
int spaceDim
Definition: mesh.hpp:49
IsoparametricTransformation EdgeTransformation
Definition: mesh.hpp:89
void AddTet(const int *vi, int attr=1)
Definition: mesh.cpp:784
MemAlloc< Tetrahedron, 1024 > TetMemory
Definition: mesh.hpp:111
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
Definition: mesh.cpp:3429
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
Definition: mesh.cpp:1046
void SetNodes(const Vector &node_coord)
Definition: mesh.cpp:4956
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Definition: mesh.cpp:3072
int NumOfVertices
Definition: mesh.hpp:51
void ScaleSubdomains(double sf)
Definition: mesh.cpp:8059
Array< int > be_to_face
Definition: mesh.hpp:79
Mesh * nc_coarse_level
Definition: mesh.hpp:100
virtual void PrintXG(std::ostream &out=std::cout) const
Print the mesh to the given stream using Netgen/Truegrid format.
Definition: mesh.cpp:6992
void Bisection(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:5942
Class for parallel meshes.
Definition: pmesh.hpp:27
Abstract data type element.
Definition: element.hpp:27
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:583
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
Definition: mesh.cpp:3906
const Table & ElementToEdgeTable() const
Definition: mesh.cpp:3855
double GetElementVolume(int i)
Definition: mesh.cpp:60
Array< int > attributes
Definition: mesh.hpp:304
void InitMesh(int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Definition: mesh.cpp:742
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL)
Definition: mesh.cpp:76
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
Definition: mesh.cpp:3475
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:428
void GeneralRefinement(Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:5894
void UpdateNodes()
Update the nodes of a curved mesh after refinement.
Definition: mesh.cpp:5003
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:159
static const int tri_orientations[6][3]
Definition: mesh.hpp:105
Class used to exrude the nodes of a mesh.
Definition: mesh.hpp:773
int NumOfFaces
Definition: mesh.hpp:52