MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mfem::NCMesh Class Reference

A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension of the Mesh class. More...

#include <ncmesh.hpp>

Inheritance diagram for mfem::NCMesh:
[legend]
Collaboration diagram for mfem::NCMesh:
[legend]

Classes

struct  Element
 
struct  Face
 
struct  GeomInfo
 This holds in one place the constants about the geometries we support. More...
 
struct  Master
 
struct  MeshId
 Identifies a vertex/edge/face in both Mesh and NCMesh. More...
 
struct  NCList
 Lists all edges/faces in the nonconforming mesh. More...
 
struct  Node
 
struct  Point
 
struct  PointMatrix
 The PointMatrix stores the coordinates of the slave face using the master face coordinate as reference. More...
 
struct  Slave
 Nonconforming edge/face within a bigger edge/face. More...
 
struct  TmpVertex
 
struct  TriFaceTraverseResults
 

Public Types

using RefCoord = std::int64_t
 

Public Member Functions

 NCMesh (const Mesh *mesh)
 
 NCMesh (std::istream &input, int version, int &curved, int &is_nc)
 
 NCMesh (const NCMesh &other)
 Deep copy of another instance.
 
NCMeshoperator= (NCMesh &)=delete
 Copy assignment not supported.
 
virtual ~NCMesh ()
 
int Dimension () const
 Return the dimension of the NCMesh.
 
int SpaceDimension () const
 Return the space dimension of the NCMesh.
 
int GetNVertices () const
 Return the number of vertices in the NCMesh.
 
int GetNEdges () const
 Return the number of edges in the NCMesh.
 
int GetNFaces () const
 Return the number of (2D) faces in the NCMesh.
 
virtual int GetNGhostElements () const
 
virtual void Refine (const Array< Refinement > &refinements)
 
virtual void LimitNCLevel (int max_nc_level)
 
const TableGetDerefinementTable ()
 
virtual void CheckDerefinementNCLevel (const Table &deref_table, Array< int > &level_ok, int max_nc_level)
 
virtual void Derefine (const Array< int > &derefs)
 
const NCListGetFaceList ()
 Return the current list of conforming and nonconforming faces.
 
const NCListGetEdgeList ()
 Return the current list of conforming and nonconforming edges.
 
const NCListGetVertexList ()
 
const NCListGetNCList (int entity)
 Return vertex/edge/face list (entity = 0/1/2, respectively).
 
void MarkCoarseLevel ()
 
const CoarseFineTransformationsGetRefinementTransforms () const
 
const CoarseFineTransformationsGetDerefinementTransforms () const
 
void ClearTransforms ()
 Free all internal data created by the above three functions.
 
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'.
 
int GetEdgeNCOrientation (const MeshId &edge_id) const
 
int GetFaceVerticesEdges (const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
 
int GetEdgeMaster (int v1, int v2) const
 
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 the specified attributes (marked in 'bdr_attr_is_ess').
 
Geometry::Type GetElementGeometry (int index) const
 Return element geometry type. index is the Mesh element number.
 
Geometry::Type GetFaceGeometry (int index) const
 Return face geometry type. index is the Mesh face number.
 
int GetNumRootElements ()
 Return the number of root elements.
 
int GetElementDepth (int i) const
 Return the distance of leaf i from the root.
 
int GetElementSizeReduction (int i) const
 
void GetElementFacesAttributes (int i, Array< int > &faces, Array< int > &fattr) const
 Return the faces and face attributes of leaf element i.
 
void SetAttribute (int i, int attr)
 Set the attribute of leaf element i, which is a Mesh element index.
 
void Print (std::ostream &out, const std::string &comments="") const
 
bool IsLegacyLoaded () const
 I/O: Return true if the mesh was loaded from the legacy v1.1 format.
 
void LegacyToNewVertexOrdering (Array< int > &order) const
 I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
 
virtual void Trim ()
 Save memory by releasing all non-essential and cached data.
 
long MemoryUsage () const
 Return total number of bytes allocated.
 
int PrintMemoryDetail () const
 
int GetNodeVertex (int node)
 Given a node index, return the vertex index associated.
 
int GetNumNodes () const
 The number of Nodes.
 
const NodeGetNode (int i) const
 Access a Node.
 
int GetNumFaces () const
 The number of faces.
 
const FaceGetFace (int i) const
 Access a Face.
 
int GetNumElements () const
 The number of elements.
 
const ElementGetElement (int i) const
 Access an Element.
 
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 parent face and replace the input nodes with the parent nodes. Additionally return the child index that the child face would be, relative to the discovered parent face.
 
std::array< int, 4 > FindFaceNodes (int face) const
 Method for finding the nodes associated to a face.
 
std::array< int, 4 > FindFaceNodes (const Face &fa) const
 
MFEM_DEPRECATED void FindFaceNodes (int face, int node[4]) const
 Backwards compatible method for finding the node associated to a face.
 
void DebugLeafOrder (std::ostream &out) const
 
void DebugDump (std::ostream &out) const
 

Static Public Member Functions

static void GridSfcOrdering2D (int width, int height, Array< int > &coords)
 
static void GridSfcOrdering3D (int width, int height, int depth, Array< int > &coords)
 

Static Public Attributes

static constexpr int MaxElemNodes
 Number of nodes an element can have.
 
static constexpr int MaxElemEdges
 Number of edges an element can have.
 
static constexpr int MaxElemFaces
 Number of faces an element can have.
 
static constexpr int MaxElemChildren
 Number of children an element can have.
 
static constexpr int MaxFaceNodes
 Number of faces an element can have.
 

Protected Types

using RefPathMap = std::map<std::string, int>
 

Protected Member Functions

 NCMesh ()=default
 
void GetMeshComponents (Mesh &mesh) const
 Fill Mesh::{vertices,elements,boundary} for the current finest level.
 
void OnMeshUpdated (Mesh *mesh)
 
void MakeTopologyOnly ()
 
virtual void Update ()
 
void UpdateLeafElements ()
 Update the leaf elements indices in leaf_elements.
 
void UpdateVertices ()
 This method assigns indices to vertices (Node::vert_index) that will be seen by the Mesh class and the rest of MFEM.
 
void CollectLeafElements (int elem, int state, Array< int > &ghosts, int &counter)
 
void InitRootState (int root_count)
 
void InitGeomFlags ()
 
bool HavePrisms () const
 Return true if the mesh contains prism elements.
 
bool HavePyramids () const
 Return true if the mesh contains pyramid elements.
 
bool HaveTets () const
 Return true if the mesh contains tetrahedral elements.
 
bool IsGhost (const Element &el) const
 Return true if the Element el is a ghost element.
 
void RefineElement (int elem, char ref_type)
 Refine the element elem with the refinement type ref_type.
 
void RefineElement (const Refinement &ref)
 Refine one element with type and scale specified by ref.
 
void DerefineElement (int elem)
 Derefine the element elem, does nothing on leaf elements.
 
void SetNodeScale (int p0, int p1, real_t scale)
 Helper function to set scale for a node with parents p0, p1.
 
int AddElement (const Element &el)
 Add an Element el to the NCMesh, optimized to reuse freed elements.
 
int AddElement (Geometry::Type geom, int attr)
 
void FreeElement (int id)
 
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)
 
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)
 
int NewTetrahedron (int n0, int n1, int n2, int n3, int attr, int fattr0, int fattr1, int fattr2, int fattr3)
 
int NewPyramid (int n0, int n1, int n2, int n3, int n4, int attr, int fattr0, int fattr1, int fattr2, int fattr3, int fattr4)
 
int NewQuadrilateral (int n0, int n1, int n2, int n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
 
int NewTriangle (int n0, int n1, int n2, int attr, int eattr0, int eattr1, int eattr2)
 
int NewSegment (int n0, int n1, int attr, int vattr1, int vattr2)
 
mfem::ElementNewMeshElement (int geom) const
 
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, and if so optionally return the mid points of those edges.
 
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 been split, and if so optionally return the midpoints.
 
bool TriFaceIsMaster (int n1, int n2, int n3) const
 Determine if a Triangle face is a master face.
 
bool QuadFaceIsMaster (int n1, int n2, int n3, int n4) const
 Determine if a Quad face is a master face.
 
void ForceRefinement (int vn1, int vn2, int vn3, int vn4)
 
void FindEdgeElements (int vn1, int vn2, int vn3, int vn4, Array< MeshId > &prisms) const
 
void CheckAnisoPrism (int vn1, int vn2, int vn3, int vn4, const Refinement *refs, int nref)
 
void CheckAnisoFace (int vn1, int vn2, int vn3, int vn4, int mid12, int mid34, int level=0)
 
void CheckIsoFace (int vn1, int vn2, int vn3, int vn4, int en1, int en2, int en3, int en4, int midf)
 
void ReparentNode (int node, int new_p1, int new_p2, real_t scale)
 
int FindMidEdgeNode (int node1, int node2) const
 
int GetMidEdgeNode (int node1, int node2)
 
int GetMidFaceNode (int en1, int en2, int en3, int en4)
 
void ReferenceElement (int elem)
 Add references to all nodes, edges and faces of the element.
 
void UnreferenceElement (int elem, Array< int > &elemFaces)
 
FaceGetFace (Element &elem, int face_no)
 
void RegisterFaces (int elem, int *fattr=NULL)
 
void DeleteUnusedFaces (const Array< int > &elemFaces)
 
void CollectDerefinements (int elem, Array< Connection > &list)
 
int RetrieveNode (const Element &el, int index)
 Return el.node[index] correctly, even if the element is refined.
 
int FindNodeExt (const Element &el, int node, bool abort=true)
 Extended version of find_node: works if 'el' is refined.
 
int ReorderFacePointMat (int v0, int v1, int v2, int v3, int elem, const PointMatrix &pm, PointMatrix &reordered) const
 
void TraverseQuadFace (int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level, Face *eface[4], MatrixMap &matrix_map)
 
TriFaceTraverseResults TraverseTriFace (int vn0, int vn1, int vn2, const PointMatrix &pm, int level, MatrixMap &matrix_map)
 
void TraverseTetEdge (int vn0, int vn1, const Point &p0, const Point &p1, MatrixMap &matrix_map)
 
void TraverseEdge (int vn0, int vn1, real_t t0, real_t t1, int flags, int level, MatrixMap &matrix_map)
 
virtual void BuildFaceList ()
 
virtual void BuildEdgeList ()
 
virtual void BuildVertexList ()
 
virtual void ElementSharesFace (int elem, int local, int face)
 
virtual void ElementSharesEdge (int elem, int local, int enode)
 
virtual void ElementSharesVertex (int elem, int local, int vnode)
 
void FindSetNeighbors (const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
 
void FindNeighbors (int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
 
void NeighborExpand (const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
 
void CollectEdgeVertices (int v0, int v1, Array< int > &indices)
 
void CollectTriFaceVertices (int v0, int v1, int v2, Array< int > &indices)
 
void CollectQuadFaceVertices (int v0, int v1, int v2, int v3, Array< int > &indices)
 
void BuildElementToVertexTable ()
 
void UpdateElementToVertexTable ()
 
int GetVertexRootCoord (int elem, RefCoord coord[3]) const
 
void CollectIncidentElements (int elem, const RefCoord coord[3], Array< int > &list) const
 
void FindVertexCousins (int elem, int local, Array< int > &cousins) const
 
void GetPointMatrix (Geometry::Type geom, const char *ref_path, DenseMatrix &matrix) const
 
void TraverseRefinements (int elem, int coarse_index, std::string &ref_path, RefPathMap &map) const
 
void InitDerefTransforms ()
 
void SetDerefMatrixCodes (int parent, Array< int > &fine_coarse)
 
const real_tCalcVertexPos (int node) const
 
int GetEdgeMaster (int node) const
 
real_t GetScale (real_t s, bool reverse) const
 Return directed scale in (0,1).
 
int EdgeSplitLevel (int vn1, int vn2) const
 Return the number of splits of this edge that have occurred in the NCMesh. If zero, this means the segment is not the master of any other segments.
 
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, this means the triangle is neither split, nor the master of a split face.
 
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. If zero, this means the quad is not the master of any other quad.
 
int QuadFaceSplitLevel (int vn1, int vn2, int vn3, int vn4) const
 Returns the total number of splits of this quad that have occurred in the NCMesh. If zero, this means the quad is not the master of any other quad.
 
void CountSplits (int elem, int splits[3]) const
 
void GetLimitRefinements (Array< Refinement > &refinements, int max_level)
 
int PrintVertexParents (std::ostream *out) const
 Print the "vertex_parents" section of the mesh file.
 
void LoadVertexParents (std::istream &input)
 Load the vertex parent hierarchy from a mesh file.
 
int PrintBoundary (std::ostream *out) const
 
void LoadBoundary (std::istream &input)
 Load the "boundary" section of the mesh file.
 
void PrintCoordinates (std::ostream &out) const
 Print the "coordinates" section of the mesh file.
 
void LoadCoordinates (std::istream &input)
 Load the "coordinates" section of the mesh file.
 
void InitRootElements ()
 Count root elements and initialize root_state.
 
int CountTopLevelNodes () const
 Return the index of the last top-level node plus one.
 
bool ZeroRootStates () const
 Return true if all root_states are zero.
 
void LoadCoarseElements (std::istream &input)
 Load the element refinement hierarchy from a legacy mesh file.
 
void CopyElements (int elem, const BlockArray< Element > &tmp_elements)
 
void LoadLegacyFormat (std::istream &input, int &curved, int &is_nc)
 Load the deprecated MFEM mesh v1.1 format for backward compatibility.
 

Static Protected Member Functions

static int find_node (const Element &el, int node)
 
static int find_element_edge (const Element &el, int vn0, int vn1, bool abort=true)
 
static int find_local_face (int geom, int a, int b, int c)
 
static const PointMatrixGetGeomIdentity (Geometry::Type geom)
 
static void CheckSupportedGeom (Geometry::Type geom)
 

Protected Attributes

int Dim
 
int spaceDim
 dimensions of the elements and the vertex coordinates
 
int MyRank
 used in parallel, or when loading a parallel file in serial
 
bool Iso
 true if the mesh only contains isotropic refinements
 
int Geoms
 bit mask of element geometries present, see InitGeomFlags()
 
bool Legacy
 true if the mesh was loaded from the legacy v1.1 format
 
HashTable< Nodenodes
 
HashTable< Facefaces
 
bool using_scaling = false
 
BlockArray< Elementelements
 
Array< int > free_element_ids
 
Array< int > root_state
 
Array< real_tcoordinates
 
int NElements
 
int NVertices
 
int NEdges
 
int NFaces
 
int NGhostElements
 
int NGhostVertices
 
int NGhostEdges
 
int NGhostFaces
 
Array< int > leaf_elements
 finest elements, in Mesh ordering (+ ghosts)
 
Array< int > leaf_sfc_index
 natural tree ordering of leaf elements
 
Array< int > vertex_nodeId
 vertex-index to node-id map, see UpdateVertices
 
NCList face_list
 lazy-initialized list of faces, see GetFaceList
 
NCList edge_list
 lazy-initialized list of edges, see GetEdgeList
 
NCList vertex_list
 lazy-initialized list of vertices, see GetVertexList
 
Array< int > boundary_faces
 subset of all faces, set by BuildFaceList
 
Array< char > face_geom
 face geometry by face index, set by OnMeshUpdated
 
Table element_vertex
 leaf-element to vertex table, see FindSetNeighbors
 
Array< Refinementref_stack
 stack of scheduled refinements (temporary)
 
HashTable< Nodeshadow
 temporary storage for reparented nodes
 
Array< Triple< int, int, int > > reparents
 scheduled node reparents (tmp)
 
Array< real_treparent_scale
 scale associated with reparents (tmp)
 
Table derefinements
 possible derefinements, see GetDerefinementTable
 
CoarseFineTransformations transforms
 storage for data returned by Get[De]RefinementTransforms()
 
Array< int > coarse_elements
 state of leaf_elements before Refine(), set by MarkCoarseLevel()
 
TmpVertextmp_vertex
 

Static Protected Attributes

static PointMatrix pm_seg_identity
 
static PointMatrix pm_tri_identity
 
static PointMatrix pm_quad_identity
 
static PointMatrix pm_tet_identity
 
static PointMatrix pm_prism_identity
 
static PointMatrix pm_pyramid_identity
 
static PointMatrix pm_hex_identity
 
static GeomInfo GI [Geometry::NumGeom]
 

Friends

class Mesh
 
class ParNCMesh
 
struct MatrixMap
 
struct PointMatrixHash
 
class NCSubMesh
 
class ParNCSubMesh
 

Detailed Description

A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension of the Mesh class.

In general, the class is used by MFEM as follows:

  1. NCMesh is constructed from elements of an existing Mesh. The elements are copied and become roots of the refinement hierarchy.
  2. Some elements are refined with the Refine() method. Both isotropic and anisotropic refinements of quads/hexes are supported.
  3. A new Mesh is created from NCMesh containing the leaf elements. This new Mesh may have non-conforming (hanging) edges and faces and is the one seen by the user.
  4. FiniteElementSpace asks NCMesh for a list of conforming, master and slave edges/faces and creates the conforming interpolation matrix P.
  5. A continuous/conforming solution is obtained by solving P'*A*P x = P'*b.
  6. Repeat from step 2.

Definition at line 139 of file ncmesh.hpp.

Member Typedef Documentation

◆ RefCoord

using mfem::NCMesh::RefCoord = std::int64_t

Definition at line 486 of file ncmesh.hpp.

◆ RefPathMap

using mfem::NCMesh::RefPathMap = std::map<std::string, int>
protected

Definition at line 1189 of file ncmesh.hpp.

Constructor & Destructor Documentation

◆ NCMesh() [1/4]

mfem::NCMesh::NCMesh ( )
protecteddefault

◆ NCMesh() [2/4]

mfem::NCMesh::NCMesh ( const Mesh * mesh)
explicit

Definition at line 127 of file ncmesh.cpp.

◆ NCMesh() [3/4]

mfem::NCMesh::NCMesh ( std::istream & input,
int version,
int & curved,
int & is_nc )

Load from a stream. The id header is assumed to have been read already from

Parameters
[in]input.
[in]versionis 10 for the v1.0 NC format, or 1 for the legacy v1.1 format.
[out]curvedis set to 1 if the curvature GridFunction follows after mesh data.
[out]is_nc(again treated as a boolean) is set to 0 if the legacy v1.1 format in fact defines a conforming mesh. See Mesh::Loader for details.

Definition at line 6388 of file ncmesh.cpp.

◆ NCMesh() [4/4]

mfem::NCMesh::NCMesh ( const NCMesh & other)

Deep copy of another instance.

Definition at line 233 of file ncmesh.cpp.

◆ ~NCMesh()

mfem::NCMesh::~NCMesh ( )
virtual

Definition at line 280 of file ncmesh.cpp.

Member Function Documentation

◆ AddElement() [1/2]

int mfem::NCMesh::AddElement ( const Element & el)
inlineprotected

Add an Element el to the NCMesh, optimized to reuse freed elements.

Definition at line 812 of file ncmesh.hpp.

◆ AddElement() [2/2]

int mfem::NCMesh::AddElement ( Geometry::Type geom,
int attr )
inlineprotected

Definition at line 823 of file ncmesh.hpp.

◆ BuildEdgeList()

void mfem::NCMesh::BuildEdgeList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 3800 of file ncmesh.cpp.

◆ BuildElementToVertexTable()

void mfem::NCMesh::BuildElementToVertexTable ( )
protected

Definition at line 4129 of file ncmesh.cpp.

◆ BuildFaceList()

void mfem::NCMesh::BuildFaceList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 3675 of file ncmesh.cpp.

◆ BuildVertexList()

void mfem::NCMesh::BuildVertexList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 3912 of file ncmesh.cpp.

◆ CalcVertexPos()

const real_t * mfem::NCMesh::CalcVertexPos ( int node) const
protected

Definition at line 2733 of file ncmesh.cpp.

◆ CheckAnisoFace()

void mfem::NCMesh::CheckAnisoFace ( int vn1,
int vn2,
int vn3,
int vn4,
int mid12,
int mid34,
int level = 0 )
protected

Definition at line 978 of file ncmesh.cpp.

◆ CheckAnisoPrism()

void mfem::NCMesh::CheckAnisoPrism ( int vn1,
int vn2,
int vn3,
int vn4,
const Refinement * refs,
int nref )
protected

Definition at line 953 of file ncmesh.cpp.

◆ CheckDerefinementNCLevel()

void mfem::NCMesh::CheckDerefinementNCLevel ( const Table & deref_table,
Array< int > & level_ok,
int max_nc_level )
virtual

Check derefinements returned by GetDerefinementTable and mark those that can be done safely so that the maximum NC level condition is not violated. On return, level_ok.Size() == deref_table.Size() and contains 0/1s.

Reimplemented in mfem::ParNCMesh.

Definition at line 2278 of file ncmesh.cpp.

◆ CheckIsoFace()

void mfem::NCMesh::CheckIsoFace ( int vn1,
int vn2,
int vn3,
int vn4,
int en1,
int en2,
int en3,
int en4,
int midf )
protected

Definition at line 1099 of file ncmesh.cpp.

◆ CheckSupportedGeom()

static void mfem::NCMesh::CheckSupportedGeom ( Geometry::Type geom)
inlinestaticprotected

Definition at line 1283 of file ncmesh.hpp.

◆ ClearTransforms()

void mfem::NCMesh::ClearTransforms ( )

Free all internal data created by the above three functions.

Definition at line 5284 of file ncmesh.cpp.

◆ CollectDerefinements()

void mfem::NCMesh::CollectDerefinements ( int elem,
Array< Connection > & list )
protected

Definition at line 2230 of file ncmesh.cpp.

◆ CollectEdgeVertices()

void mfem::NCMesh::CollectEdgeVertices ( int v0,
int v1,
Array< int > & indices )
protected

Definition at line 4058 of file ncmesh.cpp.

◆ CollectIncidentElements()

void mfem::NCMesh::CollectIncidentElements ( int elem,
const RefCoord coord[3],
Array< int > & list ) const
protected

Definition at line 4515 of file ncmesh.cpp.

◆ CollectLeafElements()

void mfem::NCMesh::CollectLeafElements ( int elem,
int state,
Array< int > & ghosts,
int & counter )
protected

Collect the leaf elements in leaf_elements, and the ghost elements in ghosts. Compute and set the element indices of elements. On quad and hex refined elements tries to order leaf elements along a space-filling curve according to the given state variable.

Definition at line 2381 of file ncmesh.cpp.

◆ CollectQuadFaceVertices()

void mfem::NCMesh::CollectQuadFaceVertices ( int v0,
int v1,
int v2,
int v3,
Array< int > & indices )
protected

Definition at line 4094 of file ncmesh.cpp.

◆ CollectTriFaceVertices()

void mfem::NCMesh::CollectTriFaceVertices ( int v0,
int v1,
int v2,
Array< int > & indices )
protected

Definition at line 4070 of file ncmesh.cpp.

◆ CopyElements()

void mfem::NCMesh::CopyElements ( int elem,
const BlockArray< Element > & tmp_elements )
protected

Definition at line 6567 of file ncmesh.cpp.

◆ CountSplits()

void mfem::NCMesh::CountSplits ( int elem,
int splits[3] ) const
protected

Definition at line 5919 of file ncmesh.cpp.

◆ CountTopLevelNodes()

int mfem::NCMesh::CountTopLevelNodes ( ) const
protected

Return the index of the last top-level node plus one.

Definition at line 6378 of file ncmesh.cpp.

◆ DebugDump()

void mfem::NCMesh::DebugDump ( std::ostream & out) const

Definition at line 6958 of file ncmesh.cpp.

◆ DebugLeafOrder()

void mfem::NCMesh::DebugLeafOrder ( std::ostream & out) const

Definition at line 6933 of file ncmesh.cpp.

◆ DeleteUnusedFaces()

void mfem::NCMesh::DeleteUnusedFaces ( const Array< int > & elemFaces)
protected

Definition at line 455 of file ncmesh.cpp.

◆ Derefine()

void mfem::NCMesh::Derefine ( const Array< int > & derefs)
virtual

Perform a subset of the possible derefinements (see GetDerefinementTable). Note that if anisotropic refinements are present in the mesh, some of the derefinements may have to be skipped to preserve mesh consistency.

Reimplemented in mfem::ParNCMesh.

Definition at line 2307 of file ncmesh.cpp.

◆ DerefineElement()

void mfem::NCMesh::DerefineElement ( int elem)
protected

Derefine the element elem, does nothing on leaf elements.

Definition at line 2036 of file ncmesh.cpp.

◆ Dimension()

int mfem::NCMesh::Dimension ( ) const
inline

Return the dimension of the NCMesh.

Definition at line 164 of file ncmesh.hpp.

◆ EdgeSplitLevel()

int mfem::NCMesh::EdgeSplitLevel ( int vn1,
int vn2 ) const
protected

Return the number of splits of this edge that have occurred in the NCMesh. If zero, this means the segment is not the master of any other segments.

Parameters
vn1The first vertex making up the segment
vn2The second vertex making up the segment
Returns
int The depth of splits of this segment that are present in the mesh.

Definition at line 5862 of file ncmesh.cpp.

◆ ElementSharesEdge()

virtual void mfem::NCMesh::ElementSharesEdge ( int elem,
int local,
int enode )
inlineprotectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 1013 of file ncmesh.hpp.

◆ ElementSharesFace()

virtual void mfem::NCMesh::ElementSharesFace ( int elem,
int local,
int face )
inlineprotectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 1012 of file ncmesh.hpp.

◆ ElementSharesVertex()

virtual void mfem::NCMesh::ElementSharesVertex ( int elem,
int local,
int vnode )
inlineprotectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 1014 of file ncmesh.hpp.

◆ find_element_edge()

int mfem::NCMesh::find_element_edge ( const Element & el,
int vn0,
int vn1,
bool abort = true )
staticprotected

Definition at line 3296 of file ncmesh.cpp.

◆ find_local_face()

int mfem::NCMesh::find_local_face ( int geom,
int a,
int b,
int c )
staticprotected

Definition at line 3314 of file ncmesh.cpp.

◆ find_node()

int mfem::NCMesh::find_node ( const Element & el,
int node )
staticprotected

Definition at line 3276 of file ncmesh.cpp.

◆ FindEdgeElements()

void mfem::NCMesh::FindEdgeElements ( int vn1,
int vn2,
int vn3,
int vn4,
Array< MeshId > & prisms ) const
protected

Definition at line 905 of file ncmesh.cpp.

◆ FindFaceNodes() [1/3]

std::array< int, 4 > mfem::NCMesh::FindFaceNodes ( const Face & fa) const

Definition at line 5753 of file ncmesh.cpp.

◆ FindFaceNodes() [2/3]

std::array< int, 4 > mfem::NCMesh::FindFaceNodes ( int face) const

Method for finding the nodes associated to a face.

Returns
Nodes making up the face

Definition at line 5748 of file ncmesh.cpp.

◆ FindFaceNodes() [3/3]

void mfem::NCMesh::FindFaceNodes ( int face,
int node[4] ) const

Backwards compatible method for finding the node associated to a face.

Definition at line 5742 of file ncmesh.cpp.

◆ FindMidEdgeNode()

int mfem::NCMesh::FindMidEdgeNode ( int node1,
int node2 ) const
protected

Definition at line 336 of file ncmesh.cpp.

◆ FindNeighbors()

void mfem::NCMesh::FindNeighbors ( int elem,
Array< int > & neighbors,
const Array< int > * search_set = NULL )
protected

Return all vertex-, edge- and face-neighbors of a single element. You can limit the number of elements being checked using 'search_set'. The complexity of the function is linear in the size of the search set.

Definition at line 4320 of file ncmesh.cpp.

◆ FindNodeExt()

int mfem::NCMesh::FindNodeExt ( const Element & el,
int node,
bool abort = true )
protected

Extended version of find_node: works if 'el' is refined.

Definition at line 3286 of file ncmesh.cpp.

◆ FindSetNeighbors()

void mfem::NCMesh::FindSetNeighbors ( const Array< char > & elem_set,
Array< int > * neighbors,
Array< char > * neighbor_set = NULL )
protected

Return all vertex-, edge- and face-neighbors of a set of elements. The neighbors are returned as a list (neighbors != NULL), as a set (neighbor_set != NULL), or both. The sizes of the set arrays must match that of leaf_elements. The function is intended to be used for large sets of elements and its complexity is linear in the number of leaf elements in the mesh.

Definition at line 4209 of file ncmesh.cpp.

◆ FindVertexCousins()

void mfem::NCMesh::FindVertexCousins ( int elem,
int local,
Array< int > & cousins ) const
protected

Return elements neighboring to a local vertex of element 'elem'. Only elements from within the same refinement tree ('cousins') are returned. Complexity is proportional to the depth of elem's refinement tree.

Definition at line 4538 of file ncmesh.cpp.

◆ ForceRefinement()

void mfem::NCMesh::ForceRefinement ( int vn1,
int vn2,
int vn3,
int vn4 )
protected

Definition at line 824 of file ncmesh.cpp.

◆ FreeElement()

void mfem::NCMesh::FreeElement ( int id)
inlineprotected

Definition at line 826 of file ncmesh.hpp.

◆ GetBoundaryClosure()

void mfem::NCMesh::GetBoundaryClosure ( const Array< int > & bdr_attr_is_ess,
Array< int > & bdr_vertices,
Array< int > & bdr_edges,
Array< int > & bdr_faces )
virtual

Get a list of vertices (2D/3D), edges (3D) and faces (3D) that coincide with boundary elements with the specified attributes (marked in 'bdr_attr_is_ess').

Get a list of vertices (2D/3D), edges (3D) and faces (3D) that coincide with boundary elements with the specified attributes (marked in 'bdr_attr_is_ess'). In 3D this function also reveals "hidden" boundary edges. In parallel it helps identifying boundary vertices/edges/faces affected by non-local boundary elements. Hidden faces can occur for an internal boundary coincident to a processor boundary.

In 3D this function also reveals "hidden" boundary edges. In parallel it helps identifying boundary vertices/edges/faces affected by non-local boundary elements. Hidden faces can occur for an internal boundary coincident to a processor boundary.

Parameters
bdr_attr_is_essIndicator if a given attribute is essential.
bdr_verticesArray of vertices that are essential.
bdr_edgesArray of edges that are essential.
bdr_facesArray of faces that are essential.

Reimplemented in mfem::ParNCMesh.

Definition at line 5776 of file ncmesh.cpp.

◆ GetDerefinementTable()

const Table & mfem::NCMesh::GetDerefinementTable ( )

Return a list of derefinement opportunities. Each row of the table contains Mesh indices of existing elements that can be derefined to form a single new coarse element. Row numbers are then passed to Derefine. This function works both in serial and parallel.

Definition at line 2263 of file ncmesh.cpp.

◆ GetDerefinementTransforms()

const CoarseFineTransformations & mfem::NCMesh::GetDerefinementTransforms ( ) const

After derefinement, calculate the relations of previous fine elements (some of which may no longer exist) to the current leaf elements. Unlike for refinement, Derefine() may only be called once before this function so there is no MarkFineLevel().

Definition at line 5204 of file ncmesh.cpp.

◆ GetEdgeList()

const NCList & mfem::NCMesh::GetEdgeList ( )
inline

Return the current list of conforming and nonconforming edges.

Definition at line 326 of file ncmesh.hpp.

◆ GetEdgeMaster() [1/2]

int mfem::NCMesh::GetEdgeMaster ( int node) const
protected

Definition at line 5657 of file ncmesh.cpp.

◆ GetEdgeMaster() [2/2]

int mfem::NCMesh::GetEdgeMaster ( int v1,
int v2 ) const

Given an edge (by its vertex indices v1 and v2) return the first (geometric) parent edge that exists in the Mesh or -1 if there is no such parent.

Definition at line 5687 of file ncmesh.cpp.

◆ GetEdgeNCOrientation()

int mfem::NCMesh::GetEdgeNCOrientation ( const MeshId & edge_id) const

Return "NC" orientation of an edge. As opposed to standard Mesh edge orientation based on vertex IDs, "NC" edge orientation follows the local edge orientation within the element 'edge_id.element' and is thus processor independent. TODO: this seems only partially true?

Definition at line 5607 of file ncmesh.cpp.

◆ GetEdgeVertices()

void mfem::NCMesh::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 at line 5588 of file ncmesh.cpp.

◆ GetElement()

const Element & mfem::NCMesh::GetElement ( int i) const
inline

Access an Element.

Parameters
iIndex of the element
Returns
const Element&

Definition at line 665 of file ncmesh.hpp.

◆ GetElementDepth()

int mfem::NCMesh::GetElementDepth ( int i) const

Return the distance of leaf i from the root.

Definition at line 5696 of file ncmesh.cpp.

◆ GetElementFacesAttributes()

void mfem::NCMesh::GetElementFacesAttributes ( int i,
Array< int > & faces,
Array< int > & fattr ) const

Return the faces and face attributes of leaf element i.

Definition at line 5722 of file ncmesh.cpp.

◆ GetElementGeometry()

Geometry::Type mfem::NCMesh::GetElementGeometry ( int index) const
inline

Return element geometry type. index is the Mesh element number.

Definition at line 442 of file ncmesh.hpp.

◆ GetElementSizeReduction()

int mfem::NCMesh::GetElementSizeReduction ( int i) const

Return the size reduction compared to the root element (ignoring local stretching and curvature).

Definition at line 5708 of file ncmesh.cpp.

◆ GetFace() [1/2]

NCMesh::Face * mfem::NCMesh::GetFace ( Element & elem,
int face_no )
protected

Definition at line 480 of file ncmesh.cpp.

◆ GetFace() [2/2]

const Face & mfem::NCMesh::GetFace ( int i) const
inline

Access a Face.

Parameters
iIndex of the face
Returns
const Face&

Definition at line 652 of file ncmesh.hpp.

◆ GetFaceGeometry()

Geometry::Type mfem::NCMesh::GetFaceGeometry ( int index) const
inline

Return face geometry type. index is the Mesh face number.

Definition at line 446 of file ncmesh.hpp.

◆ GetFaceList()

const NCList & mfem::NCMesh::GetFaceList ( )
inline

Return the current list of conforming and nonconforming faces.

Definition at line 319 of file ncmesh.hpp.

◆ GetFaceVerticesEdges()

int mfem::NCMesh::GetFaceVerticesEdges ( const MeshId & face_id,
int vert_index[4],
int edge_index[4],
int edge_orientation[4] ) const

Return Mesh vertex and edge indices of a face identified by 'face_id'. The return value is the number of face vertices.

Definition at line 5619 of file ncmesh.cpp.

◆ GetGeomIdentity()

const NCMesh::PointMatrix & mfem::NCMesh::GetGeomIdentity ( Geometry::Type geom)
staticprotected

Definition at line 4606 of file ncmesh.cpp.

◆ GetLimitRefinements()

void mfem::NCMesh::GetLimitRefinements ( Array< Refinement > & refinements,
int max_level )
protected

Definition at line 6009 of file ncmesh.cpp.

◆ GetMeshComponents()

void mfem::NCMesh::GetMeshComponents ( Mesh & mesh) const
protected

Fill Mesh::{vertices,elements,boundary} for the current finest level.

Definition at line 2758 of file ncmesh.cpp.

◆ GetMidEdgeNode()

int mfem::NCMesh::GetMidEdgeNode ( int node1,
int node2 )
protected

Definition at line 352 of file ncmesh.cpp.

◆ GetMidFaceNode()

int mfem::NCMesh::GetMidFaceNode ( int en1,
int en2,
int en3,
int en4 )
protected

Definition at line 359 of file ncmesh.cpp.

◆ GetNCList()

const NCList & mfem::NCMesh::GetNCList ( int entity)
inline

Return vertex/edge/face list (entity = 0/1/2, respectively).

Definition at line 341 of file ncmesh.hpp.

◆ GetNEdges()

int mfem::NCMesh::GetNEdges ( ) const
inline

Return the number of edges in the NCMesh.

Definition at line 171 of file ncmesh.hpp.

◆ GetNFaces()

int mfem::NCMesh::GetNFaces ( ) const
inline

Return the number of (2D) faces in the NCMesh.

Definition at line 173 of file ncmesh.hpp.

◆ GetNGhostElements()

virtual int mfem::NCMesh::GetNGhostElements ( ) const
inlinevirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 174 of file ncmesh.hpp.

◆ GetNode()

const Node & mfem::NCMesh::GetNode ( int i) const
inline

Access a Node.

Parameters
iIndex of the node
Returns
const Node&

Definition at line 639 of file ncmesh.hpp.

◆ GetNodeVertex()

int mfem::NCMesh::GetNodeVertex ( int node)
inline

Given a node index, return the vertex index associated.

Parameters
node
Returns
int

Definition at line 505 of file ncmesh.hpp.

◆ GetNumElements()

int mfem::NCMesh::GetNumElements ( ) const
inline

The number of elements.

Returns
int

Definition at line 658 of file ncmesh.hpp.

◆ GetNumFaces()

int mfem::NCMesh::GetNumFaces ( ) const
inline

The number of faces.

Returns
int

Definition at line 645 of file ncmesh.hpp.

◆ GetNumNodes()

int mfem::NCMesh::GetNumNodes ( ) const
inline

The number of Nodes.

Returns
int

Definition at line 632 of file ncmesh.hpp.

◆ GetNumRootElements()

int mfem::NCMesh::GetNumRootElements ( )
inline

Return the number of root elements.

Definition at line 450 of file ncmesh.hpp.

◆ GetNVertices()

int mfem::NCMesh::GetNVertices ( ) const
inline

Return the number of vertices in the NCMesh.

Definition at line 169 of file ncmesh.hpp.

◆ GetPointMatrix()

void mfem::NCMesh::GetPointMatrix ( Geometry::Type geom,
const char * ref_path,
DenseMatrix & matrix ) const
protected

Definition at line 4623 of file ncmesh.cpp.

◆ GetRefinementTransforms()

const CoarseFineTransformations & mfem::NCMesh::GetRefinementTransforms ( ) const

After refinement, calculate the relation of each fine element to its parent coarse element. Note that Refine() or LimitNCLevel() can be called multiple times between MarkCoarseLevel() and this function.

Definition at line 5153 of file ncmesh.cpp.

◆ GetScale()

real_t mfem::NCMesh::GetScale ( real_t s,
bool reverse ) const
inlineprotected

Return directed scale in (0,1).

Definition at line 1221 of file ncmesh.hpp.

◆ GetVertexList()

const NCList & mfem::NCMesh::GetVertexList ( )
inline

Return a list of vertices (in 'conforming'); this function is provided for uniformity/completeness. Needed in ParNCMesh/ParFESpace.

Definition at line 334 of file ncmesh.hpp.

◆ GetVertexRootCoord()

int mfem::NCMesh::GetVertexRootCoord ( int elem,
RefCoord coord[3] ) const
protected

Definition at line 4462 of file ncmesh.cpp.

◆ GridSfcOrdering2D()

void mfem::NCMesh::GridSfcOrdering2D ( int width,
int height,
Array< int > & coords )
static

Return a space filling curve for a rectangular grid of elements. Implemented is a generalized Hilbert curve for arbitrary grid dimensions. If the width is odd, height should be odd too, otherwise one diagonal (vertex-neighbor) step cannot be avoided in the curve. Even dimensions are recommended.

Definition at line 5541 of file ncmesh.cpp.

◆ GridSfcOrdering3D()

void mfem::NCMesh::GridSfcOrdering3D ( int width,
int height,
int depth,
Array< int > & coords )
static

Return a space filling curve for a 3D rectangular grid of elements. The Hilbert-curve-like algorithm works well for even dimensions. For odd width/height/depth it tends to produce some diagonal (edge-neighbor) steps. Even dimensions are recommended.

Definition at line 5556 of file ncmesh.cpp.

◆ HavePrisms()

bool mfem::NCMesh::HavePrisms ( ) const
inlineprotected

Return true if the mesh contains prism elements.

Definition at line 778 of file ncmesh.hpp.

◆ HavePyramids()

bool mfem::NCMesh::HavePyramids ( ) const
inlineprotected

Return true if the mesh contains pyramid elements.

Definition at line 781 of file ncmesh.hpp.

◆ HaveTets()

bool mfem::NCMesh::HaveTets ( ) const
inlineprotected

Return true if the mesh contains tetrahedral elements.

Definition at line 784 of file ncmesh.hpp.

◆ InitDerefTransforms()

void mfem::NCMesh::InitDerefTransforms ( )
protected

Definition at line 2343 of file ncmesh.cpp.

◆ InitGeomFlags()

void mfem::NCMesh::InitGeomFlags ( )
protected

Compute the Geometry::Type present in the root elements (coarse elements) and set Geoms bitmask accordingly.

Definition at line 259 of file ncmesh.cpp.

◆ InitRootElements()

void mfem::NCMesh::InitRootElements ( )
protected

Count root elements and initialize root_state.

Definition at line 6335 of file ncmesh.cpp.

◆ InitRootState()

void mfem::NCMesh::InitRootState ( int root_count)
protected

Try to find a space-filling curve friendly orientation of the root elements: set 'root_state' based on the ordering of coarse elements. Note that the coarse mesh itself must be ordered as an SFC by e.g. Mesh::GetGeckoElementOrdering.

Definition at line 2652 of file ncmesh.cpp.

◆ IsGhost()

bool mfem::NCMesh::IsGhost ( const Element & el) const
inlineprotected

Return true if the Element el is a ghost element.

Definition at line 787 of file ncmesh.hpp.

◆ IsLegacyLoaded()

bool mfem::NCMesh::IsLegacyLoaded ( ) const
inline

I/O: Return true if the mesh was loaded from the legacy v1.1 format.

Definition at line 473 of file ncmesh.hpp.

◆ LegacyToNewVertexOrdering()

void mfem::NCMesh::LegacyToNewVertexOrdering ( Array< int > & order) const

I/O: Return a map from old (v1.1) vertex indices to new vertex indices.

Definition at line 6814 of file ncmesh.cpp.

◆ LimitNCLevel()

void mfem::NCMesh::LimitNCLevel ( int max_nc_level)
virtual

Check the mesh and potentially refine some elements so that the maximum difference of refinement levels between adjacent elements is not greater than 'max_nc_level'.

Reimplemented in mfem::ParNCMesh.

Definition at line 6039 of file ncmesh.cpp.

◆ LoadBoundary()

void mfem::NCMesh::LoadBoundary ( std::istream & input)
protected

Load the "boundary" section of the mesh file.

Definition at line 6156 of file ncmesh.cpp.

◆ LoadCoarseElements()

void mfem::NCMesh::LoadCoarseElements ( std::istream & input)
protected

Load the element refinement hierarchy from a legacy mesh file.

Definition at line 6585 of file ncmesh.cpp.

◆ LoadCoordinates()

void mfem::NCMesh::LoadCoordinates ( std::istream & input)
protected

Load the "coordinates" section of the mesh file.

Definition at line 6214 of file ncmesh.cpp.

◆ LoadLegacyFormat()

void mfem::NCMesh::LoadLegacyFormat ( std::istream & input,
int & curved,
int & is_nc )
protected

Load the deprecated MFEM mesh v1.1 format for backward compatibility.

Definition at line 6656 of file ncmesh.cpp.

◆ LoadVertexParents()

void mfem::NCMesh::LoadVertexParents ( std::istream & input)
protected

Load the vertex parent hierarchy from a mesh file.

Definition at line 6090 of file ncmesh.cpp.

◆ MakeTopologyOnly()

void mfem::NCMesh::MakeTopologyOnly ( )
inlineprotected

Delete top-level vertex coordinates if the Mesh became curved, e.g., by calling Mesh::SetCurvature or otherwise setting the Nodes.

Definition at line 520 of file ncmesh.hpp.

◆ MarkCoarseLevel()

void mfem::NCMesh::MarkCoarseLevel ( )

Remember the current layer of leaf elements before the mesh is refined. Needed by GetRefinementTransforms(), must be called before Refine().

Definition at line 5105 of file ncmesh.cpp.

◆ MemoryUsage()

long mfem::NCMesh::MemoryUsage ( ) const

Return total number of bytes allocated.

Definition at line 6882 of file ncmesh.cpp.

◆ NeighborExpand()

void mfem::NCMesh::NeighborExpand ( const Array< int > & elems,
Array< int > & expanded,
const Array< int > * search_set = NULL )
protected

Expand a set of elements by all vertex-, edge- and face-neighbors. The output array 'expanded' will contain all items from 'elems' (provided they are in 'search_set') plus their neighbors. The neighbor search can be limited to the optional search set. The complexity is linear in the sum of the sizes of 'elems' and 'search_set'.

Definition at line 4403 of file ncmesh.cpp.

◆ NewHexahedron()

int mfem::NCMesh::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 )
protected

Definition at line 615 of file ncmesh.cpp.

◆ NewMeshElement()

mfem::Element * mfem::NCMesh::NewMeshElement ( int geom) const
protected

Definition at line 2718 of file ncmesh.cpp.

◆ NewPyramid()

int mfem::NCMesh::NewPyramid ( int n0,
int n1,
int n2,
int n3,
int n4,
int attr,
int fattr0,
int fattr1,
int fattr2,
int fattr3,
int fattr4 )
protected

Definition at line 702 of file ncmesh.cpp.

◆ NewQuadrilateral()

int mfem::NCMesh::NewQuadrilateral ( int n0,
int n1,
int n2,
int n3,
int attr,
int eattr0,
int eattr1,
int eattr2,
int eattr3 )
protected

Definition at line 732 of file ncmesh.cpp.

◆ NewSegment()

int mfem::NCMesh::NewSegment ( int n0,
int n1,
int attr,
int vattr1,
int vattr2 )
protected

Definition at line 784 of file ncmesh.cpp.

◆ NewTetrahedron()

int mfem::NCMesh::NewTetrahedron ( int n0,
int n1,
int n2,
int n3,
int attr,
int fattr0,
int fattr1,
int fattr2,
int fattr3 )
protected

Definition at line 677 of file ncmesh.cpp.

◆ NewTriangle()

int mfem::NCMesh::NewTriangle ( int n0,
int n1,
int n2,
int attr,
int eattr0,
int eattr1,
int eattr2 )
protected

Definition at line 758 of file ncmesh.cpp.

◆ NewWedge()

int mfem::NCMesh::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 )
protected

Definition at line 645 of file ncmesh.cpp.

◆ OnMeshUpdated()

void mfem::NCMesh::OnMeshUpdated ( Mesh * mesh)
protected

Get edge and face numbering from 'mesh' (i.e., set all Edge::index and Face::index) after a new mesh was created from us.

Definition at line 2885 of file ncmesh.cpp.

◆ operator=()

NCMesh & mfem::NCMesh::operator= ( NCMesh & )
delete

Copy assignment not supported.

◆ ParentFaceNodes()

int mfem::NCMesh::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 parent face and replace the input nodes with the parent nodes. Additionally return the child index that the child face would be, relative to the discovered parent face.

This method is concerned with the construction of an NCMesh structure for a d-1 manifold of an existing NCMesh. It forms a key element in a leaf -> root traversal of the parent ncmesh elements structure.

Parameters
[out]nodesThe collection of nodes whose parent we are searching for
Returns
int The child index corresponding to placing the face for the original nodes within the face defined by the returned parent nodes. If child index is -1, then the face is made up of root nodes, and nodes is unchanged.

Definition at line 3106 of file ncmesh.cpp.

◆ Print()

void mfem::NCMesh::Print ( std::ostream & out,
const std::string & comments = "" ) const

I/O: Print the mesh in "MFEM NC mesh v1.0" format. If comments is non-empty, it will be printed after the first line of the file, and each line should begin with '#'.

Definition at line 6244 of file ncmesh.cpp.

◆ PrintBoundary()

int mfem::NCMesh::PrintBoundary ( std::ostream * out) const
protected

Print the "boundary" section of the mesh file. If out == NULL, only return the number of boundary elements.

Definition at line 6117 of file ncmesh.cpp.

◆ PrintCoordinates()

void mfem::NCMesh::PrintCoordinates ( std::ostream & out) const
protected

Print the "coordinates" section of the mesh file.

Definition at line 6196 of file ncmesh.cpp.

◆ PrintMemoryDetail()

int mfem::NCMesh::PrintMemoryDetail ( ) const

Definition at line 6905 of file ncmesh.cpp.

◆ PrintVertexParents()

int mfem::NCMesh::PrintVertexParents ( std::ostream * out) const
protected

Print the "vertex_parents" section of the mesh file.

Definition at line 6055 of file ncmesh.cpp.

◆ QuadFaceIsMaster()

bool mfem::NCMesh::QuadFaceIsMaster ( int n1,
int n2,
int n3,
int n4 ) const
inlineprotected

Determine if a Quad face is a master face.

Parameters
n1The first node defining the face
n2The second node defining the face
n3The third node defining the face
n4The fourth node defining the face
Returns
true The quad face is a master face
false The quad face is not a master face

Definition at line 930 of file ncmesh.hpp.

◆ QuadFaceSplitLevel() [1/2]

int mfem::NCMesh::QuadFaceSplitLevel ( int vn1,
int vn2,
int vn3,
int vn4 ) const
protected

Returns the total number of splits of this quad that have occurred in the NCMesh. If zero, this means the quad is not the master of any other quad.

This is a convenience wrapper that sums the horizontal and vertical levels from the full method.

Parameters
vn1The first vertex making up the quad
vn2The second vertex making up the quad
vn3The third vertex making up the quad
vn4The fourth vertex making up the quad
Returns
int The depth of splits of this triangle that are present in the mesh. NB: An isotropic refinement has a level of 2, one horizontal split, followed by a vertical split.

Definition at line 5912 of file ncmesh.cpp.

◆ QuadFaceSplitLevel() [2/2]

void mfem::NCMesh::QuadFaceSplitLevel ( int vn1,
int vn2,
int vn3,
int vn4,
int & h_level,
int & v_level ) const
protected

Computes the number of horizontal and vertical splits of this quad that have occurred in the NCMesh. If zero, this means the quad is not the master of any other quad.

Parameters
vn1The first vertex making up the quad
vn2The second vertex making up the quad
vn3The third vertex making up the quad
vn4The fourth vertex making up the quad
h_levelThe number of "horizontal" splits of the quad
v_levelThe number of "vertical" splits of the quad

Definition at line 5884 of file ncmesh.cpp.

◆ QuadFaceSplitType()

int mfem::NCMesh::QuadFaceSplitType ( int n1,
int n2,
int n3,
int n4,
real_t & s,
int mid[5] = NULL ) const
protected

Given a quad face defined by four vertices, establish which edges of this face have been split, and if so optionally return the mid points of those edges.

Parameters
n1The first node defining the face
n2The second node defining the face
n3The third node defining the face
n4The fourth node defining the face
sreturns the scale of the split
midoptional return of the edge mid points.
Returns
int 0 – no split, 1 – "vertical" split, 2 – "horizontal" split

Definition at line 3028 of file ncmesh.cpp.

◆ ReferenceElement()

void mfem::NCMesh::ReferenceElement ( int elem)
protected

Add references to all nodes, edges and faces of the element.

Parameters
elemindex into elements

Definition at line 367 of file ncmesh.cpp.

◆ Refine()

void mfem::NCMesh::Refine ( const Array< Refinement > & refinements)
virtual

Perform the given batch of refinements. Please note that in the presence of anisotropic splits additional refinements may be necessary to keep the mesh consistent. However, the function always performs at least the requested refinements.

Reimplemented in mfem::ParNCMesh.

Definition at line 1945 of file ncmesh.cpp.

◆ RefineElement() [1/2]

void mfem::NCMesh::RefineElement ( const Refinement & ref)
protected

Refine one element with type and scale specified by ref.

Definition at line 1126 of file ncmesh.cpp.

◆ RefineElement() [2/2]

void mfem::NCMesh::RefineElement ( int elem,
char ref_type )
protected

Refine the element elem with the refinement type ref_type.

Definition at line 1121 of file ncmesh.cpp.

◆ RegisterFaces()

void mfem::NCMesh::RegisterFaces ( int elem,
int * fattr = NULL )
protected

Definition at line 441 of file ncmesh.cpp.

◆ ReorderFacePointMat()

int mfem::NCMesh::ReorderFacePointMat ( int v0,
int v1,
int v2,
int v3,
int elem,
const PointMatrix & pm,
PointMatrix & reordered ) const
protected

Definition at line 3421 of file ncmesh.cpp.

◆ ReparentNode()

void mfem::NCMesh::ReparentNode ( int node,
int new_p1,
int new_p2,
real_t scale )
protected

Definition at line 319 of file ncmesh.cpp.

◆ RetrieveNode()

int mfem::NCMesh::RetrieveNode ( const Element & el,
int index )
protected

Return el.node[index] correctly, even if the element is refined.

Definition at line 1994 of file ncmesh.cpp.

◆ SetAttribute()

void mfem::NCMesh::SetAttribute ( int i,
int attr )
inline

Set the attribute of leaf element i, which is a Mesh element index.

Definition at line 464 of file ncmesh.hpp.

◆ SetDerefMatrixCodes()

void mfem::NCMesh::SetDerefMatrixCodes ( int parent,
Array< int > & fine_coarse )
protected

Definition at line 2362 of file ncmesh.cpp.

◆ SetNodeScale()

void mfem::NCMesh::SetNodeScale ( int p0,
int p1,
real_t scale )
protected

Helper function to set scale for a node with parents p0, p1.

Definition at line 1115 of file ncmesh.cpp.

◆ SpaceDimension()

int mfem::NCMesh::SpaceDimension ( ) const
inline

Return the space dimension of the NCMesh.

Definition at line 166 of file ncmesh.hpp.

◆ TraverseEdge()

void mfem::NCMesh::TraverseEdge ( int vn0,
int vn1,
real_t t0,
real_t t1,
int flags,
int level,
MatrixMap & matrix_map )
protected

Definition at line 3770 of file ncmesh.cpp.

◆ TraverseQuadFace()

void mfem::NCMesh::TraverseQuadFace ( int vn0,
int vn1,
int vn2,
int vn3,
const PointMatrix & pm,
int level,
Face * eface[4],
MatrixMap & matrix_map )
protected

Definition at line 3452 of file ncmesh.cpp.

◆ TraverseRefinements()

void mfem::NCMesh::TraverseRefinements ( int elem,
int coarse_index,
std::string & ref_path,
RefPathMap & map ) const
protected

Definition at line 5119 of file ncmesh.cpp.

◆ TraverseTetEdge()

void mfem::NCMesh::TraverseTetEdge ( int vn0,
int vn1,
const Point & p0,
const Point & p1,
MatrixMap & matrix_map )
protected

Definition at line 3572 of file ncmesh.cpp.

◆ TraverseTriFace()

NCMesh::TriFaceTraverseResults mfem::NCMesh::TraverseTriFace ( int vn0,
int vn1,
int vn2,
const PointMatrix & pm,
int level,
MatrixMap & matrix_map )
protected

Definition at line 3609 of file ncmesh.cpp.

◆ TriFaceIsMaster()

bool mfem::NCMesh::TriFaceIsMaster ( int n1,
int n2,
int n3 ) const
inlineprotected

Determine if a Triangle face is a master face.

This check requires looking for the edges making up the triangle being split, if nodes exist at their midpoints, and there are vertices at them, this implies the face COULD be split. To determine if it is, we then check whether these midpoints have all been connected, this is required to discriminate between an internal master face surrounded by nonconformal refinements and a conformal boundary face surrounded by refinements.

Parameters
n1The first node defining the face
n2The second node defining the face
n3The third node defining the face
Returns
true The face is a master
false The face is not a master

Definition at line 910 of file ncmesh.hpp.

◆ TriFaceSplit()

bool mfem::NCMesh::TriFaceSplit ( int n1,
int n2,
int n3,
int mid[3] = NULL ) const
protected

Given a tri face defined by three vertices, establish whether the edges that make up this face have been split, and if so optionally return the midpoints.

This is a necessary condition for this face to have been split, but is not sufficient. Consider a triangle attached to three refined triangles, in this scenario all edges can be split but this face not be split. In this case, it is necessary to check if there is a face made up of the returned midpoint nodes.

Parameters
n1The first node defining the face
n2The second node defining the face
n3The third node defining the face
midoptional return of the edge mid points.
Returns
true Splits for all edges have been found
false

Definition at line 3082 of file ncmesh.cpp.

◆ TriFaceSplitLevel()

int mfem::NCMesh::TriFaceSplitLevel ( int vn1,
int vn2,
int vn3 ) const
protected

Return the number of splits of this triangle that have occurred in the NCMesh. If zero, this means the triangle is neither split, nor the master of a split face.

Parameters
vn1The first vertex making up the triangle
vn2The second vertex making up the triangle
vn3The third vertex making up the triangle
Returns
int The depth of splits of this triangle that are present in the mesh.

Definition at line 5869 of file ncmesh.cpp.

◆ Trim()

void mfem::NCMesh::Trim ( )
virtual

Save memory by releasing all non-essential and cached data.

Reimplemented in mfem::ParNCMesh.

Definition at line 6839 of file ncmesh.cpp.

◆ UnreferenceElement()

void mfem::NCMesh::UnreferenceElement ( int elem,
Array< int > & elemFaces )
protected

Definition at line 398 of file ncmesh.cpp.

◆ Update()

void mfem::NCMesh::Update ( )
protectedvirtual

Apart from the primary data structure, which is the element/node/face hierarchy, there is secondary data that is derived from the primary data and needs to be updated when the primary data changes. Update() takes care of that and needs to be called after each refinement and derefinement.

Reimplemented in mfem::ParNCMesh.

Definition at line 268 of file ncmesh.cpp.

◆ UpdateElementToVertexTable()

void mfem::NCMesh::UpdateElementToVertexTable ( )
inlineprotected

Definition at line 1051 of file ncmesh.hpp.

◆ UpdateLeafElements()

void mfem::NCMesh::UpdateLeafElements ( )
protected

Update the leaf elements indices in leaf_elements.

Definition at line 2449 of file ncmesh.cpp.

◆ UpdateVertices()

void mfem::NCMesh::UpdateVertices ( )
protected

This method assigns indices to vertices (Node::vert_index) that will be seen by the Mesh class and the rest of MFEM.

We must be careful to:

  1. Stay compatible with the conforming code, which expects top-level (original) vertices to be indexed first, otherwise GridFunctions defined on a conforming mesh would no longer be valid when the mesh is converted to an NC mesh.
  2. Make sure serial NCMesh is compatible with the parallel ParNCMesh, so it is possible to read parallel partial solutions in serial code (e.g., serial GLVis). This means handling ghost elements, if present.
  3. Assign vertices in a globally consistent order for parallel meshes: if two vertices i,j are shared by two ranks r1,r2, and i<j on r1, then i<j on r2 as well. This is true for top-level vertices but also for the remaining shared vertices thanks to the globally consistent SFC ordering of the leaf elements. This property reduces communication and simplifies ParNCMesh. update Vertex::index and vertex_nodeId

Definition at line 2476 of file ncmesh.cpp.

◆ ZeroRootStates()

bool mfem::NCMesh::ZeroRootStates ( ) const
protected

Return true if all root_states are zero.

Definition at line 6235 of file ncmesh.cpp.

Friends And Related Symbol Documentation

◆ MatrixMap

friend struct MatrixMap
friend

Definition at line 1349 of file ncmesh.hpp.

◆ Mesh

friend class Mesh
friend

Definition at line 509 of file ncmesh.hpp.

◆ NCSubMesh

friend class NCSubMesh
friend

Definition at line 1351 of file ncmesh.hpp.

◆ ParNCMesh

friend class ParNCMesh
friend

Definition at line 1348 of file ncmesh.hpp.

◆ ParNCSubMesh

friend class ParNCSubMesh
friend

Definition at line 1352 of file ncmesh.hpp.

◆ PointMatrixHash

friend struct PointMatrixHash
friend

Definition at line 1350 of file ncmesh.hpp.

Member Data Documentation

◆ boundary_faces

Array<int> mfem::NCMesh::boundary_faces
protected

subset of all faces, set by BuildFaceList

Definition at line 731 of file ncmesh.hpp.

◆ coarse_elements

Array<int> mfem::NCMesh::coarse_elements
protected

state of leaf_elements before Refine(), set by MarkCoarseLevel()

Definition at line 1198 of file ncmesh.hpp.

◆ coordinates

Array<real_t> mfem::NCMesh::coordinates
protected

Coordinates of top-level vertices (organized as triples). If empty, the Mesh is curved (Nodes != NULL) and NCMesh is topology-only.

Definition at line 705 of file ncmesh.hpp.

◆ derefinements

Table mfem::NCMesh::derefinements
protected

possible derefinements, see GetDerefinementTable

Definition at line 796 of file ncmesh.hpp.

◆ Dim

int mfem::NCMesh::Dim
protected

Definition at line 524 of file ncmesh.hpp.

◆ edge_list

NCList mfem::NCMesh::edge_list
protected

lazy-initialized list of edges, see GetEdgeList

Definition at line 728 of file ncmesh.hpp.

◆ element_vertex

Table mfem::NCMesh::element_vertex
protected

leaf-element to vertex table, see FindSetNeighbors

Definition at line 734 of file ncmesh.hpp.

◆ elements

BlockArray<Element> mfem::NCMesh::elements
protected

Definition at line 624 of file ncmesh.hpp.

◆ face_geom

Array<char> mfem::NCMesh::face_geom
protected

face geometry by face index, set by OnMeshUpdated

Definition at line 732 of file ncmesh.hpp.

◆ face_list

NCList mfem::NCMesh::face_list
protected

lazy-initialized list of faces, see GetFaceList

Definition at line 727 of file ncmesh.hpp.

◆ faces

HashTable<Face> mfem::NCMesh::faces
protected

Definition at line 620 of file ncmesh.hpp.

◆ free_element_ids

Array<int> mfem::NCMesh::free_element_ids
protected

Definition at line 625 of file ncmesh.hpp.

◆ Geoms

int mfem::NCMesh::Geoms
protected

bit mask of element geometries present, see InitGeomFlags()

Definition at line 527 of file ncmesh.hpp.

◆ GI

NCMesh::GeomInfo mfem::NCMesh::GI
staticprotected

Definition at line 1340 of file ncmesh.hpp.

◆ Iso

bool mfem::NCMesh::Iso
protected

true if the mesh only contains isotropic refinements

Definition at line 526 of file ncmesh.hpp.

◆ leaf_elements

Array<int> mfem::NCMesh::leaf_elements
protected

finest elements, in Mesh ordering (+ ghosts)

Definition at line 723 of file ncmesh.hpp.

◆ leaf_sfc_index

Array<int> mfem::NCMesh::leaf_sfc_index
protected

natural tree ordering of leaf elements

Definition at line 724 of file ncmesh.hpp.

◆ Legacy

bool mfem::NCMesh::Legacy
protected

true if the mesh was loaded from the legacy v1.1 format

Definition at line 528 of file ncmesh.hpp.

◆ MaxElemChildren

int mfem::NCMesh::MaxElemChildren
staticconstexpr
Initial value:
=
10

Number of children an element can have.

Definition at line 494 of file ncmesh.hpp.

◆ MaxElemEdges

int mfem::NCMesh::MaxElemEdges
staticconstexpr
Initial value:
=
12

Number of edges an element can have.

Definition at line 490 of file ncmesh.hpp.

◆ MaxElemFaces

int mfem::NCMesh::MaxElemFaces
staticconstexpr
Initial value:
=
6

Number of faces an element can have.

Definition at line 492 of file ncmesh.hpp.

◆ MaxElemNodes

int mfem::NCMesh::MaxElemNodes
staticconstexpr
Initial value:
=
8

Number of nodes an element can have.

Definition at line 488 of file ncmesh.hpp.

◆ MaxFaceNodes

int mfem::NCMesh::MaxFaceNodes
staticconstexpr
Initial value:
=
4

Number of faces an element can have.

Definition at line 496 of file ncmesh.hpp.

◆ MyRank

int mfem::NCMesh::MyRank
protected

used in parallel, or when loading a parallel file in serial

Definition at line 525 of file ncmesh.hpp.

◆ NEdges

int mfem::NCMesh::NEdges
protected

Definition at line 716 of file ncmesh.hpp.

◆ NElements

int mfem::NCMesh::NElements
protected

Definition at line 716 of file ncmesh.hpp.

◆ NFaces

int mfem::NCMesh::NFaces
protected

Definition at line 716 of file ncmesh.hpp.

◆ NGhostEdges

int mfem::NCMesh::NGhostEdges
protected

Definition at line 721 of file ncmesh.hpp.

◆ NGhostElements

int mfem::NCMesh::NGhostElements
protected

Definition at line 721 of file ncmesh.hpp.

◆ NGhostFaces

int mfem::NCMesh::NGhostFaces
protected

Definition at line 721 of file ncmesh.hpp.

◆ NGhostVertices

int mfem::NCMesh::NGhostVertices
protected

Definition at line 721 of file ncmesh.hpp.

◆ nodes

HashTable<Node> mfem::NCMesh::nodes
protected

Definition at line 619 of file ncmesh.hpp.

◆ NVertices

int mfem::NCMesh::NVertices
protected

Definition at line 716 of file ncmesh.hpp.

◆ pm_hex_identity

NCMesh::PointMatrix mfem::NCMesh::pm_hex_identity
staticprotected

Definition at line 1182 of file ncmesh.hpp.

◆ pm_prism_identity

NCMesh::PointMatrix mfem::NCMesh::pm_prism_identity
staticprotected

Definition at line 1180 of file ncmesh.hpp.

◆ pm_pyramid_identity

NCMesh::PointMatrix mfem::NCMesh::pm_pyramid_identity
staticprotected

Definition at line 1181 of file ncmesh.hpp.

◆ pm_quad_identity

NCMesh::PointMatrix mfem::NCMesh::pm_quad_identity
staticprotected

Definition at line 1178 of file ncmesh.hpp.

◆ pm_seg_identity

NCMesh::PointMatrix mfem::NCMesh::pm_seg_identity
staticprotected

Definition at line 1176 of file ncmesh.hpp.

◆ pm_tet_identity

NCMesh::PointMatrix mfem::NCMesh::pm_tet_identity
staticprotected

Definition at line 1179 of file ncmesh.hpp.

◆ pm_tri_identity

NCMesh::PointMatrix mfem::NCMesh::pm_tri_identity
staticprotected

Definition at line 1177 of file ncmesh.hpp.

◆ ref_stack

Array<Refinement> mfem::NCMesh::ref_stack
protected

stack of scheduled refinements (temporary)

Definition at line 791 of file ncmesh.hpp.

◆ reparent_scale

Array<real_t> mfem::NCMesh::reparent_scale
protected

scale associated with reparents (tmp)

Definition at line 794 of file ncmesh.hpp.

◆ reparents

Array<Triple<int, int, int> > mfem::NCMesh::reparents
protected

scheduled node reparents (tmp)

Definition at line 793 of file ncmesh.hpp.

◆ root_state

Array<int> mfem::NCMesh::root_state
protected

Initial traversal state (~ element orientation) for each root element NOTE: M = root_state.Size() is the number of root elements. NOTE: the first M items of 'elements' is the coarse mesh.

Definition at line 701 of file ncmesh.hpp.

◆ shadow

HashTable<Node> mfem::NCMesh::shadow
protected

temporary storage for reparented nodes

Definition at line 792 of file ncmesh.hpp.

◆ spaceDim

int mfem::NCMesh::spaceDim
protected

dimensions of the elements and the vertex coordinates

Definition at line 524 of file ncmesh.hpp.

◆ tmp_vertex

TmpVertex* mfem::NCMesh::tmp_vertex
mutableprotected

Definition at line 1211 of file ncmesh.hpp.

◆ transforms

CoarseFineTransformations mfem::NCMesh::transforms
mutableprotected

storage for data returned by Get[De]RefinementTransforms()

Definition at line 1195 of file ncmesh.hpp.

◆ using_scaling

bool mfem::NCMesh::using_scaling = false
protected

Definition at line 622 of file ncmesh.hpp.

◆ vertex_list

NCList mfem::NCMesh::vertex_list
protected

lazy-initialized list of vertices, see GetVertexList

Definition at line 729 of file ncmesh.hpp.

◆ vertex_nodeId

Array<int> mfem::NCMesh::vertex_nodeId
protected

vertex-index to node-id map, see UpdateVertices

Definition at line 725 of file ncmesh.hpp.


The documentation for this class was generated from the following files: