MFEM  v4.3.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-2021, 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,
72  bool get_coarse_to_fine_only = false) const;
73 
74  void GetCoarseToFineMap(const Mesh &fine_mesh,
75  Table &coarse_to_fine) const;
76 
77  void Clear();
78  bool IsInitialized() const;
79  long MemoryUsage() const;
80 };
81 
83 
84 struct MatrixMap; // for internal use
85 
86 
87 /** \brief A class for non-conforming AMR. The class is not used directly
88  * by the user, rather it is an extension of the Mesh class.
89  *
90  * In general, the class is used by MFEM as follows:
91  *
92  * 1. NCMesh is constructed from elements of an existing Mesh. The elements
93  * are copied and become roots of the refinement hierarchy.
94  *
95  * 2. Some elements are refined with the Refine() method. Both isotropic and
96  * anisotropic refinements of quads/hexes are supported.
97  *
98  * 3. A new Mesh is created from NCMesh containing the leaf elements.
99  * This new Mesh may have non-conforming (hanging) edges and faces and
100  * is the one seen by the user.
101  *
102  * 4. FiniteElementSpace asks NCMesh for a list of conforming, master and
103  * slave edges/faces and creates the conforming interpolation matrix P.
104  *
105  * 5. A continuous/conforming solution is obtained by solving P'*A*P x = P'*b.
106  *
107  * 6. Repeat from step 2.
108  */
109 class NCMesh
110 {
111 public:
112  //// Initialize with elements from an existing 'mesh'.
113  explicit NCMesh(const Mesh *mesh);
114 
115  /** Load from a stream. The id header is assumed to have been read already
116  from \param[in] input . \param[in] version is 10 for the v1.0 NC format,
117  or 1 for the legacy v1.1 format. \param[out] curved is set to 1 if the
118  curvature GridFunction follows after mesh data. \param[out] is_nc (again
119  treated as a boolean) is set to 0 if the legacy v1.1 format in fact
120  defines a conforming mesh. See Mesh::Loader for details. */
121  NCMesh(std::istream &input, int version, int &curved, int &is_nc);
122 
123  /// Deep copy of another instance.
124  NCMesh(const NCMesh &other);
125 
126  virtual ~NCMesh();
127 
128  int Dimension() const { return Dim; }
129  int SpaceDimension() const { return spaceDim; }
130 
131  int GetNVertices() const { return NVertices; }
132  int GetNEdges() const { return NEdges; }
133  int GetNFaces() const { return NFaces; }
134 
135  /** Perform the given batch of refinements. Please note that in the presence
136  of anisotropic splits additional refinements may be necessary to keep
137  the mesh consistent. However, the function always performs at least the
138  requested refinements. */
139  virtual void Refine(const Array<Refinement> &refinements);
140 
141  /** Check the mesh and potentially refine some elements so that the maximum
142  difference of refinement levels between adjacent elements is not greater
143  than 'max_nc_level'. */
144  virtual void LimitNCLevel(int max_nc_level);
145 
146  /** Return a list of derefinement opportunities. Each row of the table
147  contains Mesh indices of existing elements that can be derefined to form
148  a single new coarse element. Row numbers are then passed to Derefine.
149  This function works both in serial and parallel. */
150  const Table &GetDerefinementTable();
151 
152  /** Check derefinements returned by GetDerefinementTable and mark those that
153  can be done safely so that the maximum NC level condition is not violated.
154  On return, level_ok.Size() == deref_table.Size() and contains 0/1s. */
155  virtual void CheckDerefinementNCLevel(const Table &deref_table,
156  Array<int> &level_ok, int max_nc_level);
157 
158  /** Perform a subset of the possible derefinements (see GetDerefinementTable).
159  Note that if anisotropic refinements are present in the mesh, some of the
160  derefinements may have to be skipped to preserve mesh consistency. */
161  virtual void Derefine(const Array<int> &derefs);
162 
163 
164  // master/slave lists
165 
166  /// Identifies a vertex/edge/face in both Mesh and NCMesh.
167  struct MeshId
168  {
169  int index; ///< Mesh number
170  int element; ///< NCMesh::Element containing this vertex/edge/face
171  signed char local; ///< local number within 'element'
172  signed char geom; ///< Geometry::Type (faces only) (char to save RAM)
173 
174  Geometry::Type Geom() const { return Geometry::Type(geom); }
175 
176  MeshId() = default;
177  MeshId(int index, int element, int local, int geom = -1)
178  : index(index), element(element), local(local), geom(geom) {}
179  };
180 
181  /** Nonconforming edge/face that has more than one neighbor. The neighbors
182  are stored in NCList::slaves[i], slaves_begin <= i < slaves_end. */
183  struct Master : public MeshId
184  {
185  int slaves_begin, slaves_end; ///< slave faces
186 
187  Master() = default;
188  Master(int index, int element, int local, int geom, int sb, int se)
189  : MeshId(index, element, local, geom)
190  , slaves_begin(sb), slaves_end(se) {}
191  };
192 
193  /// Nonconforming edge/face within a bigger edge/face.
194  struct Slave : public MeshId
195  {
196  int master; ///< master number (in Mesh numbering)
197  unsigned matrix : 24; ///< index into NCList::point_matrices[geom]
198  unsigned edge_flags : 8; ///< orientation flags, see OrientedPointMatrix
199 
200  Slave() = default;
201  Slave(int index, int element, int local, int geom)
202  : MeshId(index, element, local, geom)
203  , master(-1), matrix(0), edge_flags(0) {}
204  };
205 
206  /// Lists all edges/faces in the nonconforming mesh.
207  struct NCList
208  {
212 
213  /// List of unique point matrices for each slave geometry.
215 
216  /// Return the point matrix oriented according to the master and slave edges
217  void OrientedPointMatrix(const Slave &slave,
218  DenseMatrix &oriented_matrix) const;
219 
220  void Clear();
221  bool Empty() const { return !conforming.Size() && !masters.Size(); }
222  long TotalSize() const;
223  long MemoryUsage() const;
224 
225  const MeshId& LookUp(int index, int *type = NULL) const;
226 
227  ~NCList() { Clear(); }
228  private:
229  mutable Array<int> inv_index;
230  };
231 
232  /// Return the current list of conforming and nonconforming faces.
234  {
235  if (face_list.Empty()) { BuildFaceList(); }
236  return face_list;
237  }
238 
239  /// Return the current list of conforming and nonconforming edges.
241  {
242  if (edge_list.Empty()) { BuildEdgeList(); }
243  return edge_list;
244  }
245 
246  /** Return a list of vertices (in 'conforming'); this function is provided
247  for uniformity/completeness. Needed in ParNCMesh/ParFESpace. */
249  {
250  if (vertex_list.Empty()) { BuildVertexList(); }
251  return vertex_list;
252  }
253 
254  /// Return vertex/edge/face list (entity = 0/1/2, respectively).
255  const NCList& GetNCList(int entity)
256  {
257  switch (entity)
258  {
259  case 0: return GetVertexList();
260  case 1: return GetEdgeList();
261  default: return GetFaceList();
262  }
263  }
264 
265 
266  // coarse/fine transforms
267 
268  /** Remember the current layer of leaf elements before the mesh is refined.
269  Needed by GetRefinementTransforms(), must be called before Refine(). */
270  void MarkCoarseLevel();
271 
272  /** After refinement, calculate the relation of each fine element to its
273  parent coarse element. Note that Refine() or LimitNCLevel() can be called
274  multiple times between MarkCoarseLevel() and this function. */
276 
277  /** After derefinement, calculate the relations of previous fine elements
278  (some of which may no longer exist) to the current leaf elements.
279  Unlike for refinement, Derefine() may only be called once before this
280  function so there is no MarkFineLevel(). */
282 
283  /// Free all internal data created by the above three functions.
284  void ClearTransforms();
285 
286 
287  // grid ordering
288 
289  /** Return a space filling curve for a rectangular grid of elements.
290  Implemented is a generalized Hilbert curve for arbitrary grid dimensions.
291  If the width is odd, height should be odd too, otherwise one diagonal
292  (vertex-neighbor) step cannot be avoided in the curve. Even dimensions
293  are recommended. */
294  static void GridSfcOrdering2D(int width, int height,
295  Array<int> &coords);
296 
297  /** Return a space filling curve for a 3D rectangular grid of elements.
298  The Hilbert-curve-like algorithm works well for even dimensions. For odd
299  width/height/depth it tends to produce some diagonal (edge-neighbor)
300  steps. Even dimensions are recommended. */
301  static void GridSfcOrdering3D(int width, int height, int depth,
302  Array<int> &coords);
303 
304 
305  // utility
306 
307  /// Return Mesh vertex indices of an edge identified by 'edge_id'.
308  void GetEdgeVertices(const MeshId &edge_id, int vert_index[2],
309  bool oriented = true) const;
310 
311  /** Return "NC" orientation of an edge. As opposed to standard Mesh edge
312  orientation based on vertex IDs, "NC" edge orientation follows the local
313  edge orientation within the element 'edge_id.element' and is thus
314  processor independent. TODO: this seems only partially true? */
315  int GetEdgeNCOrientation(const MeshId &edge_id) const;
316 
317  /** Return Mesh vertex and edge indices of a face identified by 'face_id'.
318  The return value is the number of face vertices. */
319  int GetFaceVerticesEdges(const MeshId &face_id,
320  int vert_index[4], int edge_index[4],
321  int edge_orientation[4]) const;
322 
323  /** Given an edge (by its vertex indices v1 and v2) return the first
324  (geometric) parent edge that exists in the Mesh or -1 if there is no such
325  parent. */
326  int GetEdgeMaster(int v1, int v2) const;
327 
328  /** Get a list of vertices (2D/3D) and edges (3D) that coincide with boundary
329  elements with the specified attributes (marked in 'bdr_attr_is_ess').
330  In 3D this function also reveals "hidden" boundary edges. In parallel it
331  helps identifying boundary vertices/edges affected by non-local boundary
332  elements. */
333  virtual void GetBoundaryClosure(const Array<int> &bdr_attr_is_ess,
334  Array<int> &bdr_vertices,
335  Array<int> &bdr_edges);
336 
337  /// Return element geometry type. @a index is the Mesh element number.
339  { return elements[leaf_elements[index]].Geom(); }
340 
341  /// Return face geometry type. @a index is the Mesh face number.
343  { return Geometry::Type(face_geom[index]); }
344 
345  /// Return the number of root elements.
346  int GetNumRootElements() { return root_state.Size(); }
347 
348  /// Return the distance of leaf 'i' from the root.
349  int GetElementDepth(int i) const;
350 
351  /** Return the size reduction compared to the root element (ignoring local
352  stretching and curvature). */
353  int GetElementSizeReduction(int i) const;
354 
355  /// Return the faces and face attributes of leaf element 'i'.
357  Array<int> &fattr) const;
358 
359 
360  /// I/O: Print the mesh in "MFEM NC mesh v1.0" format.
361  void Print(std::ostream &out) const;
362 
363  /// I/O: Return true if the mesh was loaded from the legacy v1.1 format.
364  bool IsLegacyLoaded() const { return Legacy; }
365 
366  /// I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
367  void LegacyToNewVertexOrdering(Array<int> &order) const;
368 
369  /// Save memory by releasing all non-essential and cached data.
370  virtual void Trim();
371 
372  /// Return total number of bytes allocated.
373  long MemoryUsage() const;
374 
375  int PrintMemoryDetail() const;
376 
377  typedef std::int64_t RefCoord;
378 
379 
380 protected: // non-public interface for the Mesh class
381 
382  friend class Mesh;
383 
384  /// Fill Mesh::{vertices,elements,boundary} for the current finest level.
385  void GetMeshComponents(Mesh &mesh) const;
386 
387  /** Get edge and face numbering from 'mesh' (i.e., set all Edge::index and
388  Face::index) after a new mesh was created from us. */
389  void OnMeshUpdated(Mesh *mesh);
390 
391  /** Delete top-level vertex coordinates if the Mesh became curved, e.g.,
392  by calling Mesh::SetCurvature or otherwise setting the Nodes. */
394 
395 
396 protected: // implementation
397 
398  int Dim, spaceDim; ///< dimensions of the elements and the vertex coordinates
399  int MyRank; ///< used in parallel, or when loading a parallel file in serial
400  bool Iso; ///< true if the mesh only contains isotropic refinements
401  int Geoms; ///< bit mask of element geometries present, see InitGeomFlags()
402  bool Legacy; ///< true if the mesh was loaded from the legacy v1.1 format
403 
404  /** A Node can hold a vertex, an edge, or both. Elements directly point to
405  their corner nodes, but edge nodes also exist and can be accessed using
406  a hash-table given their two end-point node IDs. All nodes can be
407  accessed in this way, with the exception of top-level vertex nodes.
408  When an element is being refined, the mid-edge nodes are readily
409  available with this mechanism. The new elements "sign in" to the nodes
410  by increasing the reference counts of their vertices and edges. The
411  parent element "signs off" its nodes by decrementing the ref counts. */
412  struct Node : public Hashed2
413  {
416 
417  Node() : vert_refc(0), edge_refc(0), vert_index(-1), edge_index(-1) {}
418  ~Node();
419 
420  bool HasVertex() const { return vert_refc > 0; }
421  bool HasEdge() const { return edge_refc > 0; }
422 
423  // decrease vertex/edge ref count, return false if Node should be deleted
424  bool UnrefVertex() { --vert_refc; return vert_refc || edge_refc; }
425  bool UnrefEdge() { --edge_refc; return vert_refc || edge_refc; }
426  };
427 
428  /** Similarly to nodes, faces can be accessed by hashing their four vertex
429  node IDs. A face knows about the one or two elements that are using it.
430  A face that is not on the boundary and only has one element referencing
431  it is either a master or a slave face. */
432  struct Face : public Hashed4
433  {
434  int attribute; ///< boundary element attribute, -1 if internal face
435  int index; ///< face number in the Mesh
436  int elem[2]; ///< up to 2 elements sharing the face
437 
438  Face() : attribute(-1), index(-1) { elem[0] = elem[1] = -1; }
439 
440  bool Boundary() const { return attribute >= 0; }
441  bool Unused() const { return elem[0] < 0 && elem[1] < 0; }
442 
443  // add or remove an element from the 'elem[2]' array
444  void RegisterElement(int e);
445  void ForgetElement(int e);
446 
447  /// Return one of elem[0] or elem[1] and make sure the other is -1.
448  int GetSingleElement() const;
449  };
450 
451  /** This is an element in the refinement hierarchy. Each element has
452  either been refined and points to its children, or is a leaf and points
453  to its vertex nodes. */
454  struct Element
455  {
456  char geom; ///< Geometry::Type of the element (char for storage only)
457  char ref_type; ///< bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
458  char tet_type; ///< tetrahedron split type, currently always 0
459  char flag; ///< generic flag/marker, can be used by algorithms
460  int index; ///< element number in the Mesh, -1 if refined
461  int rank; ///< processor number (ParNCMesh), -1 if undefined/unknown
463  union
464  {
465  int node[8]; ///< element corners (if ref_type == 0)
466  int child[8]; ///< 2-8 children (if ref_type != 0)
467  };
468  int parent; ///< parent element, -1 if this is a root element, -2 if free'd
469 
470  Element(Geometry::Type geom, int attr);
471 
472  Geometry::Type Geom() const { return Geometry::Type(geom); }
473  bool IsLeaf() const { return !ref_type && (parent != -2); }
474  };
475 
476 
477  // primary data
478 
479  HashTable<Node> nodes; // associative container holding all Nodes
480  HashTable<Face> faces; // associative container holding all Faces
481 
482  BlockArray<Element> elements; // storage for all Elements
483  Array<int> free_element_ids; // unused element ids - indices into 'elements'
484 
485  /** Initial traversal state (~ element orientation) for each root element
486  NOTE: M = root_state.Size() is the number of root elements.
487  NOTE: the first M items of 'elements' is the coarse mesh. */
489 
490  /** Coordinates of top-level vertices (organized as triples). If empty,
491  the Mesh is curved (Nodes != NULL) and NCMesh is topology-only. */
493 
494 
495  // secondary data
496 
497  /** Apart from the primary data structure, which is the element/node/face
498  hierarchy, there is secondary data that is derived from the primary
499  data and needs to be updated when the primary data changes. Update()
500  takes care of that and needs to be called after each refinement and
501  derefinement. */
502  virtual void Update();
503 
504  // set by UpdateLeafElements, UpdateVertices and OnMeshUpdated
506 
507  // NOTE: the serial code understands the bare minimum about ghost elements and
508  // other ghost entities in order to be able to load parallel partial meshes
510 
511  Array<int> leaf_elements; ///< finest elements, in Mesh ordering (+ ghosts)
512  Array<int> leaf_sfc_index; ///< natural tree ordering of leaf elements
513  Array<int> vertex_nodeId; ///< vertex-index to node-id map, see UpdateVertices
514 
515  NCList face_list; ///< lazy-initialized list of faces, see GetFaceList
516  NCList edge_list; ///< lazy-initialized list of edges, see GetEdgeList
517  NCList vertex_list; ///< lazy-initialized list of vertices, see GetVertexList
518 
519  Array<int> boundary_faces; ///< subset of all faces, set by BuildFaceList
520  Array<char> face_geom; ///< face geometry by face index, set by OnMeshUpdated
521 
522  Table element_vertex; ///< leaf-element to vertex table, see FindSetNeighbors
523 
524 
525  void UpdateLeafElements();
526  void UpdateVertices(); ///< update Vertex::index and vertex_nodeId
527  void CollectLeafElements(int elem, int state, Array<int> &ghosts,
528  int &counter);
529 
530  /** Try to find a space-filling curve friendly orientation of the root
531  elements: set 'root_state' based on the ordering of coarse elements.
532  Note that the coarse mesh itself must be ordered as an SFC by e.g.
533  Mesh::GetGeckoElementOrdering. */
534  void InitRootState(int root_count);
535 
536  void InitGeomFlags();
537 
538  bool HavePrisms() const { return Geoms & (1 << Geometry::PRISM); }
539  bool HaveTets() const { return Geoms & (1 << Geometry::TETRAHEDRON); }
540 
541  bool IsGhost(const Element &el) const { return el.rank != MyRank; }
542 
543 
544  // refinement/derefinement
545 
546  Array<Refinement> ref_stack; ///< stack of scheduled refinements (temporary)
547  HashTable<Node> shadow; ///< temporary storage for reparented nodes
548  Array<Triple<int, int, int> > reparents; ///< scheduled node reparents (tmp)
549 
550  Table derefinements; ///< possible derefinements, see GetDerefinementTable
551 
552  void RefineElement(int elem, char ref_type);
553  void DerefineElement(int elem);
554 
555  int AddElement(const Element &el)
556  {
557  if (free_element_ids.Size())
558  {
559  int idx = free_element_ids.Last();
561  elements[idx] = el;
562  return idx;
563  }
564  return elements.Append(el);
565  }
566  void FreeElement(int id)
567  {
569  elements[id].ref_type = 0;
570  elements[id].parent = -2; // mark the element as free
571  }
572 
573  int NewHexahedron(int n0, int n1, int n2, int n3,
574  int n4, int n5, int n6, int n7, int attr,
575  int fattr0, int fattr1, int fattr2,
576  int fattr3, int fattr4, int fattr5);
577 
578  int NewWedge(int n0, int n1, int n2,
579  int n3, int n4, int n5, int attr,
580  int fattr0, int fattr1,
581  int fattr2, int fattr3, int fattr4);
582 
583  int NewTetrahedron(int n0, int n1, int n2, int n3, int attr,
584  int fattr0, int fattr1, int fattr2, int fattr3);
585 
586  int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr,
587  int eattr0, int eattr1, int eattr2, int eattr3);
588 
589  int NewTriangle(int n0, int n1, int n2,
590  int attr, int eattr0, int eattr1, int eattr2);
591 
592  int NewSegment(int n0, int n1, int attr, int vattr1, int vattr2);
593 
594  mfem::Element* NewMeshElement(int geom) const;
595 
596  int QuadFaceSplitType(int v1, int v2, int v3, int v4, int mid[5]
597  = NULL /*optional output of mid-edge nodes*/) const;
598 
599  bool TriFaceSplit(int v1, int v2, int v3, int mid[3] = NULL) const;
600 
601  void ForceRefinement(int vn1, int vn2, int vn3, int vn4);
602 
603  void FindEdgeElements(int vn1, int vn2, int vn3, int vn4,
604  Array<MeshId> &prisms) const;
605 
606  void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4,
607  const Refinement *refs, int nref);
608 
609  void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4,
610  int mid12, int mid34, int level = 0);
611 
612  void CheckIsoFace(int vn1, int vn2, int vn3, int vn4,
613  int en1, int en2, int en3, int en4, int midf);
614 
615  void ReparentNode(int node, int new_p1, int new_p2);
616 
617  int FindMidEdgeNode(int node1, int node2) const;
618  int GetMidEdgeNode(int node1, int node2);
619 
620  int GetMidFaceNode(int en1, int en2, int en3, int en4);
621 
622  void ReferenceElement(int elem);
623  void UnreferenceElement(int elem, Array<int> &elemFaces);
624 
625  Face* GetFace(Element &elem, int face_no);
626  void RegisterFaces(int elem, int *fattr = NULL);
627  void DeleteUnusedFaces(const Array<int> &elemFaces);
628 
629  void CollectDerefinements(int elem, Array<Connection> &list);
630 
631  /// Return el.node[index] correctly, even if the element is refined.
632  int RetrieveNode(const Element &el, int index);
633 
634  /// Extended version of find_node: works if 'el' is refined.
635  int FindNodeExt(const Element &el, int node, bool abort = true);
636 
637 
638  // face/edge lists
639 
640  static int find_node(const Element &el, int node);
641  static int find_element_edge(const Element &el, int vn0, int vn1,
642  bool abort = true);
643  static int find_local_face(int geom, int a, int b, int c);
644 
645  struct Point;
646  struct PointMatrix;
647 
648  int ReorderFacePointMat(int v0, int v1, int v2, int v3,
649  int elem, const PointMatrix &pm,
650  PointMatrix &reordered) const;
651 
652  void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3,
653  const PointMatrix& pm, int level, Face* eface[4],
654  MatrixMap &matrix_map);
655  bool TraverseTriFace(int vn0, int vn1, int vn2,
656  const PointMatrix& pm, int level,
657  MatrixMap &matrix_map);
658  void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1,
659  MatrixMap &matrix_map);
660  void TraverseEdge(int vn0, int vn1, double t0, double t1, int flags,
661  int level, MatrixMap &matrix_map);
662 
663  virtual void BuildFaceList();
664  virtual void BuildEdgeList();
665  virtual void BuildVertexList();
666 
667  virtual void ElementSharesFace(int elem, int local, int face) {} // ParNCMesh
668  virtual void ElementSharesEdge(int elem, int local, int enode) {} // ParNCMesh
669  virtual void ElementSharesVertex(int elem, int local, int vnode) {} // ParNCMesh
670 
671 
672  // neighbors / element_vertex table
673 
674  /** Return all vertex-, edge- and face-neighbors of a set of elements.
675  The neighbors are returned as a list (neighbors != NULL), as a set
676  (neighbor_set != NULL), or both. The sizes of the set arrays must match
677  that of leaf_elements. The function is intended to be used for large
678  sets of elements and its complexity is linear in the number of leaf
679  elements in the mesh. */
680  void FindSetNeighbors(const Array<char> &elem_set,
681  Array<int> *neighbors, /* append */
682  Array<char> *neighbor_set = NULL);
683 
684  /** Return all vertex-, edge- and face-neighbors of a single element.
685  You can limit the number of elements being checked using 'search_set'.
686  The complexity of the function is linear in the size of the search set.*/
687  void FindNeighbors(int elem,
688  Array<int> &neighbors, /* append */
689  const Array<int> *search_set = NULL);
690 
691  /** Expand a set of elements by all vertex-, edge- and face-neighbors.
692  The output array 'expanded' will contain all items from 'elems'
693  (provided they are in 'search_set') plus their neighbors. The neighbor
694  search can be limited to the optional search set. The complexity is
695  linear in the sum of the sizes of 'elems' and 'search_set'. */
696  void NeighborExpand(const Array<int> &elems,
697  Array<int> &expanded,
698  const Array<int> *search_set = NULL);
699 
700 
701  void CollectEdgeVertices(int v0, int v1, Array<int> &indices);
702  void CollectTriFaceVertices(int v0, int v1, int v2, Array<int> &indices);
703  void CollectQuadFaceVertices(int v0, int v1, int v2, int v3,
704  Array<int> &indices);
706 
708  {
710  }
711 
712  int GetVertexRootCoord(int elem, RefCoord coord[3]) const;
713  void CollectIncidentElements(int elem, const RefCoord coord[3],
714  Array<int> &list) const;
715 
716  /** Return elements neighboring to a local vertex of element 'elem'. Only
717  elements from within the same refinement tree ('cousins') are returned.
718  Complexity is proportional to the depth of elem's refinement tree. */
719  void FindVertexCousins(int elem, int local, Array<int> &cousins) const;
720 
721 
722  // coarse/fine transformations
723 
724  struct Point
725  {
726  int dim;
727  double coord[3];
728 
729  Point() { dim = 0; }
730 
731  Point(double x)
732  { dim = 1; coord[0] = x; }
733 
734  Point(double x, double y)
735  { dim = 2; coord[0] = x; coord[1] = y; }
736 
737  Point(double x, double y, double z)
738  { dim = 3; coord[0] = x; coord[1] = y; coord[2] = z; }
739 
740  Point(const Point& p0, const Point& p1)
741  {
742  dim = p0.dim;
743  for (int i = 0; i < dim; i++)
744  {
745  coord[i] = (p0.coord[i] + p1.coord[i]) * 0.5;
746  }
747  }
748 
749  Point(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
750  {
751  dim = p0.dim;
752  MFEM_ASSERT(p1.dim == dim && p2.dim == dim && p3.dim == dim, "");
753  for (int i = 0; i < dim; i++)
754  {
755  coord[i] = (p0.coord[i] + p1.coord[i] + p2.coord[i] + p3.coord[i])
756  * 0.25;
757  }
758  }
759 
760  Point& operator=(const Point& src)
761  {
762  dim = src.dim;
763  for (int i = 0; i < dim; i++) { coord[i] = src.coord[i]; }
764  return *this;
765  }
766  };
767 
768  struct PointMatrix
769  {
770  int np;
772 
773  PointMatrix() : np(0) {}
774 
775  PointMatrix(const Point& p0, const Point& p1)
776  { np = 2; points[0] = p0; points[1] = p1; }
777 
778  PointMatrix(const Point& p0, const Point& p1, const Point& p2)
779  { np = 3; points[0] = p0; points[1] = p1; points[2] = p2; }
780 
781  PointMatrix(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
782  { np = 4; points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3; }
783 
784  PointMatrix(const Point& p0, const Point& p1, const Point& p2,
785  const Point& p3, const Point& p4, const Point& p5)
786  {
787  np = 6;
788  points[0] = p0; points[1] = p1; points[2] = p2;
789  points[3] = p3; points[4] = p4; points[5] = p5;
790  }
791  PointMatrix(const Point& p0, const Point& p1, const Point& p2,
792  const Point& p3, const Point& p4, const Point& p5,
793  const Point& p6, const Point& p7)
794  {
795  np = 8;
796  points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3;
797  points[4] = p4; points[5] = p5; points[6] = p6; points[7] = p7;
798  }
799 
800  Point& operator()(int i) { return points[i]; }
801  const Point& operator()(int i) const { return points[i]; }
802 
803  bool operator==(const PointMatrix &pm) const;
804 
805  void GetMatrix(DenseMatrix& point_matrix) const;
806  };
807 
814 
816 
817  void GetPointMatrix(Geometry::Type geom, const char* ref_path,
818  DenseMatrix& matrix);
819 
820  typedef std::map<std::string, int> RefPathMap;
821 
822  void TraverseRefinements(int elem, int coarse_index,
823  std::string &ref_path, RefPathMap &map);
824 
825  /// storage for data returned by Get[De]RefinementTransforms()
827 
828  /// state of leaf_elements before Refine(), set by MarkCoarseLevel()
830 
831  void InitDerefTransforms();
832  void SetDerefMatrixCodes(int parent, Array<int> &fine_coarse);
833 
834 
835  // vertex temporary data, used by GetMeshComponents
836 
837  struct TmpVertex
838  {
839  bool valid, visited;
840  double pos[3];
841  TmpVertex() : valid(false), visited(false) {}
842  };
843 
845 
846  const double *CalcVertexPos(int node) const;
847 
848 
849  // utility
850 
851  int GetEdgeMaster(int node) const;
852 
853  void FindFaceNodes(int face, int node[4]);
854 
855  int EdgeSplitLevel(int vn1, int vn2) const;
856  int TriFaceSplitLevel(int vn1, int vn2, int vn3) const;
857  void QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4,
858  int& h_level, int& v_level) const;
859 
860  void CountSplits(int elem, int splits[3]) const;
861  void GetLimitRefinements(Array<Refinement> &refinements, int max_level);
862 
863 
864  // I/O
865 
866  /// Print the "vertex_parents" section of the mesh file.
867  int PrintVertexParents(std::ostream *out) const;
868  /// Load the vertex parent hierarchy from a mesh file.
869  void LoadVertexParents(std::istream &input);
870 
871  /** Print the "boundary" section of the mesh file.
872  If out == NULL, only return the number of boundary elements. */
873  int PrintBoundary(std::ostream *out) const;
874  /// Load the "boundary" section of the mesh file.
875  void LoadBoundary(std::istream &input);
876 
877  /// Print the "coordinates" section of the mesh file.
878  void PrintCoordinates(std::ostream &out) const;
879  /// Load the "coordinates" section of the mesh file.
880  void LoadCoordinates(std::istream &input);
881 
882  /// Count root elements and initialize root_state.
883  void InitRootElements();
884  /// Return the index of the last top-level node plus one.
885  int CountTopLevelNodes() const;
886  /// Return true if all root_states are zero.
887  bool ZeroRootStates() const;
888 
889  /// Load the element refinement hierarchy from a legacy mesh file.
890  void LoadCoarseElements(std::istream &input);
891  void CopyElements(int elem, const BlockArray<Element> &tmp_elements);
892  /// Load the deprecated MFEM mesh v1.1 format for backward compatibility.
893  void LoadLegacyFormat(std::istream &input, int &curved, int &is_nc);
894 
895 
896  // geometry
897 
898  /// This holds in one place the constants about the geometries we support
899  struct GeomInfo
900  {
901  int nv, ne, nf; // number of: vertices, edges, faces
902  int edges[12][2]; // edge vertices (up to 12 edges)
903  int faces[6][4]; // face vertices (up to 6 faces)
904  int nfv[6]; // number of face vertices
905 
907  GeomInfo() : initialized(false) {}
909  };
910 
912 
913 #ifdef MFEM_DEBUG
914 public:
915  void DebugLeafOrder(std::ostream &out) const;
916  void DebugDump(std::ostream &out) const;
917 #endif
918 
919  friend class ParNCMesh; // for ParNCMesh::ElementSet
920  friend struct MatrixMap;
921  friend struct PointMatrixHash;
922 };
923 
924 }
925 
926 #endif
NCList face_list
lazy-initialized list of faces, see GetFaceList
Definition: ncmesh.hpp:515
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
Definition: ncmesh.hpp:516
void LegacyToNewVertexOrdering(Array< int > &order) const
I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
Definition: ncmesh.cpp:6026
void LoadLegacyFormat(std::istream &input, int &curved, int &is_nc)
Load the deprecated MFEM mesh v1.1 format for backward compatibility.
Definition: ncmesh.cpp:5868
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
static void GridSfcOrdering3D(int width, int height, int depth, Array< int > &coords)
Definition: ncmesh.cpp:4829
void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4, const Refinement *refs, int nref)
Definition: ncmesh.cpp:750
CoarseFineTransformations transforms
storage for data returned by Get[De]RefinementTransforms()
Definition: ncmesh.hpp:826
static void GridSfcOrdering2D(int width, int height, Array< int > &coords)
Definition: ncmesh.cpp:4814
int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
Definition: ncmesh.cpp:553
static PointMatrix pm_seg_identity
Definition: ncmesh.hpp:808
int elem[2]
up to 2 elements sharing the face
Definition: ncmesh.hpp:436
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:4892
TmpVertex * tmp_vertex
Definition: ncmesh.hpp:844
bool TriFaceSplit(int v1, int v2, int v3, int mid[3]=NULL) const
Definition: ncmesh.cpp:2556
void SetDerefMatrixCodes(int parent, Array< int > &fine_coarse)
Definition: ncmesh.cpp:1873
const CoarseFineTransformations & GetDerefinementTransforms()
Definition: ncmesh.cpp:4366
char flag
generic flag/marker, can be used by algorithms
Definition: ncmesh.hpp:459
bool HasEdge() const
Definition: ncmesh.hpp:421
PointMatrix(const Point &p0, const Point &p1, const Point &p2)
Definition: ncmesh.hpp:778
This holds in one place the constants about the geometries we support.
Definition: ncmesh.hpp:899
Array< Triple< int, int, int > > reparents
scheduled node reparents (tmp)
Definition: ncmesh.hpp:548
Point(const Point &p0, const Point &p1)
Definition: ncmesh.hpp:740
static const int NumGeom
Definition: geom.hpp:41
Array< Slave > slaves
Definition: ncmesh.hpp:211
void CollectTriFaceVertices(int v0, int v1, int v2, Array< int > &indices)
Definition: ncmesh.cpp:3315
char tet_type
tetrahedron split type, currently always 0
Definition: ncmesh.hpp:458
int NewSegment(int n0, int n1, int attr, int vattr1, int vattr2)
Definition: ncmesh.cpp:605
virtual void LimitNCLevel(int max_nc_level)
Definition: ncmesh.cpp:5279
virtual void Trim()
Save memory by releasing all non-essential and cached data.
Definition: ncmesh.cpp:6050
int index
element number in the Mesh, -1 if refined
Definition: ncmesh.hpp:460
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:781
const Geometry::Type geom
Definition: ex1.cpp:40
std::int64_t RefCoord
Definition: ncmesh.hpp:377
Embedding()=default
void ForceRefinement(int vn1, int vn2, int vn3, int vn4)
Definition: ncmesh.cpp:645
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:791
Lists all edges/faces in the nonconforming mesh.
Definition: ncmesh.hpp:207
bool Legacy
true if the mesh was loaded from the legacy v1.1 format
Definition: ncmesh.hpp:402
unsigned matrix
index into NCList::point_matrices[geom]
Definition: ncmesh.hpp:197
int NVertices
Definition: ncmesh.hpp:505
int GetNVertices() const
Definition: ncmesh.hpp:131
void GetMeshComponents(Mesh &mesh) const
Fill Mesh::{vertices,elements,boundary} for the current finest level.
Definition: ncmesh.cpp:2269
void ClearTransforms()
Free all internal data created by the above three functions.
Definition: ncmesh.cpp:4557
Table derefinements
possible derefinements, see GetDerefinementTable
Definition: ncmesh.hpp:550
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Geometry::Type Geom() const
Definition: ncmesh.hpp:174
NCList vertex_list
lazy-initialized list of vertices, see GetVertexList
Definition: ncmesh.hpp:517
void InitDerefTransforms()
Definition: ncmesh.cpp:1858
static PointMatrix pm_prism_identity
Definition: ncmesh.hpp:812
void Print(std::ostream &out) const
I/O: Print the mesh in &quot;MFEM NC mesh v1.0&quot; format.
Definition: ncmesh.cpp:5476
char geom
Geometry::Type of the element (char for storage only)
Definition: ncmesh.hpp:456
void InitRootState(int root_count)
Definition: ncmesh.cpp:2167
long TotalSize() const
Definition: ncmesh.cpp:3233
int GetEdgeNCOrientation(const MeshId &edge_id) const
Definition: ncmesh.cpp:4880
int NewTetrahedron(int n0, int n1, int n2, int n3, int attr, int fattr0, int fattr1, int fattr2, int fattr3)
Definition: ncmesh.cpp:527
void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4, int mid12, int mid34, int level=0)
Definition: ncmesh.cpp:775
bool ZeroRootStates() const
Return true if all root_states are zero.
Definition: ncmesh.cpp:5467
static PointMatrix pm_tri_identity
Definition: ncmesh.hpp:809
void FreeElement(int id)
Definition: ncmesh.hpp:566
void CountSplits(int elem, int splits[3]) const
Definition: ncmesh.cpp:5166
Geometry::Type Geom() const
Definition: ncmesh.hpp:472
void CollectDerefinements(int elem, Array< Connection > &list)
Definition: ncmesh.cpp:1745
const Point & operator()(int i) const
Definition: ncmesh.hpp:801
Point & operator()(int i)
Definition: ncmesh.hpp:800
void GetLimitRefinements(Array< Refinement > &refinements, int max_level)
Definition: ncmesh.cpp:5249
int PrintMemoryDetail() const
Definition: ncmesh.cpp:6116
void TraverseEdge(int vn0, int vn1, double t0, double t1, int flags, int level, MatrixMap &matrix_map)
Definition: ncmesh.cpp:3027
double coord[3]
Definition: ncmesh.hpp:727
bool UnrefVertex()
Definition: ncmesh.hpp:424
int spaceDim
dimensions of the elements and the vertex coordinates
Definition: ncmesh.hpp:398
bool HavePrisms() const
Definition: ncmesh.hpp:538
Array< int > vertex_nodeId
vertex-index to node-id map, see UpdateVertices
Definition: ncmesh.hpp:513
int attribute
boundary element attribute, -1 if internal face
Definition: ncmesh.hpp:434
void DebugDump(std::ostream &out) const
Definition: ncmesh.cpp:6169
void FindFaceNodes(int face, int node[4])
Definition: ncmesh.cpp:5018
Array< int > root_state
Definition: ncmesh.hpp:488
static const PointMatrix & GetGeomIdentity(Geometry::Type geom)
Definition: ncmesh.cpp:3833
void DeleteAll()
Delete the whole array.
Definition: array.hpp:841
static int find_element_edge(const Element &el, int vn0, int vn1, bool abort=true)
Definition: ncmesh.cpp:2593
int index
Mesh element number.
Definition: ncmesh.hpp:37
int master
master number (in Mesh numbering)
Definition: ncmesh.hpp:196
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition: ncmesh.hpp:233
const MeshId & LookUp(int index, int *type=NULL) const
Definition: ncmesh.cpp:3238
void OnMeshUpdated(Mesh *mesh)
Definition: ncmesh.cpp:2365
void DebugLeafOrder(std::ostream &out) const
Definition: ncmesh.cpp:6144
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:495
void CollectEdgeVertices(int v0, int v1, Array< int > &indices)
Definition: ncmesh.cpp:3303
void CollectLeafElements(int elem, int state, Array< int > &ghosts, int &counter)
Definition: ncmesh.cpp:1892
static int find_local_face(int geom, int a, int b, int c)
Definition: ncmesh.cpp:2611
Array< Refinement > ref_stack
stack of scheduled refinements (temporary)
Definition: ncmesh.hpp:546
bool Empty() const
Definition: ncmesh.hpp:221
Point(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition: ncmesh.hpp:749
void CollectIncidentElements(int elem, const RefCoord coord[3], Array< int > &list) const
Definition: ncmesh.cpp:3746
void ReferenceElement(int elem)
Definition: ncmesh.cpp:315
void GetPointMatrix(Geometry::Type geom, const char *ref_path, DenseMatrix &matrix)
Definition: ncmesh.cpp:3849
void CheckIsoFace(int vn1, int vn2, int vn3, int vn4, int en1, int en2, int en3, int en4, int midf)
Definition: ncmesh.cpp:870
void InitGeomFlags()
Definition: ncmesh.cpp:222
Face * GetFace(Element &elem, int face_no)
Definition: ncmesh.cpp:428
void TraverseRefinements(int elem, int coarse_index, std::string &ref_path, RefPathMap &map)
Definition: ncmesh.cpp:4283
int PrintBoundary(std::ostream *out) const
Definition: ncmesh.cpp:5349
int GetMidEdgeNode(int node1, int node2)
Definition: ncmesh.cpp:300
Element(Geometry::Type geom, int attr)
Definition: ncmesh.cpp:453
Refinement(int index, int type=7)
Definition: ncmesh.hpp:42
int index
face number in the Mesh
Definition: ncmesh.hpp:435
const Table & GetDerefinementTable()
Definition: ncmesh.cpp:1778
Table element_vertex
leaf-element to vertex table, see FindSetNeighbors
Definition: ncmesh.hpp:522
DenseTensor point_matrices[Geometry::NumGeom]
Matrices for IsoparametricTransformation organized by Geometry::Type.
Definition: ncmesh.hpp:63
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:746
int CountTopLevelNodes() const
Return the index of the last top-level node plus one.
Definition: ncmesh.cpp:5595
bool operator==(const PointMatrix &pm) const
Definition: ncmesh.cpp:3786
void GetMatrix(DenseMatrix &point_matrix) const
Definition: ncmesh.cpp:3800
int AddElement(const Element &el)
Definition: ncmesh.hpp:555
virtual void Derefine(const Array< int > &derefs)
Definition: ncmesh.cpp:1822
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:5139
void FindNeighbors(int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
Definition: ncmesh.cpp:3556
Geometry::Type GetElementGeometry(int index) const
Return element geometry type. index is the Mesh element number.
Definition: ncmesh.hpp:338
void InitGeom(Geometry::Type geom)
Definition: ncmesh.cpp:27
friend struct MatrixMap
Definition: ncmesh.hpp:920
Array< int > coarse_elements
state of leaf_elements before Refine(), set by MarkCoarseLevel()
Definition: ncmesh.hpp:829
void UpdateVertices()
update Vertex::index and vertex_nodeId
Definition: ncmesh.cpp:1987
int NGhostVertices
Definition: ncmesh.hpp:509
int ReorderFacePointMat(int v0, int v1, int v2, int v3, int elem, const PointMatrix &pm, PointMatrix &reordered) const
Definition: ncmesh.cpp:2689
Master(int index, int element, int local, int geom, int sb, int se)
Definition: ncmesh.hpp:188
Array< DenseMatrix * > point_matrices[Geometry::NumGeom]
List of unique point matrices for each slave geometry.
Definition: ncmesh.hpp:214
void RegisterFaces(int elem, int *fattr=NULL)
Definition: ncmesh.cpp:389
Array< Embedding > embeddings
Fine element positions in their parents.
Definition: ncmesh.hpp:65
virtual void BuildEdgeList()
Definition: ncmesh.cpp:3055
const CoarseFineTransformations & GetRefinementTransforms()
Definition: ncmesh.cpp:4315
void LoadCoarseElements(std::istream &input)
Load the element refinement hierarchy from a legacy mesh file.
Definition: ncmesh.cpp:5797
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition: ncmesh.hpp:167
int FindMidEdgeNode(int node1, int node2) const
Definition: ncmesh.cpp:284
int parent
parent element, -1 if this is a root element, -2 if free&#39;d
Definition: ncmesh.hpp:468
void MarkCoarseLevel()
Definition: ncmesh.cpp:4269
signed char local
local number within &#39;element&#39;
Definition: ncmesh.hpp:171
virtual void ElementSharesFace(int elem, int local, int face)
Definition: ncmesh.hpp:667
int GetNumRootElements()
Return the number of root elements.
Definition: ncmesh.hpp:346
void LoadCoordinates(std::istream &input)
Load the &quot;coordinates&quot; section of the mesh file.
Definition: ncmesh.cpp:5446
virtual void ElementSharesVertex(int elem, int local, int vnode)
Definition: ncmesh.hpp:669
int MyRank
used in parallel, or when loading a parallel file in serial
Definition: ncmesh.hpp:399
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
Definition: ncmesh.hpp:109
int GetElementSizeReduction(int i) const
Definition: ncmesh.cpp:4983
virtual ~NCMesh()
Definition: ncmesh.cpp:243
friend struct PointMatrixHash
Definition: ncmesh.hpp:921
int GetNFaces() const
Definition: ncmesh.hpp:133
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
Definition: ncmesh.cpp:5042
int element
NCMesh::Element containing this vertex/edge/face.
Definition: ncmesh.hpp:170
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:342
static GeomInfo GI[Geometry::NumGeom]
Definition: ncmesh.hpp:911
int slaves_end
slave faces
Definition: ncmesh.hpp:185
void MakeTopologyOnly()
Definition: ncmesh.hpp:393
void ReparentNode(int node, int new_p1, int new_p2)
Definition: ncmesh.cpp:268
int NewTriangle(int n0, int n1, int n2, int attr, int eattr0, int eattr1, int eattr2)
Definition: ncmesh.cpp:579
void DeleteUnusedFaces(const Array< int > &elemFaces)
Definition: ncmesh.cpp:403
bool Iso
true if the mesh only contains isotropic refinements
Definition: ncmesh.hpp:400
bool IsLegacyLoaded() const
I/O: Return true if the mesh was loaded from the legacy v1.1 format.
Definition: ncmesh.hpp:364
Point(double x, double y, double z)
Definition: ncmesh.hpp:737
bool HasVertex() const
Definition: ncmesh.hpp:420
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:625
std::map< std::string, int > RefPathMap
Definition: ncmesh.hpp:820
int SpaceDimension() const
Definition: ncmesh.hpp:129
virtual void BuildFaceList()
Definition: ncmesh.cpp:2935
Nonconforming edge/face within a bigger edge/face.
Definition: ncmesh.hpp:194
void DerefineElement(int elem)
Definition: ncmesh.cpp:1610
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:4997
const NCList & GetVertexList()
Definition: ncmesh.hpp:248
bool Boundary() const
Definition: ncmesh.hpp:440
bool IsLeaf() const
Definition: ncmesh.hpp:473
Point(double x)
Definition: ncmesh.hpp:731
void RegisterElement(int e)
Definition: ncmesh.cpp:414
mfem::Element * NewMeshElement(int geom) const
Definition: ncmesh.cpp:2230
Array< int > leaf_sfc_index
natural tree ordering of leaf elements
Definition: ncmesh.hpp:512
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4, const Point &p5)
Definition: ncmesh.hpp:784
void DeleteLast()
Delete the last entry of the array.
Definition: array.hpp:189
void UpdateLeafElements()
Definition: ncmesh.cpp:1960
int PrintVertexParents(std::ostream *out) const
Print the &quot;vertex_parents&quot; section of the mesh file.
Definition: ncmesh.cpp:5297
int TriFaceSplitLevel(int vn1, int vn2, int vn3) const
Definition: ncmesh.cpp:5122
int GetNEdges() const
Definition: ncmesh.hpp:132
int RetrieveNode(const Element &el, int index)
Return el.node[index] correctly, even if the element is refined.
Definition: ncmesh.cpp:1574
int Geoms
bit mask of element geometries present, see InitGeomFlags()
Definition: ncmesh.hpp:401
PointMatrix(const Point &p0, const Point &p1)
Definition: ncmesh.hpp:775
int NGhostElements
Definition: ncmesh.hpp:509
long MemoryUsage() const
Definition: ncmesh.cpp:6065
void BuildElementToVertexTable()
Definition: ncmesh.cpp:3373
HashTable< Node > nodes
Definition: ncmesh.hpp:479
double a
Definition: lissajous.cpp:41
void ForgetElement(int e)
Definition: ncmesh.cpp:421
static PointMatrix pm_tet_identity
Definition: ncmesh.hpp:811
int GetMidFaceNode(int en1, int en2, int en3, int en4)
Definition: ncmesh.cpp:307
void RefineElement(int elem, char ref_type)
Definition: ncmesh.cpp:887
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:465
void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1, MatrixMap &matrix_map)
Definition: ncmesh.cpp:2839
Array< int > boundary_faces
subset of all faces, set by BuildFaceList
Definition: ncmesh.hpp:519
int child[8]
2-8 children (if ref_type != 0)
Definition: ncmesh.hpp:466
int NGhostFaces
Definition: ncmesh.hpp:509
Array< double > coordinates
Definition: ncmesh.hpp:492
const NCList & GetNCList(int entity)
Return vertex/edge/face list (entity = 0/1/2, respectively).
Definition: ncmesh.hpp:255
Array< int > leaf_elements
finest elements, in Mesh ordering (+ ghosts)
Definition: ncmesh.hpp:511
static PointMatrix pm_quad_identity
Definition: ncmesh.hpp:810
Array< Master > masters
Definition: ncmesh.hpp:210
int EdgeSplitLevel(int vn1, int vn2) const
Definition: ncmesh.cpp:5115
Array< MeshId > conforming
Definition: ncmesh.hpp:209
unsigned edge_flags
orientation flags, see OrientedPointMatrix
Definition: ncmesh.hpp:198
void FindVertexCousins(int elem, int local, Array< int > &cousins) const
Definition: ncmesh.cpp:3769
signed char geom
Geometry::Type (faces only) (char to save RAM)
Definition: ncmesh.hpp:172
int QuadFaceSplitType(int v1, int v2, int v3, int v4, int mid[5]=NULL) const
Definition: ncmesh.cpp:2509
void OrientedPointMatrix(const Slave &slave, DenseMatrix &oriented_matrix) const
Return the point matrix oriented according to the master and slave edges.
Definition: ncmesh.cpp:3193
Array< int > free_element_ids
Definition: ncmesh.hpp:483
void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level, Face *eface[4], MatrixMap &matrix_map)
Definition: ncmesh.cpp:2720
int index(int i, int j, int nx, int ny)
Definition: life.cpp:237
int NGhostEdges
Definition: ncmesh.hpp:509
NCMesh(const Mesh *mesh)
Definition: ncmesh.cpp:104
int parent
Element index in the coarse mesh.
Definition: ncmesh.hpp:49
Point(double x, double y)
Definition: ncmesh.hpp:734
static PointMatrix pm_hex_identity
Definition: ncmesh.hpp:813
int Dimension() const
Definition: ncmesh.hpp:128
void FindEdgeElements(int vn1, int vn2, int vn3, int vn4, Array< MeshId > &prisms) const
Definition: ncmesh.cpp:703
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
Definition: ncmesh.cpp:1793
T & Last()
Return the last element in the array.
Definition: array.hpp:779
void NeighborExpand(const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
Definition: ncmesh.cpp:3639
BlockArray< Element > elements
Definition: ncmesh.hpp:482
void FindSetNeighbors(const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
Definition: ncmesh.cpp:3453
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:4861
virtual void Update()
Definition: ncmesh.cpp:231
Slave(int index, int element, int local, int geom)
Definition: ncmesh.hpp:201
bool HaveTets() const
Definition: ncmesh.hpp:539
Data type point element.
Definition: point.hpp:22
void InitRootElements()
Count root elements and initialize root_state.
Definition: ncmesh.cpp:5555
long MemoryUsage() const
Return total number of bytes allocated.
Definition: ncmesh.cpp:6093
HashTable< Node > shadow
temporary storage for reparented nodes
Definition: ncmesh.hpp:547
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:60
virtual void BuildVertexList()
Definition: ncmesh.cpp:3158
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition: ncmesh.hpp:240
MeshId(int index, int element, int local, int geom=-1)
Definition: ncmesh.hpp:177
void UpdateElementToVertexTable()
Definition: ncmesh.hpp:707
bool IsGhost(const Element &el) const
Definition: ncmesh.hpp:541
void LoadBoundary(std::istream &input)
Load the &quot;boundary&quot; section of the mesh file.
Definition: ncmesh.cpp:5388
void CopyElements(int elem, const BlockArray< Element > &tmp_elements)
Definition: ncmesh.cpp:5779
virtual void ElementSharesEdge(int elem, int local, int enode)
Definition: ncmesh.hpp:668
int GetElementDepth(int i) const
Return the distance of leaf &#39;i&#39; from the root.
Definition: ncmesh.cpp:4971
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:2583
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
Definition: ncmesh.hpp:457
int NElements
Definition: ncmesh.hpp:505
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:3339
bool Unused() const
Definition: ncmesh.hpp:441
static int find_node(const Element &el, int node)
Definition: ncmesh.cpp:2573
void UnreferenceElement(int elem, Array< int > &elemFaces)
Definition: ncmesh.cpp:346
int index
Mesh number.
Definition: ncmesh.hpp:169
void LoadVertexParents(std::istream &input)
Load the vertex parent hierarchy from a mesh file.
Definition: ncmesh.cpp:5326
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, bool get_coarse_to_fine_only=false) const
Definition: ncmesh.cpp:4455
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:745
Abstract data type element.
Definition: element.hpp:28
Point & operator=(const Point &src)
Definition: ncmesh.hpp:760
Refinement()=default
int GetSingleElement() const
Return one of elem[0] or elem[1] and make sure the other is -1.
Definition: ncmesh.cpp:436
void PrintCoordinates(std::ostream &out) const
Print the &quot;coordinates&quot; section of the mesh file.
Definition: ncmesh.cpp:5428
const double * CalcVertexPos(int node) const
Definition: ncmesh.cpp:2244
int rank
processor number (ParNCMesh), -1 if undefined/unknown
Definition: ncmesh.hpp:461
int node[8]
element corners (if ref_type == 0)
Definition: ncmesh.hpp:465
int GetVertexRootCoord(int elem, RefCoord coord[3]) const
Definition: ncmesh.cpp:3698
bool TraverseTriFace(int vn0, int vn1, int vn2, const PointMatrix &pm, int level, MatrixMap &matrix_map)
Definition: ncmesh.cpp:2876
HashTable< Face > faces
Definition: ncmesh.hpp:480
Defines the position of a fine element within a coarse element.
Definition: ncmesh.hpp:46
virtual void Refine(const Array< Refinement > &refinements)
Definition: ncmesh.cpp:1528
int matrix
Index into the DenseTensor corresponding to the parent Geometry::Type stored in CoarseFineTransformat...
Definition: ncmesh.hpp:52
Array< char > face_geom
face geometry by face index, set by OnMeshUpdated
Definition: ncmesh.hpp:520
int GetEdgeMaster(int v1, int v2) const
Definition: ncmesh.cpp:4962