MFEM  v3.4
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 "../config/config.hpp"
16 #include "../general/hash.hpp"
17 #include "../general/globals.hpp"
18 #include "../linalg/densemat.hpp"
19 #include "element.hpp"
20 #include "vertex.hpp"
21 #include "../fem/geom.hpp"
22 
23 #include <vector>
24 #include <map>
25 #include <iostream>
26 
27 namespace mfem
28 {
29 
30 /** Represents the index of an element to refine, plus a refinement type.
31  The refinement type is needed for anisotropic refinement of quads and hexes.
32  Bits 0,1 and 2 of 'ref_type' specify whether the element should be split
33  in the X, Y and Z directions, respectively (Z is ignored for quads). */
34 struct Refinement
35 {
36  int index; ///< Mesh element number
37  char ref_type; ///< refinement XYZ bit mask (7 = full isotropic)
38 
39  Refinement(int index, int type = 7) : index(index), ref_type(type) {}
40 };
41 
42 /// Defines the position of a fine element within a coarse element.
43 struct Embedding
44 {
45  int parent; ///< element index in the coarse mesh
46  int matrix; ///< index into CoarseFineTransformations::point_matrices
47 
48  Embedding(int elem, int matrix = 0) : parent(elem), matrix(matrix) {}
49 };
50 
51 /// Defines the coarse-fine transformations of all fine elements.
53 {
54  DenseTensor point_matrices; ///< matrices for IsoparametricTransformation
55  Array<Embedding> embeddings; ///< fine element positions in their parents
56 
57  void Clear() { point_matrices.Clear(); embeddings.DeleteAll(); }
58  long MemoryUsage() const;
59 };
60 
61 
62 /** \brief A class for non-conforming AMR on higher-order hexahedral,
63  * quadrilateral or triangular meshes.
64  *
65  * The class is used as follows:
66  *
67  * 1. NCMesh is constructed from elements of an existing Mesh. The elements
68  * are copied and become roots of the refinement hierarchy.
69  *
70  * 2. Some elements are refined with the Refine() method. Both isotropic and
71  * anisotropic refinements of quads/hexes are supported.
72  *
73  * 3. A new Mesh is created from NCMesh containing the leaf elements.
74  * This new mesh may have non-conforming (hanging) edges and faces.
75  *
76  * 4. FiniteElementSpace asks NCMesh for a list of conforming, master and
77  * slave edges/faces and creates the conforming interpolation matrix P.
78  *
79  * 5. A continuous/conforming solution is obtained by solving P'*A*P x = P'*b.
80  *
81  * 6. Repeat from step 2.
82  */
83 class NCMesh
84 {
85 public:
86  /** Initialize with elements from 'mesh'. If an already nonconforming mesh
87  is being loaded, 'vertex_parents' must point to a stream at the appropriate
88  section of the mesh file which contains the vertex hierarchy. */
89  explicit NCMesh(const Mesh *mesh, std::istream *vertex_parents = NULL);
90 
91  NCMesh(const NCMesh &other); // deep copy
92 
93  virtual ~NCMesh();
94 
95  int Dimension() const { return Dim; }
96  int SpaceDimension() const { return spaceDim; }
97 
98  int GetNVertices() const { return NVertices; }
99  int GetNEdges() const { return NEdges; }
100  int GetNFaces() const { return NFaces; }
101 
102  /** Perform the given batch of refinements. Please note that in the presence
103  of anisotropic splits additional refinements may be necessary to keep
104  the mesh consistent. However, the function always performs at least the
105  requested refinements. */
106  virtual void Refine(const Array<Refinement> &refinements);
107 
108  /** Check the mesh and potentially refine some elements so that the maximum
109  difference of refinement levels between adjacent elements is not greater
110  than 'max_nc_level'. */
111  virtual void LimitNCLevel(int max_nc_level);
112 
113  /** Return a list of derefinement opportunities. Each row of the table
114  contains Mesh indices of existing elements that can be derefined to form
115  a single new coarse element. Row numbers are then passed to Derefine.
116  This function works both in serial and parallel. */
117  const Table &GetDerefinementTable();
118 
119  /** Check derefinements returned by GetDerefinementTable and mark those that
120  can be done safely so that the maximum NC level condition is not violated.
121  On return, level_ok.Size() == deref_table.Size() and contains 0/1s. */
122  virtual void CheckDerefinementNCLevel(const Table &deref_table,
123  Array<int> &level_ok, int max_nc_level);
124 
125  /** Perform a subset of the possible derefinements (see GetDerefinementTable).
126  Note that if anisotropic refinements are present in the mesh, some of the
127  derefinements may have to be skipped to preserve mesh consistency. */
128  virtual void Derefine(const Array<int> &derefs);
129 
130 
131  // master/slave lists
132 
133  /// Identifies a vertex/edge/face in both Mesh and NCMesh.
134  struct MeshId
135  {
136  int index; ///< Mesh number
137  int element; ///< NCMesh::Element containing this vertex/edge/face
138  int local; ///< local number within 'element'
139 
140  MeshId(int index = -1, int element = -1, int local = -1)
142  };
143 
144  /** Nonconforming edge/face that has more than one neighbor. The neighbors
145  are stored in NCList::slaves[i], slaves_begin <= i < slaves_end. */
146  struct Master : public MeshId
147  {
148  int slaves_begin, slaves_end; ///< slave faces
149 
150  Master(int index, int element, int local, int sb, int se)
151  : MeshId(index, element, local), slaves_begin(sb), slaves_end(se) {}
152  };
153 
154  /// Nonconforming edge/face within a bigger edge/face.
155  struct Slave : public MeshId
156  {
157  int master; ///< master number (in Mesh numbering)
158  int edge_flags; ///< edge orientation flags
159  DenseMatrix point_matrix; ///< position within the master edge/face
160 
161  Slave(int index, int element, int local)
162  : MeshId(index, element, local), master(-1), edge_flags(0) {}
163 
164  /// Return the point matrix oriented according to the master and slave edges
165  void OrientedPointMatrix(DenseMatrix &oriented_matrix) const;
166  };
167 
168  /// Lists all edges/faces in the nonconforming mesh.
169  struct NCList
170  {
171  std::vector<MeshId> conforming;
172  std::vector<Master> masters;
173  std::vector<Slave> slaves;
174  // TODO: switch to Arrays when fixed for non-POD types
175  // TODO: make a list of unique slave matrices to save memory (+ time later)
176 
177  void Clear(bool hard = false);
178  bool Empty() const { return !conforming.size() && !masters.size(); }
179  long TotalSize() const;
180  long MemoryUsage() const;
181 
182  const MeshId& LookUp(int index, int *type = NULL) const;
183  private:
184  mutable Array<int> inv_index;
185  };
186 
187  /// Return the current list of conforming and nonconforming faces.
189  {
190  if (face_list.Empty()) { BuildFaceList(); }
191  return face_list;
192  }
193 
194  /// Return the current list of conforming and nonconforming edges.
196  {
197  if (edge_list.Empty()) { BuildEdgeList(); }
198  return edge_list;
199  }
200 
201  /** Return a list of vertices (in 'conforming'); this function is provided
202  for uniformity/completeness. Needed in ParNCMesh/ParFESpace. */
204  {
205  if (vertex_list.Empty()) { BuildVertexList(); }
206  return vertex_list;
207  }
208 
209  /// Return vertex/edge/face list (entity = 0/1/2, respectively).
210  const NCList& GetNCList(int entity)
211  {
212  switch (entity)
213  {
214  case 0: return GetVertexList();
215  case 1: return GetEdgeList();
216  default: return GetFaceList();
217  }
218  }
219 
220 
221  // coarse/fine transforms
222 
223  /** Remember the current layer of leaf elements before the mesh is refined.
224  Needed by GetRefinementTransforms(), must be called before Refine(). */
225  void MarkCoarseLevel();
226 
227  /** After refinement, calculate the relation of each fine element to its
228  parent coarse element. Note that Refine() or LimitNCLevel() can be called
229  multiple times between MarkCoarseLevel() and this function. */
231 
232  /** After derefinement, calculate the relations of previous fine elements
233  (some of which may no longer exist) to the current leaf elements.
234  Unlike for refinement, Derefine() may only be called once before this
235  function so there is no MarkFineLevel(). */
237 
238  /// Free all internal data created by the above three functions.
239  void ClearTransforms();
240 
241 
242  // utility
243 
244  /// Return Mesh vertex indices of an edge identified by 'edge_id'.
245  void GetEdgeVertices(const MeshId &edge_id, int vert_index[2]) const;
246 
247  /** Return "NC" orientation of an edge. As opposed standard Mesh edge
248  orientation based on vertex IDs, "NC" edge orientation follows the local
249  edge orientation within the element 'edge_id.element' and is thus
250  processor independent. */
251  int GetEdgeNCOrientation(const MeshId &edge_id) const;
252 
253  /// Return Mesh vertex and edge indices of a face identified by 'face_id'.
254  void GetFaceVerticesEdges(const MeshId &face_id,
255  int vert_index[4], int edge_index[4],
256  int edge_orientation[4]) const;
257 
258  /** Given an edge (by its vertex indices v1 and v2) return the first
259  (geometric) parent edge that exists in the Mesh or -1 if there is no such
260  parent. */
261  int GetEdgeMaster(int v1, int v2) const;
262 
263  /** Get a list of vertices (2D/3D) and edges (3D) that coincide with boundary
264  elements with the specified attributes (marked in 'bdr_attr_is_ess').
265  In 3D this function also reveals "hidden" boundary edges. In parallel it
266  helps identifying boundary vertices/edges affected by non-local boundary
267  elements. */
268  virtual void GetBoundaryClosure(const Array<int> &bdr_attr_is_ess,
269  Array<int> &bdr_vertices,
270  Array<int> &bdr_edges);
271 
272  /// Return the type of elements in the mesh.
273  int GetElementGeometry() const { return elements[0].geom; }
274 
275  int GetFaceGeometry() const { return Geometry::SQUARE; }
276 
277  /// Return the distance of leaf 'i' from the root.
278  int GetElementDepth(int i) const;
279 
280  /// I/O: Print the "vertex_parents" section of the mesh file (ver. >= 1.1).
281  void PrintVertexParents(std::ostream &out) const;
282 
283  /// I/O: Print the "coarse_elements" section of the mesh file (ver. >= 1.1).
284  void PrintCoarseElements(std::ostream &out) const;
285 
286  /** I/O: Load the vertex parent hierarchy from a mesh file. NOTE: called
287  indirectly through the constructor. */
288  void LoadVertexParents(std::istream &input);
289 
290  /// I/O: Load the element refinement hierarchy from a mesh file.
291  void LoadCoarseElements(std::istream &input);
292 
293  /// I/O: Set positions of all vertices (used by mesh loader).
294  void SetVertexPositions(const Array<mfem::Vertex> &vertices);
295 
296  /// Save memory by releasing all non-essential and cached data.
297  virtual void Trim();
298 
299  /// Return total number of bytes allocated.
300  long MemoryUsage() const;
301 
302  int PrintMemoryDetail() const;
303 
304  void PrintStats(std::ostream &out = mfem::out) const;
305 
306 
307 protected: // interface for Mesh to be able to construct itself from NCMesh
308 
309  friend class Mesh;
310 
311  /// Return the basic Mesh arrays for the current finest level.
312  void GetMeshComponents(Array<mfem::Vertex>& mvertices,
313  Array<mfem::Element*>& melements,
314  Array<mfem::Element*>& mboundary) const;
315 
316  /** Get edge and face numbering from 'mesh' (i.e., set all Edge::index and
317  Face::index) after a new mesh was created from us. */
318  virtual void OnMeshUpdated(Mesh *mesh);
319 
320 
321 protected: // implementation
322 
323  int Dim, spaceDim; ///< dimensions of the elements and the vertex coordinates
324  bool Iso; ///< true if the mesh only contains isotropic refinements
325 
326  /** A Node can hold a vertex, an edge, or both. Elements directly point to
327  their corner nodes, but edge nodes also exist and can be accessed using
328  a hash-table given their two end-point node IDs. All nodes can be
329  accessed in this way, with the exception of top-level vertex nodes.
330  When an element is being refined, the mid-edge nodes are readily
331  available with this mechanism. The new elements "sign in" to the nodes
332  by increasing the reference counts of their vertices and edges. The
333  parent element "signs off" its nodes by decrementing the ref counts. */
334  struct Node : public Hashed2
335  {
338 
339  Node() : vert_refc(0), edge_refc(0), vert_index(-1), edge_index(-1) {}
340  ~Node();
341 
342  bool HasVertex() const { return vert_refc > 0; }
343  bool HasEdge() const { return edge_refc > 0; }
344 
345  // decrease vertex/edge ref count, return false if Node should be deleted
346  bool UnrefVertex() { --vert_refc; return vert_refc || edge_refc; }
347  bool UnrefEdge() { --edge_refc; return vert_refc || edge_refc; }
348  };
349 
350  /** Similarly to nodes, faces can be accessed by hashing their four vertex
351  node IDs. A face knows about the one or two elements that are using it.
352  A face that is not on the boundary and only has one element referencing
353  it is either a master or a slave face. */
354  struct Face : public Hashed4
355  {
356  int attribute; ///< boundary element attribute, -1 if internal face
357  int index; ///< face number in the Mesh
358  int elem[2]; ///< up to 2 elements sharing the face
359 
360  Face() : attribute(-1), index(-1) { elem[0] = elem[1] = -1; }
361 
362  bool Boundary() const { return attribute >= 0; }
363  bool Unused() const { return elem[0] < 0 && elem[1] < 0;}
364 
365  // add or remove an element from the 'elem[2]' array
366  void RegisterElement(int e);
367  void ForgetElement(int e);
368 
369  /// Return one of elem[0] or elem[1] and make sure the other is -1.
370  int GetSingleElement() const;
371  };
372 
373  /** This is an element in the refinement hierarchy. Each element has
374  either been refined and points to its children, or is a leaf and points
375  to its vertex nodes. */
376  struct Element
377  {
378  char geom; ///< Geometry::Type of the element
379  char ref_type; ///< bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
380  char flag; ///< generic flag/marker, can be used by algorithms
381  int index; ///< element number in the Mesh, -1 if refined
382  int rank; ///< processor number (ParNCMesh), -1 if undefined/unknown
384  union
385  {
386  int node[8]; ///< element corners (if ref_type == 0)
387  int child[8]; ///< 2-8 children (if ref_type != 0)
388  };
389  int parent; ///< parent element, -1 if this is a root element, -2 if free
390 
391  Element(int geom, int attr);
392  };
393 
394  // primary data
395 
396  HashTable<Node> nodes; // associative container holding all Nodes
397  HashTable<Face> faces; // associative container holding all Faces
398 
399  BlockArray<Element> elements; // storage for all Elements
400  Array<int> free_element_ids; // unused element ids - indices into 'elements'
401 
402  // the first 'root_count' entries of 'elements' is the coarse mesh
404 
405  // coordinates of top-level vertices (organized as triples)
407 
413 
414 
415  // secondary data
416 
417  /** Apart from the primary data structure, which is the element/node/face
418  hierarchy, there is secondary data that is derived from the primary
419  data and needs to be updated when the primary data changes. Update()
420  takes care of that and needs to be called after refinement and
421  derefinement. */
422  virtual void Update();
423 
424  int NVertices; // set by UpdateVertices
425  int NEdges, NFaces; // set by OnMeshUpdated
426 
427  Array<int> leaf_elements; // finest level, calculated by UpdateLeafElements
428  Array<int> vertex_nodeId; // vertex-index to node-id map, see UpdateVertices
429 
430  NCList face_list; ///< lazy-initialized list of faces, see GetFaceList
431  NCList edge_list; ///< lazy-initialized list of edges, see GetEdgeList
432  NCList vertex_list; ///< lazy-initialized list of vertices, see GetVertexList
433 
434  Array<int> boundary_faces; ///< subset of all faces, set by BuildFaceList
435 
436  Table element_vertex; ///< leaf-element to vertex table, see FindSetNeighbors
437 
438 
439  virtual void UpdateVertices(); ///< update Vertex::index and vertex_nodeId
440 
441  void CollectLeafElements(int elem, int state);
442  void UpdateLeafElements();
443 
444  virtual void AssignLeafIndices();
445 
446  virtual bool IsGhost(const Element &el) const { return false; }
447  virtual int GetNumGhostElements() const { return 0; }
448  virtual int GetNumGhostVertices() const { return 0; }
449 
450 
451  // refinement/derefinement
452 
453  Array<Refinement> ref_stack; ///< stack of scheduled refinements (temporary)
454 
455  Table derefinements; ///< possible derefinements, see GetDerefinementTable
456 
457  void RefineElement(int elem, char ref_type);
458  void DerefineElement(int elem);
459 
460  int AddElement(const Element &el)
461  {
462  if (free_element_ids.Size())
463  {
464  int idx = free_element_ids.Last();
466  elements[idx] = el;
467  return idx;
468  }
469  return elements.Append(el);
470  }
471  void FreeElement(int id)
472  {
474  elements[id].parent = -2; // mark the element as free
475  }
476 
477  int NewHexahedron(int n0, int n1, int n2, int n3,
478  int n4, int n5, int n6, int n7,
479  int attr,
480  int fattr0, int fattr1, int fattr2,
481  int fattr3, int fattr4, int fattr5);
482 
483  int NewQuadrilateral(int n0, int n1, int n2, int n3,
484  int attr,
485  int eattr0, int eattr1, int eattr2, int eattr3);
486 
487  int NewTriangle(int n0, int n1, int n2,
488  int attr, int eattr0, int eattr1, int eattr2);
489 
490  mfem::Element* NewMeshElement(int geom) const;
491 
492  int GetMidEdgeNode(int vn1, int vn2);
493  int GetMidFaceNode(int en1, int en2, int en3, int en4);
494 
495  int FaceSplitType(int v1, int v2, int v3, int v4, int mid[4]
496  = NULL /*optional output of mid-edge nodes*/) const;
497 
498  void ForceRefinement(int vn1, int vn2, int vn3, int vn4);
499 
500  void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4,
501  int mid12, int mid34, int level = 0);
502 
503  void CheckIsoFace(int vn1, int vn2, int vn3, int vn4,
504  int en1, int en2, int en3, int en4, int midf);
505 
506  void RefElement(int elem);
507  void UnrefElement(int elem, Array<int> &elemFaces);
508 
509  Face* GetFace(Element &elem, int face_no);
510  void RegisterFaces(int elem, int *fattr = NULL);
511  void DeleteUnusedFaces(const Array<int> &elemFaces);
512 
513  int FindAltParents(int node1, int node2);
514 
515  bool NodeSetX1(int node, int* n);
516  bool NodeSetX2(int node, int* n);
517  bool NodeSetY1(int node, int* n);
518  bool NodeSetY2(int node, int* n);
519  bool NodeSetZ1(int node, int* n);
520  bool NodeSetZ2(int node, int* n);
521 
522  void CollectDerefinements(int elem, Array<Connection> &list);
523 
524 
525  // face/edge lists
526 
527  static int find_node(const Element &el, int node);
528  static int find_element_edge(const Element &el, int vn0, int vn1);
529  static int find_hex_face(int a, int b, int c);
530 
531  int ReorderFacePointMat(int v0, int v1, int v2, int v3,
532  int elem, DenseMatrix& mat) const;
533  struct PointMatrix;
534  void TraverseFace(int vn0, int vn1, int vn2, int vn3,
535  const PointMatrix& pm, int level);
536 
537  void TraverseEdge(int vn0, int vn1, double t0, double t1, int flags,
538  int level);
539 
540  virtual void BuildFaceList();
541  virtual void BuildEdgeList();
542  virtual void BuildVertexList();
543 
544  virtual void ElementSharesFace(int elem, int face) {} // ParNCMesh
545  virtual void ElementSharesEdge(int elem, int enode) {} // ParNCMesh
546  virtual void ElementSharesVertex(int elem, int vnode) {} // ParNCMesh
547 
548 
549  // neighbors / element_vertex table
550 
551  /** Return all vertex-, edge- and face-neighbors of a set of elements.
552  The neighbors are returned as a list (neighbors != NULL), as a set
553  (neighbor_set != NULL), or both. The sizes of the set arrays must match
554  that of leaf_elements. The function is intended to be used for large
555  sets of elements and its complexity is linear in the number of leaf
556  elements in the mesh. */
557  void FindSetNeighbors(const Array<char> &elem_set,
558  Array<int> *neighbors, /* append */
559  Array<char> *neighbor_set = NULL);
560 
561  /** Return all vertex-, edge- and face-neighbors of a single element.
562  You can limit the number of elements being checked using 'search_set'.
563  The complexity of the function is linear in the size of the search set.*/
564  void FindNeighbors(int elem,
565  Array<int> &neighbors, /* append */
566  const Array<int> *search_set = NULL);
567 
568  /** Expand a set of elements by all vertex-, edge- and face-neighbors.
569  The output array 'expanded' will contain all items from 'elems'
570  (provided they are in 'search_set') plus their neighbors. The neighbor
571  search can be limited to the optional search set. The complexity is
572  linear in the sum of the sizes of 'elems' and 'search_set'. */
573  void NeighborExpand(const Array<int> &elems,
574  Array<int> &expanded,
575  const Array<int> *search_set = NULL);
576 
577 
578  void CollectEdgeVertices(int v0, int v1, Array<int> &indices);
579  void CollectFaceVertices(int v0, int v1, int v2, int v3,
580  Array<int> &indices);
582 
584  {
586  }
587 
588 
589  // coarse/fine transformations
590 
591  struct Point
592  {
593  int dim;
594  double coord[3];
595 
596  Point() { dim = 0; }
597 
598  Point(double x, double y)
599  { dim = 2; coord[0] = x; coord[1] = y; }
600 
601  Point(double x, double y, double z)
602  { dim = 3; coord[0] = x; coord[1] = y; coord[2] = z; }
603 
604  Point(const Point& p0, const Point& p1)
605  {
606  dim = p0.dim;
607  for (int i = 0; i < dim; i++)
608  {
609  coord[i] = (p0.coord[i] + p1.coord[i]) * 0.5;
610  }
611  }
612 
613  Point(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
614  {
615  dim = p0.dim;
616  MFEM_ASSERT(p1.dim == dim && p2.dim == dim && p3.dim == dim, "");
617  for (int i = 0; i < dim; i++)
618  {
619  coord[i] = (p0.coord[i] + p1.coord[i] + p2.coord[i] + p3.coord[i])
620  * 0.25;
621  }
622  }
623 
624  Point& operator=(const Point& src)
625  {
626  dim = src.dim;
627  for (int i = 0; i < dim; i++) { coord[i] = src.coord[i]; }
628  return *this;
629  }
630  };
631 
632  struct PointMatrix
633  {
634  int np;
636 
637  PointMatrix(const Point& p0, const Point& p1, const Point& p2)
638  { np = 3; points[0] = p0; points[1] = p1; points[2] = p2; }
639 
640  PointMatrix(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
641  { np = 4; points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3; }
642 
643  PointMatrix(const Point& p0, const Point& p1, const Point& p2,
644  const Point& p3, const Point& p4, const Point& p5,
645  const Point& p6, const Point& p7)
646  {
647  np = 8;
648  points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3;
649  points[4] = p4; points[5] = p5; points[6] = p6; points[7] = p7;
650  }
651 
652  Point& operator()(int i) { return points[i]; }
653  const Point& operator()(int i) const { return points[i]; }
654 
655  void GetMatrix(DenseMatrix& point_matrix) const;
656  };
657 
661 
662  static const PointMatrix& GetGeomIdentity(int geom);
663 
664  void GetPointMatrix(int geom, const char* ref_path, DenseMatrix& matrix);
665 
666  typedef std::map<std::string, int> RefPathMap;
667 
668  void TraverseRefinements(int elem, int coarse_index,
669  std::string &ref_path, RefPathMap &map);
670 
671  /// storage for data returned by Get[De]RefinementTransforms()
673 
674  /// state of leaf_elements before Refine(), set by MarkCoarseLevel()
676 
677  void InitDerefTransforms();
678  void SetDerefMatrixCodes(int parent, Array<int> &fine_coarse);
679 
680 
681  // vertex temporary data, used by GetMeshComponents
682 
683  struct TmpVertex
684  {
685  bool valid, visited;
686  double pos[3];
687  TmpVertex() : valid(false), visited(false) {}
688  };
689 
691 
692  const double *CalcVertexPos(int node) const;
693 
694 
695  // utility
696 
697  int GetEdgeMaster(int node) const;
698 
699  void FindFaceNodes(int face, int node[4]);
700 
701  int EdgeSplitLevel(int vn1, int vn2) const;
702  void FaceSplitLevel(int vn1, int vn2, int vn3, int vn4,
703  int& h_level, int& v_level) const;
704 
705  void CountSplits(int elem, int splits[3]) const;
706  void GetLimitRefinements(Array<Refinement> &refinements, int max_level);
707 
708  int PrintElements(std::ostream &out, int elem, int &coarse_id) const;
709  void CopyElements(int elem, const BlockArray<Element> &tmp_elements,
710  Array<int> &index_map);
711 
712  // geometry
713 
714  /** This holds in one place the constants about the geometries we support
715  (triangles, quads, cubes) */
716  struct GeomInfo
717  {
718  int nv, ne, nf, nfv; // number of: vertices, edges, faces, face vertices
719  int edges[12][2]; // edge vertices (up to 12 edges)
720  int faces[6][4]; // face vertices (up to 6 faces)
721 
723  GeomInfo() : initialized(false) {}
724  void Initialize(const mfem::Element* elem);
725  };
726 
728 
730 
731 #ifdef MFEM_DEBUG
732 public:
733  void DebugDump(std::ostream &out) const;
734 #endif
735 
736  friend class ParNCMesh; // for ParNCMesh::ElementSet
737  friend struct CompareRanks;
738 };
739 
740 }
741 
742 #endif
NCList face_list
lazy-initialized list of faces, see GetFaceList
Definition: ncmesh.hpp:430
void CollectLeafElements(int elem, int state)
Definition: ncmesh.cpp:1423
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
Definition: ncmesh.hpp:431
NCMesh(const Mesh *mesh, std::istream *vertex_parents=NULL)
Definition: ncmesh.cpp:67
int Size() const
Logical size of the array.
Definition: array.hpp:133
CoarseFineTransformations transforms
storage for data returned by Get[De]RefinementTransforms()
Definition: ncmesh.hpp:672
int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
Definition: ncmesh.cpp:464
int elem[2]
up to 2 elements sharing the face
Definition: ncmesh.hpp:358
char ref_type
refinement XYZ bit mask (7 = full isotropic)
Definition: ncmesh.hpp:37
TmpVertex * tmp_vertex
Definition: ncmesh.hpp:690
Array< double > top_vertex_pos
Definition: ncmesh.hpp:406
void SetDerefMatrixCodes(int parent, Array< int > &fine_coarse)
Definition: ncmesh.cpp:1353
const CoarseFineTransformations & GetDerefinementTransforms()
Definition: ncmesh.cpp:2915
char flag
generic flag/marker, can be used by algorithms
Definition: ncmesh.hpp:380
bool HasEdge() const
Definition: ncmesh.hpp:343
PointMatrix(const Point &p0, const Point &p1, const Point &p2)
Definition: ncmesh.hpp:637
Point(const Point &p0, const Point &p1)
Definition: ncmesh.hpp:604
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:3292
void GetPointMatrix(int geom, const char *ref_path, DenseMatrix &matrix)
Definition: ncmesh.cpp:2537
virtual void LimitNCLevel(int max_nc_level)
Definition: ncmesh.cpp:3277
void GetFaceVerticesEdges(const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
Return Mesh vertex and edge indices of a face identified by &#39;face_id&#39;.
Definition: ncmesh.cpp:2998
virtual void Trim()
Save memory by releasing all non-essential and cached data.
Definition: ncmesh.cpp:3514
int index
element number in the Mesh, -1 if refined
Definition: ncmesh.hpp:381
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:640
const Geometry::Type geom
Definition: ex1.cpp:40
void CopyElements(int elem, const BlockArray< Element > &tmp_elements, Array< int > &index_map)
Definition: ncmesh.cpp:3404
void ForceRefinement(int vn1, int vn2, int vn3, int vn4)
Definition: ncmesh.cpp:549
virtual int GetNumGhostElements() const
Definition: ncmesh.hpp:447
int GetElementGeometry() const
Return the type of elements in the mesh.
Definition: ncmesh.hpp:273
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:643
Lists all edges/faces in the nonconforming mesh.
Definition: ncmesh.hpp:169
void GetMeshComponents(Array< mfem::Vertex > &mvertices, Array< mfem::Element * > &melements, Array< mfem::Element * > &mboundary) const
Return the basic Mesh arrays for the current finest level.
Definition: ncmesh.cpp:1528
void TraverseEdge(int vn0, int vn1, double t0, double t1, int flags, int level)
Definition: ncmesh.cpp:1874
int NVertices
Definition: ncmesh.hpp:424
int GetNVertices() const
Definition: ncmesh.hpp:98
void ClearTransforms()
Free all internal data created by the above three functions.
Definition: ncmesh.cpp:2958
Table derefinements
possible derefinements, see GetDerefinementTable
Definition: ncmesh.hpp:455
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void CollectFaceVertices(int v0, int v1, int v2, int v3, Array< int > &indices)
Definition: ncmesh.cpp:2157
NCList vertex_list
lazy-initialized list of vertices, see GetVertexList
Definition: ncmesh.hpp:432
void InitDerefTransforms()
Definition: ncmesh.cpp:1338
char geom
Geometry::Type of the element.
Definition: ncmesh.hpp:378
int GetFaceGeometry() const
Definition: ncmesh.hpp:275
BlockArray< Element >::iterator elem_iterator
Definition: ncmesh.hpp:412
long TotalSize() const
Definition: ncmesh.cpp:2087
void GetEdgeVertices(const MeshId &edge_id, int vert_index[2]) const
Return Mesh vertex indices of an edge identified by &#39;edge_id&#39;.
Definition: ncmesh.cpp:2968
int GetEdgeNCOrientation(const MeshId &edge_id) const
Definition: ncmesh.cpp:2986
void OrientedPointMatrix(DenseMatrix &oriented_matrix) const
Return the point matrix oriented according to the master and slave edges.
Definition: ncmesh.cpp:2048
static GeomInfo & gi_tri
Definition: ncmesh.hpp:729
void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4, int mid12, int mid34, int level=0)
Definition: ncmesh.cpp:583
virtual int GetNumGhostVertices() const
Definition: ncmesh.hpp:448
static PointMatrix pm_tri_identity
Definition: ncmesh.hpp:658
void FreeElement(int id)
Definition: ncmesh.hpp:471
void CountSplits(int elem, int splits[3]) const
Definition: ncmesh.cpp:3199
void CollectDerefinements(int elem, Array< Connection > &list)
Definition: ncmesh.cpp:1225
const Point & operator()(int i) const
Definition: ncmesh.hpp:653
std::vector< MeshId > conforming
Definition: ncmesh.hpp:171
Point & operator()(int i)
Definition: ncmesh.hpp:652
void GetLimitRefinements(Array< Refinement > &refinements, int max_level)
Definition: ncmesh.cpp:3247
int PrintMemoryDetail() const
Definition: ncmesh.cpp:3566
double coord[3]
Definition: ncmesh.hpp:594
bool UnrefVertex()
Definition: ncmesh.hpp:346
int spaceDim
dimensions of the elements and the vertex coordinates
Definition: ncmesh.hpp:323
Array< int > vertex_nodeId
Definition: ncmesh.hpp:428
int attribute
boundary element attribute, -1 if internal face
Definition: ncmesh.hpp:356
void DebugDump(std::ostream &out) const
Definition: ncmesh.cpp:3644
void FindFaceNodes(int face, int node[4])
Definition: ncmesh.cpp:3077
int index
Mesh element number.
Definition: ncmesh.hpp:36
int master
master number (in Mesh numbering)
Definition: ncmesh.hpp:157
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition: ncmesh.hpp:188
const MeshId & LookUp(int index, int *type=NULL) const
Definition: ncmesh.cpp:2092
virtual void ElementSharesEdge(int elem, int enode)
Definition: ncmesh.hpp:545
virtual void OnMeshUpdated(Mesh *mesh)
Definition: ncmesh.cpp:1603
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:61
void CollectEdgeVertices(int v0, int v1, Array< int > &indices)
Definition: ncmesh.cpp:2145
Array< Refinement > ref_stack
stack of scheduled refinements (temporary)
Definition: ncmesh.hpp:453
bool Empty() const
Definition: ncmesh.hpp:178
Point(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:613
MeshId(int index=-1, int element=-1, int local=-1)
Definition: ncmesh.hpp:140
virtual bool IsGhost(const Element &el) const
Definition: ncmesh.hpp:446
void CheckIsoFace(int vn1, int vn2, int vn3, int vn4, int en1, int en2, int en3, int en4, int midf)
Definition: ncmesh.cpp:649
int root_count
Definition: ncmesh.hpp:403
Face * GetFace(Element &elem, int face_no)
Definition: ncmesh.cpp:344
void TraverseRefinements(int elem, int coarse_index, std::string &ref_path, RefPathMap &map)
Definition: ncmesh.cpp:2848
int PrintElements(std::ostream &out, int elem, int &coarse_id) const
Definition: ncmesh.cpp:3364
Refinement(int index, int type=7)
Definition: ncmesh.hpp:39
int index
face number in the Mesh
Definition: ncmesh.hpp:357
const Table & GetDerefinementTable()
Definition: ncmesh.cpp:1258
Table element_vertex
leaf-element to vertex table, see FindSetNeighbors
Definition: ncmesh.hpp:436
DenseTensor point_matrices
matrices for IsoparametricTransformation
Definition: ncmesh.hpp:54
int local
local number within &#39;element&#39;
Definition: ncmesh.hpp:138
int FindAltParents(int node1, int node2)
Definition: ncmesh.cpp:366
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:614
void GetMatrix(DenseMatrix &point_matrix) const
Definition: ncmesh.cpp:2501
int AddElement(const Element &el)
Definition: ncmesh.hpp:460
bool NodeSetY1(int node, int *n)
Definition: ncmesh.cpp:536
virtual void Derefine(const Array< int > &derefs)
Definition: ncmesh.cpp:1302
void FindNeighbors(int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
Definition: ncmesh.cpp:2356
std::vector< Master > masters
Definition: ncmesh.hpp:172
Array< int > coarse_elements
state of leaf_elements before Refine(), set by MarkCoarseLevel()
Definition: ncmesh.hpp:675
virtual void UpdateVertices()
update Vertex::index and vertex_nodeId
Definition: ncmesh.cpp:1372
void RegisterFaces(int elem, int *fattr=NULL)
Definition: ncmesh.cpp:305
Element(int geom, int attr)
Definition: ncmesh.cpp:423
Array< Embedding > embeddings
fine element positions in their parents
Definition: ncmesh.hpp:55
virtual void BuildEdgeList()
Definition: ncmesh.cpp:1913
const CoarseFineTransformations & GetRefinementTransforms()
Definition: ncmesh.cpp:2878
void LoadCoarseElements(std::istream &input)
I/O: Load the element refinement hierarchy from a mesh file.
Definition: ncmesh.cpp:3424
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:134
int parent
parent element, -1 if this is a root element, -2 if free
Definition: ncmesh.hpp:389
void MarkCoarseLevel()
Definition: ncmesh.cpp:2834
Slave(int index, int element, int local)
Definition: ncmesh.hpp:161
Master(int index, int element, int local, int sb, int se)
Definition: ncmesh.hpp:150
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Definition: ncmesh.hpp:83
virtual ~NCMesh()
Definition: ncmesh.cpp:200
int GetNFaces() const
Definition: ncmesh.hpp:100
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
Definition: ncmesh.cpp:3100
int element
NCMesh::Element containing this vertex/edge/face.
Definition: ncmesh.hpp:137
static GeomInfo & gi_hex
Definition: ncmesh.hpp:729
Embedding(int elem, int matrix=0)
Definition: ncmesh.hpp:48
static GeomInfo GI[Geometry::NumGeom]
Definition: ncmesh.hpp:727
void Clear(bool hard=false)
Definition: ncmesh.cpp:2069
HashTable< Node >::iterator node_iterator
Definition: ncmesh.hpp:408
int slaves_end
slave faces
Definition: ncmesh.hpp:148
bool NodeSetY2(int node, int *n)
Definition: ncmesh.cpp:539
int NewTriangle(int n0, int n1, int n2, int attr, int eattr0, int eattr1, int eattr2)
Definition: ncmesh.cpp:489
void DeleteUnusedFaces(const Array< int > &elemFaces)
Definition: ncmesh.cpp:319
virtual void ElementSharesFace(int elem, int face)
Definition: ncmesh.hpp:544
bool Iso
true if the mesh only contains isotropic refinements
Definition: ncmesh.hpp:324
Point(double x, double y, double z)
Definition: ncmesh.hpp:601
bool HasVertex() const
Definition: ncmesh.hpp:342
std::map< std::string, int > RefPathMap
Definition: ncmesh.hpp:666
int FaceSplitType(int v1, int v2, int v3, int v4, int mid[4]=NULL) const
Definition: ncmesh.cpp:1652
int SpaceDimension() const
Definition: ncmesh.hpp:96
std::vector< Slave > slaves
Definition: ncmesh.hpp:173
virtual void BuildFaceList()
Definition: ncmesh.cpp:1803
Nonconforming edge/face within a bigger edge/face.
Definition: ncmesh.hpp:155
void DerefineElement(int elem)
Definition: ncmesh.cpp:1119
int ReorderFacePointMat(int v0, int v1, int v2, int v3, int elem, DenseMatrix &mat) const
Definition: ncmesh.cpp:1721
void SetVertexPositions(const Array< mfem::Vertex > &vertices)
I/O: Set positions of all vertices (used by mesh loader).
Definition: ncmesh.cpp:3341
int edge_flags
edge orientation flags
Definition: ncmesh.hpp:158
HashTable< Node >::const_iterator node_const_iterator
Definition: ncmesh.hpp:410
const NCList & GetVertexList()
Definition: ncmesh.hpp:203
bool Boundary() const
Definition: ncmesh.hpp:362
void RegisterElement(int e)
Definition: ncmesh.cpp:330
static const PointMatrix & GetGeomIdentity(int geom)
Definition: ncmesh.cpp:2524
mfem::Element * NewMeshElement(int geom) const
Definition: ncmesh.cpp:1487
void DeleteLast()
Delete the last entry.
Definition: array.hpp:181
void UpdateLeafElements()
Definition: ncmesh.cpp:1464
int GetNEdges() const
Definition: ncmesh.hpp:99
bool NodeSetZ2(int node, int *n)
Definition: ncmesh.cpp:545
static GeomInfo & gi_quad
Definition: ncmesh.hpp:729
long MemoryUsage() const
Definition: ncmesh.cpp:3526
void BuildElementToVertexTable()
Definition: ncmesh.cpp:2181
HashTable< Node > nodes
Definition: ncmesh.hpp:396
bool NodeSetZ1(int node, int *n)
Definition: ncmesh.cpp:542
void ForgetElement(int e)
Definition: ncmesh.cpp:337
int GetMidFaceNode(int en1, int en2, int en3, int en4)
Definition: ncmesh.cpp:521
void RefineElement(int elem, char ref_type)
Definition: ncmesh.cpp:666
int NewHexahedron(int n0, int n1, int n2, int n3, int n4, int n5, int n6, int n7, int attr, int fattr0, int fattr1, int fattr2, int fattr3, int fattr4, int fattr5)
Definition: ncmesh.cpp:435
Array< int > boundary_faces
subset of all faces, set by BuildFaceList
Definition: ncmesh.hpp:434
int child[8]
2-8 children (if ref_type != 0)
Definition: ncmesh.hpp:387
const NCList & GetNCList(int entity)
Return vertex/edge/face list (entity = 0/1/2, respectively).
Definition: ncmesh.hpp:210
Array< int > leaf_elements
Definition: ncmesh.hpp:427
static PointMatrix pm_quad_identity
Definition: ncmesh.hpp:659
void TraverseFace(int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level)
Definition: ncmesh.cpp:1754
int EdgeSplitLevel(int vn1, int vn2) const
Definition: ncmesh.cpp:3159
static int find_element_edge(const Element &el, int vn0, int vn1)
Definition: ncmesh.cpp:1689
void Initialize(const mfem::Element *elem)
Definition: ncmesh.cpp:28
Array< int > free_element_ids
Definition: ncmesh.hpp:400
bool NodeSetX2(int node, int *n)
Definition: ncmesh.cpp:533
int parent
element index in the coarse mesh
Definition: ncmesh.hpp:45
Point(double x, double y)
Definition: ncmesh.hpp:598
virtual void AssignLeafIndices()
Definition: ncmesh.cpp:1478
virtual void ElementSharesVertex(int elem, int vnode)
Definition: ncmesh.hpp:546
static PointMatrix pm_hex_identity
Definition: ncmesh.hpp:660
static int find_hex_face(int a, int b, int c)
Definition: ncmesh.cpp:1705
int Dimension() const
Definition: ncmesh.hpp:95
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
Definition: ncmesh.cpp:1273
T & Last()
Return the last element in the array.
Definition: array.hpp:647
void NeighborExpand(const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
Definition: ncmesh.cpp:2439
BlockArray< Element > elements
Definition: ncmesh.hpp:399
void FindSetNeighbors(const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
Definition: ncmesh.cpp:2253
virtual void Update()
Definition: ncmesh.cpp:188
long MemoryUsage() const
Return total number of bytes allocated.
Definition: ncmesh.cpp:3545
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:3390
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:52
virtual void BuildVertexList()
Definition: ncmesh.cpp:2013
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition: ncmesh.hpp:195
void UpdateElementToVertexTable()
Definition: ncmesh.hpp:583
friend struct CompareRanks
Definition: ncmesh.hpp:737
int GetElementDepth(int i) const
Return the distance of leaf &#39;i&#39; from the root.
Definition: ncmesh.cpp:3065
HashTable< Face >::iterator face_iterator
Definition: ncmesh.hpp:409
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
Definition: ncmesh.hpp:379
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
bool Unused() const
Definition: ncmesh.hpp:363
static int find_node(const Element &el, int node)
Definition: ncmesh.cpp:1679
int index
Mesh number.
Definition: ncmesh.hpp:136
void LoadVertexParents(std::istream &input)
Definition: ncmesh.cpp:3319
int GetMidEdgeNode(int vn1, int vn2)
Definition: ncmesh.cpp:513
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:638
Abstract data type element.
Definition: element.hpp:27
void UnrefElement(int elem, Array< int > &elemFaces)
Definition: ncmesh.cpp:262
Point & operator=(const Point &src)
Definition: ncmesh.hpp:624
void PrintStats(std::ostream &out=mfem::out) const
Definition: ncmesh.cpp:3591
int GetSingleElement() const
Return one of elem[0] or elem[1] and make sure the other is -1.
Definition: ncmesh.cpp:352
const double * CalcVertexPos(int node) const
Definition: ncmesh.cpp:1499
int rank
processor number (ParNCMesh), -1 if undefined/unknown
Definition: ncmesh.hpp:382
bool NodeSetX1(int node, int *n)
Definition: ncmesh.cpp:530
int node[8]
element corners (if ref_type == 0)
Definition: ncmesh.hpp:386
DenseMatrix point_matrix
position within the master edge/face
Definition: ncmesh.hpp:159
HashTable< Face > faces
Definition: ncmesh.hpp:397
Defines the position of a fine element within a coarse element.
Definition: ncmesh.hpp:43
void RefElement(int elem)
Definition: ncmesh.cpp:231
void FaceSplitLevel(int vn1, int vn2, int vn3, int vn4, int &h_level, int &v_level) const
Definition: ncmesh.cpp:3166
HashTable< Face >::const_iterator face_const_iterator
Definition: ncmesh.hpp:411
virtual void Refine(const Array< Refinement > &refinements)
Definition: ncmesh.cpp:1075
int matrix
index into CoarseFineTransformations::point_matrices
Definition: ncmesh.hpp:46
int GetEdgeMaster(int v1, int v2) const
Definition: ncmesh.cpp:3056