MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Protected Attributes | Friends | List of all members
mfem::NCMesh Class Reference

A class for non-conforming AMR on higher-order hexahedral, prismatic, quadrilateral or triangular meshes. More...

#include <ncmesh.hpp>

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

Classes

struct  Element
 
struct  Face
 
struct  GeomInfo
 
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
 
struct  Slave
 Nonconforming edge/face within a bigger edge/face. More...
 
struct  TmpVertex
 

Public Types

typedef int64_t RefCoord
 

Public Member Functions

 NCMesh (const Mesh *mesh, std::istream *vertex_parents=NULL)
 
 NCMesh (const NCMesh &other)
 
virtual ~NCMesh ()
 
int Dimension () const
 
int SpaceDimension () const
 
int GetNVertices () const
 
int GetNEdges () const
 
int GetNFaces () 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. More...
 
const NCListGetEdgeList ()
 Return the current list of conforming and nonconforming edges. More...
 
const NCListGetVertexList ()
 
const NCListGetNCList (int entity)
 Return vertex/edge/face list (entity = 0/1/2, respectively). More...
 
void MarkCoarseLevel ()
 
const CoarseFineTransformationsGetRefinementTransforms ()
 
const CoarseFineTransformationsGetDerefinementTransforms ()
 
void ClearTransforms ()
 Free all internal data created by the above three functions. More...
 
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'. More...
 
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)
 
Geometry::Type GetElementGeometry (int index) const
 Return element geometry type. index is the Mesh element number. More...
 
Geometry::Type GetFaceGeometry (int index) const
 Return face geometry type. index is the Mesh face number. More...
 
int GetNumRootElements ()
 Return the number of root elements. More...
 
int GetElementDepth (int i) const
 Return the distance of leaf 'i' from the root. More...
 
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'. More...
 
void PrintVertexParents (std::ostream &out) const
 I/O: Print the "vertex_parents" section of the mesh file (ver. >= 1.1). More...
 
void PrintCoarseElements (std::ostream &out) const
 I/O: Print the "coarse_elements" section of the mesh file (ver. >= 1.1). More...
 
void LoadVertexParents (std::istream &input)
 
void LoadCoarseElements (std::istream &input)
 I/O: Load the element refinement hierarchy from a mesh file. More...
 
void SetVertexPositions (const Array< mfem::Vertex > &vertices)
 I/O: Set positions of all vertices (used by mesh loader). More...
 
virtual void Trim ()
 Save memory by releasing all non-essential and cached data. More...
 
long MemoryUsage () const
 Return total number of bytes allocated. More...
 
int PrintMemoryDetail () const
 
void PrintStats (std::ostream &out=mfem::out) const
 
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)
 

Protected Types

typedef HashTable< Node >::iterator node_iterator
 
typedef HashTable< Face >::iterator face_iterator
 
typedef HashTable< Node >
::const_iterator 
node_const_iterator
 
typedef HashTable< Face >
::const_iterator 
face_const_iterator
 
typedef BlockArray< Element >
::iterator 
elem_iterator
 
typedef std::map< std::string,
int > 
RefPathMap
 

Protected Member Functions

void GetMeshComponents (Mesh &mesh) const
 Fill Mesh::{vertices,elements,boundary} for the current finest level. More...
 
virtual void OnMeshUpdated (Mesh *mesh)
 
virtual void Update ()
 
virtual void UpdateVertices ()
 update Vertex::index and vertex_nodeId More...
 
void CollectLeafElements (int elem, int state)
 
void UpdateLeafElements ()
 
virtual void AssignLeafIndices ()
 
void InitRootState (int root_count)
 
virtual bool IsGhost (const Element &el) const
 
virtual int GetNumGhostElements () const
 
virtual int GetNumGhostVertices () const
 
void InitGeomFlags ()
 
bool HavePrisms () const
 
bool HaveTets () const
 
void RefineElement (int elem, char ref_type)
 
void DerefineElement (int elem)
 
int AddElement (const Element &el)
 
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 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)
 
mfem::ElementNewMeshElement (int geom) const
 
int QuadFaceSplitType (int v1, int v2, int v3, int v4, int mid[5]=NULL) const
 
bool TriFaceSplit (int v1, int v2, int v3, int mid[3]=NULL) const
 
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)
 
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)
 
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. More...
 
int FindNodeExt (const Element &el, int node, bool abort=true)
 Extended version of find_node: works if 'el' is refined. More...
 
int ReorderFacePointMat (int v0, int v1, int v2, int v3, int elem, DenseMatrix &mat) const
 
void TraverseQuadFace (int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level, Face *eface[4])
 
bool TraverseTriFace (int vn0, int vn1, int vn2, const PointMatrix &pm, int level)
 
void TraverseTetEdge (int vn0, int vn1, const Point &p0, const Point &p1)
 
void TraverseEdge (int vn0, int vn1, double t0, double t1, int flags, int level)
 
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)
 
void TraverseRefinements (int elem, int coarse_index, std::string &ref_path, RefPathMap &map)
 
void InitDerefTransforms ()
 
void SetDerefMatrixCodes (int parent, Array< int > &fine_coarse)
 
const double * CalcVertexPos (int node) const
 
int GetEdgeMaster (int node) const
 
void FindFaceNodes (int face, int node[4])
 
int EdgeSplitLevel (int vn1, int vn2) const
 
int TriFaceSplitLevel (int vn1, int vn2, int vn3) const
 
void QuadFaceSplitLevel (int vn1, int vn2, int vn3, int vn4, int &h_level, int &v_level) const
 
void CountSplits (int elem, int splits[3]) const
 
void GetLimitRefinements (Array< Refinement > &refinements, int max_level)
 
int PrintElements (std::ostream &out, int elem, int &coarse_id) const
 
void CopyElements (int elem, const BlockArray< Element > &tmp_elements, Array< int > &index_map)
 

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)
 

Protected Attributes

int Dim
 
int spaceDim
 dimensions of the elements and the vertex coordinates More...
 
bool Iso
 true if the mesh only contains isotropic refinements More...
 
int Geoms
 bit mask of element geometries present, see InitGeomFlags() More...
 
HashTable< Nodenodes
 
HashTable< Facefaces
 
BlockArray< Elementelements
 
Array< int > free_element_ids
 
Array< int > root_state
 
Array< double > top_vertex_pos
 coordinates of top-level vertices (organized as triples) More...
 
int NVertices
 
int NEdges
 
int NFaces
 
Array< int > leaf_elements
 
Array< int > vertex_nodeId
 
NCList face_list
 lazy-initialized list of faces, see GetFaceList More...
 
NCList edge_list
 lazy-initialized list of edges, see GetEdgeList More...
 
NCList vertex_list
 lazy-initialized list of vertices, see GetVertexList More...
 
Array< int > boundary_faces
 subset of all faces, set by BuildFaceList More...
 
Array< char > face_geom
 face geometry by face index, set by OnMeshUpdated More...
 
Table element_vertex
 leaf-element to vertex table, see FindSetNeighbors More...
 
Array< Refinementref_stack
 stack of scheduled refinements (temporary) More...
 
HashTable< Nodeshadow
 temporary storage for reparented nodes More...
 
Array< Triple< int, int, int > > reparents
 scheduled node reparents (tmp) More...
 
Table derefinements
 possible derefinements, see GetDerefinementTable More...
 
CoarseFineTransformations transforms
 storage for data returned by Get[De]RefinementTransforms() More...
 
Array< int > coarse_elements
 state of leaf_elements before Refine(), set by MarkCoarseLevel() More...
 
TmpVertextmp_vertex
 

Static Protected Attributes

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

Friends

class Mesh
 
class ParNCMesh
 
struct CompareRanks
 

Detailed Description

A class for non-conforming AMR on higher-order hexahedral, prismatic, quadrilateral or triangular meshes.

The class is used 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.
  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 100 of file ncmesh.hpp.

Member Typedef Documentation

typedef BlockArray<Element>::iterator mfem::NCMesh::elem_iterator
protected

Definition at line 477 of file ncmesh.hpp.

typedef HashTable<Face>::const_iterator mfem::NCMesh::face_const_iterator
protected

Definition at line 476 of file ncmesh.hpp.

typedef HashTable<Face>::iterator mfem::NCMesh::face_iterator
protected

Definition at line 474 of file ncmesh.hpp.

typedef HashTable<Node>::const_iterator mfem::NCMesh::node_const_iterator
protected

Definition at line 475 of file ncmesh.hpp.

typedef HashTable<Node>::iterator mfem::NCMesh::node_iterator
protected

Definition at line 473 of file ncmesh.hpp.

typedef int64_t mfem::NCMesh::RefCoord

Definition at line 365 of file ncmesh.hpp.

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

Definition at line 784 of file ncmesh.hpp.

Constructor & Destructor Documentation

mfem::NCMesh::NCMesh ( const Mesh mesh,
std::istream *  vertex_parents = NULL 
)
explicit

Initialize with elements from 'mesh'. If an already nonconforming mesh is being loaded, 'vertex_parents' must point to a stream at the appropriate section of the mesh file which contains the vertex hierarchy.

Definition at line 72 of file ncmesh.cpp.

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

Definition at line 188 of file ncmesh.cpp.

mfem::NCMesh::~NCMesh ( )
virtual

Definition at line 225 of file ncmesh.cpp.

Member Function Documentation

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

Definition at line 539 of file ncmesh.hpp.

void mfem::NCMesh::AssignLeafIndices ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 1918 of file ncmesh.cpp.

void mfem::NCMesh::BuildEdgeList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 2656 of file ncmesh.cpp.

void mfem::NCMesh::BuildElementToVertexTable ( )
protected

Definition at line 2971 of file ncmesh.cpp.

void mfem::NCMesh::BuildFaceList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 2542 of file ncmesh.cpp.

void mfem::NCMesh::BuildVertexList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 2757 of file ncmesh.cpp.

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

Definition at line 2004 of file ncmesh.cpp.

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

Definition at line 747 of file ncmesh.cpp.

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

Definition at line 722 of file ncmesh.cpp.

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 1748 of file ncmesh.cpp.

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 842 of file ncmesh.cpp.

void mfem::NCMesh::ClearTransforms ( )

Free all internal data created by the above three functions.

Definition at line 4127 of file ncmesh.cpp.

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

Definition at line 1700 of file ncmesh.cpp.

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

Definition at line 2901 of file ncmesh.cpp.

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

Definition at line 3352 of file ncmesh.cpp.

void mfem::NCMesh::CollectLeafElements ( int  elem,
int  state 
)
protected

Definition at line 1865 of file ncmesh.cpp.

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

Definition at line 2937 of file ncmesh.cpp.

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

Definition at line 2913 of file ncmesh.cpp.

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

Definition at line 4965 of file ncmesh.cpp.

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

Definition at line 4727 of file ncmesh.cpp.

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

Definition at line 5241 of file ncmesh.cpp.

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

Definition at line 5216 of file ncmesh.cpp.

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

Definition at line 391 of file ncmesh.cpp.

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 1777 of file ncmesh.cpp.

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

Definition at line 1574 of file ncmesh.cpp.

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

Definition at line 112 of file ncmesh.hpp.

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

Definition at line 4676 of file ncmesh.cpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 643 of file ncmesh.hpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 642 of file ncmesh.hpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 644 of file ncmesh.hpp.

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

Definition at line 2270 of file ncmesh.cpp.

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

Definition at line 2288 of file ncmesh.cpp.

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

Definition at line 2250 of file ncmesh.cpp.

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

Definition at line 675 of file ncmesh.cpp.

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

Definition at line 4579 of file ncmesh.cpp.

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

Definition at line 272 of file ncmesh.cpp.

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 3154 of file ncmesh.cpp.

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 2260 of file ncmesh.cpp.

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 3051 of file ncmesh.cpp.

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 3375 of file ncmesh.cpp.

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

Definition at line 617 of file ncmesh.cpp.

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

Definition at line 550 of file ncmesh.hpp.

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

Get a list of vertices (2D/3D) and edges (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 affected by non-local boundary elements.

Reimplemented in mfem::ParNCMesh.

Definition at line 4603 of file ncmesh.cpp.

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 1733 of file ncmesh.cpp.

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

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 3954 of file ncmesh.cpp.

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

Return the current list of conforming and nonconforming edges.

Definition at line 219 of file ncmesh.hpp.

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 4523 of file ncmesh.cpp.

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

Definition at line 4491 of file ncmesh.cpp.

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 4441 of file ncmesh.cpp.

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 4422 of file ncmesh.cpp.

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

Return the distance of leaf 'i' from the root.

Definition at line 4532 of file ncmesh.cpp.

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 4558 of file ncmesh.cpp.

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

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

Definition at line 317 of file ncmesh.hpp.

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

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

Definition at line 4544 of file ncmesh.cpp.

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

Definition at line 416 of file ncmesh.cpp.

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

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

Definition at line 321 of file ncmesh.hpp.

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

Return the current list of conforming and nonconforming faces.

Definition at line 212 of file ncmesh.hpp.

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 4453 of file ncmesh.cpp.

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

Definition at line 3422 of file ncmesh.cpp.

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

Definition at line 4810 of file ncmesh.cpp.

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

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

Definition at line 2033 of file ncmesh.cpp.

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

Definition at line 288 of file ncmesh.cpp.

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

Definition at line 295 of file ncmesh.cpp.

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

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

Definition at line 234 of file ncmesh.hpp.

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

Definition at line 116 of file ncmesh.hpp.

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

Definition at line 117 of file ncmesh.hpp.

virtual int mfem::NCMesh::GetNumGhostElements ( ) const
inlineprotectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 519 of file ncmesh.hpp.

virtual int mfem::NCMesh::GetNumGhostVertices ( ) const
inlineprotectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 520 of file ncmesh.hpp.

int mfem::NCMesh::GetNumRootElements ( )
inline

Return the number of root elements.

Definition at line 325 of file ncmesh.hpp.

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

Definition at line 115 of file ncmesh.hpp.

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

Definition at line 3437 of file ncmesh.cpp.

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

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 3903 of file ncmesh.cpp.

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 227 of file ncmesh.hpp.

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

Definition at line 3304 of file ncmesh.cpp.

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 4375 of file ncmesh.cpp.

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 4390 of file ncmesh.cpp.

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

Definition at line 524 of file ncmesh.hpp.

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

Definition at line 525 of file ncmesh.hpp.

void mfem::NCMesh::InitDerefTransforms ( )
protected

Definition at line 1813 of file ncmesh.cpp.

void mfem::NCMesh::InitGeomFlags ( )
protected

Definition at line 204 of file ncmesh.cpp.

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 1927 of file ncmesh.cpp.

virtual bool mfem::NCMesh::IsGhost ( const Element el) const
inlineprotectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 518 of file ncmesh.hpp.

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 4840 of file ncmesh.cpp.

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

I/O: Load the element refinement hierarchy from a mesh file.

Definition at line 4985 of file ncmesh.cpp.

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

I/O: Load the vertex parent hierarchy from a mesh file. NOTE: called indirectly through the constructor.

Definition at line 4882 of file ncmesh.cpp.

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 3857 of file ncmesh.cpp.

long mfem::NCMesh::MemoryUsage ( ) const

Return total number of bytes allocated.

Definition at line 5114 of file ncmesh.cpp.

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 3237 of file ncmesh.cpp.

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 453 of file ncmesh.cpp.

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

Definition at line 1990 of file ncmesh.cpp.

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 541 of file ncmesh.cpp.

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 515 of file ncmesh.cpp.

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

Definition at line 567 of file ncmesh.cpp.

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 483 of file ncmesh.cpp.

void mfem::NCMesh::OnMeshUpdated ( Mesh mesh)
protectedvirtual

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

Reimplemented in mfem::ParNCMesh.

Definition at line 2120 of file ncmesh.cpp.

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

I/O: Print the "coarse_elements" section of the mesh file (ver. >= 1.1).

Definition at line 4951 of file ncmesh.cpp.

int mfem::NCMesh::PrintElements ( std::ostream &  out,
int  elem,
int &  coarse_id 
) const
protected

Definition at line 4925 of file ncmesh.cpp.

int mfem::NCMesh::PrintMemoryDetail ( ) const

Definition at line 5136 of file ncmesh.cpp.

void mfem::NCMesh::PrintStats ( std::ostream &  out = mfem::out) const

Definition at line 5162 of file ncmesh.cpp.

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

I/O: Print the "vertex_parents" section of the mesh file (ver. >= 1.1).

Definition at line 4855 of file ncmesh.cpp.

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

Definition at line 4700 of file ncmesh.cpp.

int mfem::NCMesh::QuadFaceSplitType ( int  v1,
int  v2,
int  v3,
int  v4,
int  mid[5] = NULL 
) const
protected

Definition at line 2186 of file ncmesh.cpp.

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

Definition at line 303 of file ncmesh.cpp.

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 1492 of file ncmesh.cpp.

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

Definition at line 859 of file ncmesh.cpp.

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

Definition at line 377 of file ncmesh.cpp.

int mfem::NCMesh::ReorderFacePointMat ( int  v0,
int  v1,
int  v2,
int  v3,
int  elem,
DenseMatrix mat 
) const
protected

Definition at line 2305 of file ncmesh.cpp.

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

Definition at line 256 of file ncmesh.cpp.

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

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

Definition at line 1538 of file ncmesh.cpp.

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

Definition at line 1828 of file ncmesh.cpp.

void mfem::NCMesh::SetVertexPositions ( const Array< mfem::Vertex > &  vertices)

I/O: Set positions of all vertices (used by mesh loader).

Definition at line 4904 of file ncmesh.cpp.

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

Definition at line 113 of file ncmesh.hpp.

void mfem::NCMesh::TraverseEdge ( int  vn0,
int  vn1,
double  t0,
double  t1,
int  flags,
int  level 
)
protected

Definition at line 2626 of file ncmesh.cpp.

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

Definition at line 2339 of file ncmesh.cpp.

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

Definition at line 3871 of file ncmesh.cpp.

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

Definition at line 2451 of file ncmesh.cpp.

bool mfem::NCMesh::TraverseTriFace ( int  vn0,
int  vn1,
int  vn2,
const PointMatrix pm,
int  level 
)
protected

Definition at line 2487 of file ncmesh.cpp.

bool mfem::NCMesh::TriFaceSplit ( int  v1,
int  v2,
int  v3,
int  mid[3] = NULL 
) const
protected

Definition at line 2233 of file ncmesh.cpp.

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

Definition at line 4683 of file ncmesh.cpp.

void mfem::NCMesh::Trim ( )
virtual

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

Reimplemented in mfem::ParNCMesh.

Definition at line 5078 of file ncmesh.cpp.

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

Definition at line 334 of file ncmesh.cpp.

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 refinement and derefinement.

Reimplemented in mfem::ParNCMesh.

Definition at line 213 of file ncmesh.cpp.

void mfem::NCMesh::UpdateElementToVertexTable ( )
inlineprotected

Definition at line 682 of file ncmesh.hpp.

void mfem::NCMesh::UpdateLeafElements ( )
protected

Definition at line 1907 of file ncmesh.cpp.

void mfem::NCMesh::UpdateVertices ( )
protectedvirtual

update Vertex::index and vertex_nodeId

Reimplemented in mfem::ParNCMesh.

Definition at line 1847 of file ncmesh.cpp.

Friends And Related Function Documentation

friend struct CompareRanks
friend

Definition at line 856 of file ncmesh.hpp.

friend class Mesh
friend

Definition at line 370 of file ncmesh.hpp.

friend class ParNCMesh
friend

Definition at line 855 of file ncmesh.hpp.

Member Data Documentation

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

subset of all faces, set by BuildFaceList

Definition at line 499 of file ncmesh.hpp.

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

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

Definition at line 793 of file ncmesh.hpp.

Table mfem::NCMesh::derefinements
protected

possible derefinements, see GetDerefinementTable

Definition at line 534 of file ncmesh.hpp.

int mfem::NCMesh::Dim
protected

Definition at line 382 of file ncmesh.hpp.

NCList mfem::NCMesh::edge_list
protected

lazy-initialized list of edges, see GetEdgeList

Definition at line 496 of file ncmesh.hpp.

Table mfem::NCMesh::element_vertex
protected

leaf-element to vertex table, see FindSetNeighbors

Definition at line 502 of file ncmesh.hpp.

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

Definition at line 462 of file ncmesh.hpp.

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

face geometry by face index, set by OnMeshUpdated

Definition at line 500 of file ncmesh.hpp.

NCList mfem::NCMesh::face_list
protected

lazy-initialized list of faces, see GetFaceList

Definition at line 495 of file ncmesh.hpp.

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

Definition at line 460 of file ncmesh.hpp.

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

Definition at line 463 of file ncmesh.hpp.

int mfem::NCMesh::Geoms
protected

bit mask of element geometries present, see InitGeomFlags()

Definition at line 384 of file ncmesh.hpp.

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

Definition at line 847 of file ncmesh.hpp.

bool mfem::NCMesh::Iso
protected

true if the mesh only contains isotropic refinements

Definition at line 383 of file ncmesh.hpp.

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

Definition at line 492 of file ncmesh.hpp.

int mfem::NCMesh::NEdges
protected

Definition at line 490 of file ncmesh.hpp.

int mfem::NCMesh::NFaces
protected

Definition at line 490 of file ncmesh.hpp.

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

Definition at line 459 of file ncmesh.hpp.

int mfem::NCMesh::NVertices
protected

Definition at line 489 of file ncmesh.hpp.

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

Definition at line 777 of file ncmesh.hpp.

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

Definition at line 776 of file ncmesh.hpp.

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

Definition at line 774 of file ncmesh.hpp.

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

Definition at line 775 of file ncmesh.hpp.

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

Definition at line 773 of file ncmesh.hpp.

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

stack of scheduled refinements (temporary)

Definition at line 530 of file ncmesh.hpp.

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

scheduled node reparents (tmp)

Definition at line 532 of file ncmesh.hpp.

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 468 of file ncmesh.hpp.

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

temporary storage for reparented nodes

Definition at line 531 of file ncmesh.hpp.

int mfem::NCMesh::spaceDim
protected

dimensions of the elements and the vertex coordinates

Definition at line 382 of file ncmesh.hpp.

TmpVertex* mfem::NCMesh::tmp_vertex
mutableprotected

Definition at line 808 of file ncmesh.hpp.

Array<double> mfem::NCMesh::top_vertex_pos
protected

coordinates of top-level vertices (organized as triples)

Definition at line 471 of file ncmesh.hpp.

CoarseFineTransformations mfem::NCMesh::transforms
protected

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

Definition at line 790 of file ncmesh.hpp.

NCList mfem::NCMesh::vertex_list
protected

lazy-initialized list of vertices, see GetVertexList

Definition at line 497 of file ncmesh.hpp.

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

Definition at line 493 of file ncmesh.hpp.


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