MFEM  v4.1.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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
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 "../general/sort_pairs.hpp"
19 #include "../linalg/densemat.hpp"
20 #include "element.hpp"
21 #include "vertex.hpp"
22 #include "../fem/geom.hpp"
23 
24 #include <vector>
25 #include <map>
26 #include <iostream>
27 
28 namespace mfem
29 {
30 
31 /** Represents the index of an element to refine, plus a refinement type.
32  The refinement type is needed for anisotropic refinement of quads and hexes.
33  Bits 0,1 and 2 of 'ref_type' specify whether the element should be split
34  in the X, Y and Z directions, respectively (Z is ignored for quads). */
35 struct Refinement
36 {
37  int index; ///< Mesh element number
38  char ref_type; ///< refinement XYZ bit mask (7 = full isotropic)
39 
40  Refinement() = default;
41 
42  Refinement(int index, int type = 7) : index(index), ref_type(type) {}
43 };
44 
45 /// Defines the position of a fine element within a coarse element.
46 struct Embedding
47 {
48  /// %Element index in the coarse mesh.
49  int parent;
50  /** @brief Index into the DenseTensor corresponding to the parent
51  Geometry::Type stored in CoarseFineTransformations::point_matrices. */
52  int matrix;
53 
54  Embedding() = default;
55 
56  Embedding(int elem, int matrix = 0) : parent(elem), matrix(matrix) {}
57 };
58 
59 /// Defines the coarse-fine transformations of all fine elements.
61 {
62  /// Matrices for IsoparametricTransformation organized by Geometry::Type
64  /// Fine element positions in their parents.
66 
67  void GetCoarseToFineMap(const Mesh &fine_mesh,
68  Table &coarse_to_fine,
69  Array<int> &coarse_to_ref_type,
70  Table &ref_type_to_matrix,
71  Array<Geometry::Type> &ref_type_to_geom) const;
72 
73  void Clear();
74  bool IsInitialized() const;
75  long MemoryUsage() const;
76 };
77 
78 
79 /** \brief A class for non-conforming AMR on higher-order hexahedral, prismatic,
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  signed char local; ///< local number within 'element'
156  signed char geom; /**< Geometry::Type (faces only) (char storage to save
157  RAM) */
158 
159  MeshId(int index = -1, int element = -1, signed char local = -1,
160  signed char geom = -1)
162 
163  Geometry::Type Geom() const { return Geometry::Type(geom); }
164  };
165 
166  /** Nonconforming edge/face that has more than one neighbor. The neighbors
167  are stored in NCList::slaves[i], slaves_begin <= i < slaves_end. */
168  struct Master : public MeshId
169  {
170  int slaves_begin, slaves_end; ///< slave faces
171 
172  Master(int index, int element, char local, char geom, int sb, int se)
173  : MeshId(index, element, local, geom)
174  , slaves_begin(sb), slaves_end(se) {}
175  };
176 
177  /// Nonconforming edge/face within a bigger edge/face.
178  struct Slave : public MeshId
179  {
180  int master; ///< master number (in Mesh numbering)
181  int edge_flags; ///< edge orientation flags
182  DenseMatrix point_matrix; ///< position within the master edge/face
183 
184  Slave(int index, int element, signed char local, signed char geom)
185  : MeshId(index, element, local, geom)
186  , master(-1), edge_flags(0) {}
187 
188  /// Return the point matrix oriented according to the master and slave edges
189  void OrientedPointMatrix(DenseMatrix &oriented_matrix) const;
190  };
191 
192  /// Lists all edges/faces in the nonconforming mesh.
193  struct NCList
194  {
195  std::vector<MeshId> conforming;
196  std::vector<Master> masters;
197  std::vector<Slave> slaves;
198  // TODO: switch to Arrays when fixed for non-POD types
199  // TODO: make a list of unique slave matrices to save memory (+ time later)
200 
201  void Clear(bool hard = false);
202  bool Empty() const { return !conforming.size() && !masters.size(); }
203  long TotalSize() const;
204  long MemoryUsage() const;
205 
206  const MeshId& LookUp(int index, int *type = NULL) const;
207  private:
208  mutable Array<int> inv_index;
209  };
210 
211  /// Return the current list of conforming and nonconforming faces.
213  {
214  if (face_list.Empty()) { BuildFaceList(); }
215  return face_list;
216  }
217 
218  /// Return the current list of conforming and nonconforming edges.
220  {
221  if (edge_list.Empty()) { BuildEdgeList(); }
222  return edge_list;
223  }
224 
225  /** Return a list of vertices (in 'conforming'); this function is provided
226  for uniformity/completeness. Needed in ParNCMesh/ParFESpace. */
228  {
229  if (vertex_list.Empty()) { BuildVertexList(); }
230  return vertex_list;
231  }
232 
233  /// Return vertex/edge/face list (entity = 0/1/2, respectively).
234  const NCList& GetNCList(int entity)
235  {
236  switch (entity)
237  {
238  case 0: return GetVertexList();
239  case 1: return GetEdgeList();
240  default: return GetFaceList();
241  }
242  }
243 
244 
245  // coarse/fine transforms
246 
247  /** Remember the current layer of leaf elements before the mesh is refined.
248  Needed by GetRefinementTransforms(), must be called before Refine(). */
249  void MarkCoarseLevel();
250 
251  /** After refinement, calculate the relation of each fine element to its
252  parent coarse element. Note that Refine() or LimitNCLevel() can be called
253  multiple times between MarkCoarseLevel() and this function. */
255 
256  /** After derefinement, calculate the relations of previous fine elements
257  (some of which may no longer exist) to the current leaf elements.
258  Unlike for refinement, Derefine() may only be called once before this
259  function so there is no MarkFineLevel(). */
261 
262  /// Free all internal data created by the above three functions.
263  void ClearTransforms();
264 
265 
266  // grid ordering
267 
268  /** Return a space filling curve for a rectangular grid of elements.
269  Implemented is a generalized Hilbert curve for arbitrary grid dimensions.
270  If the width is odd, height should be odd too, otherwise one diagonal
271  (vertex-neighbor) step cannot be avoided in the curve. Even dimensions
272  are recommended. */
273  static void GridSfcOrdering2D(int width, int height,
274  Array<int> &coords);
275 
276  /** Return a space filling curve for a 3D rectangular grid of elements.
277  The Hilbert-curve-like algorithm works well for even dimensions. For odd
278  width/height/depth it tends to produce some diagonal (edge-neighbor)
279  steps. Even dimensions are recommended. */
280  static void GridSfcOrdering3D(int width, int height, int depth,
281  Array<int> &coords);
282 
283 
284  // utility
285 
286  /// Return Mesh vertex indices of an edge identified by 'edge_id'.
287  void GetEdgeVertices(const MeshId &edge_id, int vert_index[2],
288  bool oriented = true) const;
289 
290  /** Return "NC" orientation of an edge. As opposed to standard Mesh edge
291  orientation based on vertex IDs, "NC" edge orientation follows the local
292  edge orientation within the element 'edge_id.element' and is thus
293  processor independent. TODO: this seems only partially true? */
294  int GetEdgeNCOrientation(const MeshId &edge_id) const;
295 
296  /** Return Mesh vertex and edge indices of a face identified by 'face_id'.
297  The return value is the number of face vertices. */
298  int GetFaceVerticesEdges(const MeshId &face_id,
299  int vert_index[4], int edge_index[4],
300  int edge_orientation[4]) const;
301 
302  /** Given an edge (by its vertex indices v1 and v2) return the first
303  (geometric) parent edge that exists in the Mesh or -1 if there is no such
304  parent. */
305  int GetEdgeMaster(int v1, int v2) const;
306 
307  /** Get a list of vertices (2D/3D) and edges (3D) that coincide with boundary
308  elements with the specified attributes (marked in 'bdr_attr_is_ess').
309  In 3D this function also reveals "hidden" boundary edges. In parallel it
310  helps identifying boundary vertices/edges affected by non-local boundary
311  elements. */
312  virtual void GetBoundaryClosure(const Array<int> &bdr_attr_is_ess,
313  Array<int> &bdr_vertices,
314  Array<int> &bdr_edges);
315 
316  /// Return element geometry type. @a index is the Mesh element number.
318  { return elements[leaf_elements[index]].Geom(); }
319 
320  /// Return face geometry type. @a index is the Mesh face number.
322  { return Geometry::Type(face_geom[index]); }
323 
324  /// Return the number of root elements.
325  int GetNumRootElements() { return root_state.Size(); }
326 
327  /// Return the distance of leaf 'i' from the root.
328  int GetElementDepth(int i) const;
329 
330  /** Return the size reduction compared to the root element (ignoring local
331  stretching and curvature). */
332  int GetElementSizeReduction(int i) const;
333 
334  /// Return the faces and face attributes of leaf element 'i'.
336  Array<int> &fattr) const;
337 
338 
339  /// I/O: Print the "vertex_parents" section of the mesh file (ver. >= 1.1).
340  void PrintVertexParents(std::ostream &out) const;
341 
342  /// I/O: Print the "coarse_elements" section of the mesh file (ver. >= 1.1).
343  void PrintCoarseElements(std::ostream &out) const;
344 
345  /** I/O: Load the vertex parent hierarchy from a mesh file. NOTE: called
346  indirectly through the constructor. */
347  void LoadVertexParents(std::istream &input);
348 
349  /// I/O: Load the element refinement hierarchy from a mesh file.
350  void LoadCoarseElements(std::istream &input);
351 
352  /// I/O: Set positions of all vertices (used by mesh loader).
353  void SetVertexPositions(const Array<mfem::Vertex> &vertices);
354 
355  /// Save memory by releasing all non-essential and cached data.
356  virtual void Trim();
357 
358  /// Return total number of bytes allocated.
359  long MemoryUsage() const;
360 
361  int PrintMemoryDetail() const;
362 
363  void PrintStats(std::ostream &out = mfem::out) const;
364 
365  typedef int64_t RefCoord;
366 
367 
368 protected: // interface for Mesh to be able to construct itself from NCMesh
369 
370  friend class Mesh;
371 
372  /// Fill Mesh::{vertices,elements,boundary} for the current finest level.
373  void GetMeshComponents(Mesh &mesh) const;
374 
375  /** Get edge and face numbering from 'mesh' (i.e., set all Edge::index and
376  Face::index) after a new mesh was created from us. */
377  virtual void OnMeshUpdated(Mesh *mesh);
378 
379 
380 protected: // implementation
381 
382  int Dim, spaceDim; ///< dimensions of the elements and the vertex coordinates
383  bool Iso; ///< true if the mesh only contains isotropic refinements
384  int Geoms; ///< bit mask of element geometries present, see InitGeomFlags()
385 
386  /** A Node can hold a vertex, an edge, or both. Elements directly point to
387  their corner nodes, but edge nodes also exist and can be accessed using
388  a hash-table given their two end-point node IDs. All nodes can be
389  accessed in this way, with the exception of top-level vertex nodes.
390  When an element is being refined, the mid-edge nodes are readily
391  available with this mechanism. The new elements "sign in" to the nodes
392  by increasing the reference counts of their vertices and edges. The
393  parent element "signs off" its nodes by decrementing the ref counts. */
394  struct Node : public Hashed2
395  {
398 
399  Node() : vert_refc(0), edge_refc(0), vert_index(-1), edge_index(-1) {}
400  ~Node();
401 
402  bool HasVertex() const { return vert_refc > 0; }
403  bool HasEdge() const { return edge_refc > 0; }
404 
405  // decrease vertex/edge ref count, return false if Node should be deleted
406  bool UnrefVertex() { --vert_refc; return vert_refc || edge_refc; }
407  bool UnrefEdge() { --edge_refc; return vert_refc || edge_refc; }
408  };
409 
410  /** Similarly to nodes, faces can be accessed by hashing their four vertex
411  node IDs. A face knows about the one or two elements that are using it.
412  A face that is not on the boundary and only has one element referencing
413  it is either a master or a slave face. */
414  struct Face : public Hashed4
415  {
416  int attribute; ///< boundary element attribute, -1 if internal face
417  int index; ///< face number in the Mesh
418  int elem[2]; ///< up to 2 elements sharing the face
419 
420  Face() : attribute(-1), index(-1) { elem[0] = elem[1] = -1; }
421 
422  bool Boundary() const { return attribute >= 0; }
423  bool Unused() const { return elem[0] < 0 && elem[1] < 0; }
424 
425  // add or remove an element from the 'elem[2]' array
426  void RegisterElement(int e);
427  void ForgetElement(int e);
428 
429  /// Return one of elem[0] or elem[1] and make sure the other is -1.
430  int GetSingleElement() const;
431  };
432 
433  /** This is an element in the refinement hierarchy. Each element has
434  either been refined and points to its children, or is a leaf and points
435  to its vertex nodes. */
436  struct Element
437  {
438  char geom; ///< Geometry::Type of the element (char for storage only)
439  char ref_type; ///< bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
440  char tet_type; ///< tetrahedron split type, currently always 0
441  char flag; ///< generic flag/marker, can be used by algorithms
442  int index; ///< element number in the Mesh, -1 if refined
443  int rank; ///< processor number (ParNCMesh), -1 if undefined/unknown
445  union
446  {
447  int node[8]; ///< element corners (if ref_type == 0)
448  int child[8]; ///< 2-8 children (if ref_type != 0)
449  };
450  int parent; ///< parent element, -1 if this is a root element, -2 if free
451 
452  Element(Geometry::Type geom, int attr);
453 
454  Geometry::Type Geom() const { return Geometry::Type(geom); }
455  };
456 
457  // primary data
458 
459  HashTable<Node> nodes; // associative container holding all Nodes
460  HashTable<Face> faces; // associative container holding all Faces
461 
462  BlockArray<Element> elements; // storage for all Elements
463  Array<int> free_element_ids; // unused element ids - indices into 'elements'
464 
465  /** Initial traversal state (~ element orientation) for each root element
466  NOTE: M = root_state.Size() is the number of root elements.
467  NOTE: the first M items of 'elements' is the coarse mesh. */
469 
470  /// coordinates of top-level vertices (organized as triples)
472 
478 
479 
480  // secondary data
481 
482  /** Apart from the primary data structure, which is the element/node/face
483  hierarchy, there is secondary data that is derived from the primary
484  data and needs to be updated when the primary data changes. Update()
485  takes care of that and needs to be called after refinement and
486  derefinement. */
487  virtual void Update();
488 
489  int NVertices; // set by UpdateVertices
490  int NEdges, NFaces; // set by OnMeshUpdated
491 
492  Array<int> leaf_elements; // finest level, calculated by UpdateLeafElements
493  Array<int> vertex_nodeId; // vertex-index to node-id map, see UpdateVertices
494 
495  NCList face_list; ///< lazy-initialized list of faces, see GetFaceList
496  NCList edge_list; ///< lazy-initialized list of edges, see GetEdgeList
497  NCList vertex_list; ///< lazy-initialized list of vertices, see GetVertexList
498 
499  Array<int> boundary_faces; ///< subset of all faces, set by BuildFaceList
500  Array<char> face_geom; ///< face geometry by face index, set by OnMeshUpdated
501 
502  Table element_vertex; ///< leaf-element to vertex table, see FindSetNeighbors
503 
504 
505  virtual void UpdateVertices(); ///< update Vertex::index and vertex_nodeId
506 
507  void CollectLeafElements(int elem, int state);
508  void UpdateLeafElements();
509 
510  virtual void AssignLeafIndices();
511 
512  /** Try to find a space-filling curve friendly orientation of the root
513  elements: set 'root_state' based on the ordering of coarse elements.
514  Note that the coarse mesh itself must be ordered as an SFC by e.g.
515  Mesh::GetGeckoElementOrdering. */
516  void InitRootState(int root_count);
517 
518  virtual bool IsGhost(const Element &el) const { return false; }
519  virtual int GetNumGhostElements() const { return 0; }
520  virtual int GetNumGhostVertices() const { return 0; }
521 
522  void InitGeomFlags();
523 
524  bool HavePrisms() const { return Geoms & (1 << Geometry::PRISM); }
525  bool HaveTets() const { return Geoms & (1 << Geometry::TETRAHEDRON); }
526 
527 
528  // refinement/derefinement
529 
530  Array<Refinement> ref_stack; ///< stack of scheduled refinements (temporary)
531  HashTable<Node> shadow; ///< temporary storage for reparented nodes
532  Array<Triple<int, int, int> > reparents; ///< scheduled node reparents (tmp)
533 
534  Table derefinements; ///< possible derefinements, see GetDerefinementTable
535 
536  void RefineElement(int elem, char ref_type);
537  void DerefineElement(int elem);
538 
539  int AddElement(const Element &el)
540  {
541  if (free_element_ids.Size())
542  {
543  int idx = free_element_ids.Last();
545  elements[idx] = el;
546  return idx;
547  }
548  return elements.Append(el);
549  }
550  void FreeElement(int id)
551  {
553  elements[id].parent = -2; // mark the element as free
554  }
555 
556  int NewHexahedron(int n0, int n1, int n2, int n3,
557  int n4, int n5, int n6, int n7, int attr,
558  int fattr0, int fattr1, int fattr2,
559  int fattr3, int fattr4, int fattr5);
560 
561  int NewWedge(int n0, int n1, int n2,
562  int n3, int n4, int n5, int attr,
563  int fattr0, int fattr1,
564  int fattr2, int fattr3, int fattr4);
565 
566  int NewTetrahedron(int n0, int n1, int n2, int n3, int attr,
567  int fattr0, int fattr1, int fattr2, int fattr3);
568 
569  int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr,
570  int eattr0, int eattr1, int eattr2, int eattr3);
571 
572  int NewTriangle(int n0, int n1, int n2,
573  int attr, int eattr0, int eattr1, int eattr2);
574 
575  mfem::Element* NewMeshElement(int geom) const;
576 
577  int QuadFaceSplitType(int v1, int v2, int v3, int v4, int mid[5]
578  = NULL /*optional output of mid-edge nodes*/) const;
579 
580  bool TriFaceSplit(int v1, int v2, int v3, int mid[3] = NULL) const;
581 
582  void ForceRefinement(int vn1, int vn2, int vn3, int vn4);
583 
584  void FindEdgeElements(int vn1, int vn2, int vn3, int vn4,
585  Array<MeshId> &prisms) const;
586 
587  void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4,
588  const Refinement *refs, int nref);
589 
590  void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4,
591  int mid12, int mid34, int level = 0);
592 
593  void CheckIsoFace(int vn1, int vn2, int vn3, int vn4,
594  int en1, int en2, int en3, int en4, int midf);
595 
596  void ReparentNode(int node, int new_p1, int new_p2);
597 
598  int FindMidEdgeNode(int node1, int node2) const;
599  int GetMidEdgeNode(int node1, int node2);
600 
601  int GetMidFaceNode(int en1, int en2, int en3, int en4);
602 
603  void ReferenceElement(int elem);
604  void UnreferenceElement(int elem, Array<int> &elemFaces);
605 
606  Face* GetFace(Element &elem, int face_no);
607  void RegisterFaces(int elem, int *fattr = NULL);
608  void DeleteUnusedFaces(const Array<int> &elemFaces);
609 
610  void CollectDerefinements(int elem, Array<Connection> &list);
611 
612  /// Return el.node[index] correctly, even if the element is refined.
613  int RetrieveNode(const Element &el, int index);
614 
615  /// Extended version of find_node: works if 'el' is refined.
616  int FindNodeExt(const Element &el, int node, bool abort = true);
617 
618 
619  // face/edge lists
620 
621  static int find_node(const Element &el, int node);
622  static int find_element_edge(const Element &el, int vn0, int vn1,
623  bool abort = true);
624  static int find_local_face(int geom, int a, int b, int c);
625 
626  int ReorderFacePointMat(int v0, int v1, int v2, int v3,
627  int elem, DenseMatrix& mat) const;
628  struct Point;
629  struct PointMatrix;
630  void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3,
631  const PointMatrix& pm, int level, Face* eface[4]);
632  bool TraverseTriFace(int vn0, int vn1, int vn2,
633  const PointMatrix& pm, int level);
634  void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1);
635  void TraverseEdge(int vn0, int vn1, double t0, double t1, int flags,
636  int level);
637 
638  virtual void BuildFaceList();
639  virtual void BuildEdgeList();
640  virtual void BuildVertexList();
641 
642  virtual void ElementSharesFace(int elem, int local, int face) {} // ParNCMesh
643  virtual void ElementSharesEdge(int elem, int local, int enode) {} // ParNCMesh
644  virtual void ElementSharesVertex(int elem, int local, int vnode) {} // ParNCMesh
645 
646 
647  // neighbors / element_vertex table
648 
649  /** Return all vertex-, edge- and face-neighbors of a set of elements.
650  The neighbors are returned as a list (neighbors != NULL), as a set
651  (neighbor_set != NULL), or both. The sizes of the set arrays must match
652  that of leaf_elements. The function is intended to be used for large
653  sets of elements and its complexity is linear in the number of leaf
654  elements in the mesh. */
655  void FindSetNeighbors(const Array<char> &elem_set,
656  Array<int> *neighbors, /* append */
657  Array<char> *neighbor_set = NULL);
658 
659  /** Return all vertex-, edge- and face-neighbors of a single element.
660  You can limit the number of elements being checked using 'search_set'.
661  The complexity of the function is linear in the size of the search set.*/
662  void FindNeighbors(int elem,
663  Array<int> &neighbors, /* append */
664  const Array<int> *search_set = NULL);
665 
666  /** Expand a set of elements by all vertex-, edge- and face-neighbors.
667  The output array 'expanded' will contain all items from 'elems'
668  (provided they are in 'search_set') plus their neighbors. The neighbor
669  search can be limited to the optional search set. The complexity is
670  linear in the sum of the sizes of 'elems' and 'search_set'. */
671  void NeighborExpand(const Array<int> &elems,
672  Array<int> &expanded,
673  const Array<int> *search_set = NULL);
674 
675 
676  void CollectEdgeVertices(int v0, int v1, Array<int> &indices);
677  void CollectTriFaceVertices(int v0, int v1, int v2, Array<int> &indices);
678  void CollectQuadFaceVertices(int v0, int v1, int v2, int v3,
679  Array<int> &indices);
681 
683  {
685  }
686 
687  int GetVertexRootCoord(int elem, RefCoord coord[3]) const;
688  void CollectIncidentElements(int elem, const RefCoord coord[3],
689  Array<int> &list) const;
690 
691  /** Return elements neighboring to a local vertex of element 'elem'. Only
692  elements from within the same refinement tree ('cousins') are returned.
693  Complexity is proportional to the depth of elem's refinement tree. */
694  void FindVertexCousins(int elem, int local, Array<int> &cousins) const;
695 
696 
697  // coarse/fine transformations
698 
699  struct Point
700  {
701  int dim;
702  double coord[3];
703 
704  Point() { dim = 0; }
705 
706  Point(double x, double y)
707  { dim = 2; coord[0] = x; coord[1] = y; }
708 
709  Point(double x, double y, double z)
710  { dim = 3; coord[0] = x; coord[1] = y; coord[2] = z; }
711 
712  Point(const Point& p0, const Point& p1)
713  {
714  dim = p0.dim;
715  for (int i = 0; i < dim; i++)
716  {
717  coord[i] = (p0.coord[i] + p1.coord[i]) * 0.5;
718  }
719  }
720 
721  Point(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
722  {
723  dim = p0.dim;
724  MFEM_ASSERT(p1.dim == dim && p2.dim == dim && p3.dim == dim, "");
725  for (int i = 0; i < dim; i++)
726  {
727  coord[i] = (p0.coord[i] + p1.coord[i] + p2.coord[i] + p3.coord[i])
728  * 0.25;
729  }
730  }
731 
732  Point& operator=(const Point& src)
733  {
734  dim = src.dim;
735  for (int i = 0; i < dim; i++) { coord[i] = src.coord[i]; }
736  return *this;
737  }
738  };
739 
740  struct PointMatrix
741  {
742  int np;
744 
745  PointMatrix(const Point& p0, const Point& p1, const Point& p2)
746  { np = 3; points[0] = p0; points[1] = p1; points[2] = p2; }
747 
748  PointMatrix(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
749  { np = 4; points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3; }
750 
751  PointMatrix(const Point& p0, const Point& p1, const Point& p2,
752  const Point& p3, const Point& p4, const Point& p5)
753  {
754  np = 6;
755  points[0] = p0; points[1] = p1; points[2] = p2;
756  points[3] = p3; points[4] = p4; points[5] = p5;
757  }
758  PointMatrix(const Point& p0, const Point& p1, const Point& p2,
759  const Point& p3, const Point& p4, const Point& p5,
760  const Point& p6, const Point& p7)
761  {
762  np = 8;
763  points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3;
764  points[4] = p4; points[5] = p5; points[6] = p6; points[7] = p7;
765  }
766 
767  Point& operator()(int i) { return points[i]; }
768  const Point& operator()(int i) const { return points[i]; }
769 
770  void GetMatrix(DenseMatrix& point_matrix) const;
771  };
772 
778 
780 
781  void GetPointMatrix(Geometry::Type geom, const char* ref_path,
782  DenseMatrix& matrix);
783 
784  typedef std::map<std::string, int> RefPathMap;
785 
786  void TraverseRefinements(int elem, int coarse_index,
787  std::string &ref_path, RefPathMap &map);
788 
789  /// storage for data returned by Get[De]RefinementTransforms()
791 
792  /// state of leaf_elements before Refine(), set by MarkCoarseLevel()
794 
795  void InitDerefTransforms();
796  void SetDerefMatrixCodes(int parent, Array<int> &fine_coarse);
797 
798 
799  // vertex temporary data, used by GetMeshComponents
800 
801  struct TmpVertex
802  {
803  bool valid, visited;
804  double pos[3];
805  TmpVertex() : valid(false), visited(false) {}
806  };
807 
809 
810  const double *CalcVertexPos(int node) const;
811 
812 
813  // utility
814 
815  int GetEdgeMaster(int node) const;
816 
817  void FindFaceNodes(int face, int node[4]);
818 
819  int EdgeSplitLevel(int vn1, int vn2) const;
820  int TriFaceSplitLevel(int vn1, int vn2, int vn3) const;
821  void QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4,
822  int& h_level, int& v_level) const;
823 
824  void CountSplits(int elem, int splits[3]) const;
825  void GetLimitRefinements(Array<Refinement> &refinements, int max_level);
826 
827  int PrintElements(std::ostream &out, int elem, int &coarse_id) const;
828  void CopyElements(int elem, const BlockArray<Element> &tmp_elements,
829  Array<int> &index_map);
830 
831  // geometry
832 
833  /** This holds in one place the constants about the geometries we support
834  (triangles, quads, cubes) */
835  struct GeomInfo
836  {
837  int nv, ne, nf; // number of: vertices, edges, faces
838  int edges[12][2]; // edge vertices (up to 12 edges)
839  int faces[6][4]; // face vertices (up to 6 faces)
840  int nfv[6]; // number of face vertices
841 
843  GeomInfo() : initialized(false) {}
844  void Initialize(const mfem::Element* elem);
845  };
846 
848 
849 #ifdef MFEM_DEBUG
850 public:
851  void DebugLeafOrder(std::ostream &out) const;
852  void DebugDump(std::ostream &out) const;
853 #endif
854 
855  friend class ParNCMesh; // for ParNCMesh::ElementSet
856  friend struct CompareRanks;
857 };
858 
859 }
860 
861 #endif
NCList face_list
lazy-initialized list of faces, see GetFaceList
Definition: ncmesh.hpp:495
void CollectLeafElements(int elem, int state)
Definition: ncmesh.cpp:1865
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
Definition: ncmesh.hpp:496
NCMesh(const Mesh *mesh, std::istream *vertex_parents=NULL)
Definition: ncmesh.cpp:72
bool TraverseTriFace(int vn0, int vn1, int vn2, const PointMatrix &pm, int level)
Definition: ncmesh.cpp:2487
int Size() const
Logical size of the array.
Definition: array.hpp:124
static void GridSfcOrdering3D(int width, int height, int depth, Array< int > &coords)
Definition: ncmesh.cpp:4390
void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4, const Refinement *refs, int nref)
Definition: ncmesh.cpp:722
CoarseFineTransformations transforms
storage for data returned by Get[De]RefinementTransforms()
Definition: ncmesh.hpp:790
static void GridSfcOrdering2D(int width, int height, Array< int > &coords)
Definition: ncmesh.cpp:4375
int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
Definition: ncmesh.cpp:541
int elem[2]
up to 2 elements sharing the face
Definition: ncmesh.hpp:418
char ref_type
refinement XYZ bit mask (7 = full isotropic)
Definition: ncmesh.hpp:38
int GetFaceVerticesEdges(const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
Definition: ncmesh.cpp:4453
TmpVertex * tmp_vertex
Definition: ncmesh.hpp:808
Array< double > top_vertex_pos
coordinates of top-level vertices (organized as triples)
Definition: ncmesh.hpp:471
bool TriFaceSplit(int v1, int v2, int v3, int mid[3]=NULL) const
Definition: ncmesh.cpp:2233
void SetDerefMatrixCodes(int parent, Array< int > &fine_coarse)
Definition: ncmesh.cpp:1828
const CoarseFineTransformations & GetDerefinementTransforms()
Definition: ncmesh.cpp:3954
char flag
generic flag/marker, can be used by algorithms
Definition: ncmesh.hpp:441
bool HasEdge() const
Definition: ncmesh.hpp:403
PointMatrix(const Point &p0, const Point &p1, const Point &p2)
Definition: ncmesh.hpp:745
Array< Triple< int, int, int > > reparents
scheduled node reparents (tmp)
Definition: ncmesh.hpp:532
Point(const Point &p0, const Point &p1)
Definition: ncmesh.hpp:712
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:4855
void CollectTriFaceVertices(int v0, int v1, int v2, Array< int > &indices)
Definition: ncmesh.cpp:2913
char tet_type
tetrahedron split type, currently always 0
Definition: ncmesh.hpp:440
virtual void LimitNCLevel(int max_nc_level)
Definition: ncmesh.cpp:4840
virtual void Trim()
Save memory by releasing all non-essential and cached data.
Definition: ncmesh.cpp:5078
int index
element number in the Mesh, -1 if refined
Definition: ncmesh.hpp:442
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:748
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:4965
void ForceRefinement(int vn1, int vn2, int vn3, int vn4)
Definition: ncmesh.cpp:617
virtual int GetNumGhostElements() const
Definition: ncmesh.hpp:519
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:758
Lists all edges/faces in the nonconforming mesh.
Definition: ncmesh.hpp:193
void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level, Face *eface[4])
Definition: ncmesh.cpp:2339
void TraverseEdge(int vn0, int vn1, double t0, double t1, int flags, int level)
Definition: ncmesh.cpp:2626
int NVertices
Definition: ncmesh.hpp:489
int GetNVertices() const
Definition: ncmesh.hpp:115
void GetMeshComponents(Mesh &mesh) const
Fill Mesh::{vertices,elements,boundary} for the current finest level.
Definition: ncmesh.cpp:2033
void ClearTransforms()
Free all internal data created by the above three functions.
Definition: ncmesh.cpp:4127
Table derefinements
possible derefinements, see GetDerefinementTable
Definition: ncmesh.hpp:534
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Geometry::Type Geom() const
Definition: ncmesh.hpp:163
NCList vertex_list
lazy-initialized list of vertices, see GetVertexList
Definition: ncmesh.hpp:497
void InitDerefTransforms()
Definition: ncmesh.cpp:1813
static PointMatrix pm_prism_identity
Definition: ncmesh.hpp:776
char geom
Geometry::Type of the element (char for storage only)
Definition: ncmesh.hpp:438
BlockArray< Element >::iterator elem_iterator
Definition: ncmesh.hpp:477
void InitRootState(int root_count)
Definition: ncmesh.cpp:1927
long TotalSize() const
Definition: ncmesh.cpp:2831
MeshId(int index=-1, int element=-1, signed char local=-1, signed char geom=-1)
Definition: ncmesh.hpp:159
int GetEdgeNCOrientation(const MeshId &edge_id) const
Definition: ncmesh.cpp:4441
void OrientedPointMatrix(DenseMatrix &oriented_matrix) const
Return the point matrix oriented according to the master and slave edges.
Definition: ncmesh.cpp:2792
int NewTetrahedron(int n0, int n1, int n2, int n3, int attr, int fattr0, int fattr1, int fattr2, int fattr3)
Definition: ncmesh.cpp:515
void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4, int mid12, int mid34, int level=0)
Definition: ncmesh.cpp:747
virtual int GetNumGhostVertices() const
Definition: ncmesh.hpp:520
static PointMatrix pm_tri_identity
Definition: ncmesh.hpp:773
void FreeElement(int id)
Definition: ncmesh.hpp:550
void CountSplits(int elem, int splits[3]) const
Definition: ncmesh.cpp:4727
Geometry::Type Geom() const
Definition: ncmesh.hpp:454
void CollectDerefinements(int elem, Array< Connection > &list)
Definition: ncmesh.cpp:1700
const Point & operator()(int i) const
Definition: ncmesh.hpp:768
std::vector< MeshId > conforming
Definition: ncmesh.hpp:195
Point & operator()(int i)
Definition: ncmesh.hpp:767
void GetLimitRefinements(Array< Refinement > &refinements, int max_level)
Definition: ncmesh.cpp:4810
int PrintMemoryDetail() const
Definition: ncmesh.cpp:5136
double coord[3]
Definition: ncmesh.hpp:702
bool UnrefVertex()
Definition: ncmesh.hpp:406
int spaceDim
dimensions of the elements and the vertex coordinates
Definition: ncmesh.hpp:382
bool HavePrisms() const
Definition: ncmesh.hpp:524
Array< int > vertex_nodeId
Definition: ncmesh.hpp:493
int attribute
boundary element attribute, -1 if internal face
Definition: ncmesh.hpp:416
void DebugDump(std::ostream &out) const
Definition: ncmesh.cpp:5241
void FindFaceNodes(int face, int node[4])
Definition: ncmesh.cpp:4579
Array< int > root_state
Definition: ncmesh.hpp:468
static const PointMatrix & GetGeomIdentity(Geometry::Type geom)
Definition: ncmesh.cpp:3422
static int find_element_edge(const Element &el, int vn0, int vn1, bool abort=true)
Definition: ncmesh.cpp:2270
int index
Mesh element number.
Definition: ncmesh.hpp:37
int master
master number (in Mesh numbering)
Definition: ncmesh.hpp:180
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition: ncmesh.hpp:212
const MeshId & LookUp(int index, int *type=NULL) const
Definition: ncmesh.cpp:2836
virtual void OnMeshUpdated(Mesh *mesh)
Definition: ncmesh.cpp:2120
void DebugLeafOrder(std::ostream &out) const
Definition: ncmesh.cpp:5216
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:64
int NewWedge(int n0, int n1, int n2, int n3, int n4, int n5, int attr, int fattr0, int fattr1, int fattr2, int fattr3, int fattr4)
Definition: ncmesh.cpp:483
void CollectEdgeVertices(int v0, int v1, Array< int > &indices)
Definition: ncmesh.cpp:2901
static int find_local_face(int geom, int a, int b, int c)
Definition: ncmesh.cpp:2288
Array< Refinement > ref_stack
stack of scheduled refinements (temporary)
Definition: ncmesh.hpp:530
bool Empty() const
Definition: ncmesh.hpp:202
Point(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:721
void CollectIncidentElements(int elem, const RefCoord coord[3], Array< int > &list) const
Definition: ncmesh.cpp:3352
void ReferenceElement(int elem)
Definition: ncmesh.cpp:303
virtual bool IsGhost(const Element &el) const
Definition: ncmesh.hpp:518
void GetPointMatrix(Geometry::Type geom, const char *ref_path, DenseMatrix &matrix)
Definition: ncmesh.cpp:3437
void CheckIsoFace(int vn1, int vn2, int vn3, int vn4, int en1, int en2, int en3, int en4, int midf)
Definition: ncmesh.cpp:842
void InitGeomFlags()
Definition: ncmesh.cpp:204
Face * GetFace(Element &elem, int face_no)
Definition: ncmesh.cpp:416
void TraverseRefinements(int elem, int coarse_index, std::string &ref_path, RefPathMap &map)
Definition: ncmesh.cpp:3871
int GetMidEdgeNode(int node1, int node2)
Definition: ncmesh.cpp:288
Element(Geometry::Type geom, int attr)
Definition: ncmesh.cpp:441
int PrintElements(std::ostream &out, int elem, int &coarse_id) const
Definition: ncmesh.cpp:4925
Refinement(int index, int type=7)
Definition: ncmesh.hpp:42
int index
face number in the Mesh
Definition: ncmesh.hpp:417
const Table & GetDerefinementTable()
Definition: ncmesh.cpp:1733
Table element_vertex
leaf-element to vertex table, see FindSetNeighbors
Definition: ncmesh.hpp:502
DenseTensor point_matrices[Geometry::NumGeom]
Matrices for IsoparametricTransformation organized by Geometry::Type.
Definition: ncmesh.hpp:63
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:707
void GetMatrix(DenseMatrix &point_matrix) const
Definition: ncmesh.cpp:3392
int AddElement(const Element &el)
Definition: ncmesh.hpp:539
virtual void Derefine(const Array< int > &derefs)
Definition: ncmesh.cpp:1777
double b
Definition: lissajous.cpp:42
void QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4, int &h_level, int &v_level) const
Definition: ncmesh.cpp:4700
void FindNeighbors(int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
Definition: ncmesh.cpp:3154
Geometry::Type GetElementGeometry(int index) const
Return element geometry type. index is the Mesh element number.
Definition: ncmesh.hpp:317
std::vector< Master > masters
Definition: ncmesh.hpp:196
Array< int > coarse_elements
state of leaf_elements before Refine(), set by MarkCoarseLevel()
Definition: ncmesh.hpp:793
virtual void UpdateVertices()
update Vertex::index and vertex_nodeId
Definition: ncmesh.cpp:1847
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:4043
void RegisterFaces(int elem, int *fattr=NULL)
Definition: ncmesh.cpp:377
Array< Embedding > embeddings
Fine element positions in their parents.
Definition: ncmesh.hpp:65
virtual void BuildEdgeList()
Definition: ncmesh.cpp:2656
const CoarseFineTransformations & GetRefinementTransforms()
Definition: ncmesh.cpp:3903
void LoadCoarseElements(std::istream &input)
I/O: Load the element refinement hierarchy from a mesh file.
Definition: ncmesh.cpp:4985
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:151
int FindMidEdgeNode(int node1, int node2) const
Definition: ncmesh.cpp:272
int parent
parent element, -1 if this is a root element, -2 if free
Definition: ncmesh.hpp:450
void MarkCoarseLevel()
Definition: ncmesh.cpp:3857
Slave(int index, int element, signed char local, signed char geom)
Definition: ncmesh.hpp:184
signed char local
local number within &#39;element&#39;
Definition: ncmesh.hpp:155
virtual void ElementSharesFace(int elem, int local, int face)
Definition: ncmesh.hpp:642
int GetNumRootElements()
Return the number of root elements.
Definition: ncmesh.hpp:325
virtual void ElementSharesVertex(int elem, int local, int vnode)
Definition: ncmesh.hpp:644
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
A class for non-conforming AMR on higher-order hexahedral, prismatic, quadrilateral or triangular mes...
Definition: ncmesh.hpp:100
int GetElementSizeReduction(int i) const
Definition: ncmesh.cpp:4544
virtual ~NCMesh()
Definition: ncmesh.cpp:225
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:4603
int element
NCMesh::Element containing this vertex/edge/face.
Definition: ncmesh.hpp:154
Embedding(int elem, int matrix=0)
Definition: ncmesh.hpp:56
Geometry::Type GetFaceGeometry(int index) const
Return face geometry type. index is the Mesh face number.
Definition: ncmesh.hpp:321
void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1)
Definition: ncmesh.cpp:2451
static GeomInfo GI[Geometry::NumGeom]
Definition: ncmesh.hpp:847
void Clear(bool hard=false)
Definition: ncmesh.cpp:2813
HashTable< Node >::iterator node_iterator
Definition: ncmesh.hpp:473
int slaves_end
slave faces
Definition: ncmesh.hpp:170
NCMesh::RefCoord RefCoord
void ReparentNode(int node, int new_p1, int new_p2)
Definition: ncmesh.cpp:256
int NewTriangle(int n0, int n1, int n2, int attr, int eattr0, int eattr1, int eattr2)
Definition: ncmesh.cpp:567
void DeleteUnusedFaces(const Array< int > &elemFaces)
Definition: ncmesh.cpp:391
bool Iso
true if the mesh only contains isotropic refinements
Definition: ncmesh.hpp:383
Point(double x, double y, double z)
Definition: ncmesh.hpp:709
bool HasVertex() const
Definition: ncmesh.hpp:402
std::map< std::string, int > RefPathMap
Definition: ncmesh.hpp:784
int SpaceDimension() const
Definition: ncmesh.hpp:113
std::vector< Slave > slaves
Definition: ncmesh.hpp:197
virtual void BuildFaceList()
Definition: ncmesh.cpp:2542
Nonconforming edge/face within a bigger edge/face.
Definition: ncmesh.hpp:178
void DerefineElement(int elem)
Definition: ncmesh.cpp:1574
void GetElementFacesAttributes(int i, Array< int > &faces, Array< int > &fattr) const
Return the faces and face attributes of leaf element &#39;i&#39;.
Definition: ncmesh.cpp:4558
int ReorderFacePointMat(int v0, int v1, int v2, int v3, int elem, DenseMatrix &mat) const
Definition: ncmesh.cpp:2305
void SetVertexPositions(const Array< mfem::Vertex > &vertices)
I/O: Set positions of all vertices (used by mesh loader).
Definition: ncmesh.cpp:4904
int edge_flags
edge orientation flags
Definition: ncmesh.hpp:181
HashTable< Node >::const_iterator node_const_iterator
Definition: ncmesh.hpp:475
const NCList & GetVertexList()
Definition: ncmesh.hpp:227
bool Boundary() const
Definition: ncmesh.hpp:422
void RegisterElement(int e)
Definition: ncmesh.cpp:402
mfem::Element * NewMeshElement(int geom) const
Definition: ncmesh.cpp:1990
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4, const Point &p5)
Definition: ncmesh.hpp:751
void DeleteLast()
Delete the last entry.
Definition: array.hpp:177
void UpdateLeafElements()
Definition: ncmesh.cpp:1907
int TriFaceSplitLevel(int vn1, int vn2, int vn3) const
Definition: ncmesh.cpp:4683
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:1538
int Geoms
bit mask of element geometries present, see InitGeomFlags()
Definition: ncmesh.hpp:384
long MemoryUsage() const
Definition: ncmesh.cpp:5090
void BuildElementToVertexTable()
Definition: ncmesh.cpp:2971
HashTable< Node > nodes
Definition: ncmesh.hpp:459
double a
Definition: lissajous.cpp:41
void ForgetElement(int e)
Definition: ncmesh.cpp:409
static PointMatrix pm_tet_identity
Definition: ncmesh.hpp:775
int GetMidFaceNode(int en1, int en2, int en3, int en4)
Definition: ncmesh.cpp:295
void RefineElement(int elem, char ref_type)
Definition: ncmesh.cpp:859
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:453
Array< int > boundary_faces
subset of all faces, set by BuildFaceList
Definition: ncmesh.hpp:499
int child[8]
2-8 children (if ref_type != 0)
Definition: ncmesh.hpp:448
const NCList & GetNCList(int entity)
Return vertex/edge/face list (entity = 0/1/2, respectively).
Definition: ncmesh.hpp:234
Array< int > leaf_elements
Definition: ncmesh.hpp:492
static PointMatrix pm_quad_identity
Definition: ncmesh.hpp:774
Master(int index, int element, char local, char geom, int sb, int se)
Definition: ncmesh.hpp:172
int EdgeSplitLevel(int vn1, int vn2) const
Definition: ncmesh.cpp:4676
void Initialize(const mfem::Element *elem)
Definition: ncmesh.cpp:30
void FindVertexCousins(int elem, int local, Array< int > &cousins) const
Definition: ncmesh.cpp:3375
signed char geom
Definition: ncmesh.hpp:156
int QuadFaceSplitType(int v1, int v2, int v3, int v4, int mid[5]=NULL) const
Definition: ncmesh.cpp:2186
Array< int > free_element_ids
Definition: ncmesh.hpp:463
int index(int i, int j, int nx, int ny)
Definition: life.cpp:241
int parent
Element index in the coarse mesh.
Definition: ncmesh.hpp:49
Point(double x, double y)
Definition: ncmesh.hpp:706
virtual void AssignLeafIndices()
Definition: ncmesh.cpp:1918
static PointMatrix pm_hex_identity
Definition: ncmesh.hpp:777
int Dimension() const
Definition: ncmesh.hpp:112
void FindEdgeElements(int vn1, int vn2, int vn3, int vn4, Array< MeshId > &prisms) const
Definition: ncmesh.cpp:675
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
Definition: ncmesh.cpp:1748
T & Last()
Return the last element in the array.
Definition: array.hpp:740
void NeighborExpand(const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
Definition: ncmesh.cpp:3237
BlockArray< Element > elements
Definition: ncmesh.hpp:462
void FindSetNeighbors(const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
Definition: ncmesh.cpp:3051
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:4422
virtual void Update()
Definition: ncmesh.cpp:213
bool HaveTets() const
Definition: ncmesh.hpp:525
Data type point element.
Definition: point.hpp:22
long MemoryUsage() const
Return total number of bytes allocated.
Definition: ncmesh.cpp:5114
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:4951
HashTable< Node > shadow
temporary storage for reparented nodes
Definition: ncmesh.hpp:531
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:60
virtual void BuildVertexList()
Definition: ncmesh.cpp:2757
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition: ncmesh.hpp:219
void UpdateElementToVertexTable()
Definition: ncmesh.hpp:682
friend struct CompareRanks
Definition: ncmesh.hpp:856
virtual void ElementSharesEdge(int elem, int local, int enode)
Definition: ncmesh.hpp:643
int GetElementDepth(int i) const
Return the distance of leaf &#39;i&#39; from the root.
Definition: ncmesh.cpp:4532
HashTable< Face >::iterator face_iterator
Definition: ncmesh.hpp:474
int FindNodeExt(const Element &el, int node, bool abort=true)
Extended version of find_node: works if &#39;el&#39; is refined.
Definition: ncmesh.cpp:2260
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
Definition: ncmesh.hpp:439
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:66
void CollectQuadFaceVertices(int v0, int v1, int v2, int v3, Array< int > &indices)
Definition: ncmesh.cpp:2937
bool Unused() const
Definition: ncmesh.hpp:423
static int find_node(const Element &el, int node)
Definition: ncmesh.cpp:2250
void UnreferenceElement(int elem, Array< int > &elemFaces)
Definition: ncmesh.cpp:334
int index
Mesh number.
Definition: ncmesh.hpp:153
void LoadVertexParents(std::istream &input)
Definition: ncmesh.cpp:4882
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:725
Abstract data type element.
Definition: element.hpp:28
Point & operator=(const Point &src)
Definition: ncmesh.hpp:732
Refinement()=default
void PrintStats(std::ostream &out=mfem::out) const
Definition: ncmesh.cpp:5162
int GetSingleElement() const
Return one of elem[0] or elem[1] and make sure the other is -1.
Definition: ncmesh.cpp:424
const double * CalcVertexPos(int node) const
Definition: ncmesh.cpp:2004
int rank
processor number (ParNCMesh), -1 if undefined/unknown
Definition: ncmesh.hpp:443
int node[8]
element corners (if ref_type == 0)
Definition: ncmesh.hpp:447
int GetVertexRootCoord(int elem, RefCoord coord[3]) const
Definition: ncmesh.cpp:3304
DenseMatrix point_matrix
position within the master edge/face
Definition: ncmesh.hpp:182
HashTable< Face > faces
Definition: ncmesh.hpp:460
Defines the position of a fine element within a coarse element.
Definition: ncmesh.hpp:46
HashTable< Face >::const_iterator face_const_iterator
Definition: ncmesh.hpp:476
virtual void Refine(const Array< Refinement > &refinements)
Definition: ncmesh.cpp:1492
int matrix
Index into the DenseTensor corresponding to the parent Geometry::Type stored in CoarseFineTransformat...
Definition: ncmesh.hpp:52
int64_t RefCoord
Definition: ncmesh.hpp:365
Array< char > face_geom
face geometry by face index, set by OnMeshUpdated
Definition: ncmesh.hpp:500
int GetEdgeMaster(int v1, int v2) const
Definition: ncmesh.cpp:4523