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