MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
ncmesh.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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). The
36 refinement spacing or scale in each direction is a number in (0,1), with the
37 default 0.5 meaning bisection. This linear scale parameter defines the
38 position of the refinement between the beginning (0) and end (1) of the
39 reference element in each direction. */
41{
42 int index; ///< Mesh element number
43 enum : char { X = 1, Y = 2, Z = 4, XY = 3, XZ = 5, YZ = 6, XYZ = 7 };
44 using ScaledType = std::pair<char, real_t>;
45 real_t s[3]; /// Refinement scale in each dimension
46 Refinement() = default;
47 /// Refinement type XYZ, with scale 0.5.
48 Refinement(int index);
49 /// Default case of empty list @a refs is XYZ with scale 0.5.
50 Refinement(int index, const std::initializer_list<ScaledType> &refs);
51 /// Refine element with a single type and scale in all dimensions.
52 Refinement(int index, char type, real_t scale = 0.5);
53 /// Return the type as char.
54 char GetType() const;
55 /// Set the element, type, and scale.
56 void Set(int element, char type,
57 real_t scale = 0.5); /// Uses @a scale in all dimensions
58 /// Set the type and scale, assuming the element is already set.
59 void SetType(char type,
60 real_t scale = 0.5); /// Uses @a scale in all dimensions
61 /// Set the scale in the directions for the currently set type.
62 void SetScaleForType(const real_t *scale);
63private :
64 void SetScale(const ScaledType &ref);
65};
66
67/// Defines the position of a fine element within a coarse element.
69{
70 /// Coarse %Element index in the coarse mesh.
71 int parent;
72
73 /** The (geom, matrix) pair determines the sub-element transformation for the
74 fine element: CoarseFineTransformations::point_matrices[geom](matrix) is
75 the point matrix of the region within the coarse element reference
76 domain.*/
77 unsigned geom : 4;
78 unsigned matrix : 27;
79
80 /// For internal use: 0 if regular fine element, 1 if parallel ghost element.
81 unsigned ghost : 1;
82
83 Embedding() = default;
84 Embedding(int elem, Geometry::Type geom, int matrix = 0, bool ghost = false)
85 : parent(elem), geom(geom), matrix(matrix), ghost(ghost) {}
86};
87
88/// Defines the coarse-fine transformations of all fine elements.
90{
91 /// Fine element positions in their parents.
93
94 /** A "dictionary" of matrices for IsoparametricTransformation. Use
95 Embedding::{geom,matrix} to access a fine element point matrix. */
97
98 /** Invert the 'embeddings' array: create a Table with coarse elements as
99 rows and fine elements as columns. If 'want_ghosts' is false, parallel
100 ghost fine elements are not included in the table. */
101 void MakeCoarseToFineTable(Table &coarse_to_fine,
102 bool want_ghosts = false) const;
103
104 void Clear();
105 bool IsInitialized() const;
106 long MemoryUsage() const;
107
108 MFEM_DEPRECATED
109 void GetCoarseToFineMap(const Mesh &fine_mesh, Table &coarse_to_fine) const
110 { MakeCoarseToFineTable(coarse_to_fine, true); (void) fine_mesh; }
111};
112
113void Swap(CoarseFineTransformations &a, CoarseFineTransformations &b);
114
115struct MatrixMap; // for internal use
116
117/** @brief For a NURBS mesh with nonconforming patch topology, this struct
118 provides a map from hanging vertices in the patch topology to the knotvector
119 of a neighboring patch. This facilitates ensuring mesh conformity.
120 */
122{
123public:
124 /// Set the spatial dimension and number of vertices.
125 void SetSize(int dimension, int numVertices);
126
127 // The following set and get functions are for a single entry in the array of
128 // data, for a hanging vertex in the patch topology, with the given 'index'.
129 // The vertex index is 'v', parent vertices are 'pv', and knot-span is 'ks'.
130
131 /// Set the data for a vertex in 2D.
132 void SetVertex2D(int index, int v, int ks,
133 const std::array<int, 2> &pv);
134
135 /// Set the data for a vertex in 3D.
136 void SetVertex3D(int index, int v, const std::array<int, 2> &ks,
137 const std::array<int, 4> &pv);
138
139 /// Set the knot-span index for a vertex in 2D.
140 void SetKnotSpan2D(int index, int ks);
141
142 /// Set the knot-span indices for a vertex in 3D.
143 void SetKnotSpans3D(int index, const std::array<int, 2> &ks);
144
145 /// Get the data for a vertex in 2D.
146 void GetVertex2D(int index, int &v, int &ks,
147 std::array<int, 2> &pv) const;
148
149 /// Get the data for a vertex in 3D.
150 void GetVertex3D(int index, int &v, std::array<int, 2> &ks,
151 std::array<int, 4> &pv) const;
152
153 /// Print all the data.
154 void Print(std::ostream &os) const;
155
156 /// Return the number of vertices.
157 int Size() const { return data.NumRows(); }
158
159 /// Return the vertex pair representing the parent edge (2D) or face (3D).
160 std::pair<int, int> GetVertexParentPair(int index) const;
161
162private:
163 int dim; /// Spatial dimension
164 Array2D<int> data; /// Row-wise data for each vertex.
165};
166
167/** @brief A class for non-conforming AMR. The class is not used directly by the
168 * user, rather it is an extension of the Mesh class.
169 *
170 * In general, the class is used by MFEM as follows:
171 *
172 * 1. NCMesh is constructed from elements of an existing Mesh. The elements are
173 * copied and become roots of the refinement hierarchy.
174 *
175 * 2. Some elements are refined with the Refine() method. Both isotropic and
176 * anisotropic refinements of quads/hexes are supported.
177 *
178 * 3. A new Mesh is created from NCMesh containing the leaf elements. This new
179 * Mesh may have non-conforming (hanging) edges and faces and is the one
180 * seen by the user.
181 *
182 * 4. FiniteElementSpace asks NCMesh for a list of conforming, master and slave
183 * edges/faces and creates the conforming interpolation matrix P.
184 *
185 * 5. A continuous/conforming solution is obtained by solving P'*A*P x = P'*b.
186 *
187 * 6. Repeat from step 2.
188 */
190{
191protected:
192 NCMesh() = default;
193public:
194 //// Initialize with elements from an existing Mesh.
195 explicit NCMesh(const Mesh *mesh);
196
197 /** Load from a stream. The id header is assumed to have been read already
198 from \param[in] input . \param[in] version is 10 for the v1.0 NC format,
199 or 1 for the legacy v1.1 format. \param[out] curved is set to 1 if the
200 curvature GridFunction follows after mesh data. \param[out] is_nc (again
201 treated as a boolean) is set to 0 if the legacy v1.1 format in fact
202 defines a conforming mesh. See Mesh::Loader for details. */
203 NCMesh(std::istream &input, int version, int &curved, int &is_nc);
204
205 /// Deep copy of another instance.
206 NCMesh(const NCMesh &other);
207
208 /// Copy assignment not supported
210
211 virtual ~NCMesh();
212
213 /// Return the dimension of the NCMesh.
214 int Dimension() const { return Dim; }
215 /// Return the space dimension of the NCMesh.
216 int SpaceDimension() const { return spaceDim; }
217
218 /// Return the number of vertices in the NCMesh.
219 int GetNVertices() const { return NVertices; }
220 /// Return the number of edges in the NCMesh.
221 int GetNEdges() const { return NEdges; }
222 /// Return the number of (2D) faces in the NCMesh.
223 int GetNFaces() const { return NFaces; }
224 virtual int GetNGhostElements() const { return 0; }
225
226 /** Perform the given batch of refinements. Please note that in the presence
227 of anisotropic splits additional refinements may be necessary to keep the
228 mesh consistent. However, the function always performs at least the
229 requested refinements. */
230 virtual void Refine(const Array<Refinement> &refinements);
231
232 /** Check the mesh and potentially refine some elements so that the maximum
233 difference of refinement levels between adjacent elements is not greater
234 than 'max_nc_level'. */
235 virtual void LimitNCLevel(int max_nc_level);
236
237 /** Return a list of derefinement opportunities. Each row of the table
238 contains Mesh indices of existing elements that can be derefined to form
239 a single new coarse element. Row numbers are then passed to Derefine.
240 This function works both in serial and parallel. */
242
243 /** Check derefinements returned by GetDerefinementTable and mark those that
244 can be done safely so that the maximum NC level condition is not
245 violated. On return, level_ok.Size() == deref_table.Size() and contains
246 0/1s. */
247 virtual void CheckDerefinementNCLevel(const Table &deref_table,
248 Array<int> &level_ok, int max_nc_level);
249
250 /** Perform a subset of the possible derefinements (see
251 GetDerefinementTable). Note that if anisotropic refinements are present
252 in the mesh, some of the derefinements may have to be skipped to preserve
253 mesh consistency. */
254 virtual void Derefine(const Array<int> &derefs);
255
256 // master/slave lists
257
258 /// Identifies a vertex/edge/face in both Mesh and NCMesh.
259 struct MeshId
260 {
261 int index; ///< Mesh number
262 int element; ///< NCMesh::Element containing this vertex/edge/face
263 signed char local; ///< local number within 'element'
264 signed char geom; ///< Geometry::Type (faces only) (char to save RAM)
265
267
268 MeshId() = default;
269 MeshId(int index, int element, int local, int geom = -1)
271 };
272
273 /** Nonconforming edge/face that has more than one neighbor. The neighbors
274 are stored in NCList::slaves[i], slaves_begin <= i < slaves_end. */
275 struct Master : public MeshId
276 {
277 int slaves_begin, slaves_end; ///< slave faces
278
279 Master() = default;
280 Master(int index, int element, int local, int geom, int sb, int se)
282 , slaves_begin(sb), slaves_end(se) {}
283 };
284
285 /// Nonconforming edge/face within a bigger edge/face.
286 struct Slave : public MeshId
287 {
288 int master; ///< master number (in Mesh numbering)
289 unsigned matrix : 24; ///< index into NCList::point_matrices[geom]
290 unsigned edge_flags : 8; ///< orientation flags, see OrientedPointMatrix
291
292 Slave() = default;
293 Slave(int index, int element, int local, int geom)
295 , master(-1), matrix(0), edge_flags(0) {}
296 };
297
298
299 /// Lists all edges/faces in the nonconforming mesh.
300 struct NCList
301 {
302 Array<MeshId> conforming; ///< All MeshIds corresponding to conformal faces
303 Array<Master> masters; ///< All MeshIds corresponding to master faces
304 Array<Slave> slaves; ///< All MeshIds corresponding to slave faces
305
306 /// List of unique point matrices for each slave geometry.
308
309 void OrientedPointMatrix(const Slave &slave,
310 DenseMatrix &oriented_matrix) const;
311
312 /// Particular MeshId type, used for allowing static casting to the
313 /// appropriate child type after searching the NCList. UNRECOGNIZED
314 /// denotes that an instance is not known within the NCList, meaning that
315 /// it does not play a part in NC mechanics. This can be because the index
316 /// did not exist in the original Mesh, or because the entry is a boundary
317 /// face, whose NC status is always conforming.
319
320 /// Helper storing a reference to a MeshId type, and the face type it can
321 /// be cast to
323 {
324 const MeshId * const id; ///< Pointer to a possible MeshId, nullptr if not found
325 /// MeshIdType corresponding to the MeshId. UNRECOGNIZED if unfound.
327 };
328 /// Return a mesh id and type for a given nc index.
330
331 /// Return a face type for a given nc index.
332 MeshIdType GetMeshIdType(int index) const;
333
334 /// Given an index, check if this is a certain face type.
335 bool CheckMeshIdType(int index, MeshIdType type) const;
336
337 /// Erase the contents of the conforming, master and slave arrays.
338 void Clear();
339 /// Whether the NCList is empty.
340 bool Empty() const
341 {
342 return conforming.Size() == 0
343 && masters.Size() == 0
344 && slaves.Size() == 0;
345 }
346 /// The total size of the component arrays in the NCList.
347 long TotalSize() const
348 {
349 return conforming.Size() + masters.Size() + slaves.Size();
350 }
351 /// The memory usage of the three public arrays. Does not account for the
352 /// inverse index.
353 long MemoryUsage() const;
354 ~NCList() { Clear(); }
355 private:
356 // Check for existence or construct the inv_index list map if necessary.
357 // const because only modifies the mutable member inv_index.
358 void BuildIndex() const;
359
360 /// A lazily constructed map from index to MeshId. Built whenever
361 /// GetMeshIdAndType, GetMeshIdType or CheckMeshIdType is called for the
362 /// first time. The MeshIdType is stored with, to enable casting to Slave
363 /// or Master elements appropriately.
364 mutable std::unordered_map<int, std::pair<MeshIdType, int>> inv_index;
365 };
366
367
368 /// Return the current list of conforming and nonconforming faces.
370 {
371 if (face_list.Empty()) { BuildFaceList(); }
372 return face_list;
373 }
374
375 /// Return the current list of conforming and nonconforming edges.
377 {
378 if (edge_list.Empty()) { BuildEdgeList(); }
379 return edge_list;
380 }
381
382 /** Return a list of vertices (in 'conforming'); this function is provided
383 for uniformity/completeness. Needed in ParNCMesh/ParFESpace. */
385 {
386 if (vertex_list.Empty()) { BuildVertexList(); }
387 return vertex_list;
388 }
389
390 /// Return vertex/edge/face list (entity = 0/1/2, respectively).
391 const NCList& GetNCList(int entity)
392 {
393 switch (entity)
394 {
395 case 0: return GetVertexList();
396 case 1: return GetEdgeList();
397 default: return GetFaceList();
398 }
399 }
400
402 {
403 return vertex_to_knotspan;
404 }
405
406 /// Remap knot-span indices @a vertex_to_knotspan after refinement.
407 void RefineVertexToKnotSpan(const std::vector<Array<int>> &kvf,
408 const Array<KnotVector*> &kvext,
409 std::map<std::pair<int, int>,
410 std::array<int, 2>> &parentToKV);
411
412 // coarse/fine transforms
413
414 /** Remember the current layer of leaf elements before the mesh is refined.
415 Needed by GetRefinementTransforms(), must be called before Refine(). */
416 void MarkCoarseLevel();
417
418 /** After refinement, calculate the relation of each fine element to its
419 parent coarse element. Note that Refine() or LimitNCLevel() can be called
420 multiple times between MarkCoarseLevel() and this function. */
422
423 /** After derefinement, calculate the relations of previous fine elements
424 (some of which may no longer exist) to the current leaf elements. Unlike
425 for refinement, Derefine() may only be called once before this function
426 so there is no MarkFineLevel(). */
428
429 /// Free all internal data created by the above three functions.
430 void ClearTransforms();
431
432
433 // grid ordering
434
435 /** Return a space filling curve for a rectangular grid of elements.
436 Implemented is a generalized Hilbert curve for arbitrary grid dimensions.
437 If the width is odd, height should be odd too, otherwise one diagonal
438 (vertex-neighbor) step cannot be avoided in the curve. Even dimensions
439 are recommended. */
440 static void GridSfcOrdering2D(int width, int height,
441 Array<int> &coords);
442
443 /** Return a space filling curve for a 3D rectangular grid of elements. The
444 Hilbert-curve-like algorithm works well for even dimensions. For odd
445 width/height/depth it tends to produce some diagonal (edge-neighbor)
446 steps. Even dimensions are recommended. */
447 static void GridSfcOrdering3D(int width, int height, int depth,
448 Array<int> &coords);
449
450
451 // utility
452
453 /// Return Mesh vertex indices of an edge identified by 'edge_id'.
454 void GetEdgeVertices(const MeshId &edge_id, int vert_index[2],
455 bool oriented = true) const;
456
457 /** Return "NC" orientation of an edge. As opposed to standard Mesh edge
458 orientation based on vertex IDs, "NC" edge orientation follows the local
459 edge orientation within the element 'edge_id.element' and is thus
460 processor independent. TODO: this seems only partially true? */
461 int GetEdgeNCOrientation(const MeshId &edge_id) const;
462
463 /** Return Mesh vertex and edge indices of a face identified by 'face_id'.
464 The return value is the number of face vertices. */
465 int GetFaceVerticesEdges(const MeshId &face_id,
466 int vert_index[4], int edge_index[4],
467 int edge_orientation[4]) const;
468
469 /** Given an edge (by its vertex indices v1 and v2) return the first
470 (geometric) parent edge that exists in the Mesh or -1 if there is no such
471 parent. */
472 int GetEdgeMaster(int v1, int v2) const;
473
474 /** Get a list of vertices (2D/3D), edges (3D) and faces (3D) that coincide
475 with boundary elements with the specified attributes (marked in
476 'bdr_attr_is_ess'). In 3D this function also reveals "hidden" boundary
477 edges. In parallel it helps identifying boundary vertices/edges/faces
478 affected by non-local boundary elements. Hidden faces can occur for an
479 internal boundary coincident to a processor boundary.
480 */
481
482 /**
483 * @brief Get a list of vertices (2D/3D), edges (3D) and faces (3D) that
484 * coincide with boundary elements with the specified attributes (marked in
485 * 'bdr_attr_is_ess').
486 *
487 * @details In 3D this function also reveals "hidden" boundary edges. In
488 * parallel it helps identifying boundary vertices/edges/faces affected by
489 * non-local boundary elements. Hidden faces can occur for an internal
490 * boundary coincident to a processor boundary.
491 *
492 * @param bdr_attr_is_ess Indicator if a given attribute is essential.
493 * @param bdr_vertices Array of vertices that are essential.
494 * @param bdr_edges Array of edges that are essential.
495 * @param bdr_faces Array of faces that are essential.
496 */
497 virtual void GetBoundaryClosure(const Array<int> &bdr_attr_is_ess,
498 Array<int> &bdr_vertices,
499 Array<int> &bdr_edges, Array<int> &bdr_faces);
500
501 /// Return element geometry type. @a index is the Mesh element number.
504
505 /// Return face geometry type. @a index is the Mesh face number.
508
509 /// Return the number of root elements.
511
512 /// Return the distance of leaf @a i from the root.
513 int GetElementDepth(int i) const;
514
515 /** Return the size reduction compared to the root element (ignoring local
516 stretching and curvature). */
517 int GetElementSizeReduction(int i) const;
518
519 /// Return the faces and face attributes of leaf element @a i.
521 Array<int> &fattr) const;
522
523 /// Set the attribute of leaf element @a i, which is a Mesh element index.
524 void SetAttribute(int i, int attr)
525 { elements[leaf_elements[i]].attribute = attr; }
526
527 /** I/O: Print the mesh in "MFEM NC mesh v1.0" format. If @a comments is
528 non-empty, it will be printed after the first line of the file, and each
529 line should begin with '#'. */
530 void Print(std::ostream &out, const std::string &comments = "",
531 bool nurbs=false) const;
532
533 /// I/O: Return true if the mesh was loaded from the legacy v1.1 format.
534 bool IsLegacyLoaded() const { return Legacy; }
535
536 /// I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
537 void LegacyToNewVertexOrdering(Array<int> &order) const;
538
539 /// Save memory by releasing all non-essential and cached data.
540 virtual void Trim();
541
542 /// Return total number of bytes allocated.
543 long MemoryUsage() const;
544
545 int PrintMemoryDetail() const;
546
547 /// Return true for ParNCMesh with more than one MPI process.
548 virtual bool IsParallel() const { return false; }
549
550 using RefCoord = std::int64_t;
551
552 static constexpr int MaxElemNodes =
553 8; ///< Number of nodes an element can have
554 static constexpr int MaxElemEdges =
555 12; ///< Number of edges an element can have
556 static constexpr int MaxElemFaces =
557 6; ///< Number of faces an element can have
558 static constexpr int MaxElemChildren =
559 10; ///< Number of children an element can have
560 static constexpr int MaxFaceNodes =
561 4; ///< Number of faces an element can have
562
563 /**
564 * @brief Given a node index, return the vertex index associated
565 *
566 * @param node
567 * @return int
568 */
569 int GetNodeVertex(int node) { return nodes[node].vert_index; }
570
571protected: // non-public interface for the Mesh class
572
573 friend class Mesh;
574
575 /// Fill Mesh::{vertices,elements,boundary} for the current finest level.
576 void GetMeshComponents(Mesh &mesh) const;
577
578 /** Get edge and face numbering from 'mesh' (i.e., set all Edge::index and
579 Face::index) after a new mesh was created from us. */
580 void OnMeshUpdated(Mesh *mesh);
581
582 /** Delete top-level vertex coordinates if the Mesh became curved, e.g., by
583 calling Mesh::SetCurvature or otherwise setting the Nodes. */
585
586protected: // implementation
587
588 int Dim, spaceDim; ///< dimensions of the elements and the vertex coordinates
589 int MyRank; ///< used in parallel, or when loading a parallel file in serial
590 bool Iso; ///< true if the mesh only contains isotropic refinements
591 int Geoms; ///< bit mask of element geometries present, see InitGeomFlags()
592 bool Legacy; ///< true if the mesh was loaded from the legacy v1.1 format
593
594
595 /** A Node can hold a vertex, an edge, or both. Elements directly point to
596 their corner nodes, but edge nodes also exist and can be accessed using a
597 hash-table given their two end-point node IDs. All nodes can be accessed
598 in this way, with the exception of top-level vertex nodes. When an
599 element is being refined, the mid-edge nodes are readily available with
600 this mechanism. The new elements "sign in" to the nodes by increasing the
601 reference counts of their vertices and edges. The parent element "signs
602 off" its nodes by decrementing the ref counts. */
603 struct Node : public Hashed2
604 {
607
609 scale(0.5), scaleSet(false) {}
610 ~Node();
611
612 bool HasVertex() const { return vert_refc > 0; }
613 bool HasEdge() const { return edge_refc > 0; }
614
615 // decrease vertex/edge ref count, return false if Node should be deleted
616 bool UnrefVertex() { --vert_refc; return vert_refc || edge_refc; }
617 bool UnrefEdge() { --edge_refc; return vert_refc || edge_refc; }
618
619 real_t GetScale() const { return scale; }
620 void SetScale(real_t s, bool overwrite = false);
621
622 private:
623 real_t scale; ///< Scale from struct Refinement, default 0.5
624 bool scaleSet; ///< Indicates whether scale is set and cannot be changed
625#ifdef MFEM_USE_DOUBLE
626 static constexpr real_t scaleTol = 1.0e-8; ///< Scale comparison tolerance
627#else
628 static constexpr real_t scaleTol = 1.0e-5; ///< Scale comparison tolerance
629#endif
630 };
631
632 /** Similarly to nodes, faces can be accessed by hashing their four vertex
633 node IDs. A face knows about the one or two elements that are using it. A
634 face that is not on the boundary and only has one element referencing it
635 is either a master or a slave face. */
636 struct Face : public Hashed4
637 {
638 int attribute; ///< boundary element attribute, -1 if internal face
639 int index; ///< face number in the Mesh
640 int elem[2]; ///< up to 2 elements sharing the face
641
642 Face() : attribute(-1), index(-1) { elem[0] = elem[1] = -1; }
643
644 bool Boundary() const { return attribute >= 0; }
645 bool Unused() const { return elem[0] < 0 && elem[1] < 0; }
646
647 // add or remove an element from the 'elem[2]' array
648 void RegisterElement(int e);
649 void ForgetElement(int e);
650
651 /// Return one of elem[0] or elem[1] and make sure the other is -1.
652 int GetSingleElement() const;
653 int GetAttribute() const { return attribute; }
654 };
655
656 /** This is an element in the refinement hierarchy. Each element has either
657 been refined and points to its children, or is a leaf and points to its
658 vertex nodes. */
659 struct Element
660 {
661 char geom; ///< Geometry::Type of the element (char for storage only)
662 char ref_type; ///< bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
663 char tet_type; ///< tetrahedron split type, currently always 0
664 char flag; ///< generic flag/marker, can be used by algorithms
665 int index; ///< element number in the Mesh, -1 if refined
666 int rank; ///< processor number (ParNCMesh), -1 if undefined/unknown
668 union
669 {
670 int node[MaxElemNodes]; ///< element corners (if ref_type == 0)
671 int child[MaxElemChildren]; ///< 2-10 children (if ref_type != 0)
672 };
673 int parent; ///< parent element, -1 if this is a root element, -2 if free'd
674 Element(Geometry::Type geom, int attr);
675
677 bool IsLeaf() const { return !ref_type && (parent != -2); }
678 int GetAttribute() const { return attribute; }
679 };
680
681
682 // primary data
683 HashTable<Node> nodes; // associative container holding all Nodes
684 HashTable<Face> faces; // associative container holding all Faces
685
686 bool using_scaling = false; // Whether Node::scale is being used
687
688 BlockArray<Element> elements; // storage for all Elements
689 Array<int> free_element_ids; // unused element ids - indices into 'elements'
690public:
691 /**
692 * @brief The number of Nodes.
693 *
694 * @return int
695 */
696 int GetNumNodes() const { return nodes.Size(); }
697 /**
698 * @brief Access a Node
699 *
700 * @param i Index of the node
701 * @return const Node&
702 */
703 const Node& GetNode(int i) const {return nodes[i]; }
704 /**
705 * @brief The number of faces
706 *
707 * @return int
708 */
709 int GetNumFaces() const { return faces.Size(); }
710 /**
711 * @brief Access a Face
712 *
713 * @param i Index of the face
714 * @return const Face&
715 */
716 const Face& GetFace(int i) const {return faces[i]; }
717 /**
718 * @brief The number of elements
719 *
720 * @return int
721 */
722 int GetNumElements() const { return elements.Size(); }
723 /**
724 * @brief Access an Element
725 *
726 * @param i Index of the element
727 * @return const Element&
728 */
729 const Element& GetElement(int i) const { return elements[i]; }
730
731 /**
732 * @brief Given a set of nodes defining a face, traverse the nodes structure
733 * to find the nodes that make up the parent face and replace the input nodes
734 * with the parent nodes. Additionally return the child index that the child
735 * face would be, relative to the discovered parent face.
736 * @details This method is concerned with the construction of an NCMesh
737 * structure for a d-1 manifold of an existing NCMesh. It forms a key element
738 * in a leaf -> root traversal of the parent ncmesh elements structure.
739 *
740 * @param[out] nodes The collection of nodes whose parent we are searching
741 * for
742 * @return int The child index corresponding to placing the face for the
743 * original nodes within the face defined by the returned parent nodes. If
744 * child index is -1, then the face is made up of root nodes, and nodes is
745 * unchanged.
746 */
747 int ParentFaceNodes(std::array<int, 4> &nodes) const;
748
749 /**
750 * @brief Method for finding the nodes associated to a @a face
751 * @return Nodes making up the face
752 */
753 std::array<int, 4> FindFaceNodes(int face) const;
754 std::array<int, 4> FindFaceNodes(const Face &fa) const;
755 /**
756 * @brief Backwards compatible method for finding the @a node associated to a
757 * @a face
758 */
759 MFEM_DEPRECATED void FindFaceNodes(int face, int node[4]) const;
760protected:
761
762 /** Initial traversal state (~ element orientation) for each root element
763 NOTE: M = root_state.Size() is the number of root elements. NOTE: the
764 first M items of 'elements' is the coarse mesh. */
766
767 /** Coordinates of top-level vertices (organized as triples). If empty, the
768 Mesh is curved (Nodes != NULL) and NCMesh is topology-only. */
770
771 // secondary data
772 /** Apart from the primary data structure, which is the element/node/face
773 hierarchy, there is secondary data that is derived from the primary data
774 and needs to be updated when the primary data changes. Update() takes
775 care of that and needs to be called after each refinement and
776 derefinement. */
777 virtual void Update();
778
779 // set by UpdateLeafElements, UpdateVertices and OnMeshUpdated
781
782 // NOTE: the serial code understands the bare minimum about ghost elements
783 // and other ghost entities in order to be able to load parallel partial
784 // meshes
786
787 Array<int> leaf_elements; ///< finest elements, in Mesh ordering (+ ghosts)
788 Array<int> leaf_sfc_index; ///< natural tree ordering of leaf elements
789 Array<int> vertex_nodeId; ///< vertex-index to node-id map, see UpdateVertices
790
791 NCList face_list; ///< lazy-initialized list of faces, see GetFaceList
792 NCList edge_list; ///< lazy-initialized list of edges, see GetEdgeList
793 NCList vertex_list; ///< lazy-initialized list of vertices, see GetVertexList
794
795 Array<int> boundary_faces; ///< subset of all faces, set by BuildFaceList
796 Array<char> face_geom; ///< face geometry by face index, set by OnMeshUpdated
797
798 Table element_vertex; ///< leaf-element to vertex table, see FindSetNeighbors
799
800 /// Update the leaf elements indices in leaf_elements
801 void UpdateLeafElements();
802
803 /** @brief This method assigns indices to vertices (Node::vert_index) that
804 will be seen by the Mesh class and the rest of MFEM.
805
806 We must be careful to:
807 1. Stay compatible with the conforming code, which expects top-level
808 (original) vertices to be indexed first, otherwise GridFunctions
809 defined on a conforming mesh would no longer be valid when the mesh is
810 converted to an NC mesh.
811
812 2. Make sure serial NCMesh is compatible with the parallel ParNCMesh, so
813 it is possible to read parallel partial solutions in serial code
814 (e.g., serial GLVis). This means handling ghost elements, if present.
815
816 3. Assign vertices in a globally consistent order for parallel meshes: if
817 two vertices i,j are shared by two ranks r1,r2, and i<j on r1, then
818 i<j on r2 as well. This is true for top-level vertices but also for
819 the remaining shared vertices thanks to the globally consistent SFC
820 ordering of the leaf elements. This property reduces communication and
821 simplifies ParNCMesh. */
822 void UpdateVertices(); ///< update Vertex::index and vertex_nodeId
823
824 /** Collect the leaf elements in leaf_elements, and the ghost elements in
825 ghosts. Compute and set the element indices of @a elements. On quad and
826 hex refined elements tries to order leaf elements along a space-filling
827 curve according to the given @a state variable. */
828 void CollectLeafElements(int elem, int state, Array<int> &ghosts,
829 int &counter);
830
831 /** Try to find a space-filling curve friendly orientation of the root
832 elements: set 'root_state' based on the ordering of coarse elements. Note
833 that the coarse mesh itself must be ordered as an SFC by e.g.
834 Mesh::GetGeckoElementOrdering. */
835 void InitRootState(int root_count);
836
837 /** Compute the Geometry::Type present in the root elements (coarse elements)
838 and set @a Geoms bitmask accordingly. */
839 void InitGeomFlags();
840
841 /// Return true if the mesh contains prism elements.
842 bool HavePrisms() const { return Geoms & (1 << Geometry::PRISM); }
843
844 /// Return true if the mesh contains pyramid elements.
845 bool HavePyramids() const { return Geoms & (1 << Geometry::PYRAMID); }
846
847 /// Return true if the mesh contains tetrahedral elements.
848 bool HaveTets() const { return Geoms & (1 << Geometry::TETRAHEDRON); }
849
850 /// Return true if the Element @a el is a ghost element.
851 bool IsGhost(const Element &el) const { return el.rank != MyRank; }
852
853 // refinement/derefinement
854
855 Array<Refinement> ref_stack; ///< stack of scheduled refinements (temporary)
856 HashTable<Node> shadow; ///< temporary storage for reparented nodes
857 Array<Triple<int, int, int> > reparents; ///< scheduled node reparents (tmp)
858 Array<real_t> reparent_scale; ///< scale associated with reparents (tmp)
859
860 Table derefinements; ///< possible derefinements, see GetDerefinementTable
861
862
863 /// Refine the element @a elem with the refinement type @a ref_type.
864 void RefineElement(int elem, char ref_type);
865
866 /// Refine one element with type and scale specified by @a ref.
867 void RefineElement(const Refinement & ref);
868
869 /// Derefine the element @a elem, does nothing on leaf elements.
870 void DerefineElement(int elem);
871
872 /// Helper function to set scale for a node with parents @a p0, @a p1.
873 void SetNodeScale(int p0, int p1, real_t scale);
874
875 /// Add an Element @a el to the NCMesh, optimized to reuse freed elements.
876 int AddElement(const Element &el)
877 {
879 {
880 int idx = free_element_ids.Last();
882 elements[idx] = el;
883 return idx;
884 }
885 return elements.Append(el);
886 }
887 int AddElement(Geometry::Type geom, int attr) { return AddElement(Element(geom,attr)); }
888
889 // Free the element with index @a id.
890 void FreeElement(int id)
891 {
893 elements[id].ref_type = 0;
894 elements[id].parent = -2; // mark the element as free
895 }
896
897 int NewHexahedron(int n0, int n1, int n2, int n3,
898 int n4, int n5, int n6, int n7, int attr,
899 int fattr0, int fattr1, int fattr2,
900 int fattr3, int fattr4, int fattr5);
901
902 int NewWedge(int n0, int n1, int n2,
903 int n3, int n4, int n5, int attr,
904 int fattr0, int fattr1,
905 int fattr2, int fattr3, int fattr4);
906
907 int NewTetrahedron(int n0, int n1, int n2, int n3, int attr,
908 int fattr0, int fattr1, int fattr2, int fattr3);
909
910 int NewPyramid(int n0, int n1, int n2, int n3, int n4, int attr,
911 int fattr0, int fattr1, int fattr2, int fattr3,
912 int fattr4);
913
914 int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr,
915 int eattr0, int eattr1, int eattr2, int eattr3);
916
917 int NewTriangle(int n0, int n1, int n2,
918 int attr, int eattr0, int eattr1, int eattr2);
919
920 int NewSegment(int n0, int n1, int attr, int vattr1, int vattr2);
921
922 mfem::Element* NewMeshElement(int geom) const;
923
924 /**
925 * @brief Given a quad face defined by four vertices, establish which edges
926 * of this face have been split, and if so optionally return the mid points
927 * of those edges.
928 *
929 * @param n1 The first node defining the face
930 * @param n2 The second node defining the face
931 * @param n3 The third node defining the face
932 * @param n4 The fourth node defining the face
933 * @param s returns the scale of the split
934 * @param mid optional return of the edge mid points.
935 * @return int 0 -- no split, 1 -- "vertical" split, 2 -- "horizontal" split
936 */
937 int QuadFaceSplitType(int n1, int n2, int n3, int n4, real_t & s,
938 int mid[5] = NULL /*optional output of mid-edge nodes*/) const;
939
940 /**
941 * @brief Given a tri face defined by three vertices, establish whether the
942 * edges that make up this face have been split, and if so optionally return
943 * the midpoints.
944 * @details This is a necessary condition for this face to have been split,
945 * but is not sufficient. Consider a triangle attached to three refined
946 * triangles, in this scenario all edges can be split but this face not be
947 * split. In this case, it is necessary to check if there is a face made up
948 * of the returned midpoint nodes.
949 *
950 * @param n1 The first node defining the face
951 * @param n2 The second node defining the face
952 * @param n3 The third node defining the face
953 * @param mid optional return of the edge mid points.
954 * @return true Splits for all edges have been found
955 * @return false
956 */
957 bool TriFaceSplit(int n1, int n2, int n3, int mid[3] = NULL) const;
958
959 /**
960 * @brief Determine if a Triangle face is a master face
961 * @details This check requires looking for the edges making up the triangle
962 * being split, if nodes exist at their midpoints, and there are vertices at
963 * them, this implies the face COULD be split. To determine if it is, we then
964 * check whether these midpoints have all been connected, this is required to
965 * discriminate between an internal master face surrounded by nonconformal
966 * refinements and a conformal boundary face surrounded by refinements.
967 *
968 * @param n1 The first node defining the face
969 * @param n2 The second node defining the face
970 * @param n3 The third node defining the face
971 * @return true The face is a master
972 * @return false The face is not a master
973 */
974 inline bool TriFaceIsMaster(int n1, int n2, int n3) const
975 {
976 int mid[3];
977 return !(!TriFaceSplit(n1, n2, n3, mid) // The edges aren't split
978 // OR none of the midpoints are connected.
979 || (nodes.FindId(mid[0], mid[1]) < 0 &&
980 nodes.FindId(mid[0], mid[2]) < 0 &&
981 nodes.FindId(mid[1], mid[2]) < 0));
982 }
983
984 /**
985 * @brief Determine if a Quad face is a master face
986 *
987 * @param n1 The first node defining the face
988 * @param n2 The second node defining the face
989 * @param n3 The third node defining the face
990 * @param n4 The fourth node defining the face
991 * @return true The quad face is a master face
992 * @return false The quad face is not a master face
993 */
994 inline bool QuadFaceIsMaster(int n1, int n2, int n3, int n4) const
995 {
996 real_t s;
997 return QuadFaceSplitType(n1, n2, n3, n4, s) != 0;
998 }
999
1000 void ForceRefinement(int vn1, int vn2, int vn3, int vn4);
1001
1002 void FindEdgeElements(int vn1, int vn2, int vn3, int vn4,
1003 Array<MeshId> &prisms) const;
1004
1005 void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4,
1006 const Refinement *refs, int nref);
1007
1008 void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4,
1009 int mid12, int mid34, int level = 0);
1010
1011 void CheckIsoFace(int vn1, int vn2, int vn3, int vn4,
1012 int en1, int en2, int en3, int en4, int midf);
1013
1014 void ReparentNode(int node, int new_p1, int new_p2, real_t scale);
1015
1016 int FindMidEdgeNode(int node1, int node2) const;
1017 int GetMidEdgeNode(int node1, int node2);
1018
1019 int GetMidFaceNode(int en1, int en2, int en3, int en4);
1020
1021 /**
1022 * @brief Add references to all nodes, edges and faces of the element
1023 *
1024 * @param elem index into elements
1025 */
1026 void ReferenceElement(int elem);
1027 void UnreferenceElement(int elem, Array<int> &elemFaces);
1028
1029 Face* GetFace(Element &elem, int face_no);
1030 void RegisterFaces(int elem, int *fattr = NULL);
1031 void DeleteUnusedFaces(const Array<int> &elemFaces);
1032
1033 void CollectDerefinements(int elem, Array<Connection> &list);
1034
1035 /// Return el.node[index] correctly, even if the element is refined.
1036 int RetrieveNode(const Element &el, int index);
1037
1038 /// Extended version of find_node: works if 'el' is refined.
1039 int FindNodeExt(const Element &el, int node, bool abort = true);
1040
1041
1042 // face/edge lists
1043
1044 static int find_node(const Element &el, int node);
1045 static int find_element_edge(const Element &el, int vn0, int vn1,
1046 bool abort = true);
1047 static int find_local_face(int geom, int a, int b, int c);
1048
1049 struct Point;
1050 struct PointMatrix;
1051
1052 int ReorderFacePointMat(int v0, int v1, int v2, int v3,
1053 int elem, const PointMatrix &pm,
1054 PointMatrix &reordered) const;
1055
1056 void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3,
1057 const PointMatrix& pm, int level, Face* eface[4],
1058 MatrixMap &matrix_map);
1060 {
1061 bool unsplit; ///< Whether this face has no further splits.
1062 bool ghost_neighbor; ///< Whether the face neighbor is a ghost.
1063 };
1064 TriFaceTraverseResults TraverseTriFace(int vn0, int vn1, int vn2,
1065 const PointMatrix& pm, int level,
1066 MatrixMap &matrix_map);
1067 void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1,
1068 MatrixMap &matrix_map);
1069 void TraverseEdge(int vn0, int vn1, real_t t0, real_t t1, int flags,
1070 int level, MatrixMap &matrix_map);
1071
1072 virtual void BuildFaceList();
1073 virtual void BuildEdgeList();
1074 virtual void BuildVertexList();
1075
1076 virtual void ElementSharesFace(int elem, int local, int face) {} // ParNCMesh
1077 virtual void ElementSharesEdge(int elem, int local, int enode) {} // ParNCMesh
1078 virtual void ElementSharesVertex(int elem, int local, int vnode) {} // ParNCMesh
1079
1080 // neighbors / element_vertex table
1081
1082 /** Return all vertex-, edge- and face-neighbors of a set of elements. The
1083 neighbors are returned as a list (neighbors != NULL), as a set
1084 (neighbor_set != NULL), or both. The sizes of the set arrays must match
1085 that of leaf_elements. The function is intended to be used for large sets
1086 of elements and its complexity is linear in the number of leaf elements
1087 in the mesh. */
1088 void FindSetNeighbors(const Array<char> &elem_set,
1089 Array<int> *neighbors, /* append */
1090 Array<char> *neighbor_set = NULL);
1091
1092 /** Return all vertex-, edge- and face-neighbors of a single element. You can
1093 limit the number of elements being checked using 'search_set'. The
1094 complexity of the function is linear in the size of the search set.*/
1095 void FindNeighbors(int elem,
1096 Array<int> &neighbors, /* append */
1097 const Array<int> *search_set = NULL);
1098
1099 /** Expand a set of elements by all vertex-, edge- and face-neighbors. The
1100 output array 'expanded' will contain all items from 'elems' (provided
1101 they are in 'search_set') plus their neighbors. The neighbor search can
1102 be limited to the optional search set. The complexity is linear in the
1103 sum of the sizes of 'elems' and 'search_set'. */
1104 void NeighborExpand(const Array<int> &elems,
1105 Array<int> &expanded,
1106 const Array<int> *search_set = NULL);
1107
1108
1109 void CollectEdgeVertices(int v0, int v1, Array<int> &indices);
1110 void CollectTriFaceVertices(int v0, int v1, int v2, Array<int> &indices);
1111 void CollectQuadFaceVertices(int v0, int v1, int v2, int v3,
1112 Array<int> &indices);
1114
1119
1120 int GetVertexRootCoord(int elem, RefCoord coord[3]) const;
1121 void CollectIncidentElements(int elem, const RefCoord coord[3],
1122 Array<int> &list) const;
1123
1124 /** Return elements neighboring to a local vertex of element 'elem'. Only
1125 elements from within the same refinement tree ('cousins') are returned.
1126 Complexity is proportional to the depth of elem's refinement tree. */
1127 void FindVertexCousins(int elem, int local, Array<int> &cousins) const;
1128
1129
1130 // coarse/fine transformations
1131
1132 struct Point
1133 {
1134 int dim;
1136
1137 Point() { dim = 0; }
1138
1139 Point(const Point &) = default;
1140
1142 { dim = 1; coord[0] = x; }
1143
1145 { dim = 2; coord[0] = x; coord[1] = y; }
1146
1148 { dim = 3; coord[0] = x; coord[1] = y; coord[2] = z; }
1149
1150 Point(const Point& p0, const Point& p1, real_t s = 0.5)
1151 {
1152 dim = p0.dim;
1153 for (int i = 0; i < dim; i++)
1154 {
1155 coord[i] = ((1.0 - s) * p0.coord[i]) + (s * p1.coord[i]);
1156 }
1157 }
1158
1159 Point(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
1160 {
1161 dim = p0.dim;
1162 MFEM_ASSERT(p1.dim == dim && p2.dim == dim && p3.dim == dim, "");
1163 for (int i = 0; i < dim; i++)
1164 {
1165 coord[i] = (p0.coord[i] + p1.coord[i] + p2.coord[i] + p3.coord[i])
1166 * 0.25;
1167 }
1168 }
1169
1170 Point& operator=(const Point& src)
1171 {
1172 dim = src.dim;
1173 for (int i = 0; i < dim; i++) { coord[i] = src.coord[i]; }
1174 return *this;
1175 }
1176 };
1177
1178 /** @brief The PointMatrix stores the coordinates of the slave face using the
1179 master face coordinate as reference.
1180
1181 In 2D, the point matrix has the orientation of the parent edge, so its
1182 columns need to be flipped when applying it, see
1183 ApplyLocalSlaveTransformation.
1184
1185 In 3D, the orientation part of Elem2Inf is encoded in the point matrix.
1186
1187 The following transformation gives the relation between the reference
1188 quad face coordinates (xi, eta) in [0,1]^2, and the fine quad face
1189 coordinates (x, y):
1190 x = a0*(1-xi)*(1-eta) + a1*xi*(1-eta) + a2*xi*eta + a3*(1-xi)*eta
1191 y = b0*(1-xi)*(1-eta) + b1*xi*(1-eta) + b2*xi*eta + b3*(1-xi)*eta
1192 */
1194 {
1195 int np;
1197
1198 PointMatrix() : np(0) {}
1199
1200 PointMatrix(const Point& p0, const Point& p1)
1201 { np = 2; points[0] = p0; points[1] = p1; }
1202
1203 PointMatrix(const Point& p0, const Point& p1, const Point& p2)
1204 { np = 3; points[0] = p0; points[1] = p1; points[2] = p2; }
1205
1206 PointMatrix(const Point& p0, const Point& p1, const Point& p2, const Point& p3)
1207 { np = 4; points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3; }
1208
1209 PointMatrix(const Point& p0, const Point& p1, const Point& p2,
1210 const Point& p3, const Point& p4)
1211 {
1212 np = 5;
1213 points[0] = p0; points[1] = p1; points[2] = p2;
1214 points[3] = p3; points[4] = p4;
1215 }
1216 PointMatrix(const Point& p0, const Point& p1, const Point& p2,
1217 const Point& p3, const Point& p4, const Point& p5)
1218 {
1219 np = 6;
1220 points[0] = p0; points[1] = p1; points[2] = p2;
1221 points[3] = p3; points[4] = p4; points[5] = p5;
1222 }
1223 PointMatrix(const Point& p0, const Point& p1, const Point& p2,
1224 const Point& p3, const Point& p4, const Point& p5,
1225 const Point& p6, const Point& p7)
1226 {
1227 np = 8;
1228 points[0] = p0; points[1] = p1; points[2] = p2; points[3] = p3;
1229 points[4] = p4; points[5] = p5; points[6] = p6; points[7] = p7;
1230 }
1231
1232 Point& operator()(int i) { return points[i]; }
1233 const Point& operator()(int i) const { return points[i]; }
1234
1235 bool operator==(const PointMatrix &pm) const;
1236
1237 void GetMatrix(DenseMatrix& point_matrix) const;
1238 };
1239
1247
1248 static const PointMatrix& GetGeomIdentity(Geometry::Type geom);
1249
1250 void GetPointMatrix(Geometry::Type geom, const char* ref_path,
1251 DenseMatrix& matrix) const;
1252
1253 using RefPathMap = std::map<std::string, int>;
1254
1255 void TraverseRefinements(int elem, int coarse_index,
1256 std::string &ref_path, RefPathMap &map) const;
1257
1258 /// storage for data returned by Get[De]RefinementTransforms()
1260
1261 /// state of leaf_elements before Refine(), set by MarkCoarseLevel()
1263
1264 void InitDerefTransforms();
1265 void SetDerefMatrixCodes(int parent, Array<int> &fine_coarse);
1266
1267 // vertex temporary data, used by GetMeshComponents
1269 {
1272 TmpVertex() : valid(false), visited(false) {}
1273 };
1274
1276
1277 const real_t *CalcVertexPos(int node) const;
1278
1279
1280 // utility
1281
1282 int GetEdgeMaster(int node) const;
1283
1284 /// Return directed scale in (0,1).
1285 inline real_t GetScale(real_t s, bool reverse) const
1286 { return reverse ? 1.0 - s : s; }
1287
1288 /**
1289 * @brief Return the number of splits of this edge that have occurred in the
1290 * NCMesh. If zero, this means the segment is not the master of any other
1291 * segments.
1292 *
1293 * @param vn1 The first vertex making up the segment
1294 * @param vn2 The second vertex making up the segment
1295 * @return int The depth of splits of this segment that are present in the
1296 * mesh.
1297 */
1298 int EdgeSplitLevel(int vn1, int vn2) const;
1299 /**
1300 * @brief Return the number of splits of this triangle that have occurred in
1301 * the NCMesh. If zero, this means the triangle is neither split, nor the
1302 * master of a split face.
1303 *
1304 * @param vn1 The first vertex making up the triangle
1305 * @param vn2 The second vertex making up the triangle
1306 * @param vn3 The third vertex making up the triangle
1307 * @return int The depth of splits of this triangle that are present in the
1308 * mesh.
1309 */
1310 int TriFaceSplitLevel(int vn1, int vn2, int vn3) const;
1311 /**
1312 * @brief Computes the number of horizontal and vertical splits of this quad
1313 * that have occurred in the NCMesh. If zero, this means the quad is not the
1314 * master of any other quad.
1315 *
1316 * @param vn1 The first vertex making up the quad
1317 * @param vn2 The second vertex making up the quad
1318 * @param vn3 The third vertex making up the quad
1319 * @param vn4 The fourth vertex making up the quad
1320 * @param h_level The number of "horizontal" splits of the quad
1321 * @param v_level The number of "vertical" splits of the quad
1322 */
1323 void QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4,
1324 int& h_level, int& v_level) const;
1325 /**
1326 * @brief Returns the total number of splits of this quad that have occurred
1327 * in the NCMesh. If zero, this means the quad is not the master of any other
1328 * quad.
1329 * @details This is a convenience wrapper that sums the horizontal and
1330 * vertical levels from the full method.
1331 *
1332 * @param vn1 The first vertex making up the quad
1333 * @param vn2 The second vertex making up the quad
1334 * @param vn3 The third vertex making up the quad
1335 * @param vn4 The fourth vertex making up the quad
1336 * @return int The depth of splits of this triangle that are present in the
1337 * mesh. NB: An isotropic refinement has a level of 2, one horizontal split,
1338 * followed by a vertical split.
1339 */
1340 int QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4) const;
1341
1342 void CountSplits(int elem, int splits[3]) const;
1343 void GetLimitRefinements(Array<Refinement> &refinements, int max_level);
1344
1345 // Checker helpers
1346
1348 {
1349 MFEM_VERIFY(geom == Geometry::SEGMENT ||
1350 geom == Geometry::TRIANGLE || geom == Geometry::SQUARE ||
1351 geom == Geometry::CUBE || geom == Geometry::PRISM ||
1352 geom == Geometry::PYRAMID || geom == Geometry::TETRAHEDRON,
1353 "Element type " << geom << " is not supported by NCMesh.");
1354 }
1355
1356
1357 // I/O
1358
1359 /// Print the "vertex_parents" section of the mesh file.
1360 int PrintVertexParents(std::ostream *out) const;
1361 /// Load the vertex parent hierarchy from a mesh file.
1362 void LoadVertexParents(std::istream &input);
1363
1364 /// Load VertexToKnotSpan data for the NC patch topology mesh of a 2D or 3D
1365 /// MFEM NURBS NC-patch mesh.
1366 void LoadVertexToKnotSpan(std::istream &input);
1367 void LoadVertexToKnotSpan2D(std::istream &input);
1368 void LoadVertexToKnotSpan3D(std::istream &input);
1369
1370 /** Print the "boundary" section of the mesh file. If out == NULL, only
1371 return the number of boundary elements. */
1372 int PrintBoundary(std::ostream *out) const;
1373 /// Load the "boundary" section of the mesh file.
1374 void LoadBoundary(std::istream &input);
1375
1376 /// Print the "coordinates" section of the mesh file.
1377 void PrintCoordinates(std::ostream &out) const;
1378 /// Load the "coordinates" section of the mesh file.
1379 void LoadCoordinates(std::istream &input);
1380
1381 /// Count root elements and initialize root_state.
1382 void InitRootElements();
1383 /// Return the index of the last top-level node plus one.
1384 int CountTopLevelNodes() const;
1385 /// Return true if all root_states are zero.
1386 bool ZeroRootStates() const;
1387
1388 /// Load the element refinement hierarchy from a legacy mesh file.
1389 void LoadCoarseElements(std::istream &input);
1390 void CopyElements(int elem, const BlockArray<Element> &tmp_elements);
1391 /// Load the deprecated MFEM mesh v1.1 format for backward compatibility.
1392 void LoadLegacyFormat(std::istream &input, int &curved, int &is_nc);
1393
1394 // geometry
1395
1396 /// This holds in one place the constants about the geometries we support
1398 {
1399 int nv, ne, nf; // number of: vertices, edges, faces
1400 int edges[MaxElemEdges][2]; // edge vertices (up to 12 edges)
1401 int faces[MaxElemFaces][4]; // face vertices (up to 6 faces)
1402 int nfv[MaxElemFaces]; // number of face vertices
1403
1405 GeomInfo() : initialized(false) {}
1407 void InitGeom(Geometry::Type geom);
1408 };
1409
1411
1412 /// This is used for a NURBS mesh with this NCMesh as its patch topology.
1414
1415#ifdef MFEM_DEBUG
1416public:
1417 void DebugLeafOrder(std::ostream &out) const;
1418 void DebugDump(std::ostream &out) const;
1419#endif
1420
1421 friend class ParNCMesh; // for ParNCMesh::ElementSet
1422 friend struct MatrixMap;
1423 friend struct PointMatrixHash;
1424 friend class NCSubMesh; // for faces, nodes
1425 friend class ParNCSubMesh; // for faces, nodes
1426};
1427
1428}
1429
1430#endif
Dynamic 2D array using row-major layout.
Definition array.hpp:430
int NumRows() const
Definition array.hpp:447
int Size() const
Return the logical size of the array.
Definition array.hpp:166
void DeleteAll()
Delete the whole array.
Definition array.hpp:1033
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:912
void DeleteLast()
Delete the last entry of the array.
Definition array.hpp:224
T & Last()
Return the last element in the array.
Definition array.hpp:945
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:46
Mesh data type.
Definition mesh.hpp:65
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
Definition ncmesh.hpp:190
int GetEdgeMaster(int v1, int v2) const
Definition ncmesh.cpp:5738
static constexpr int MaxElemEdges
Number of edges an element can have.
Definition ncmesh.hpp:554
static PointMatrix pm_tet_identity
Definition ncmesh.hpp:1243
void OnMeshUpdated(Mesh *mesh)
Definition ncmesh.cpp:2887
int GetNumFaces() const
The number of faces.
Definition ncmesh.hpp:709
int GetNodeVertex(int node)
Given a node index, return the vertex index associated.
Definition ncmesh.hpp:569
void GetElementFacesAttributes(int i, Array< int > &faces, Array< int > &fattr) const
Return the faces and face attributes of leaf element i.
Definition ncmesh.cpp:5773
static GeomInfo GI[Geometry::NumGeom]
Definition ncmesh.hpp:1410
virtual void Update()
Definition ncmesh.cpp:268
void FindNeighbors(int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
Definition ncmesh.cpp:4322
void FreeElement(int id)
Definition ncmesh.hpp:890
virtual void Trim()
Save memory by releasing all non-essential and cached data.
Definition ncmesh.cpp:6964
virtual void ElementSharesVertex(int elem, int local, int vnode)
Definition ncmesh.hpp:1078
HashTable< Node > shadow
temporary storage for reparented nodes
Definition ncmesh.hpp:856
void SetDerefMatrixCodes(int parent, Array< int > &fine_coarse)
Definition ncmesh.cpp:2364
void CopyElements(int elem, const BlockArray< Element > &tmp_elements)
Definition ncmesh.cpp:6692
Array< Refinement > ref_stack
stack of scheduled refinements (temporary)
Definition ncmesh.hpp:855
std::int64_t RefCoord
Definition ncmesh.hpp:550
const Face & GetFace(int i) const
Access a Face.
Definition ncmesh.hpp:716
const VertexToKnotSpan & GetVertexToKnotSpan() const
Definition ncmesh.hpp:401
mfem::Element * NewMeshElement(int geom) const
Definition ncmesh.cpp:2720
void SetNodeScale(int p0, int p1, real_t scale)
Helper function to set scale for a node with parents p0, p1.
Definition ncmesh.cpp:1117
int NewTriangle(int n0, int n1, int n2, int attr, int eattr0, int eattr1, int eattr2)
Definition ncmesh.cpp:758
void MakeTopologyOnly()
Definition ncmesh.hpp:584
void GetMeshComponents(Mesh &mesh) const
Fill Mesh::{vertices,elements,boundary} for the current finest level.
Definition ncmesh.cpp:2760
void LoadVertexToKnotSpan3D(std::istream &input)
Definition ncmesh.cpp:6194
NCMesh()=default
int NGhostElements
Definition ncmesh.hpp:785
void LoadCoarseElements(std::istream &input)
Load the element refinement hierarchy from a legacy mesh file.
Definition ncmesh.cpp:6710
void NeighborExpand(const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
Definition ncmesh.cpp:4405
HashTable< Node > nodes
Definition ncmesh.hpp:683
int PrintVertexParents(std::ostream *out) const
Print the "vertex_parents" section of the mesh file.
Definition ncmesh.cpp:6106
static PointMatrix pm_prism_identity
Definition ncmesh.hpp:1244
int NewSegment(int n0, int n1, int attr, int vattr1, int vattr2)
Definition ncmesh.cpp:784
virtual bool IsParallel() const
Return true for ParNCMesh with more than one MPI process.
Definition ncmesh.hpp:548
static int find_node(const Element &el, int node)
Definition ncmesh.cpp:3278
void DeleteUnusedFaces(const Array< int > &elemFaces)
Definition ncmesh.cpp:455
const CoarseFineTransformations & GetRefinementTransforms() const
Definition ncmesh.cpp:5204
void BuildElementToVertexTable()
Definition ncmesh.cpp:4131
virtual int GetNGhostElements() const
Definition ncmesh.hpp:224
void TraverseEdge(int vn0, int vn1, real_t t0, real_t t1, int flags, int level, MatrixMap &matrix_map)
Definition ncmesh.cpp:3772
bool using_scaling
Definition ncmesh.hpp:686
static int find_element_edge(const Element &el, int vn0, int vn1, bool abort=true)
Definition ncmesh.cpp:3298
Array< int > free_element_ids
Definition ncmesh.hpp:689
void LoadBoundary(std::istream &input)
Load the "boundary" section of the mesh file.
Definition ncmesh.cpp:6261
int Dimension() const
Return the dimension of the NCMesh.
Definition ncmesh.hpp:214
virtual void ElementSharesFace(int elem, int local, int face)
Definition ncmesh.hpp:1076
int PrintMemoryDetail() const
Definition ncmesh.cpp:7030
Array< real_t > reparent_scale
scale associated with reparents (tmp)
Definition ncmesh.hpp:858
static constexpr int MaxElemNodes
Number of nodes an element can have.
Definition ncmesh.hpp:552
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:615
HashTable< Face > faces
Definition ncmesh.hpp:684
bool HaveTets() const
Return true if the mesh contains tetrahedral elements.
Definition ncmesh.hpp:848
int AddElement(Geometry::Type geom, int attr)
Definition ncmesh.hpp:887
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:5639
static PointMatrix pm_seg_identity
Definition ncmesh.hpp:1240
static constexpr int MaxFaceNodes
Number of faces an element can have.
Definition ncmesh.hpp:560
const NCList & GetNCList(int entity)
Return vertex/edge/face list (entity = 0/1/2, respectively).
Definition ncmesh.hpp:391
bool QuadFaceIsMaster(int n1, int n2, int n3, int n4) const
Determine if a Quad face is a master face.
Definition ncmesh.hpp:994
Array< int > boundary_faces
subset of all faces, set by BuildFaceList
Definition ncmesh.hpp:795
void Print(std::ostream &out, const std::string &comments="", bool nurbs=false) const
Definition ncmesh.cpp:6349
std::array< int, 4 > FindFaceNodes(int face) const
Method for finding the nodes associated to a face.
Definition ncmesh.cpp:5799
static PointMatrix pm_quad_identity
Definition ncmesh.hpp:1242
void ReferenceElement(int elem)
Add references to all nodes, edges and faces of the element.
Definition ncmesh.cpp:367
void UpdateElementToVertexTable()
Definition ncmesh.hpp:1115
void LoadVertexToKnotSpan2D(std::istream &input)
Definition ncmesh.cpp:6174
BlockArray< Element > elements
Definition ncmesh.hpp:688
Array< int > coarse_elements
state of leaf_elements before Refine(), set by MarkCoarseLevel()
Definition ncmesh.hpp:1262
const CoarseFineTransformations & GetDerefinementTransforms() const
Definition ncmesh.cpp:5255
Table derefinements
possible derefinements, see GetDerefinementTable
Definition ncmesh.hpp:860
int FindNodeExt(const Element &el, int node, bool abort=true)
Extended version of find_node: works if 'el' is refined.
Definition ncmesh.cpp:3288
int ReorderFacePointMat(int v0, int v1, int v2, int v3, int elem, const PointMatrix &pm, PointMatrix &reordered) const
Definition ncmesh.cpp:3423
int FindMidEdgeNode(int node1, int node2) const
Definition ncmesh.cpp:336
Array< char > face_geom
face geometry by face index, set by OnMeshUpdated
Definition ncmesh.hpp:796
int GetEdgeNCOrientation(const MeshId &edge_id) const
Definition ncmesh.cpp:5658
Table element_vertex
leaf-element to vertex table, see FindSetNeighbors
Definition ncmesh.hpp:798
Array< Triple< int, int, int > > reparents
scheduled node reparents (tmp)
Definition ncmesh.hpp:857
virtual void BuildFaceList()
Definition ncmesh.cpp:3677
void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4, const Refinement *refs, int nref)
Definition ncmesh.cpp:955
void LoadVertexParents(std::istream &input)
Load the vertex parent hierarchy from a mesh file.
Definition ncmesh.cpp:6141
TmpVertex * tmp_vertex
Definition ncmesh.hpp:1275
static void GridSfcOrdering3D(int width, int height, int depth, Array< int > &coords)
Definition ncmesh.cpp:5607
Array< int > leaf_elements
finest elements, in Mesh ordering (+ ghosts)
Definition ncmesh.hpp:787
const real_t * CalcVertexPos(int node) const
Definition ncmesh.cpp:2735
void CollectLeafElements(int elem, int state, Array< int > &ghosts, int &counter)
Definition ncmesh.cpp:2383
const Element & GetElement(int i) const
Access an Element.
Definition ncmesh.hpp:729
int CountTopLevelNodes() const
Return the index of the last top-level node plus one.
Definition ncmesh.cpp:6494
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:5935
int Geoms
bit mask of element geometries present, see InitGeomFlags()
Definition ncmesh.hpp:591
void InitDerefTransforms()
Definition ncmesh.cpp:2345
int GetElementSizeReduction(int i) const
Definition ncmesh.cpp:5759
virtual void LimitNCLevel(int max_nc_level)
Definition ncmesh.cpp:6090
int NGhostVertices
Definition ncmesh.hpp:785
const NCList & GetVertexList()
Definition ncmesh.hpp:384
bool ZeroRootStates() const
Return true if all root_states are zero.
Definition ncmesh.cpp:6340
Array< int > vertex_nodeId
vertex-index to node-id map, see UpdateVertices
Definition ncmesh.hpp:789
int ParentFaceNodes(std::array< int, 4 > &nodes) const
Given a set of nodes defining a face, traverse the nodes structure to find the nodes that make up the...
Definition ncmesh.cpp:3108
int GetNumRootElements()
Return the number of root elements.
Definition ncmesh.hpp:510
int RetrieveNode(const Element &el, int index)
Return el.node[index] correctly, even if the element is refined.
Definition ncmesh.cpp:1996
void FindEdgeElements(int vn1, int vn2, int vn3, int vn4, Array< MeshId > &prisms) const
Definition ncmesh.cpp:907
virtual ~NCMesh()
Definition ncmesh.cpp:280
void CollectDerefinements(int elem, Array< Connection > &list)
Definition ncmesh.cpp:2232
bool IsLegacyLoaded() const
I/O: Return true if the mesh was loaded from the legacy v1.1 format.
Definition ncmesh.hpp:534
const Node & GetNode(int i) const
Access a Node.
Definition ncmesh.hpp:703
virtual void ElementSharesEdge(int elem, int local, int enode)
Definition ncmesh.hpp:1077
VertexToKnotSpan vertex_to_knotspan
This is used for a NURBS mesh with this NCMesh as its patch topology.
Definition ncmesh.hpp:1413
int GetFaceVerticesEdges(const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
Definition ncmesh.cpp:5670
void RefineElement(int elem, char ref_type)
Refine the element elem with the refinement type ref_type.
Definition ncmesh.cpp:1123
void InitRootState(int root_count)
Definition ncmesh.cpp:2654
Geometry::Type GetElementGeometry(int index) const
Return element geometry type. index is the Mesh element number.
Definition ncmesh.hpp:502
int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
Definition ncmesh.cpp:732
Array< int > leaf_sfc_index
natural tree ordering of leaf elements
Definition ncmesh.hpp:788
void ForceRefinement(int vn1, int vn2, int vn3, int vn4)
Definition ncmesh.cpp:824
void CollectEdgeVertices(int v0, int v1, Array< int > &indices)
Definition ncmesh.cpp:4060
void RefineVertexToKnotSpan(const std::vector< Array< int > > &kvf, const Array< KnotVector * > &kvext, std::map< std::pair< int, int >, std::array< int, 2 > > &parentToKV)
Remap knot-span indices vertex_to_knotspan after refinement.
Definition ncmesh.cpp:5110
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
Definition ncmesh.hpp:792
static PointMatrix pm_pyramid_identity
Definition ncmesh.hpp:1245
Array< int > root_state
Definition ncmesh.hpp:765
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition ncmesh.hpp:369
real_t GetScale(real_t s, bool reverse) const
Return directed scale in (0,1).
Definition ncmesh.hpp:1285
Geometry::Type GetFaceGeometry(int index) const
Return face geometry type. index is the Mesh face number.
Definition ncmesh.hpp:506
void CollectTriFaceVertices(int v0, int v1, int v2, Array< int > &indices)
Definition ncmesh.cpp:4072
void InitRootElements()
Count root elements and initialize root_state.
Definition ncmesh.cpp:6451
int GetMidEdgeNode(int node1, int node2)
Definition ncmesh.cpp:352
virtual void BuildEdgeList()
Definition ncmesh.cpp:3802
void DebugDump(std::ostream &out) const
Definition ncmesh.cpp:7083
void DebugLeafOrder(std::ostream &out) const
Definition ncmesh.cpp:7058
std::map< std::string, int > RefPathMap
Definition ncmesh.hpp:1253
int GetElementDepth(int i) const
Return the distance of leaf i from the root.
Definition ncmesh.cpp:5747
void LoadVertexToKnotSpan(std::istream &input)
Definition ncmesh.cpp:6168
void TraverseRefinements(int elem, int coarse_index, std::string &ref_path, RefPathMap &map) const
Definition ncmesh.cpp:5170
void CheckIsoFace(int vn1, int vn2, int vn3, int vn4, int en1, int en2, int en3, int en4, int midf)
Definition ncmesh.cpp:1101
void GetLimitRefinements(Array< Refinement > &refinements, int max_level)
Definition ncmesh.cpp:6060
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:702
virtual void Derefine(const Array< int > &derefs)
Definition ncmesh.cpp:2309
bool Iso
true if the mesh only contains isotropic refinements
Definition ncmesh.hpp:590
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:5827
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:3084
bool Legacy
true if the mesh was loaded from the legacy v1.1 format
Definition ncmesh.hpp:592
virtual void BuildVertexList()
Definition ncmesh.cpp:3914
static const PointMatrix & GetGeomIdentity(Geometry::Type geom)
Definition ncmesh.cpp:4608
int GetNVertices() const
Return the number of vertices in the NCMesh.
Definition ncmesh.hpp:219
static constexpr int MaxElemFaces
Number of faces an element can have.
Definition ncmesh.hpp:556
int GetNumElements() const
The number of elements.
Definition ncmesh.hpp:722
void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4, int mid12, int mid34, int level=0)
Definition ncmesh.cpp:980
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:645
int GetVertexRootCoord(int elem, RefCoord coord[3]) const
Definition ncmesh.cpp:4464
static void CheckSupportedGeom(Geometry::Type geom)
Definition ncmesh.hpp:1347
Array< real_t > coordinates
Definition ncmesh.hpp:769
void FindVertexCousins(int elem, int local, Array< int > &cousins) const
Definition ncmesh.cpp:4540
bool IsGhost(const Element &el) const
Return true if the Element el is a ghost element.
Definition ncmesh.hpp:851
static constexpr int MaxElemChildren
Number of children an element can have.
Definition ncmesh.hpp:558
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition ncmesh.hpp:376
bool TriFaceIsMaster(int n1, int n2, int n3) const
Determine if a Triangle face is a master face.
Definition ncmesh.hpp:974
int MyRank
used in parallel, or when loading a parallel file in serial
Definition ncmesh.hpp:589
int AddElement(const Element &el)
Add an Element el to the NCMesh, optimized to reuse freed elements.
Definition ncmesh.hpp:876
void CountSplits(int elem, int splits[3]) const
Definition ncmesh.cpp:5970
int spaceDim
dimensions of the elements and the vertex coordinates
Definition ncmesh.hpp:588
void UnreferenceElement(int elem, Array< int > &elemFaces)
Definition ncmesh.cpp:398
void LegacyToNewVertexOrdering(Array< int > &order) const
I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
Definition ncmesh.cpp:6939
int NGhostEdges
Definition ncmesh.hpp:785
void CollectIncidentElements(int elem, const RefCoord coord[3], Array< int > &list) const
Definition ncmesh.cpp:4517
void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1, MatrixMap &matrix_map)
Definition ncmesh.cpp:3574
int QuadFaceSplitType(int n1, int n2, int n3, int n4, real_t &s, 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:3030
void RegisterFaces(int elem, int *fattr=NULL)
Definition ncmesh.cpp:441
int PrintBoundary(std::ostream *out) const
Definition ncmesh.cpp:6222
CoarseFineTransformations transforms
storage for data returned by Get[De]RefinementTransforms()
Definition ncmesh.hpp:1259
void MarkCoarseLevel()
Definition ncmesh.cpp:5156
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:5920
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
Definition ncmesh.cpp:2280
void UpdateVertices()
This method assigns indices to vertices (Node::vert_index) that will be seen by the Mesh class and th...
Definition ncmesh.cpp:2478
int GetNFaces() const
Return the number of (2D) faces in the NCMesh.
Definition ncmesh.hpp:223
NCList vertex_list
lazy-initialized list of vertices, see GetVertexList
Definition ncmesh.hpp:793
bool HavePyramids() const
Return true if the mesh contains pyramid elements.
Definition ncmesh.hpp:845
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:677
void LoadCoordinates(std::istream &input)
Load the "coordinates" section of the mesh file.
Definition ncmesh.cpp:6319
void GetPointMatrix(Geometry::Type geom, const char *ref_path, DenseMatrix &matrix) const
Definition ncmesh.cpp:4625
void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level, Face *eface[4], MatrixMap &matrix_map)
Definition ncmesh.cpp:3454
const Table & GetDerefinementTable()
Definition ncmesh.cpp:2265
void ClearTransforms()
Free all internal data created by the above three functions.
Definition ncmesh.cpp:5335
void SetAttribute(int i, int attr)
Set the attribute of leaf element i, which is a Mesh element index.
Definition ncmesh.hpp:524
int SpaceDimension() const
Return the space dimension of the NCMesh.
Definition ncmesh.hpp:216
void FindSetNeighbors(const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
Definition ncmesh.cpp:4211
int GetNumNodes() const
The number of Nodes.
Definition ncmesh.hpp:696
NCMesh & operator=(NCMesh &)=delete
Copy assignment not supported.
void CollectQuadFaceVertices(int v0, int v1, int v2, int v3, Array< int > &indices)
Definition ncmesh.cpp:4096
void PrintCoordinates(std::ostream &out) const
Print the "coordinates" section of the mesh file.
Definition ncmesh.cpp:6301
TriFaceTraverseResults TraverseTriFace(int vn0, int vn1, int vn2, const PointMatrix &pm, int level, MatrixMap &matrix_map)
Definition ncmesh.cpp:3611
int GetMidFaceNode(int en1, int en2, int en3, int en4)
Definition ncmesh.cpp:359
bool HavePrisms() const
Return true if the mesh contains prism elements.
Definition ncmesh.hpp:842
long MemoryUsage() const
Return total number of bytes allocated.
Definition ncmesh.cpp:7007
void DerefineElement(int elem)
Derefine the element elem, does nothing on leaf elements.
Definition ncmesh.cpp:2038
NCList face_list
lazy-initialized list of faces, see GetFaceList
Definition ncmesh.hpp:791
virtual void Refine(const Array< Refinement > &refinements)
Definition ncmesh.cpp:1947
void ReparentNode(int node, int new_p1, int new_p2, real_t scale)
Definition ncmesh.cpp:319
static PointMatrix pm_tri_identity
Definition ncmesh.hpp:1241
int NGhostFaces
Definition ncmesh.hpp:785
static PointMatrix pm_hex_identity
Definition ncmesh.hpp:1246
void UpdateLeafElements()
Update the leaf elements indices in leaf_elements.
Definition ncmesh.cpp:2451
void LoadLegacyFormat(std::istream &input, int &curved, int &is_nc)
Load the deprecated MFEM mesh v1.1 format for backward compatibility.
Definition ncmesh.cpp:6781
int GetNEdges() const
Return the number of edges in the NCMesh.
Definition ncmesh.hpp:221
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:5913
static void GridSfcOrdering2D(int width, int height, Array< int > &coords)
Definition ncmesh.cpp:5592
static int find_local_face(int geom, int a, int b, int c)
Definition ncmesh.cpp:3316
Class representing a Nonconformal SubMesh. This is only used by SubMesh.
Definition ncsubmesh.hpp:28
A parallel extension of the NCMesh class.
Definition pncmesh.hpp:63
Class representing a Parallel Nonconformal SubMesh. This is only used by ParSubMesh.
Data type point element.
Definition point.hpp:23
Table stores the connectivity of elements of TYPE I to elements of TYPE II. For example,...
Definition table.hpp:43
int Size() const
Returns the number of TYPE I elements.
Definition table.hpp:103
For a NURBS mesh with nonconforming patch topology, this struct provides a map from hanging vertices ...
Definition ncmesh.hpp:122
void GetVertex3D(int index, int &v, std::array< int, 2 > &ks, std::array< int, 4 > &pv) const
Get the data for a vertex in 3D.
Definition ncnurbs.cpp:3293
void SetVertex2D(int index, int v, int ks, const std::array< int, 2 > &pv)
Set the data for a vertex in 2D.
Definition ncnurbs.cpp:3251
void SetVertex3D(int index, int v, const std::array< int, 2 > &ks, const std::array< int, 4 > &pv)
Set the data for a vertex in 3D.
Definition ncnurbs.cpp:3260
void SetSize(int dimension, int numVertices)
Set the spatial dimension and number of vertices.
Definition ncnurbs.cpp:3244
void GetVertex2D(int index, int &v, int &ks, std::array< int, 2 > &pv) const
Get the data for a vertex in 2D.
Definition ncnurbs.cpp:3284
void Print(std::ostream &os) const
Print all the data.
Definition ncnurbs.cpp:3305
void SetKnotSpans3D(int index, const std::array< int, 2 > &ks)
Set the knot-span indices for a vertex in 3D.
Definition ncnurbs.cpp:3278
std::pair< int, int > GetVertexParentPair(int index) const
Return the vertex pair representing the parent edge (2D) or face (3D).
Definition ncnurbs.cpp:3321
int Size() const
Return the number of vertices.
Definition ncmesh.hpp:157
void SetKnotSpan2D(int index, int ks)
Set the knot-span index for a vertex in 2D.
Definition ncnurbs.cpp:3273
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
Definition hooke.cpp:45
int index(int i, int j, int nx, int ny)
Definition life.cpp:236
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
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 Swap(T &a, T &b)
Swap objects of type T. The operation is performed using the most specialized swap function from the ...
Definition array.hpp:738
NCMesh::RefCoord RefCoord
float real_t
Definition config.hpp:46
Defines the coarse-fine transformations of all fine elements.
Definition ncmesh.hpp:90
MFEM_DEPRECATED void GetCoarseToFineMap(const Mesh &fine_mesh, Table &coarse_to_fine) const
Definition ncmesh.hpp:109
Array< Embedding > embeddings
Fine element positions in their parents.
Definition ncmesh.hpp:92
void MakeCoarseToFineTable(Table &coarse_to_fine, bool want_ghosts=false) const
Definition ncmesh.cpp:5313
DenseTensor point_matrices[Geometry::NumGeom]
Definition ncmesh.hpp:96
Defines the position of a fine element within a coarse element.
Definition ncmesh.hpp:69
unsigned geom
Definition ncmesh.hpp:77
Embedding()=default
unsigned ghost
For internal use: 0 if regular fine element, 1 if parallel ghost element.
Definition ncmesh.hpp:81
Embedding(int elem, Geometry::Type geom, int matrix=0, bool ghost=false)
Definition ncmesh.hpp:84
int parent
Coarse Element index in the coarse mesh.
Definition ncmesh.hpp:71
unsigned matrix
Definition ncmesh.hpp:78
int rank
processor number (ParNCMesh), -1 if undefined/unknown
Definition ncmesh.hpp:666
bool IsLeaf() const
Definition ncmesh.hpp:677
int child[MaxElemChildren]
2-10 children (if ref_type != 0)
Definition ncmesh.hpp:671
char tet_type
tetrahedron split type, currently always 0
Definition ncmesh.hpp:663
Element(Geometry::Type geom, int attr)
Definition ncmesh.cpp:602
int GetAttribute() const
Definition ncmesh.hpp:678
char flag
generic flag/marker, can be used by algorithms
Definition ncmesh.hpp:664
int node[MaxElemNodes]
element corners (if ref_type == 0)
Definition ncmesh.hpp:670
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
Definition ncmesh.hpp:662
char geom
Geometry::Type of the element (char for storage only)
Definition ncmesh.hpp:661
int index
element number in the Mesh, -1 if refined
Definition ncmesh.hpp:665
int parent
parent element, -1 if this is a root element, -2 if free'd
Definition ncmesh.hpp:673
Geometry::Type Geom() const
Definition ncmesh.hpp:676
int GetAttribute() const
Definition ncmesh.hpp:653
void ForgetElement(int e)
Definition ncmesh.cpp:473
int elem[2]
up to 2 elements sharing the face
Definition ncmesh.hpp:640
bool Unused() const
Definition ncmesh.hpp:645
void RegisterElement(int e)
Definition ncmesh.cpp:466
bool Boundary() const
Definition ncmesh.hpp:644
int index
face number in the Mesh
Definition ncmesh.hpp:639
int GetSingleElement() const
Return one of elem[0] or elem[1] and make sure the other is -1.
Definition ncmesh.cpp:488
int attribute
boundary element attribute, -1 if internal face
Definition ncmesh.hpp:638
This holds in one place the constants about the geometries we support.
Definition ncmesh.hpp:1398
int nfv[MaxElemFaces]
Definition ncmesh.hpp:1402
int faces[MaxElemFaces][4]
Definition ncmesh.hpp:1401
int edges[MaxElemEdges][2]
Definition ncmesh.hpp:1400
void InitGeom(Geometry::Type geom)
Definition ncmesh.cpp:57
GeomInfo(Geometry::Type geom)
Definition ncmesh.hpp:1406
Master(int index, int element, int local, int geom, int sb, int se)
Definition ncmesh.hpp:280
int slaves_end
slave faces
Definition ncmesh.hpp:277
Identifies a vertex/edge/face in both Mesh and NCMesh.
Definition ncmesh.hpp:260
int element
NCMesh::Element containing this vertex/edge/face.
Definition ncmesh.hpp:262
int index
Mesh number.
Definition ncmesh.hpp:261
MeshId(int index, int element, int local, int geom=-1)
Definition ncmesh.hpp:269
Geometry::Type Geom() const
Definition ncmesh.hpp:266
signed char geom
Geometry::Type (faces only) (char to save RAM)
Definition ncmesh.hpp:264
signed char local
local number within 'element'
Definition ncmesh.hpp:263
const MeshIdType type
MeshIdType corresponding to the MeshId. UNRECOGNIZED if unfound.
Definition ncmesh.hpp:326
Lists all edges/faces in the nonconforming mesh.
Definition ncmesh.hpp:301
MeshIdType GetMeshIdType(int index) const
Return a face type for a given nc index.
Definition ncmesh.cpp:4010
Array< MeshId > conforming
All MeshIds corresponding to conformal faces.
Definition ncmesh.hpp:302
long MemoryUsage() const
Definition ncmesh.cpp:6979
Array< Slave > slaves
All MeshIds corresponding to slave faces.
Definition ncmesh.hpp:304
bool CheckMeshIdType(int index, MeshIdType type) const
Given an index, check if this is a certain face type.
Definition ncmesh.cpp:4018
bool Empty() const
Whether the NCList is empty.
Definition ncmesh.hpp:340
void Clear()
Erase the contents of the conforming, master and slave arrays.
Definition ncmesh.cpp:3971
void OrientedPointMatrix(const Slave &slave, DenseMatrix &oriented_matrix) const
Definition ncmesh.cpp:3949
Array< Master > masters
All MeshIds corresponding to master faces.
Definition ncmesh.hpp:303
long TotalSize() const
The total size of the component arrays in the NCList.
Definition ncmesh.hpp:347
Array< DenseMatrix * > point_matrices[Geometry::NumGeom]
List of unique point matrices for each slave geometry.
Definition ncmesh.hpp:307
MeshIdAndType GetMeshIdAndType(int index) const
Return a mesh id and type for a given nc index.
Definition ncmesh.cpp:3990
bool HasVertex() const
Definition ncmesh.hpp:612
bool HasEdge() const
Definition ncmesh.hpp:613
real_t GetScale() const
Definition ncmesh.hpp:619
void SetScale(real_t s, bool overwrite=false)
Definition ncmesh.cpp:298
The PointMatrix stores the coordinates of the slave face using the master face coordinate as referenc...
Definition ncmesh.hpp:1194
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4)
Definition ncmesh.hpp:1209
Point & operator()(int i)
Definition ncmesh.hpp:1232
void GetMatrix(DenseMatrix &point_matrix) const
Definition ncmesh.cpp:4571
PointMatrix(const Point &p0, const Point &p1, const Point &p2)
Definition ncmesh.hpp:1203
Point points[MaxElemNodes]
Definition ncmesh.hpp:1196
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition ncmesh.hpp:1206
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4, const Point &p5)
Definition ncmesh.hpp:1216
const Point & operator()(int i) const
Definition ncmesh.hpp:1233
bool operator==(const PointMatrix &pm) const
Definition ncmesh.cpp:4557
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:1223
PointMatrix(const Point &p0, const Point &p1)
Definition ncmesh.hpp:1200
Point(const Point &p0, const Point &p1, real_t s=0.5)
Definition ncmesh.hpp:1150
Point(const Point &)=default
Point(real_t x, real_t y, real_t z)
Definition ncmesh.hpp:1147
Point(real_t x, real_t y)
Definition ncmesh.hpp:1144
Point & operator=(const Point &src)
Definition ncmesh.hpp:1170
Point(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Definition ncmesh.hpp:1159
Nonconforming edge/face within a bigger edge/face.
Definition ncmesh.hpp:287
unsigned edge_flags
orientation flags, see OrientedPointMatrix
Definition ncmesh.hpp:290
unsigned matrix
index into NCList::point_matrices[geom]
Definition ncmesh.hpp:289
int master
master number (in Mesh numbering)
Definition ncmesh.hpp:288
Slave(int index, int element, int local, int geom)
Definition ncmesh.hpp:293
bool unsplit
Whether this face has no further splits.
Definition ncmesh.hpp:1061
bool ghost_neighbor
Whether the face neighbor is a ghost.
Definition ncmesh.hpp:1062
Refinement()=default
Refinement scale in each dimension.
real_t s[3]
Definition ncmesh.hpp:45
void SetScaleForType(const real_t *scale)
Set the scale in the directions for the currently set type.
Definition ncmesh.cpp:540
std::pair< char, real_t > ScaledType
Definition ncmesh.hpp:44
int index
Mesh element number.
Definition ncmesh.hpp:42
void Set(int element, char type, real_t scale=0.5)
Set the element, type, and scale.
Definition ncmesh.cpp:589
char GetType() const
Return the type as char.
Definition ncmesh.cpp:580
void SetType(char type, real_t scale=0.5)
Set the type and scale, assuming the element is already set.
Definition ncmesh.cpp:596