MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ncmesh.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_NCMESH
13 #define MFEM_NCMESH
14 
15 #include <vector>
16 #include <iostream>
17 
18 #include "../config/config.hpp"
19 #include "../general/hash.hpp"
20 #include "../linalg/densemat.hpp"
21 #include "element.hpp"
22 #include "vertex.hpp"
23 #include "../fem/geom.hpp"
24 
25 namespace mfem
26 {
27 
32 struct Refinement
33 {
34  int index;
35  char ref_type;
36 
37  Refinement(int index, int type = 7)
38  : index(index), ref_type(type) {}
39 };
40 
41 
63 class NCMesh
64 {
65 protected:
66  struct Element; // forward
67 
68 public:
69 
73  NCMesh(const Mesh *mesh, std::istream *vertex_parents = NULL);
74 
76  NCMesh(const NCMesh &other);
77 
78  int Dimension() const { return Dim; }
79  int SpaceDimension() const { return spaceDim; }
80 
85  virtual void Refine(const Array<Refinement> &refinements);
86 
90  virtual void LimitNCLevel(int max_level);
91 
93  struct MeshId
94  {
95  int index;
96  int local;
98 
99  MeshId(int index = -1, Element* element = NULL, int local = -1)
101  };
102 
105  struct Master : public MeshId
106  {
108 
109  Master(int index, Element* element, int local, int sb, int se)
110  : MeshId(index, element, local), slaves_begin(sb), slaves_end(se) {}
111  };
112 
115  struct Slave : public MeshId
116  {
117  int master;
119 
120  Slave(int index) : MeshId(index), master(-1) {}
121  };
122 
124  struct NCList
125  {
126  std::vector<MeshId> conforming;
127  std::vector<Master> masters;
128  std::vector<Slave> slaves;
129  // TODO: switch to Arrays when fixed for non-POD types
130 
131  void Clear() { conforming.clear(); masters.clear(); slaves.clear(); }
132  bool Empty() const { return !conforming.size() && !masters.size(); }
133  long MemoryUsage() const;
134  };
135 
138  {
139  if (face_list.Empty()) { BuildFaceList(); }
140  return face_list;
141  }
142 
145  {
146  if (edge_list.Empty()) { BuildEdgeList(); }
147  return edge_list;
148  }
149 
155  {
158 
160  bool IsIdentity() const { return !point_matrix.Data(); }
161  };
162 
166  void MarkCoarseLevel();
167 
169  void ClearCoarseLevel() { coarse_elements.DeleteAll(); }
170 
176  FineTransform* GetFineTransforms();
177 
181  int GetEdgeMaster(int v1, int v2) const;
182 
188  virtual void GetBoundaryClosure(const Array<int> &bdr_attr_is_ess,
189  Array<int> &bdr_vertices,
190  Array<int> &bdr_edges);
191 
193  void PrintVertexParents(std::ostream &out) const;
194 
196  void PrintCoarseElements(std::ostream &out) const;
197 
200  void LoadVertexParents(std::istream &input);
201 
203  void LoadCoarseElements(std::istream &input);
204 
206  void SetVertexPositions(const Array<mfem::Vertex> &vertices);
207 
209  long MemoryUsage() const;
210 
211  void PrintMemoryDetail() const;
212 
213  virtual ~NCMesh();
214 
215 
216 protected: // interface for Mesh to be able to construct itself from NCMesh
217 
218  friend class Mesh;
219 
221  void GetMeshComponents(Array<mfem::Vertex>& vertices,
222  Array<mfem::Element*>& elements,
223  Array<mfem::Element*>& boundary) const;
224 
227  virtual void OnMeshUpdated(Mesh *mesh);
228 
229 
230 protected: // implementation
231 
232  int Dim, spaceDim;
233  bool Iso;
234 
235  Element* CopyHierarchy(Element* elem);
236  void DeleteHierarchy(Element* elem);
237 
238 
239  // primary data
240 
243  struct RefCount
244  {
246 
247  RefCount() : ref_count(0) {}
248 
249  int Ref()
250  {
251  return ++ref_count;
252  }
253  int Unref()
254  {
255  int ret = --ref_count;
256  if (!ret) { delete this; }
257  return ret;
258  }
259  };
260 
263  struct Vertex : public RefCount
264  {
265  int index;
266  double pos[3];
267 
268  Vertex() {}
269  Vertex(double x, double y, double z) : index(-1)
270  { pos[0] = x, pos[1] = y, pos[2] = z; }
271  Vertex(const Vertex& other) { std::memcpy(this, &other, sizeof(*this)); }
272  };
273 
275  struct Edge : public RefCount
276  {
277  int attribute;
278  int index;
279 
280  Edge() : attribute(-1), index(-1) {}
281  Edge(const Edge &other) { std::memcpy(this, &other, sizeof(*this)); }
282  bool Boundary() const { return attribute >= 0; }
283  };
284 
295  struct Node : public Hashed2<Node>
296  {
299 
300  Node(int id) : Hashed2<Node>(id), vertex(NULL), edge(NULL) {}
301  Node(const Node &other);
302  ~Node();
303 
304  // Bump ref count on a vertex or an edge, or create them. Used when an
305  // element starts using a vertex or an edge.
306  void RefVertex();
307  void RefEdge();
308 
309  // Decrement ref on vertex or edge when an element is not using them
310  // anymore. The vertex, edge or the whole Node can autodestruct.
311  // (The hash-table pointer needs to be known then to remove the node.)
314  };
315 
320  struct Face : public RefCount, public Hashed4<Face>
321  {
322  int attribute;
323  int index;
324  Element* elem[2];
325 
326  Face(int id);
327  Face(const Face& other);
328 
329  bool Boundary() const { return attribute >= 0; }
330 
331  // add or remove an element from the 'elem[2]' array
332  void RegisterElement(Element* e);
333  void ForgetElement(Element* e);
334 
335  // return one of elem[0] or elem[1] and make sure the other is NULL
336  Element* GetSingleElement() const;
337 
338  // overloaded Unref without auto-destruction
339  int Unref() { return --ref_count; }
340  };
341 
345  struct Element
346  {
347  char geom;
348  char ref_type;
349  int index;
350  int rank;
352  union
353  {
354  Node* node[8];
356  };
358 
359  Element(int geom, int attr);
360  Element(const Element& other) { std::memcpy(this, &other, sizeof(*this)); }
361  };
362 
363  Array<Element*> root_elements; // coarsest mesh, initialized by constructor
364 
365  HashTable<Node> nodes; // associative container holding all Nodes
366  HashTable<Face> faces; // associative container holding all Faces
367 
368 
369  // secondary data
370 
377  virtual void Update();
378 
379  Array<Element*> leaf_elements; // finest level, updated by UpdateLeafElements
380 
381  Array<int> vertex_nodeId; // vertex-index to node-id map, see UpdateVertices
382 
385 
388 
391 
392  virtual void UpdateVertices();
393 
394  void CollectLeafElements(Element* elem);
395  void UpdateLeafElements();
396 
397  virtual void AssignLeafIndices();
398 
399  virtual bool IsGhost(const Element* elem) const { return false; }
400  virtual int GetNumGhosts() const { return 0; }
401 
402 
403  // refinement
404 
405  struct ElemRefType
406  {
408  int ref_type;
409 
411  : elem(elem), ref_type(type) {}
412  };
413 
415 
416  void RefineElement(Element* elem, char ref_type);
417  void DerefineElement(Element* elem);
418 
419  Element* NewHexahedron(Node* n0, Node* n1, Node* n2, Node* n3,
420  Node* n4, Node* n5, Node* n6, Node* n7,
421  int attr,
422  int fattr0, int fattr1, int fattr2,
423  int fattr3, int fattr4, int fattr5);
424 
425  Element* NewQuadrilateral(Node* n0, Node* n1, Node* n2, Node* n3,
426  int attr,
427  int eattr0, int eattr1, int eattr2, int eattr3);
428 
429  Element* NewTriangle(Node* n0, Node* n1, Node* n2,
430  int attr, int eattr0, int eattr1, int eattr2);
431 
432  Vertex* NewVertex(Node* v1, Node* v2);
433 
434  Node* GetMidEdgeVertex(Node* v1, Node* v2);
436  Node* GetMidFaceVertex(Node* e1, Node* e2, Node* e3, Node* e4);
437 
438  int FaceSplitType(Node* v1, Node* v2, Node* v3, Node* v4,
439  Node* mid[4] = NULL /* optional output of mid-edge nodes*/) const;
440 
441  void ForceRefinement(Node* v1, Node* v2, Node* v3, Node* v4);
442 
443  void CheckAnisoFace(Node* v1, Node* v2, Node* v3, Node* v4,
444  Node* mid12, Node* mid34, int level = 0);
445 
446  void CheckIsoFace(Node* v1, Node* v2, Node* v3, Node* v4,
447  Node* e1, Node* e2, Node* e3, Node* e4, Node* midf);
448 
449  void RefElementNodes(Element *elem);
450  void UnrefElementNodes(Element *elem);
451  void RegisterFaces(Element* elem, int *fattr = NULL);
452 
453  Node* PeekAltParents(Node* v1, Node* v2);
454 
455  bool NodeSetX1(Node* node, Node** n);
456  bool NodeSetX2(Node* node, Node** n);
457  bool NodeSetY1(Node* node, Node** n);
458  bool NodeSetY2(Node* node, Node** n);
459  bool NodeSetZ1(Node* node, Node** n);
460  bool NodeSetZ2(Node* node, Node** n);
461 
462 
463  // face/edge lists
464 
465  static int find_node(Element* elem, Node* node);
466  static int find_node(Element* elem, int node_id);
467  static int find_hex_face(int a, int b, int c);
468 
469  void ReorderFacePointMat(Node* v0, Node* v1, Node* v2, Node* v3,
470  Element* elem, DenseMatrix& mat) const;
471  struct PointMatrix;
472 
473  void TraverseFace(Node* v0, Node* v1, Node* v2, Node* v3,
474  const PointMatrix& pm, int level);
475 
476  void TraverseEdge(Node* v0, Node* v1, double t0, double t1, int level);
477 
478  virtual void BuildFaceList();
479  virtual void BuildEdgeList();
480 
481  virtual void ElementSharesEdge(Element* elem, Edge* edge) {} // ParNCMesh
482  virtual void ElementSharesFace(Element* elem, Face* face) {} // ParNCMesh
483 
484 
485  // neighbors / element_vertex table
486 
493  void FindSetNeighbors(const Array<char> &elem_set,
494  Array<Element*> *neighbors,
495  Array<char> *neighbor_set = NULL);
496 
500  void FindNeighbors(const Element* elem,
501  Array<Element*> &neighbors,
502  const Array<Element*> *search_set = NULL);
503 
504  void CollectEdgeVertices(Node *v0, Node *v1, Array<int> &indices);
505  void CollectFaceVertices(Node* v0, Node* v1, Node* v2, Node* v3,
506  Array<int> &indices);
508 
510  {
512  }
513 
514 
515  // coarse to fine transformations
516 
517  struct Point
518  {
519  int dim;
520  double coord[3];
521 
522  Point() { dim = 0; }
523 
524  Point(double x, double y)
525  {
526  dim = 2; coord[0] = x; coord[1] = y;
527  }
528 
529  Point(double x, double y, double z)
530  {
531  dim = 3; coord[0] = x; coord[1] = y; coord[2] = z;
532  }
533 
534  Point(const Point& p0, const Point& p1)
535  {
536  dim = p0.dim;
537  for (int i = 0; i < dim; i++)
538  {
539  coord[i] = (p0.coord[i] + p1.coord[i]) * 0.5;
540  }
541  }
542 
543  Point(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
544  {
545  dim = p0.dim;
546  MFEM_ASSERT(p1.dim == dim && p2.dim == dim && p3.dim == dim, "");
547  for (int i = 0; i < dim; i++)
548  {
549  coord[i] = (p0.coord[i] + p1.coord[i] + p2.coord[i] + p3.coord[i])
550  * 0.25;
551  }
552  }
553 
554  Point& operator=(const Point& src)
555  {
556  dim = src.dim;
557  for (int i = 0; i < dim; i++)
558  {
559  coord[i] = src.coord[i];
560  }
561  return *this;
562  }
563  };
564 
565  struct PointMatrix
566  {
567  int np;
569 
570  PointMatrix(const Point& p0, const Point& p1, const Point& p2)
571  { np = 3; points[0] = p0; points[1] = p1; points[2] = p2; }
572 
573  PointMatrix(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
574  { np = 4; points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3; }
575 
576  PointMatrix(const Point& p0, const Point& p1, const Point& p2,
577  const Point& p3, const Point& p4, const Point& p5,
578  const Point& p6, const Point& p7)
579  {
580  np = 8;
581  points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3;
582  points[4] = p4; points[5] = p5; points[6] = p6; points[7] = p7;
583  }
584 
585  Point& operator()(int i) { return points[i]; }
586  const Point& operator()(int i) const { return points[i]; }
587 
588  void GetMatrix(DenseMatrix& point_matrix) const;
589  };
590 
593 
594  void GetFineTransforms(Element* elem, int coarse_index,
595  FineTransform *transforms, const PointMatrix &pm);
596 
597 
598  // utility
599 
600  Node* GetEdgeMaster(Node* node) const;
601 
602  static void find_face_nodes(const Face *face, Node* node[4]);
603 
604  int EdgeSplitLevel(Node* v1, Node* v2) const;
605  void FaceSplitLevel(Node* v1, Node* v2, Node* v3, Node* v4,
606  int& h_level, int& v_level) const;
607 
608  void CountSplits(Element* elem, int splits[3]) const;
609 
610  int CountElements(Element* elem) const;
611 
612  int PrintElements(std::ostream &out, Element* elem, int &coarse_id) const;
613 
614  void CountObjects(int &nelem, int &nvert, int &nedges) const;
615 
616 
617 public: // TODO: maybe make this part of mfem::Geometry?
618 
621  struct GeomInfo
622  {
623  int nv, ne, nf, nfv; // number of: vertices, edges, faces, face vertices
624  int edges[12][2]; // edge vertices (up to 12 edges)
625  int faces[6][4]; // face vertices (up to 6 faces)
626 
628  GeomInfo() : initialized(false) {}
629  void Initialize(const mfem::Element* elem);
630  };
631 
633 
634 #ifdef MFEM_DEBUG
635 public:
636  void DebugNeighbors(Array<char> &elem_set);
637 #endif
638 
639  friend class ParNCMesh;
640 };
641 
642 }
643 
644 #endif
NCList face_list
lazy-initialized list of faces, see GetFaceList
Definition: ncmesh.hpp:383
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
Definition: ncmesh.hpp:384
NCMesh(const Mesh *mesh, std::istream *vertex_parents=NULL)
Definition: ncmesh.cpp:55
virtual void LimitNCLevel(int max_level)
Definition: ncmesh.cpp:2530
Element * NewHexahedron(Node *n0, Node *n1, Node *n2, Node *n3, Node *n4, Node *n5, Node *n6, Node *n7, int attr, int fattr0, int fattr1, int fattr2, int fattr3, int fattr4, int fattr5)
Definition: ncmesh.cpp:472
void FindNeighbors(const Element *elem, Array< Element * > &neighbors, const Array< Element * > *search_set=NULL)
Definition: ncmesh.cpp:1965
char ref_type
refinement XYZ bit mask (7 = full isotropic)
Definition: ncmesh.hpp:35
int index
edge number in the Mesh
Definition: ncmesh.hpp:278
Array< Element * > coarse_elements
state of leaf_elements before Refine(), set by MarkCoarseLevel()
Definition: ncmesh.hpp:592
int PrintElements(std::ostream &out, Element *elem, int &coarse_id) const
Definition: ncmesh.cpp:2630
PointMatrix(const Point &p0, const Point &p1, const Point &p2)
Definition: ncmesh.hpp:570
FineTransform * GetFineTransforms()
Definition: ncmesh.cpp:2280
Point(const Point &p0, const Point &p1)
Definition: ncmesh.hpp:534
static const int NumGeom
Definition: geom.hpp:34
void PrintVertexParents(std::ostream &out) const
I/O: Print the &quot;vertex_parents&quot; section of the mesh file (ver. &gt;= 1.1).
Definition: ncmesh.cpp:2568
Node * PeekAltParents(Node *v1, Node *v2)
Definition: ncmesh.cpp:406
void PrintMemoryDetail() const
Definition: ncmesh.cpp:2806
int index
element number in the Mesh, -1 if refined
Definition: ncmesh.hpp:349
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:573
virtual int GetNumGhosts() const
Definition: ncmesh.hpp:400
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4, const Point &p5, const Point &p6, const Point &p7)
Definition: ncmesh.hpp:576
Lists all edges/faces in the nonconforming mesh.
Definition: ncmesh.hpp:124
int attribute
boundary element attribute, -1 if internal edge (2D)
Definition: ncmesh.hpp:277
bool IsIdentity() const
As an optimization, identity transform is &quot;stored&quot; as empty matrix.
Definition: ncmesh.hpp:160
bool Boundary() const
Definition: ncmesh.hpp:282
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
void RefElementNodes(Element *elem)
Definition: ncmesh.cpp:293
char geom
Geometry::Type of the element.
Definition: ncmesh.hpp:347
static int find_node(Element *elem, Node *node)
Definition: ncmesh.cpp:1479
int FaceSplitType(Node *v1, Node *v2, Node *v3, Node *v4, Node *mid[4]=NULL) const
Definition: ncmesh.cpp:1453
virtual void ElementSharesFace(Element *elem, Face *face)
Definition: ncmesh.hpp:482
int index
vertex number in the Mesh
Definition: ncmesh.hpp:265
MeshId(int index=-1, Element *element=NULL, int local=-1)
Definition: ncmesh.hpp:99
void CountObjects(int &nelem, int &nvert, int &nedges) const
Definition: ncmesh.cpp:2769
const Point & operator()(int i) const
Definition: ncmesh.hpp:586
std::vector< MeshId > conforming
Definition: ncmesh.hpp:126
Point & operator()(int i)
Definition: ncmesh.hpp:585
Element * parent
parent element, NULL if this is a root element
Definition: ncmesh.hpp:357
Slave(int index)
Definition: ncmesh.hpp:120
double coord[3]
Definition: ncmesh.hpp:520
Vertex * NewVertex(Node *v1, Node *v2)
Definition: ncmesh.cpp:551
bool NodeSetY1(Node *node, Node **n)
Definition: ncmesh.cpp:607
void ClearCoarseLevel()
Free the internally stored array of coarse leaf elements.
Definition: ncmesh.hpp:169
int spaceDim
dimensions of the elements and the vertex coordinates
Definition: ncmesh.hpp:232
Array< int > vertex_nodeId
Definition: ncmesh.hpp:381
int attribute
boundary element attribute, -1 if internal face
Definition: ncmesh.hpp:322
int index
Mesh element number.
Definition: ncmesh.hpp:34
int master
master number (in Mesh numbering)
Definition: ncmesh.hpp:117
Element * element
NCMesh::Element containing this vertex/edge/face.
Definition: ncmesh.hpp:97
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition: ncmesh.hpp:137
virtual void OnMeshUpdated(Mesh *mesh)
Definition: ncmesh.cpp:1426
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:66
Element * NewQuadrilateral(Node *n0, Node *n1, Node *n2, Node *n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
Definition: ncmesh.cpp:500
void DebugNeighbors(Array< char > &elem_set)
Definition: ncmesh.cpp:1995
int coarse_index
coarse Mesh element index
Definition: ncmesh.hpp:156
void CollectEdgeVertices(Node *v0, Node *v1, Array< int > &indices)
Definition: ncmesh.cpp:1771
virtual bool IsGhost(const Element *elem) const
Definition: ncmesh.hpp:399
bool Empty() const
Definition: ncmesh.hpp:132
Point(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:543
Vertex(double x, double y, double z)
Definition: ncmesh.hpp:269
int num_vertices
width of the table
Definition: ncmesh.hpp:390
Refinement(int index, int type=7)
Definition: ncmesh.hpp:37
int index
face number in the Mesh
Definition: ncmesh.hpp:323
Table element_vertex
leaf-element to vertex table, see FindSetNeighbors
Definition: ncmesh.hpp:389
int local
local number within &#39;element&#39;
Definition: ncmesh.hpp:96
Element * elem[2]
up to 2 elements sharing the face
Definition: ncmesh.hpp:324
void GetMatrix(DenseMatrix &point_matrix) const
Definition: ncmesh.cpp:2010
void CheckAnisoFace(Node *v1, Node *v2, Node *v3, Node *v4, Node *mid12, Node *mid34, int level=0)
Definition: ncmesh.cpp:654
void RegisterFaces(Element *elem, int *fattr=NULL)
Definition: ncmesh.cpp:378
std::vector< Master > masters
Definition: ncmesh.hpp:127
virtual void UpdateVertices()
update Vertex::index and vertex_nodeId
Definition: ncmesh.cpp:1291
void RegisterElement(Element *e)
Definition: ncmesh.cpp:364
DenseMatrix point_matrix
for use in IsoparametricTransformation
Definition: ncmesh.hpp:157
Element * CopyHierarchy(Element *elem)
Definition: ncmesh.cpp:185
void ReorderFacePointMat(Node *v0, Node *v1, Node *v2, Node *v3, Element *elem, DenseMatrix &mat) const
Definition: ncmesh.cpp:1513
void RefineElement(Element *elem, char ref_type)
Definition: ncmesh.cpp:725
void FaceSplitLevel(Node *v1, Node *v2, Node *v3, Node *v4, int &h_level, int &v_level) const
Definition: ncmesh.cpp:2450
Element(int geom, int attr)
Definition: ncmesh.cpp:460
virtual void BuildEdgeList()
Definition: ncmesh.cpp:1693
void LoadCoarseElements(std::istream &input)
I/O: Load the element refinement hierarchy from a mesh file.
Definition: ncmesh.cpp:2677
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:93
void CollectLeafElements(Element *elem)
Definition: ncmesh.cpp:1309
void MarkCoarseLevel()
Definition: ncmesh.cpp:2268
virtual void ElementSharesEdge(Element *elem, Edge *edge)
Definition: ncmesh.hpp:481
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:86
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Definition: ncmesh.hpp:63
bool NodeSetX1(Node *node, Node **n)
Definition: ncmesh.cpp:601
virtual ~NCMesh()
Definition: ncmesh.cpp:241
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
Definition: ncmesh.cpp:2387
static GeomInfo GI[Geometry::NumGeom]
Definition: ncmesh.hpp:632
void DeleteHierarchy(Element *elem)
Definition: ncmesh.cpp:211
int slaves_end
slave faces
Definition: ncmesh.hpp:107
Vertex * vertex
Definition: ncmesh.hpp:297
void TraverseEdge(Node *v0, Node *v1, double t0, double t1, int level)
Definition: ncmesh.cpp:1666
bool Iso
true if the mesh only contains isotropic refinements
Definition: ncmesh.hpp:233
Point(double x, double y, double z)
Definition: ncmesh.hpp:529
double * Data() const
Returns vector of the elements.
Definition: densemat.hpp:77
void DerefineElement(Element *elem)
Definition: ncmesh.cpp:1179
Array< Node * > boundary_edges
subset of all edges, set by BuildEdgeList
Definition: ncmesh.hpp:387
void UnrefEdge(HashTable< Node > &nodes)
Definition: ncmesh.cpp:271
void FindSetNeighbors(const Array< char > &elem_set, Array< Element * > *neighbors, Array< char > *neighbor_set=NULL)
Definition: ncmesh.cpp:1881
int SpaceDimension() const
Definition: ncmesh.hpp:79
std::vector< Slave > slaves
Definition: ncmesh.hpp:128
virtual void BuildFaceList()
Definition: ncmesh.cpp:1591
void ForceRefinement(Node *v1, Node *v2, Node *v3, Node *v4)
Definition: ncmesh.cpp:620
void SetVertexPositions(const Array< mfem::Vertex > &vertices)
I/O: Set positions of all vertices (used by mesh loader).
Definition: ncmesh.cpp:2616
bool Boundary() const
Definition: ncmesh.hpp:329
double pos[3]
3D position
Definition: ncmesh.hpp:266
Element * child[8]
2-8 children (if ref_type != 0)
Definition: ncmesh.hpp:355
Array< Element * > leaf_elements
Definition: ncmesh.hpp:379
void UpdateLeafElements()
Definition: ncmesh.cpp:1324
Vertex(const Vertex &other)
Definition: ncmesh.hpp:271
bool NodeSetZ1(Node *node, Node **n)
Definition: ncmesh.cpp:613
Node * GetMidEdgeVertexSimple(Node *v1, Node *v2)
Definition: ncmesh.cpp:574
long MemoryUsage() const
Definition: ncmesh.cpp:2754
void BuildElementToVertexTable()
Definition: ncmesh.cpp:1807
HashTable< Node > nodes
Definition: ncmesh.hpp:365
Edge(const Edge &other)
Definition: ncmesh.hpp:281
void UnrefElementNodes(Element *elem)
Definition: ncmesh.cpp:322
Master(int index, Element *element, int local, int sb, int se)
Definition: ncmesh.hpp:109
void CheckIsoFace(Node *v1, Node *v2, Node *v3, Node *v4, Node *e1, Node *e2, Node *e3, Node *e4, Node *midf)
Definition: ncmesh.cpp:707
Node * GetMidFaceVertex(Node *e1, Node *e2, Node *e3, Node *e4)
Definition: ncmesh.cpp:583
Element(const Element &other)
Definition: ncmesh.hpp:360
void CollectFaceVertices(Node *v0, Node *v1, Node *v2, Node *v3, Array< int > &indices)
Definition: ncmesh.cpp:1783
bool NodeSetZ2(Node *node, Node **n)
Definition: ncmesh.cpp:616
void Initialize(const mfem::Element *elem)
Definition: ncmesh.cpp:28
bool NodeSetX2(Node *node, Node **n)
Definition: ncmesh.cpp:604
void GetMeshComponents(Array< mfem::Vertex > &vertices, Array< mfem::Element * > &elements, Array< mfem::Element * > &boundary) const
Return the basic Mesh arrays for the current finest level.
Definition: ncmesh.cpp:1344
Node * node[8]
element corners (if ref_type == 0)
Definition: ncmesh.hpp:354
Element * NewTriangle(Node *n0, Node *n1, Node *n2, int attr, int eattr0, int eattr1, int eattr2)
Definition: ncmesh.cpp:527
Point(double x, double y)
Definition: ncmesh.hpp:524
virtual void AssignLeafIndices()
Definition: ncmesh.cpp:1335
static int find_hex_face(int a, int b, int c)
Definition: ncmesh.cpp:1497
int Dimension() const
Definition: ncmesh.hpp:78
Array< ElemRefType > ref_stack
stack of scheduled refinements (temporary)
Definition: ncmesh.hpp:414
static void find_face_nodes(const Face *face, Node *node[4])
Definition: ncmesh.cpp:2367
Array< Element * > root_elements
Definition: ncmesh.hpp:363
virtual void Update()
Definition: ncmesh.cpp:174
long MemoryUsage() const
Return total number of bytes allocated.
Definition: ncmesh.cpp:2783
Node * GetMidEdgeVertex(Node *v1, Node *v2)
Definition: ncmesh.cpp:565
int EdgeSplitLevel(Node *v1, Node *v2) const
Definition: ncmesh.cpp:2443
Element * GetSingleElement() const
Definition: ncmesh.cpp:392
void PrintCoarseElements(std::ostream &out) const
I/O: Print the &quot;coarse_elements&quot; section of the mesh file (ver. &gt;= 1.1).
Definition: ncmesh.cpp:2659
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition: ncmesh.hpp:144
void TraverseFace(Node *v0, Node *v1, Node *v2, Node *v3, const PointMatrix &pm, int level)
Definition: ncmesh.cpp:1544
void UpdateElementToVertexTable()
Definition: ncmesh.hpp:509
void CountSplits(Element *elem, int splits[3]) const
Definition: ncmesh.cpp:2483
void ForgetElement(Element *e)
Definition: ncmesh.cpp:371
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
Definition: ncmesh.hpp:348
ElemRefType(Element *elem, int type)
Definition: ncmesh.hpp:410
bool NodeSetY2(Node *node, Node **n)
Definition: ncmesh.cpp:610
int index
Mesh number.
Definition: ncmesh.hpp:95
void LoadVertexParents(std::istream &input)
Definition: ncmesh.cpp:2596
int CountElements(Element *elem) const
Definition: ncmesh.cpp:2741
Abstract data type element.
Definition: element.hpp:27
Point & operator=(const Point &src)
Definition: ncmesh.hpp:554
Array< Face * > boundary_faces
subset of all faces, set by BuildFaceList
Definition: ncmesh.hpp:386
int rank
processor number (ParNCMesh)
Definition: ncmesh.hpp:350
DenseMatrix point_matrix
position within the master
Definition: ncmesh.hpp:118
HashTable< Face > faces
Definition: ncmesh.hpp:366
virtual void Refine(const Array< Refinement > &refinements)
Definition: ncmesh.cpp:1139
int GetEdgeMaster(int v1, int v2) const
Definition: ncmesh.cpp:2358
void UnrefVertex(HashTable< Node > &nodes)
Definition: ncmesh.cpp:264