MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ncmesh.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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"
20#include "element.hpp"
21#include "vertex.hpp"
22#include "../fem/geom.hpp"
23
24#include <vector>
25#include <map>
26#include <iostream>
27#include <unordered_map>
28
29namespace mfem
30{
31
32/** Represents the index of an element to refine, plus a refinement type.
33 The refinement type is needed for anisotropic refinement of quads and hexes.
34 Bits 0,1 and 2 of 'ref_type' specify whether the element should be split
35 in the X, Y and Z directions, respectively (Z is ignored for quads). */
37{
38 enum : char { X = 1, Y = 2, Z = 4, XY = 3, XZ = 5, YZ = 6, XYZ = 7 };
39 int index; ///< Mesh element number
40 char ref_type; ///< refinement XYZ bit mask (7 = full isotropic)
41
42 Refinement() = default;
43
45 : index(index), ref_type(type) {}
46};
47
48
49/// Defines the position of a fine element within a coarse element.
51{
52 /// Coarse %Element index in the coarse mesh.
53 int parent;
54
55 /** The (geom, matrix) pair determines the sub-element transformation for the
56 fine element: CoarseFineTransformations::point_matrices[geom](matrix) is
57 the point matrix of the region within the coarse element reference domain.*/
58 unsigned geom : 4;
59 unsigned matrix : 27;
60
61 /// For internal use: 0 if regular fine element, 1 if parallel ghost element.
62 unsigned ghost : 1;
63
64 Embedding() = default;
65 Embedding(int elem, Geometry::Type geom, int matrix = 0, bool ghost = false)
66 : parent(elem), geom(geom), matrix(matrix), ghost(ghost) {}
67};
68
69
70/// Defines the coarse-fine transformations of all fine elements.
72{
73 /// Fine element positions in their parents.
75
76 /** A "dictionary" of matrices for IsoparametricTransformation. Use
77 Embedding::{geom,matrix} to access a fine element point matrix. */
79
80 /** Invert the 'embeddings' array: create a Table with coarse elements as
81 rows and fine elements as columns. If 'want_ghosts' is false, parallel
82 ghost fine elements are not included in the table. */
83 void MakeCoarseToFineTable(Table &coarse_to_fine,
84 bool want_ghosts = false) const;
85
86 void Clear();
87 bool IsInitialized() const;
88 long MemoryUsage() const;
89
90 MFEM_DEPRECATED
91 void GetCoarseToFineMap(const Mesh &fine_mesh, Table &coarse_to_fine) const
92 { MakeCoarseToFineTable(coarse_to_fine, true); (void) fine_mesh; }
93};
94
95void Swap(CoarseFineTransformations &a, CoarseFineTransformations &b);
96
97struct MatrixMap; // for internal use
98
99
100/** \brief A class for non-conforming AMR. The class is not used directly
101 * by the user, rather it is an extension of the Mesh class.
102 *
103 * In general, the class is used by MFEM as follows:
104 *
105 * 1. NCMesh is constructed from elements of an existing Mesh. The elements
106 * are copied and become roots of the refinement hierarchy.
107 *
108 * 2. Some elements are refined with the Refine() method. Both isotropic and
109 * anisotropic refinements of quads/hexes are supported.
110 *
111 * 3. A new Mesh is created from NCMesh containing the leaf elements.
112 * This new Mesh may have non-conforming (hanging) edges and faces and
113 * is the one seen by the user.
114 *
115 * 4. FiniteElementSpace asks NCMesh for a list of conforming, master and
116 * slave edges/faces and creates the conforming interpolation matrix P.
117 *
118 * 5. A continuous/conforming solution is obtained by solving P'*A*P x = P'*b.
119 *
120 * 6. Repeat from step 2.
121 */
123{
124public:
125 //// Initialize with elements from an existing 'mesh'.
126 explicit NCMesh(const Mesh *mesh);
127
128 /** Load from a stream. The id header is assumed to have been read already
129 from \param[in] input . \param[in] version is 10 for the v1.0 NC format,
130 or 1 for the legacy v1.1 format. \param[out] curved is set to 1 if the
131 curvature GridFunction follows after mesh data. \param[out] is_nc (again
132 treated as a boolean) is set to 0 if the legacy v1.1 format in fact
133 defines a conforming mesh. See Mesh::Loader for details. */
134 NCMesh(std::istream &input, int version, int &curved, int &is_nc);
135
136 /// Deep copy of another instance.
137 NCMesh(const NCMesh &other);
138
139 /// Copy assignment not supported
141
142 virtual ~NCMesh();
143
144 /// Return the dimension of the NCMesh.
145 int Dimension() const { return Dim; }
146 /// Return the space dimension of the NCMesh.
147 int SpaceDimension() const { return spaceDim; }
148
149 /// Return the number of vertices in the NCMesh.
150 int GetNVertices() const { return NVertices; }
151 /// Return the number of edges in the NCMesh.
152 int GetNEdges() const { return NEdges; }
153 /// Return the number of (2D) faces in the NCMesh.
154 int GetNFaces() const { return NFaces; }
155 virtual int GetNGhostElements() const { return 0; }
156
157 /** Perform the given batch of refinements. Please note that in the presence
158 of anisotropic splits additional refinements may be necessary to keep
159 the mesh consistent. However, the function always performs at least the
160 requested refinements. */
161 virtual void Refine(const Array<Refinement> &refinements);
162
163 /** Check the mesh and potentially refine some elements so that the maximum
164 difference of refinement levels between adjacent elements is not greater
165 than 'max_nc_level'. */
166 virtual void LimitNCLevel(int max_nc_level);
167
168 /** Return a list of derefinement opportunities. Each row of the table
169 contains Mesh indices of existing elements that can be derefined to form
170 a single new coarse element. Row numbers are then passed to Derefine.
171 This function works both in serial and parallel. */
173
174 /** Check derefinements returned by GetDerefinementTable and mark those that
175 can be done safely so that the maximum NC level condition is not violated.
176 On return, level_ok.Size() == deref_table.Size() and contains 0/1s. */
177 virtual void CheckDerefinementNCLevel(const Table &deref_table,
178 Array<int> &level_ok, int max_nc_level);
179
180 /** Perform a subset of the possible derefinements (see GetDerefinementTable).
181 Note that if anisotropic refinements are present in the mesh, some of the
182 derefinements may have to be skipped to preserve mesh consistency. */
183 virtual void Derefine(const Array<int> &derefs);
184
185 // master/slave lists
186
187 /// Identifies a vertex/edge/face in both Mesh and NCMesh.
188 struct MeshId
189 {
190 int index; ///< Mesh number
191 int element; ///< NCMesh::Element containing this vertex/edge/face
192 signed char local; ///< local number within 'element'
193 signed char geom; ///< Geometry::Type (faces only) (char to save RAM)
194
196
197 MeshId() = default;
198 MeshId(int index, int element, int local, int geom = -1)
200 };
201
202 /** Nonconforming edge/face that has more than one neighbor. The neighbors
203 are stored in NCList::slaves[i], slaves_begin <= i < slaves_end. */
204 struct Master : public MeshId
205 {
206 int slaves_begin, slaves_end; ///< slave faces
207
208 Master() = default;
209 Master(int index, int element, int local, int geom, int sb, int se)
211 , slaves_begin(sb), slaves_end(se) {}
212 };
213
214 /// Nonconforming edge/face within a bigger edge/face.
215 struct Slave : public MeshId
216 {
217 int master; ///< master number (in Mesh numbering)
218 unsigned matrix : 24; ///< index into NCList::point_matrices[geom]
219 unsigned edge_flags : 8; ///< orientation flags, see OrientedPointMatrix
220
221 Slave() = default;
222 Slave(int index, int element, int local, int geom)
224 , master(-1), matrix(0), edge_flags(0) {}
225 };
226
227
228 /// Lists all edges/faces in the nonconforming mesh.
229 struct NCList
230 {
231 Array<MeshId> conforming; ///< All MeshIds corresponding to conformal faces
232 Array<Master> masters; ///< All MeshIds corresponding to master faces
233 Array<Slave> slaves; ///< All MeshIds corresponding to slave faces
234
235 /// List of unique point matrices for each slave geometry.
237
238 void OrientedPointMatrix(const Slave &slave,
239 DenseMatrix &oriented_matrix) const;
240
241 /// Particular MeshId type, used for allowing static casting to the
242 /// appropriate child type after searching the NCList. UNRECOGNIZED
243 /// denotes that an instance is not known within the NCList, meaning that
244 /// it does not play a part in NC mechanics. This can be because the index
245 /// did not exist in the original Mesh, or because the entry is a boundary
246 /// face, whose NC status is always conforming.
248
249 /// Helper storing a reference to a MeshId type, and the face type it can
250 /// be cast to
252 {
253 const MeshId * const id; ///< Pointer to a possible MeshId, nullptr if not found
254 /// MeshIdType corresponding to the MeshId. UNRECOGNIZED if unfound.
256 };
257 /// Return a mesh id and type for a given nc index.
259
260 /// Return a face type for a given nc index.
261 MeshIdType GetMeshIdType(int index) const;
262
263 /// Given an index, check if this is a certain face type.
264 bool CheckMeshIdType(int index, MeshIdType type) const;
265
266 /// Erase the contents of the conforming, master and slave arrays.
267 void Clear();
268 /// Whether the NCList is empty.
269 bool Empty() const
270 {
271 return conforming.Size() == 0
272 && masters.Size() == 0
273 && slaves.Size() == 0;
274 }
275 /// The total size of the component arrays in the NCList.
276 long TotalSize() const
277 {
278 return conforming.Size() + masters.Size() + slaves.Size();
279 }
280 /// The memory usage of the three public arrays. Does not account for the
281 /// inverse index.
282 long MemoryUsage() const;
283 ~NCList() { Clear(); }
284 private:
285 // Check for existence or construct the inv_index list map if necessary.
286 // const because only modifies the mutable member inv_index.
287 void BuildIndex() const;
288
289 /// A lazily constructed map from index to MeshId. Built whenever
290 /// GetMeshIdAndType, GetMeshIdType or CheckMeshIdType is called for the
291 /// first time. The MeshIdType is stored with, to enable casting to Slave
292 /// or Master elements appropriately.
293 mutable std::unordered_map<int, std::pair<MeshIdType, int>> inv_index;
294 };
295
296
297 /// Return the current list of conforming and nonconforming faces.
299 {
300 if (face_list.Empty()) { BuildFaceList(); }
301 return face_list;
302 }
303
304 /// Return the current list of conforming and nonconforming edges.
306 {
307 if (edge_list.Empty()) { BuildEdgeList(); }
308 return edge_list;
309 }
310
311 /** Return a list of vertices (in 'conforming'); this function is provided
312 for uniformity/completeness. Needed in ParNCMesh/ParFESpace. */
314 {
315 if (vertex_list.Empty()) { BuildVertexList(); }
316 return vertex_list;
317 }
318
319 /// Return vertex/edge/face list (entity = 0/1/2, respectively).
320 const NCList& GetNCList(int entity)
321 {
322 switch (entity)
323 {
324 case 0: return GetVertexList();
325 case 1: return GetEdgeList();
326 default: return GetFaceList();
327 }
328 }
329
330
331 // coarse/fine transforms
332
333 /** Remember the current layer of leaf elements before the mesh is refined.
334 Needed by GetRefinementTransforms(), must be called before Refine(). */
335 void MarkCoarseLevel();
336
337 /** After refinement, calculate the relation of each fine element to its
338 parent coarse element. Note that Refine() or LimitNCLevel() can be called
339 multiple times between MarkCoarseLevel() and this function. */
341
342 /** After derefinement, calculate the relations of previous fine elements
343 (some of which may no longer exist) to the current leaf elements.
344 Unlike for refinement, Derefine() may only be called once before this
345 function so there is no MarkFineLevel(). */
347
348 /// Free all internal data created by the above three functions.
349 void ClearTransforms();
350
351
352 // grid ordering
353
354 /** Return a space filling curve for a rectangular grid of elements.
355 Implemented is a generalized Hilbert curve for arbitrary grid dimensions.
356 If the width is odd, height should be odd too, otherwise one diagonal
357 (vertex-neighbor) step cannot be avoided in the curve. Even dimensions
358 are recommended. */
359 static void GridSfcOrdering2D(int width, int height,
360 Array<int> &coords);
361
362 /** Return a space filling curve for a 3D rectangular grid of elements.
363 The Hilbert-curve-like algorithm works well for even dimensions. For odd
364 width/height/depth it tends to produce some diagonal (edge-neighbor)
365 steps. Even dimensions are recommended. */
366 static void GridSfcOrdering3D(int width, int height, int depth,
367 Array<int> &coords);
368
369
370 // utility
371
372 /// Return Mesh vertex indices of an edge identified by 'edge_id'.
373 void GetEdgeVertices(const MeshId &edge_id, int vert_index[2],
374 bool oriented = true) const;
375
376 /** Return "NC" orientation of an edge. As opposed to standard Mesh edge
377 orientation based on vertex IDs, "NC" edge orientation follows the local
378 edge orientation within the element 'edge_id.element' and is thus
379 processor independent. TODO: this seems only partially true? */
380 int GetEdgeNCOrientation(const MeshId &edge_id) const;
381
382 /** Return Mesh vertex and edge indices of a face identified by 'face_id'.
383 The return value is the number of face vertices. */
384 int GetFaceVerticesEdges(const MeshId &face_id,
385 int vert_index[4], int edge_index[4],
386 int edge_orientation[4]) const;
387
388 /** Given an edge (by its vertex indices v1 and v2) return the first
389 (geometric) parent edge that exists in the Mesh or -1 if there is no such
390 parent. */
391 int GetEdgeMaster(int v1, int v2) const;
392
393 /** Get a list of vertices (2D/3D), edges (3D) and faces (3D) that coincide
394 with boundary elements with the specified attributes (marked in
395 'bdr_attr_is_ess'). In 3D this function also reveals "hidden" boundary
396 edges. In parallel it helps identifying boundary vertices/edges/faces
397 affected by non-local boundary elements. Hidden faces can occur for an
398 internal boundary coincident to a processor boundary.
399 */
400
401 /**
402 * @brief Get a list of vertices (2D/3D), edges (3D) and faces (3D) that
403 * coincide with boundary elements with the specified attributes (marked in
404 * 'bdr_attr_is_ess').
405 *
406 * @details In 3D this function also reveals "hidden" boundary edges. In
407 * parallel it helps identifying boundary vertices/edges/faces affected by
408 * non-local boundary elements. Hidden faces can occur for an internal
409 * boundary coincident to a processor boundary.
410 *
411 * @param bdr_attr_is_ess Indicator if a given attribute is essential.
412 * @param bdr_vertices Array of vertices that are essential.
413 * @param bdr_edges Array of edges that are essential.
414 * @param bdr_faces Array of faces that are essential.
415 */
416 virtual void GetBoundaryClosure(const Array<int> &bdr_attr_is_ess,
417 Array<int> &bdr_vertices,
418 Array<int> &bdr_edges, Array<int> &bdr_faces);
419
420 /// Return element geometry type. @a index is the Mesh element number.
423
424 /// Return face geometry type. @a index is the Mesh face number.
427
428 /// Return the number of root elements.
430
431 /// Return the distance of leaf 'i' from the root.
432 int GetElementDepth(int i) const;
433
434 /** Return the size reduction compared to the root element (ignoring local
435 stretching and curvature). */
436 int GetElementSizeReduction(int i) const;
437
438 /// Return the faces and face attributes of leaf element 'i'.
440 Array<int> &fattr) const;
441
442
443 /** I/O: Print the mesh in "MFEM NC mesh v1.0" format. If @a comments is
444 non-empty, it will be printed after the first line of the file, and each
445 line should begin with '#'. */
446 void Print(std::ostream &out, const std::string &comments = "") const;
447
448 /// I/O: Return true if the mesh was loaded from the legacy v1.1 format.
449 bool IsLegacyLoaded() const { return Legacy; }
450
451 /// I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
452 void LegacyToNewVertexOrdering(Array<int> &order) const;
453
454 /// Save memory by releasing all non-essential and cached data.
455 virtual void Trim();
456
457 /// Return total number of bytes allocated.
458 long MemoryUsage() const;
459
460 int PrintMemoryDetail() const;
461
462 typedef std::int64_t RefCoord;
463
464
465protected: // non-public interface for the Mesh class
466
467 friend class Mesh;
468
469 /// Fill Mesh::{vertices,elements,boundary} for the current finest level.
470 void GetMeshComponents(Mesh &mesh) const;
471
472 /** Get edge and face numbering from 'mesh' (i.e., set all Edge::index and
473 Face::index) after a new mesh was created from us. */
474 void OnMeshUpdated(Mesh *mesh);
475
476 /** Delete top-level vertex coordinates if the Mesh became curved, e.g.,
477 by calling Mesh::SetCurvature or otherwise setting the Nodes. */
479
480protected: // implementation
481
482 int Dim, spaceDim; ///< dimensions of the elements and the vertex coordinates
483 int MyRank; ///< used in parallel, or when loading a parallel file in serial
484 bool Iso; ///< true if the mesh only contains isotropic refinements
485 int Geoms; ///< bit mask of element geometries present, see InitGeomFlags()
486 bool Legacy; ///< true if the mesh was loaded from the legacy v1.1 format
487
488 static const int MaxElemNodes =
489 8; ///< Number of nodes of an element can have
490 static const int MaxElemEdges =
491 12; ///< Number of edges of an element can have
492 static const int MaxElemFaces =
493 6; ///< Number of faces of an element can have
494 static const int MaxElemChildren =
495 10; ///< Number of children of an element can have
496
497 /** A Node can hold a vertex, an edge, or both. Elements directly point to
498 their corner nodes, but edge nodes also exist and can be accessed using
499 a hash-table given their two end-point node IDs. All nodes can be
500 accessed in this way, with the exception of top-level vertex nodes.
501 When an element is being refined, the mid-edge nodes are readily
502 available with this mechanism. The new elements "sign in" to the nodes
503 by increasing the reference counts of their vertices and edges. The
504 parent element "signs off" its nodes by decrementing the ref counts. */
505 struct Node : public Hashed2
506 {
509
511 ~Node();
512
513 bool HasVertex() const { return vert_refc > 0; }
514 bool HasEdge() const { return edge_refc > 0; }
515
516 // decrease vertex/edge ref count, return false if Node should be deleted
517 bool UnrefVertex() { --vert_refc; return vert_refc || edge_refc; }
518 bool UnrefEdge() { --edge_refc; return vert_refc || edge_refc; }
519 };
520
521 /** Similarly to nodes, faces can be accessed by hashing their four vertex
522 node IDs. A face knows about the one or two elements that are using it.
523 A face that is not on the boundary and only has one element referencing
524 it is either a master or a slave face. */
525 struct Face : public Hashed4
526 {
527 int attribute; ///< boundary element attribute, -1 if internal face
528 int index; ///< face number in the Mesh
529 int elem[2]; ///< up to 2 elements sharing the face
530
531 Face() : attribute(-1), index(-1) { elem[0] = elem[1] = -1; }
532
533 bool Boundary() const { return attribute >= 0; }
534 bool Unused() const { return elem[0] < 0 && elem[1] < 0; }
535
536 // add or remove an element from the 'elem[2]' array
537 void RegisterElement(int e);
538 void ForgetElement(int e);
539
540 /// Return one of elem[0] or elem[1] and make sure the other is -1.
541 int GetSingleElement() const;
542 };
543
544 /** This is an element in the refinement hierarchy. Each element has
545 either been refined and points to its children, or is a leaf and points
546 to its vertex nodes. */
547 struct Element
548 {
549 char geom; ///< Geometry::Type of the element (char for storage only)
550 char ref_type; ///< bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
551 char tet_type; ///< tetrahedron split type, currently always 0
552 char flag; ///< generic flag/marker, can be used by algorithms
553 int index; ///< element number in the Mesh, -1 if refined
554 int rank; ///< processor number (ParNCMesh), -1 if undefined/unknown
556 union
557 {
558 int node[MaxElemNodes]; ///< element corners (if ref_type == 0)
559 int child[MaxElemChildren]; ///< 2-10 children (if ref_type != 0)
560 };
561 int parent; ///< parent element, -1 if this is a root element, -2 if free'd
562
563 Element(Geometry::Type geom, int attr);
564
566 bool IsLeaf() const { return !ref_type && (parent != -2); }
567 };
568
569
570 // primary data
571
572 HashTable<Node> nodes; // associative container holding all Nodes
573 HashTable<Face> faces; // associative container holding all Faces
574
575 BlockArray<Element> elements; // storage for all Elements
576 Array<int> free_element_ids; // unused element ids - indices into 'elements'
577
578 /** Initial traversal state (~ element orientation) for each root element
579 NOTE: M = root_state.Size() is the number of root elements.
580 NOTE: the first M items of 'elements' is the coarse mesh. */
582
583 /** Coordinates of top-level vertices (organized as triples). If empty,
584 the Mesh is curved (Nodes != NULL) and NCMesh is topology-only. */
586
587
588 // secondary data
589
590 /** Apart from the primary data structure, which is the element/node/face
591 hierarchy, there is secondary data that is derived from the primary
592 data and needs to be updated when the primary data changes. Update()
593 takes care of that and needs to be called after each refinement and
594 derefinement. */
595 virtual void Update();
596
597 // set by UpdateLeafElements, UpdateVertices and OnMeshUpdated
599
600 // NOTE: the serial code understands the bare minimum about ghost elements and
601 // other ghost entities in order to be able to load parallel partial meshes
603
604 Array<int> leaf_elements; ///< finest elements, in Mesh ordering (+ ghosts)
605 Array<int> leaf_sfc_index; ///< natural tree ordering of leaf elements
606 Array<int> vertex_nodeId; ///< vertex-index to node-id map, see UpdateVertices
607
608 NCList face_list; ///< lazy-initialized list of faces, see GetFaceList
609 NCList edge_list; ///< lazy-initialized list of edges, see GetEdgeList
610 NCList vertex_list; ///< lazy-initialized list of vertices, see GetVertexList
611
612 Array<int> boundary_faces; ///< subset of all faces, set by BuildFaceList
613 Array<char> face_geom; ///< face geometry by face index, set by OnMeshUpdated
614
615 Table element_vertex; ///< leaf-element to vertex table, see FindSetNeighbors
616
617 /// Update the leaf elements indices in leaf_elements
618 void UpdateLeafElements();
619
620 /** @brief This method assigns indices to vertices (Node::vert_index) that
621 will be seen by the Mesh class and the rest of MFEM.
622
623 We must be careful to:
624 1. Stay compatible with the conforming code, which expects top-level
625 (original) vertices to be indexed first, otherwise GridFunctions
626 defined on a conforming mesh would no longer be valid when the
627 mesh is converted to an NC mesh.
628
629 2. Make sure serial NCMesh is compatible with the parallel ParNCMesh,
630 so it is possible to read parallel partial solutions in serial code
631 (e.g., serial GLVis). This means handling ghost elements, if present.
632
633 3. Assign vertices in a globally consistent order for parallel meshes:
634 if two vertices i,j are shared by two ranks r1,r2, and i<j on r1,
635 then i<j on r2 as well. This is true for top-level vertices but also
636 for the remaining shared vertices thanks to the globally consistent
637 SFC ordering of the leaf elements. This property reduces communication
638 and simplifies ParNCMesh. */
639 void UpdateVertices(); ///< update Vertex::index and vertex_nodeId
640
641 /** Collect the leaf elements in leaf_elements, and the ghost elements in
642 ghosts. Compute and set the element indices of @a elements. On quad and
643 hex refined elements tries to order leaf elements along a space-filling
644 curve according to the given @a state variable. */
645 void CollectLeafElements(int elem, int state, Array<int> &ghosts,
646 int &counter);
647
648 /** Try to find a space-filling curve friendly orientation of the root
649 elements: set 'root_state' based on the ordering of coarse elements.
650 Note that the coarse mesh itself must be ordered as an SFC by e.g.
651 Mesh::GetGeckoElementOrdering. */
652 void InitRootState(int root_count);
653
654 /** Compute the Geometry::Type present in the root elements (coarse elements)
655 and set @a Geoms bitmask accordingly. */
656 void InitGeomFlags();
657
658 /// Return true if the mesh contains prism elements.
659 bool HavePrisms() const { return Geoms & (1 << Geometry::PRISM); }
660
661 /// Return true if the mesh contains pyramid elements.
662 bool HavePyramids() const { return Geoms & (1 << Geometry::PYRAMID); }
663
664 /// Return true if the mesh contains tetrahedral elements.
665 bool HaveTets() const { return Geoms & (1 << Geometry::TETRAHEDRON); }
666
667 /// Return true if the Element @a el is a ghost element.
668 bool IsGhost(const Element &el) const { return el.rank != MyRank; }
669
670
671 // refinement/derefinement
672
673 Array<Refinement> ref_stack; ///< stack of scheduled refinements (temporary)
674 HashTable<Node> shadow; ///< temporary storage for reparented nodes
675 Array<Triple<int, int, int> > reparents; ///< scheduled node reparents (tmp)
676
677 Table derefinements; ///< possible derefinements, see GetDerefinementTable
678
679 /** Refine the element @a elem with the refinement @a ref_type
680 (c.f. Refinement::enum) */
681 void RefineElement(int elem, char ref_type);
682
683 /// Derefine the element @a elem, does nothing on leaf elements.
684 void DerefineElement(int elem);
685
686 // Add an Element @a el to the NCMesh, optimized to reuse freed elements.
687 int AddElement(const Element &el)
688 {
690 {
691 int idx = free_element_ids.Last();
693 elements[idx] = el;
694 return idx;
695 }
696 return elements.Append(el);
697 }
698
699 // Free the element with index @a id.
700 void FreeElement(int id)
701 {
703 elements[id].ref_type = 0;
704 elements[id].parent = -2; // mark the element as free
705 }
706
707 int NewHexahedron(int n0, int n1, int n2, int n3,
708 int n4, int n5, int n6, int n7, int attr,
709 int fattr0, int fattr1, int fattr2,
710 int fattr3, int fattr4, int fattr5);
711
712 int NewWedge(int n0, int n1, int n2,
713 int n3, int n4, int n5, int attr,
714 int fattr0, int fattr1,
715 int fattr2, int fattr3, int fattr4);
716
717 int NewTetrahedron(int n0, int n1, int n2, int n3, int attr,
718 int fattr0, int fattr1, int fattr2, int fattr3);
719
720 int NewPyramid(int n0, int n1, int n2, int n3, int n4, int attr,
721 int fattr0, int fattr1, int fattr2, int fattr3,
722 int fattr4);
723
724 int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr,
725 int eattr0, int eattr1, int eattr2, int eattr3);
726
727 int NewTriangle(int n0, int n1, int n2,
728 int attr, int eattr0, int eattr1, int eattr2);
729
730 int NewSegment(int n0, int n1, int attr, int vattr1, int vattr2);
731
732 mfem::Element* NewMeshElement(int geom) const;
733
734 /**
735 * @brief Given a quad face defined by four vertices, establish which edges
736 * of this face have been split, and if so optionally return the mid points
737 * of those edges.
738 *
739 * @param n1 The first node defining the face
740 * @param n2 The second node defining the face
741 * @param n3 The third node defining the face
742 * @param n4 The fourth node defining the face
743 * @param mid optional return of the edge mid points.
744 * @return int 0 -- no split, 1 -- "vertical" split, 2 -- "horizontal" split
745 */
746 int QuadFaceSplitType(int n1, int n2, int n3, int n4, int mid[5]
747 = NULL /*optional output of mid-edge nodes*/) const;
748
749 /**
750 * @brief Given a tri face defined by three vertices, establish whether the
751 * edges that make up this face have been split, and if so optionally return
752 * the midpoints.
753 * @details This is a necessary condition for this face to have been split,
754 * but is not sufficient. Consider a triangle attached to three refined
755 * triangles, in this scenario all edges can be split but this face not be
756 * split. In this case, it is necessary to check if there is a face made up
757 * of the returned midpoint nodes.
758 *
759 * @param n1 The first node defining the face
760 * @param n2 The second node defining the face
761 * @param n3 The third node defining the face
762 * @param mid optional return of the edge mid points.
763 * @return true Splits for all edges have been found
764 * @return false
765 */
766 bool TriFaceSplit(int n1, int n2, int n3, int mid[3] = NULL) const;
767
768 /**
769 * @brief Determine if a Triangle face is a master face
770 * @details This check requires looking for the edges making up the triangle
771 * being split, if nodes exist at their midpoints, and there are vertices at
772 * them, this implies the face COULD be split. To determine if it is, we then
773 * check whether these midpoints have all been connected, this is required to
774 * discriminate between an internal master face surrounded by nonconformal
775 * refinements and a conformal boundary face surrounded by refinements.
776 *
777 * @param n1 The first node defining the face
778 * @param n2 The second node defining the face
779 * @param n3 The third node defining the face
780 * @return true The face is a master
781 * @return false The face is not a master
782 */
783 inline bool TriFaceIsMaster(int n1, int n2, int n3) const
784 {
785 int mid[3];
786 return !(!TriFaceSplit(n1, n2, n3, mid) // The edges aren't split
787 // OR none of the midpoints are connected.
788 || (nodes.FindId(mid[0], mid[1]) < 0 &&
789 nodes.FindId(mid[0], mid[2]) < 0 &&
790 nodes.FindId(mid[1], mid[2]) < 0));
791 }
792
793 /**
794 * @brief Determine if a Quad face is a master face
795 *
796 * @param n1 The first node defining the face
797 * @param n2 The second node defining the face
798 * @param n3 The third node defining the face
799 * @param n4 The fourth node defining the face
800 * @return true The quad face is a master face
801 * @return false The quad face is not a master face
802 */
803 inline bool QuadFaceIsMaster(int n1, int n2, int n3, int n4) const
804 {
805 return QuadFaceSplitType(n1, n2, n3, n4) != 0;
806 }
807
808 void ForceRefinement(int vn1, int vn2, int vn3, int vn4);
809
810 void FindEdgeElements(int vn1, int vn2, int vn3, int vn4,
811 Array<MeshId> &prisms) const;
812
813 void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4,
814 const Refinement *refs, int nref);
815
816 void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4,
817 int mid12, int mid34, int level = 0);
818
819 void CheckIsoFace(int vn1, int vn2, int vn3, int vn4,
820 int en1, int en2, int en3, int en4, int midf);
821
822 void ReparentNode(int node, int new_p1, int new_p2);
823
824 int FindMidEdgeNode(int node1, int node2) const;
825 int GetMidEdgeNode(int node1, int node2);
826
827 int GetMidFaceNode(int en1, int en2, int en3, int en4);
828
829 void ReferenceElement(int elem);
830 void UnreferenceElement(int elem, Array<int> &elemFaces);
831
832 Face* GetFace(Element &elem, int face_no);
833 void RegisterFaces(int elem, int *fattr = NULL);
834 void DeleteUnusedFaces(const Array<int> &elemFaces);
835
836 void CollectDerefinements(int elem, Array<Connection> &list);
837
838 /// Return el.node[index] correctly, even if the element is refined.
839 int RetrieveNode(const Element &el, int index);
840
841 /// Extended version of find_node: works if 'el' is refined.
842 int FindNodeExt(const Element &el, int node, bool abort = true);
843
844
845 // face/edge lists
846
847 static int find_node(const Element &el, int node);
848 static int find_element_edge(const Element &el, int vn0, int vn1,
849 bool abort = true);
850 static int find_local_face(int geom, int a, int b, int c);
851
852 struct Point;
853 struct PointMatrix;
854
855 int ReorderFacePointMat(int v0, int v1, int v2, int v3,
856 int elem, const PointMatrix &pm,
857 PointMatrix &reordered) const;
858
859 void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3,
860 const PointMatrix& pm, int level, Face* eface[4],
861 MatrixMap &matrix_map);
863 {
864 bool unsplit; ///< Whether this face has no further splits.
865 bool ghost_neighbor; ///< Whether the face neighbor is a ghost.
866 };
867 TriFaceTraverseResults TraverseTriFace(int vn0, int vn1, int vn2,
868 const PointMatrix& pm, int level,
869 MatrixMap &matrix_map);
870 void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1,
871 MatrixMap &matrix_map);
872 void TraverseEdge(int vn0, int vn1, real_t t0, real_t t1, int flags,
873 int level, MatrixMap &matrix_map);
874
875 virtual void BuildFaceList();
876 virtual void BuildEdgeList();
877 virtual void BuildVertexList();
878
879 virtual void ElementSharesFace(int elem, int local, int face) {} // ParNCMesh
880 virtual void ElementSharesEdge(int elem, int local, int enode) {} // ParNCMesh
881 virtual void ElementSharesVertex(int elem, int local, int vnode) {} // ParNCMesh
882
883 // neighbors / element_vertex table
884
885 /** Return all vertex-, edge- and face-neighbors of a set of elements.
886 The neighbors are returned as a list (neighbors != NULL), as a set
887 (neighbor_set != NULL), or both. The sizes of the set arrays must match
888 that of leaf_elements. The function is intended to be used for large
889 sets of elements and its complexity is linear in the number of leaf
890 elements in the mesh. */
891 void FindSetNeighbors(const Array<char> &elem_set,
892 Array<int> *neighbors, /* append */
893 Array<char> *neighbor_set = NULL);
894
895 /** Return all vertex-, edge- and face-neighbors of a single element.
896 You can limit the number of elements being checked using 'search_set'.
897 The complexity of the function is linear in the size of the search set.*/
898 void FindNeighbors(int elem,
899 Array<int> &neighbors, /* append */
900 const Array<int> *search_set = NULL);
901
902 /** Expand a set of elements by all vertex-, edge- and face-neighbors.
903 The output array 'expanded' will contain all items from 'elems'
904 (provided they are in 'search_set') plus their neighbors. The neighbor
905 search can be limited to the optional search set. The complexity is
906 linear in the sum of the sizes of 'elems' and 'search_set'. */
907 void NeighborExpand(const Array<int> &elems,
908 Array<int> &expanded,
909 const Array<int> *search_set = NULL);
910
911
912 void CollectEdgeVertices(int v0, int v1, Array<int> &indices);
913 void CollectTriFaceVertices(int v0, int v1, int v2, Array<int> &indices);
914 void CollectQuadFaceVertices(int v0, int v1, int v2, int v3,
915 Array<int> &indices);
917
922
923 int GetVertexRootCoord(int elem, RefCoord coord[3]) const;
924 void CollectIncidentElements(int elem, const RefCoord coord[3],
925 Array<int> &list) const;
926
927 /** Return elements neighboring to a local vertex of element 'elem'. Only
928 elements from within the same refinement tree ('cousins') are returned.
929 Complexity is proportional to the depth of elem's refinement tree. */
930 void FindVertexCousins(int elem, int local, Array<int> &cousins) const;
931
932
933 // coarse/fine transformations
934
935 struct Point
936 {
937 int dim;
939
940 Point() { dim = 0; }
941
942 Point(const Point &) = default;
943
945 { dim = 1; coord[0] = x; }
946
948 { dim = 2; coord[0] = x; coord[1] = y; }
949
951 { dim = 3; coord[0] = x; coord[1] = y; coord[2] = z; }
952
953 Point(const Point& p0, const Point& p1)
954 {
955 dim = p0.dim;
956 for (int i = 0; i < dim; i++)
957 {
958 coord[i] = (p0.coord[i] + p1.coord[i]) * 0.5;
959 }
960 }
961
962 Point(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
963 {
964 dim = p0.dim;
965 MFEM_ASSERT(p1.dim == dim && p2.dim == dim && p3.dim == dim, "");
966 for (int i = 0; i < dim; i++)
967 {
968 coord[i] = (p0.coord[i] + p1.coord[i] + p2.coord[i] + p3.coord[i])
969 * 0.25;
970 }
971 }
972
973 Point& operator=(const Point& src)
974 {
975 dim = src.dim;
976 for (int i = 0; i < dim; i++) { coord[i] = src.coord[i]; }
977 return *this;
978 }
979 };
980
981 /** @brief The PointMatrix stores the coordinates of the slave face using the
982 master face coordinate as reference.
983
984 In 2D, the point matrix has the orientation of the parent
985 edge, so its columns need to be flipped when applying it, see
986 ApplyLocalSlaveTransformation.
987
988 In 3D, the orientation part of Elem2Inf is encoded in the point
989 matrix.
990
991 The following transformation gives the relation between the
992 reference quad face coordinates (xi, eta) in [0,1]^2, and the fine quad
993 face coordinates (x, y):
994 x = a0*(1-xi)*(1-eta) + a1*xi*(1-eta) + a2*xi*eta + a3*(1-xi)*eta
995 y = b0*(1-xi)*(1-eta) + b1*xi*(1-eta) + b2*xi*eta + b3*(1-xi)*eta
996 */
998 {
999 int np;
1001
1002 PointMatrix() : np(0) {}
1003
1004 PointMatrix(const Point& p0, const Point& p1)
1005 { np = 2; points[0] = p0; points[1] = p1; }
1006
1007 PointMatrix(const Point& p0, const Point& p1, const Point& p2)
1008 { np = 3; points[0] = p0; points[1] = p1; points[2] = p2; }
1009
1010 PointMatrix(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
1011 { np = 4; points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3; }
1012
1013 PointMatrix(const Point& p0, const Point& p1, const Point& p2,
1014 const Point& p3, const Point& p4)
1015 {
1016 np = 5;
1017 points[0] = p0; points[1] = p1; points[2] = p2;
1018 points[3] = p3; points[4] = p4;
1019 }
1020 PointMatrix(const Point& p0, const Point& p1, const Point& p2,
1021 const Point& p3, const Point& p4, const Point& p5)
1022 {
1023 np = 6;
1024 points[0] = p0; points[1] = p1; points[2] = p2;
1025 points[3] = p3; points[4] = p4; points[5] = p5;
1026 }
1027 PointMatrix(const Point& p0, const Point& p1, const Point& p2,
1028 const Point& p3, const Point& p4, const Point& p5,
1029 const Point& p6, const Point& p7)
1030 {
1031 np = 8;
1032 points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3;
1033 points[4] = p4; points[5] = p5; points[6] = p6; points[7] = p7;
1034 }
1035
1036 Point& operator()(int i) { return points[i]; }
1037 const Point& operator()(int i) const { return points[i]; }
1038
1039 bool operator==(const PointMatrix &pm) const;
1040
1041 void GetMatrix(DenseMatrix& point_matrix) const;
1042 };
1043
1051
1052 static const PointMatrix& GetGeomIdentity(Geometry::Type geom);
1053
1054 void GetPointMatrix(Geometry::Type geom, const char* ref_path,
1055 DenseMatrix& matrix) const;
1056
1057 typedef std::map<std::string, int> RefPathMap;
1058
1059 void TraverseRefinements(int elem, int coarse_index,
1060 std::string &ref_path, RefPathMap &map) const;
1061
1062 /// storage for data returned by Get[De]RefinementTransforms()
1064
1065 /// state of leaf_elements before Refine(), set by MarkCoarseLevel()
1067
1068 void InitDerefTransforms();
1069 void SetDerefMatrixCodes(int parent, Array<int> &fine_coarse);
1070
1071 // vertex temporary data, used by GetMeshComponents
1073 {
1076 TmpVertex() : valid(false), visited(false) {}
1077 };
1078
1080
1081 const real_t *CalcVertexPos(int node) const;
1082
1083
1084 // utility
1085
1086 int GetEdgeMaster(int node) const;
1087
1088 void FindFaceNodes(int face, int node[4]) const;
1089
1090 /**
1091 * @brief Return the number of splits of this edge that have occurred in the
1092 * NCMesh. If zero, this means the segment is not the master of any other segments.
1093 *
1094 * @param vn1 The first vertex making up the segment
1095 * @param vn2 The second vertex making up the segment
1096 * @return int The depth of splits of this segment that are present in the mesh.
1097 */
1098 int EdgeSplitLevel(int vn1, int vn2) const;
1099 /**
1100 * @brief Return the number of splits of this triangle that have occurred in
1101 * the NCMesh. If zero, this means the triangle is neither split, nor the
1102 * master of a split face.
1103 *
1104 * @param vn1 The first vertex making up the triangle
1105 * @param vn2 The second vertex making up the triangle
1106 * @param vn3 The third vertex making up the triangle
1107 * @return int The depth of splits of this triangle that are present in the mesh.
1108 */
1109 int TriFaceSplitLevel(int vn1, int vn2, int vn3) const;
1110 /**
1111 * @brief Computes the number of horizontal and vertical splits of this quad
1112 * that have occurred in the NCMesh. If zero, this means the quad is not
1113 * the master of any other quad.
1114 *
1115 * @param vn1 The first vertex making up the quad
1116 * @param vn2 The second vertex making up the quad
1117 * @param vn3 The third vertex making up the quad
1118 * @param vn4 The fourth vertex making up the quad
1119 * @param h_level The number of "horizontal" splits of the quad
1120 * @param v_level The number of "vertical" splits of the quad
1121 */
1122 void QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4,
1123 int& h_level, int& v_level) const;
1124 /**
1125 * @brief Returns the total number of splits of this quad that have occurred
1126 * in the NCMesh. If zero, this means the quad is not
1127 * the master of any other quad.
1128 * @details This is a convenience wrapper that sums the horizontal and
1129 * vertical levels from the full method.
1130 *
1131 * @param vn1 The first vertex making up the quad
1132 * @param vn2 The second vertex making up the quad
1133 * @param vn3 The third vertex making up the quad
1134 * @param vn4 The fourth vertex making up the quad
1135 * @return int The depth of splits of this triangle that are present in the
1136 * mesh. NB: An isotropic refinement has a level of 2, one horizontal split,
1137 * followed by a vertical split.
1138 */
1139 int QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4) const;
1140
1141 void CountSplits(int elem, int splits[3]) const;
1142 void GetLimitRefinements(Array<Refinement> &refinements, int max_level);
1143
1144
1145 // I/O
1146
1147 /// Print the "vertex_parents" section of the mesh file.
1148 int PrintVertexParents(std::ostream *out) const;
1149 /// Load the vertex parent hierarchy from a mesh file.
1150 void LoadVertexParents(std::istream &input);
1151
1152 /** Print the "boundary" section of the mesh file.
1153 If out == NULL, only return the number of boundary elements. */
1154 int PrintBoundary(std::ostream *out) const;
1155 /// Load the "boundary" section of the mesh file.
1156 void LoadBoundary(std::istream &input);
1157
1158 /// Print the "coordinates" section of the mesh file.
1159 void PrintCoordinates(std::ostream &out) const;
1160 /// Load the "coordinates" section of the mesh file.
1161 void LoadCoordinates(std::istream &input);
1162
1163 /// Count root elements and initialize root_state.
1164 void InitRootElements();
1165 /// Return the index of the last top-level node plus one.
1166 int CountTopLevelNodes() const;
1167 /// Return true if all root_states are zero.
1168 bool ZeroRootStates() const;
1169
1170 /// Load the element refinement hierarchy from a legacy mesh file.
1171 void LoadCoarseElements(std::istream &input);
1172 void CopyElements(int elem, const BlockArray<Element> &tmp_elements);
1173 /// Load the deprecated MFEM mesh v1.1 format for backward compatibility.
1174 void LoadLegacyFormat(std::istream &input, int &curved, int &is_nc);
1175
1176 // geometry
1177
1178 /// This holds in one place the constants about the geometries we support
1180 {
1181 int nv, ne, nf; // number of: vertices, edges, faces
1182 int edges[MaxElemEdges][2]; // edge vertices (up to 12 edges)
1183 int faces[MaxElemFaces][4]; // face vertices (up to 6 faces)
1184 int nfv[MaxElemFaces]; // number of face vertices
1185
1187 GeomInfo() : initialized(false) {}
1188 void InitGeom(Geometry::Type geom);
1189 };
1190
1192
1193#ifdef MFEM_DEBUG
1194public:
1195 void DebugLeafOrder(std::ostream &out) const;
1196 void DebugDump(std::ostream &out) const;
1197#endif
1198
1199 friend class ParNCMesh; // for ParNCMesh::ElementSet
1200 friend struct MatrixMap;
1201 friend struct PointMatrixHash;
1202};
1203
1204}
1205
1206#endif
int Size() const
Return the logical size of the array.
Definition array.hpp:144
void DeleteAll()
Delete the whole array.
Definition array.hpp:864
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
void DeleteLast()
Delete the last entry of the array.
Definition array.hpp:199
T & Last()
Return the last element in the array.
Definition array.hpp:802
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Rank 3 tensor (array of matrices)
Abstract data type element.
Definition element.hpp:29
static const int NumGeom
Definition geom.hpp:42
Mesh data type.
Definition mesh.hpp:56
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
Definition ncmesh.hpp:123
int GetEdgeMaster(int v1, int v2) const
Definition ncmesh.cpp:5208
static PointMatrix pm_tet_identity
Definition ncmesh.hpp:1047
void OnMeshUpdated(Mesh *mesh)
Definition ncmesh.cpp:2590
void GetElementFacesAttributes(int i, Array< int > &faces, Array< int > &fattr) const
Return the faces and face attributes of leaf element 'i'.
Definition ncmesh.cpp:5243
static GeomInfo GI[Geometry::NumGeom]
Definition ncmesh.hpp:1191
virtual void Update()
Definition ncmesh.cpp:268
void FindNeighbors(int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
Definition ncmesh.cpp:3841
void FreeElement(int id)
Definition ncmesh.hpp:700
virtual void Trim()
Save memory by releasing all non-essential and cached data.
Definition ncmesh.cpp:6326
virtual void ElementSharesVertex(int elem, int local, int vnode)
Definition ncmesh.hpp:881
HashTable< Node > shadow
temporary storage for reparented nodes
Definition ncmesh.hpp:674
void SetDerefMatrixCodes(int parent, Array< int > &fine_coarse)
Definition ncmesh.cpp:2064
void CopyElements(int elem, const BlockArray< Element > &tmp_elements)
Definition ncmesh.cpp:6054
Array< Refinement > ref_stack
stack of scheduled refinements (temporary)
Definition ncmesh.hpp:673
int QuadFaceSplitType(int n1, int n2, int n3, int n4, int mid[5]=NULL) const
Given a quad face defined by four vertices, establish which edges of this face have been split,...
Definition ncmesh.cpp:2735
mfem::Element * NewMeshElement(int geom) const
Definition ncmesh.cpp:2423
int NewTriangle(int n0, int n1, int n2, int attr, int eattr0, int eattr1, int eattr2)
Definition ncmesh.cpp:646
void MakeTopologyOnly()
Definition ncmesh.hpp:478
void GetMeshComponents(Mesh &mesh) const
Fill Mesh::{vertices,elements,boundary} for the current finest level.
Definition ncmesh.cpp:2463
int NGhostElements
Definition ncmesh.hpp:602
void LoadCoarseElements(std::istream &input)
Load the element refinement hierarchy from a legacy mesh file.
Definition ncmesh.cpp:6072
void NeighborExpand(const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
Definition ncmesh.cpp:3924
HashTable< Node > nodes
Definition ncmesh.hpp:572
int PrintVertexParents(std::ostream *out) const
Print the "vertex_parents" section of the mesh file.
Definition ncmesh.cpp:5567
static PointMatrix pm_prism_identity
Definition ncmesh.hpp:1048
int NewSegment(int n0, int n1, int attr, int vattr1, int vattr2)
Definition ncmesh.cpp:672
static int find_node(const Element &el, int node)
Definition ncmesh.cpp:2800
void DeleteUnusedFaces(const Array< int > &elemFaces)
Definition ncmesh.cpp:440
const CoarseFineTransformations & GetRefinementTransforms() const
Definition ncmesh.cpp:4674
void BuildElementToVertexTable()
Definition ncmesh.cpp:3650
virtual int GetNGhostElements() const
Definition ncmesh.hpp:155
void TraverseEdge(int vn0, int vn1, real_t t0, real_t t1, int flags, int level, MatrixMap &matrix_map)
Definition ncmesh.cpp:3293
static int find_element_edge(const Element &el, int vn0, int vn1, bool abort=true)
Definition ncmesh.cpp:2820
Array< int > free_element_ids
Definition ncmesh.hpp:576
void LoadBoundary(std::istream &input)
Load the "boundary" section of the mesh file.
Definition ncmesh.cpp:5658
int Dimension() const
Return the dimension of the NCMesh.
Definition ncmesh.hpp:145
virtual void ElementSharesFace(int elem, int local, int face)
Definition ncmesh.hpp:879
int PrintMemoryDetail() const
Definition ncmesh.cpp:6392
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:503
HashTable< Face > faces
Definition ncmesh.hpp:573
bool HaveTets() const
Return true if the mesh contains tetrahedral elements.
Definition ncmesh.hpp:665
static const int MaxElemChildren
Number of children of an element can have.
Definition ncmesh.hpp:494
void GetEdgeVertices(const MeshId &edge_id, int vert_index[2], bool oriented=true) const
Return Mesh vertex indices of an edge identified by 'edge_id'.
Definition ncmesh.cpp:5109
static PointMatrix pm_seg_identity
Definition ncmesh.hpp:1044
static const int MaxElemEdges
Number of edges of an element can have.
Definition ncmesh.hpp:490
const NCList & GetNCList(int entity)
Return vertex/edge/face list (entity = 0/1/2, respectively).
Definition ncmesh.hpp:320
bool QuadFaceIsMaster(int n1, int n2, int n3, int n4) const
Determine if a Quad face is a master face.
Definition ncmesh.hpp:803
Array< int > boundary_faces
subset of all faces, set by BuildFaceList
Definition ncmesh.hpp:612
static PointMatrix pm_quad_identity
Definition ncmesh.hpp:1046
void ReferenceElement(int elem)
Definition ncmesh.cpp:352
void UpdateElementToVertexTable()
Definition ncmesh.hpp:918
void ReparentNode(int node, int new_p1, int new_p2)
Definition ncmesh.cpp:305
BlockArray< Element > elements
Definition ncmesh.hpp:575
Array< int > coarse_elements
state of leaf_elements before Refine(), set by MarkCoarseLevel()
Definition ncmesh.hpp:1066
const CoarseFineTransformations & GetDerefinementTransforms() const
Definition ncmesh.cpp:4725
Table derefinements
possible derefinements, see GetDerefinementTable
Definition ncmesh.hpp:677
int FindNodeExt(const Element &el, int node, bool abort=true)
Extended version of find_node: works if 'el' is refined.
Definition ncmesh.cpp:2810
int ReorderFacePointMat(int v0, int v1, int v2, int v3, int elem, const PointMatrix &pm, PointMatrix &reordered) const
Definition ncmesh.cpp:2945
int FindMidEdgeNode(int node1, int node2) const
Definition ncmesh.cpp:321
Array< char > face_geom
face geometry by face index, set by OnMeshUpdated
Definition ncmesh.hpp:613
int GetEdgeNCOrientation(const MeshId &edge_id) const
Definition ncmesh.cpp:5128
Table element_vertex
leaf-element to vertex table, see FindSetNeighbors
Definition ncmesh.hpp:615
Array< Triple< int, int, int > > reparents
scheduled node reparents (tmp)
Definition ncmesh.hpp:675
virtual void BuildFaceList()
Definition ncmesh.cpp:3198
void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4, const Refinement *refs, int nref)
Definition ncmesh.cpp:817
void LoadVertexParents(std::istream &input)
Load the vertex parent hierarchy from a mesh file.
Definition ncmesh.cpp:5596
TmpVertex * tmp_vertex
Definition ncmesh.hpp:1079
static void GridSfcOrdering3D(int width, int height, int depth, Array< int > &coords)
Definition ncmesh.cpp:5077
Array< int > leaf_elements
finest elements, in Mesh ordering (+ ghosts)
Definition ncmesh.hpp:604
const real_t * CalcVertexPos(int node) const
Definition ncmesh.cpp:2438
void CollectLeafElements(int elem, int state, Array< int > &ghosts, int &counter)
Definition ncmesh.cpp:2083
int CountTopLevelNodes() const
Return the index of the last top-level node plus one.
Definition ncmesh.cpp:5870
void QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4, int &h_level, int &v_level) const
Computes the number of horizontal and vertical splits of this quad that have occurred in the NCMesh....
Definition ncmesh.cpp:5395
int Geoms
bit mask of element geometries present, see InitGeomFlags()
Definition ncmesh.hpp:485
void InitDerefTransforms()
Definition ncmesh.cpp:2045
int GetElementSizeReduction(int i) const
Definition ncmesh.cpp:5229
virtual void LimitNCLevel(int max_nc_level)
Definition ncmesh.cpp:5549
int NGhostVertices
Definition ncmesh.hpp:602
const NCList & GetVertexList()
Definition ncmesh.hpp:313
bool ZeroRootStates() const
Return true if all root_states are zero.
Definition ncmesh.cpp:5737
NCMesh(const Mesh *mesh)
Definition ncmesh.cpp:135
Array< int > vertex_nodeId
vertex-index to node-id map, see UpdateVertices
Definition ncmesh.hpp:606
int GetNumRootElements()
Return the number of root elements.
Definition ncmesh.hpp:429
int RetrieveNode(const Element &el, int index)
Return el.node[index] correctly, even if the element is refined.
Definition ncmesh.cpp:1696
void FindEdgeElements(int vn1, int vn2, int vn3, int vn4, Array< MeshId > &prisms) const
Definition ncmesh.cpp:770
virtual ~NCMesh()
Definition ncmesh.cpp:280
void CollectDerefinements(int elem, Array< Connection > &list)
Definition ncmesh.cpp:1932
bool IsLegacyLoaded() const
I/O: Return true if the mesh was loaded from the legacy v1.1 format.
Definition ncmesh.hpp:449
static const int MaxElemFaces
Number of faces of an element can have.
Definition ncmesh.hpp:492
virtual void ElementSharesEdge(int elem, int local, int enode)
Definition ncmesh.hpp:880
int GetFaceVerticesEdges(const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
Definition ncmesh.cpp:5140
void RefineElement(int elem, char ref_type)
Definition ncmesh.cpp:954
void InitRootState(int root_count)
Definition ncmesh.cpp:2359
Geometry::Type GetElementGeometry(int index) const
Return element geometry type. index is the Mesh element number.
Definition ncmesh.hpp:421
int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
Definition ncmesh.cpp:620
Array< int > leaf_sfc_index
natural tree ordering of leaf elements
Definition ncmesh.hpp:605
void ForceRefinement(int vn1, int vn2, int vn3, int vn4)
Definition ncmesh.cpp:712
void CollectEdgeVertices(int v0, int v1, Array< int > &indices)
Definition ncmesh.cpp:3580
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
Definition ncmesh.hpp:609
static PointMatrix pm_pyramid_identity
Definition ncmesh.hpp:1049
Array< int > root_state
Definition ncmesh.hpp:581
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition ncmesh.hpp:298
Geometry::Type GetFaceGeometry(int index) const
Return face geometry type. index is the Mesh face number.
Definition ncmesh.hpp:425
void CollectTriFaceVertices(int v0, int v1, int v2, Array< int > &indices)
Definition ncmesh.cpp:3592
void InitRootElements()
Count root elements and initialize root_state.
Definition ncmesh.cpp:5830
int GetMidEdgeNode(int node1, int node2)
Definition ncmesh.cpp:337
virtual void BuildEdgeList()
Definition ncmesh.cpp:3321
void DebugDump(std::ostream &out) const
Definition ncmesh.cpp:6445
void DebugLeafOrder(std::ostream &out) const
Definition ncmesh.cpp:6420
int GetElementDepth(int i) const
Return the distance of leaf 'i' from the root.
Definition ncmesh.cpp:5217
void TraverseRefinements(int elem, int coarse_index, std::string &ref_path, RefPathMap &map) const
Definition ncmesh.cpp:4640
void CheckIsoFace(int vn1, int vn2, int vn3, int vn4, int en1, int en2, int en3, int en4, int midf)
Definition ncmesh.cpp:937
void GetLimitRefinements(Array< Refinement > &refinements, int max_level)
Definition ncmesh.cpp:5519
int NewPyramid(int n0, int n1, int n2, int n3, int n4, int attr, int fattr0, int fattr1, int fattr2, int fattr3, int fattr4)
Definition ncmesh.cpp:590
virtual void Derefine(const Array< int > &derefs)
Definition ncmesh.cpp:2009
bool Iso
true if the mesh only contains isotropic refinements
Definition ncmesh.hpp:484
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges, Array< int > &bdr_faces)
Get a list of vertices (2D/3D), edges (3D) and faces (3D) that coincide with boundary elements with t...
Definition ncmesh.cpp:5286
bool TriFaceSplit(int n1, int n2, int n3, int mid[3]=NULL) const
Given a tri face defined by three vertices, establish whether the edges that make up this face have b...
Definition ncmesh.cpp:2782
bool Legacy
true if the mesh was loaded from the legacy v1.1 format
Definition ncmesh.hpp:486
virtual void BuildVertexList()
Definition ncmesh.cpp:3433
static const PointMatrix & GetGeomIdentity(Geometry::Type geom)
Definition ncmesh.cpp:4127
int GetNVertices() const
Return the number of vertices in the NCMesh.
Definition ncmesh.hpp:150
void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4, int mid12, int mid34, int level=0)
Definition ncmesh.cpp:842
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:533
int GetVertexRootCoord(int elem, RefCoord coord[3]) const
Definition ncmesh.cpp:3983
Array< real_t > coordinates
Definition ncmesh.hpp:585
void FindVertexCousins(int elem, int local, Array< int > &cousins) const
Definition ncmesh.cpp:4059
Face * GetFace(Element &elem, int face_no)
Definition ncmesh.cpp:465
std::map< std::string, int > RefPathMap
Definition ncmesh.hpp:1057
bool IsGhost(const Element &el) const
Return true if the Element el is a ghost element.
Definition ncmesh.hpp:668
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition ncmesh.hpp:305
bool TriFaceIsMaster(int n1, int n2, int n3) const
Determine if a Triangle face is a master face.
Definition ncmesh.hpp:783
int MyRank
used in parallel, or when loading a parallel file in serial
Definition ncmesh.hpp:483
void Print(std::ostream &out, const std::string &comments="") const
Definition ncmesh.cpp:5746
int AddElement(const Element &el)
Definition ncmesh.hpp:687
void CountSplits(int elem, int splits[3]) const
Definition ncmesh.cpp:5429
int spaceDim
dimensions of the elements and the vertex coordinates
Definition ncmesh.hpp:482
void UnreferenceElement(int elem, Array< int > &elemFaces)
Definition ncmesh.cpp:383
void LegacyToNewVertexOrdering(Array< int > &order) const
I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
Definition ncmesh.cpp:6301
int NGhostEdges
Definition ncmesh.hpp:602
void CollectIncidentElements(int elem, const RefCoord coord[3], Array< int > &list) const
Definition ncmesh.cpp:4036
void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1, MatrixMap &matrix_map)
Definition ncmesh.cpp:3095
void RegisterFaces(int elem, int *fattr=NULL)
Definition ncmesh.cpp:426
int PrintBoundary(std::ostream *out) const
Definition ncmesh.cpp:5619
static const int MaxElemNodes
Number of nodes of an element can have.
Definition ncmesh.hpp:488
CoarseFineTransformations transforms
storage for data returned by Get[De]RefinementTransforms()
Definition ncmesh.hpp:1063
void MarkCoarseLevel()
Definition ncmesh.cpp:4626
int TriFaceSplitLevel(int vn1, int vn2, int vn3) const
Return the number of splits of this triangle that have occurred in the NCMesh. If zero,...
Definition ncmesh.cpp:5380
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
Definition ncmesh.cpp:1980
void UpdateVertices()
This method assigns indices to vertices (Node::vert_index) that will be seen by the Mesh class and th...
Definition ncmesh.cpp:2179
int GetNFaces() const
Return the number of (2D) faces in the NCMesh.
Definition ncmesh.hpp:154
NCList vertex_list
lazy-initialized list of vertices, see GetVertexList
Definition ncmesh.hpp:610
bool HavePyramids() const
Return true if the mesh contains pyramid elements.
Definition ncmesh.hpp:662
void InitGeomFlags()
Definition ncmesh.cpp:259
int NewTetrahedron(int n0, int n1, int n2, int n3, int attr, int fattr0, int fattr1, int fattr2, int fattr3)
Definition ncmesh.cpp:565
void LoadCoordinates(std::istream &input)
Load the "coordinates" section of the mesh file.
Definition ncmesh.cpp:5716
void GetPointMatrix(Geometry::Type geom, const char *ref_path, DenseMatrix &matrix) const
Definition ncmesh.cpp:4144
void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level, Face *eface[4], MatrixMap &matrix_map)
Definition ncmesh.cpp:2976
const Table & GetDerefinementTable()
Definition ncmesh.cpp:1965
void FindFaceNodes(int face, int node[4]) const
Definition ncmesh.cpp:5264
void ClearTransforms()
Free all internal data created by the above three functions.
Definition ncmesh.cpp:4805
int SpaceDimension() const
Return the space dimension of the NCMesh.
Definition ncmesh.hpp:147
void FindSetNeighbors(const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
Definition ncmesh.cpp:3730
NCMesh & operator=(NCMesh &)=delete
Copy assignment not supported.
void CollectQuadFaceVertices(int v0, int v1, int v2, int v3, Array< int > &indices)
Definition ncmesh.cpp:3616
void PrintCoordinates(std::ostream &out) const
Print the "coordinates" section of the mesh file.
Definition ncmesh.cpp:5698
TriFaceTraverseResults TraverseTriFace(int vn0, int vn1, int vn2, const PointMatrix &pm, int level, MatrixMap &matrix_map)
Definition ncmesh.cpp:3132
int GetMidFaceNode(int en1, int en2, int en3, int en4)
Definition ncmesh.cpp:344
bool HavePrisms() const
Return true if the mesh contains prism elements.
Definition ncmesh.hpp:659
long MemoryUsage() const
Return total number of bytes allocated.
Definition ncmesh.cpp:6369
void DerefineElement(int elem)
Derefine the element elem, does nothing on leaf elements.
Definition ncmesh.cpp:1738
NCList face_list
lazy-initialized list of faces, see GetFaceList
Definition ncmesh.hpp:608
virtual void Refine(const Array< Refinement > &refinements)
Definition ncmesh.cpp:1648
std::int64_t RefCoord
Definition ncmesh.hpp:462
static PointMatrix pm_tri_identity
Definition ncmesh.hpp:1045
int NGhostFaces
Definition ncmesh.hpp:602
static PointMatrix pm_hex_identity
Definition ncmesh.hpp:1050
void UpdateLeafElements()
Update the leaf elements indices in leaf_elements.
Definition ncmesh.cpp:2151
void LoadLegacyFormat(std::istream &input, int &curved, int &is_nc)
Load the deprecated MFEM mesh v1.1 format for backward compatibility.
Definition ncmesh.cpp:6143
int GetNEdges() const
Return the number of edges in the NCMesh.
Definition ncmesh.hpp:152
int EdgeSplitLevel(int vn1, int vn2) const
Return the number of splits of this edge that have occurred in the NCMesh. If zero,...
Definition ncmesh.cpp:5373
static void GridSfcOrdering2D(int width, int height, Array< int > &coords)
Definition ncmesh.cpp:5062
static int find_local_face(int geom, int a, int b, int c)
Definition ncmesh.cpp:2838
A parallel extension of the NCMesh class.
Definition pncmesh.hpp:65
Data type point element.
Definition point.hpp:23
int Size() const
Returns the number of TYPE I elements.
Definition table.hpp:92
int index(int i, int j, int nx, int ny)
Definition life.cpp:235
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
void Swap(Array< T > &, Array< T > &)
Definition array.hpp:648
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
float real_t
Definition config.hpp:43
Defines the coarse-fine transformations of all fine elements.
Definition ncmesh.hpp:72
MFEM_DEPRECATED void GetCoarseToFineMap(const Mesh &fine_mesh, Table &coarse_to_fine) const
Definition ncmesh.hpp:91
Array< Embedding > embeddings
Fine element positions in their parents.
Definition ncmesh.hpp:74
void MakeCoarseToFineTable(Table &coarse_to_fine, bool want_ghosts=false) const
Definition ncmesh.cpp:4783
DenseTensor point_matrices[Geometry::NumGeom]
Definition ncmesh.hpp:78
Defines the position of a fine element within a coarse element.
Definition ncmesh.hpp:51
unsigned geom
Definition ncmesh.hpp:58
Embedding()=default
unsigned ghost
For internal use: 0 if regular fine element, 1 if parallel ghost element.
Definition ncmesh.hpp:62
Embedding(int elem, Geometry::Type geom, int matrix=0, bool ghost=false)
Definition ncmesh.hpp:65
int parent
Coarse Element index in the coarse mesh.
Definition ncmesh.hpp:53
unsigned matrix
Definition ncmesh.hpp:59
int rank
processor number (ParNCMesh), -1 if undefined/unknown
Definition ncmesh.hpp:554
bool IsLeaf() const
Definition ncmesh.hpp:566
int child[MaxElemChildren]
2-10 children (if ref_type != 0)
Definition ncmesh.hpp:559
char tet_type
tetrahedron split type, currently always 0
Definition ncmesh.hpp:551
Element(Geometry::Type geom, int attr)
Definition ncmesh.cpp:490
char flag
generic flag/marker, can be used by algorithms
Definition ncmesh.hpp:552
int node[MaxElemNodes]
element corners (if ref_type == 0)
Definition ncmesh.hpp:558
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
Definition ncmesh.hpp:550
char geom
Geometry::Type of the element (char for storage only)
Definition ncmesh.hpp:549
int index
element number in the Mesh, -1 if refined
Definition ncmesh.hpp:553
int parent
parent element, -1 if this is a root element, -2 if free'd
Definition ncmesh.hpp:561
Geometry::Type Geom() const
Definition ncmesh.hpp:565
void ForgetElement(int e)
Definition ncmesh.cpp:458
int elem[2]
up to 2 elements sharing the face
Definition ncmesh.hpp:529
bool Unused() const
Definition ncmesh.hpp:534
void RegisterElement(int e)
Definition ncmesh.cpp:451
bool Boundary() const
Definition ncmesh.hpp:533
int index
face number in the Mesh
Definition ncmesh.hpp:528
int GetSingleElement() const
Return one of elem[0] or elem[1] and make sure the other is -1.
Definition ncmesh.cpp:473
int attribute
boundary element attribute, -1 if internal face
Definition ncmesh.hpp:527
This holds in one place the constants about the geometries we support.
Definition ncmesh.hpp:1180
int nfv[MaxElemFaces]
Definition ncmesh.hpp:1184
int faces[MaxElemFaces][4]
Definition ncmesh.hpp:1183
int edges[MaxElemEdges][2]
Definition ncmesh.hpp:1182
void InitGeom(Geometry::Type geom)
Definition ncmesh.cpp:57
Master(int index, int element, int local, int geom, int sb, int se)
Definition ncmesh.hpp:209
int slaves_end
slave faces
Definition ncmesh.hpp:206
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition ncmesh.hpp:189
int element
NCMesh::Element containing this vertex/edge/face.
Definition ncmesh.hpp:191
int index
Mesh number.
Definition ncmesh.hpp:190
MeshId(int index, int element, int local, int geom=-1)
Definition ncmesh.hpp:198
Geometry::Type Geom() const
Definition ncmesh.hpp:195
signed char geom
Geometry::Type (faces only) (char to save RAM)
Definition ncmesh.hpp:193
signed char local
local number within 'element'
Definition ncmesh.hpp:192
const MeshIdType type
MeshIdType corresponding to the MeshId. UNRECOGNIZED if unfound.
Definition ncmesh.hpp:255
Lists all edges/faces in the nonconforming mesh.
Definition ncmesh.hpp:230
MeshIdType GetMeshIdType(int index) const
Return a face type for a given nc index.
Definition ncmesh.cpp:3529
Array< MeshId > conforming
All MeshIds corresponding to conformal faces.
Definition ncmesh.hpp:231
long MemoryUsage() const
Definition ncmesh.cpp:6341
Array< Slave > slaves
All MeshIds corresponding to slave faces.
Definition ncmesh.hpp:233
bool CheckMeshIdType(int index, MeshIdType type) const
Given an index, check if this is a certain face type.
Definition ncmesh.cpp:3537
bool Empty() const
Whether the NCList is empty.
Definition ncmesh.hpp:269
void Clear()
Erase the contents of the conforming, master and slave arrays.
Definition ncmesh.cpp:3490
void OrientedPointMatrix(const Slave &slave, DenseMatrix &oriented_matrix) const
Definition ncmesh.cpp:3468
Array< Master > masters
All MeshIds corresponding to master faces.
Definition ncmesh.hpp:232
long TotalSize() const
The total size of the component arrays in the NCList.
Definition ncmesh.hpp:276
Array< DenseMatrix * > point_matrices[Geometry::NumGeom]
List of unique point matrices for each slave geometry.
Definition ncmesh.hpp:236
MeshIdAndType GetMeshIdAndType(int index) const
Return a mesh id and type for a given nc index.
Definition ncmesh.cpp:3509
bool HasVertex() const
Definition ncmesh.hpp:513
bool HasEdge() const
Definition ncmesh.hpp:514
The PointMatrix stores the coordinates of the slave face using the master face coordinate as referenc...
Definition ncmesh.hpp:998
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4)
Definition ncmesh.hpp:1013
Point & operator()(int i)
Definition ncmesh.hpp:1036
void GetMatrix(DenseMatrix &point_matrix) const
Definition ncmesh.cpp:4090
PointMatrix(const Point &p0, const Point &p1, const Point &p2)
Definition ncmesh.hpp:1007
Point points[MaxElemNodes]
Definition ncmesh.hpp:1000
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition ncmesh.hpp:1010
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4, const Point &p5)
Definition ncmesh.hpp:1020
const Point & operator()(int i) const
Definition ncmesh.hpp:1037
bool operator==(const PointMatrix &pm) const
Definition ncmesh.cpp:4076
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:1027
PointMatrix(const Point &p0, const Point &p1)
Definition ncmesh.hpp:1004
Point(const Point &p0, const Point &p1)
Definition ncmesh.hpp:953
Point(const Point &)=default
Point(real_t x, real_t y, real_t z)
Definition ncmesh.hpp:950
Point(real_t x, real_t y)
Definition ncmesh.hpp:947
Point & operator=(const Point &src)
Definition ncmesh.hpp:973
Point(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition ncmesh.hpp:962
Nonconforming edge/face within a bigger edge/face.
Definition ncmesh.hpp:216
unsigned edge_flags
orientation flags, see OrientedPointMatrix
Definition ncmesh.hpp:219
unsigned matrix
index into NCList::point_matrices[geom]
Definition ncmesh.hpp:218
int master
master number (in Mesh numbering)
Definition ncmesh.hpp:217
Slave(int index, int element, int local, int geom)
Definition ncmesh.hpp:222
bool unsplit
Whether this face has no further splits.
Definition ncmesh.hpp:864
bool ghost_neighbor
Whether the face neighbor is a ghost.
Definition ncmesh.hpp:865
Refinement()=default
char ref_type
refinement XYZ bit mask (7 = full isotropic)
Definition ncmesh.hpp:40
int index
Mesh element number.
Definition ncmesh.hpp:39
Refinement(int index, int type=Refinement::XYZ)
Definition ncmesh.hpp:44