MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Classes | 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, 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 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]) const
 Return Mesh vertex indices of an edge identified by 'edge_id'. More...
 
int GetEdgeNCOrientation (const MeshId &edge_id) const
 
void 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'. More...
 
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)
 
int GetElementGeometry () const
 Return the type of elements in the mesh. More...
 
int GetFaceGeometry () const
 
int GetElementDepth (int i) const
 Return the distance of leaf 'i' from the root. 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 DebugDump (std::ostream &out) const
 

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 (Array< mfem::Vertex > &mvertices, Array< mfem::Element * > &melements, Array< mfem::Element * > &mboundary) const
 Return the basic Mesh arrays 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 ()
 
virtual bool IsGhost (const Element &el) const
 
virtual int GetNumGhostElements () const
 
virtual int GetNumGhostVertices () 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 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 GetMidEdgeNode (int vn1, int vn2)
 
int GetMidFaceNode (int en1, int en2, int en3, int en4)
 
int FaceSplitType (int v1, int v2, int v3, int v4, int mid[4]=NULL) const
 
void ForceRefinement (int vn1, int vn2, int vn3, int vn4)
 
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 RefElement (int elem)
 
void UnrefElement (int elem, Array< int > &elemFaces)
 
FaceGetFace (Element &elem, int face_no)
 
void RegisterFaces (int elem, int *fattr=NULL)
 
void DeleteUnusedFaces (const Array< int > &elemFaces)
 
int FindAltParents (int node1, int node2)
 
bool NodeSetX1 (int node, int *n)
 
bool NodeSetX2 (int node, int *n)
 
bool NodeSetY1 (int node, int *n)
 
bool NodeSetY2 (int node, int *n)
 
bool NodeSetZ1 (int node, int *n)
 
bool NodeSetZ2 (int node, int *n)
 
void CollectDerefinements (int elem, Array< Connection > &list)
 
int ReorderFacePointMat (int v0, int v1, int v2, int v3, int elem, DenseMatrix &mat) const
 
void TraverseFace (int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level)
 
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 face)
 
virtual void ElementSharesEdge (int elem, int enode)
 
virtual void ElementSharesVertex (int elem, 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 CollectFaceVertices (int v0, int v1, int v2, int v3, Array< int > &indices)
 
void BuildElementToVertexTable ()
 
void UpdateElementToVertexTable ()
 
void GetPointMatrix (int 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
 
void FaceSplitLevel (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)
 
static int find_hex_face (int a, int b, int c)
 
static const PointMatrixGetGeomIdentity (int 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...
 
HashTable< Nodenodes
 
HashTable< Facefaces
 
BlockArray< Elementelements
 
Array< int > free_element_ids
 
int root_count
 
Array< double > top_vertex_pos
 
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...
 
Table element_vertex
 leaf-element to vertex table, see FindSetNeighbors More...
 
Array< Refinementref_stack
 stack of scheduled refinements (temporary) 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_hex_identity
 
static GeomInfo GI [Geometry::NumGeom]
 
static GeomInfogi_hex = NCMesh::GI[Geometry::CUBE]
 
static GeomInfogi_quad = NCMesh::GI[Geometry::SQUARE]
 
static GeomInfogi_tri = NCMesh::GI[Geometry::TRIANGLE]
 

Friends

class Mesh
 
class ParNCMesh
 
struct CompareRanks
 

Detailed Description

A class for non-conforming AMR on higher-order hexahedral, 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 83 of file ncmesh.hpp.

Member Typedef Documentation

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

Definition at line 412 of file ncmesh.hpp.

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

Definition at line 411 of file ncmesh.hpp.

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

Definition at line 409 of file ncmesh.hpp.

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

Definition at line 410 of file ncmesh.hpp.

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

Definition at line 408 of file ncmesh.hpp.

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

Definition at line 666 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 67 of file ncmesh.cpp.

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

Definition at line 174 of file ncmesh.cpp.

mfem::NCMesh::~NCMesh ( )
virtual

Definition at line 200 of file ncmesh.cpp.

Member Function Documentation

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

Definition at line 460 of file ncmesh.hpp.

void mfem::NCMesh::AssignLeafIndices ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 1478 of file ncmesh.cpp.

void mfem::NCMesh::BuildEdgeList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 1913 of file ncmesh.cpp.

void mfem::NCMesh::BuildElementToVertexTable ( )
protected

Definition at line 2181 of file ncmesh.cpp.

void mfem::NCMesh::BuildFaceList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 1803 of file ncmesh.cpp.

void mfem::NCMesh::BuildVertexList ( )
protectedvirtual

Reimplemented in mfem::ParNCMesh.

Definition at line 2013 of file ncmesh.cpp.

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

Definition at line 1499 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 583 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 1273 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 649 of file ncmesh.cpp.

void mfem::NCMesh::ClearTransforms ( )

Free all internal data created by the above three functions.

Definition at line 2958 of file ncmesh.cpp.

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

Definition at line 1225 of file ncmesh.cpp.

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

Definition at line 2145 of file ncmesh.cpp.

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

Definition at line 2157 of file ncmesh.cpp.

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

Definition at line 1423 of file ncmesh.cpp.

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

Definition at line 3404 of file ncmesh.cpp.

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

Definition at line 3199 of file ncmesh.cpp.

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

Definition at line 3644 of file ncmesh.cpp.

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

Definition at line 319 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 1302 of file ncmesh.cpp.

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

Definition at line 1119 of file ncmesh.cpp.

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

Definition at line 95 of file ncmesh.hpp.

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

Definition at line 3159 of file ncmesh.cpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 545 of file ncmesh.hpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 544 of file ncmesh.hpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 546 of file ncmesh.hpp.

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

Definition at line 3166 of file ncmesh.cpp.

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

Definition at line 1652 of file ncmesh.cpp.

int mfem::NCMesh::find_element_edge ( const Element el,
int  vn0,
int  vn1 
)
staticprotected

Definition at line 1689 of file ncmesh.cpp.

int mfem::NCMesh::find_hex_face ( int  a,
int  b,
int  c 
)
staticprotected

Definition at line 1705 of file ncmesh.cpp.

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

Definition at line 1679 of file ncmesh.cpp.

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

Definition at line 366 of file ncmesh.cpp.

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

Definition at line 3077 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 2356 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 2253 of file ncmesh.cpp.

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

Definition at line 549 of file ncmesh.cpp.

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

Definition at line 471 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 3100 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 1258 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 2915 of file ncmesh.cpp.

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

Return the current list of conforming and nonconforming edges.

Definition at line 195 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 3056 of file ncmesh.cpp.

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

Definition at line 3024 of file ncmesh.cpp.

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

Return "NC" orientation of an edge. As opposed 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.

Definition at line 2986 of file ncmesh.cpp.

void mfem::NCMesh::GetEdgeVertices ( const MeshId edge_id,
int  vert_index[2] 
) const

Return Mesh vertex indices of an edge identified by 'edge_id'.

Definition at line 2968 of file ncmesh.cpp.

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

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

Definition at line 3065 of file ncmesh.cpp.

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

Return the type of elements in the mesh.

Definition at line 273 of file ncmesh.hpp.

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

Definition at line 344 of file ncmesh.cpp.

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

Definition at line 275 of file ncmesh.hpp.

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

Return the current list of conforming and nonconforming faces.

Definition at line 188 of file ncmesh.hpp.

void 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'.

Definition at line 2998 of file ncmesh.cpp.

const NCMesh::PointMatrix & mfem::NCMesh::GetGeomIdentity ( int  geom)
staticprotected

Definition at line 2524 of file ncmesh.cpp.

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

Definition at line 3247 of file ncmesh.cpp.

void mfem::NCMesh::GetMeshComponents ( Array< mfem::Vertex > &  mvertices,
Array< mfem::Element * > &  melements,
Array< mfem::Element * > &  mboundary 
) const
protected

Return the basic Mesh arrays for the current finest level.

Definition at line 1528 of file ncmesh.cpp.

int mfem::NCMesh::GetMidEdgeNode ( int  vn1,
int  vn2 
)
protected

Definition at line 513 of file ncmesh.cpp.

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

Definition at line 521 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 210 of file ncmesh.hpp.

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

Definition at line 99 of file ncmesh.hpp.

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

Definition at line 100 of file ncmesh.hpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 447 of file ncmesh.hpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 448 of file ncmesh.hpp.

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

Definition at line 98 of file ncmesh.hpp.

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

Definition at line 2537 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 2878 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 203 of file ncmesh.hpp.

void mfem::NCMesh::InitDerefTransforms ( )
protected

Definition at line 1338 of file ncmesh.cpp.

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

Reimplemented in mfem::ParNCMesh.

Definition at line 446 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 3277 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 3424 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 3319 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 2834 of file ncmesh.cpp.

long mfem::NCMesh::MemoryUsage ( ) const

Return total number of bytes allocated.

Definition at line 3545 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 2439 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 435 of file ncmesh.cpp.

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

Definition at line 1487 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 464 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 489 of file ncmesh.cpp.

bool mfem::NCMesh::NodeSetX1 ( int  node,
int *  n 
)
inlineprotected

Definition at line 530 of file ncmesh.cpp.

bool mfem::NCMesh::NodeSetX2 ( int  node,
int *  n 
)
inlineprotected

Definition at line 533 of file ncmesh.cpp.

bool mfem::NCMesh::NodeSetY1 ( int  node,
int *  n 
)
inlineprotected

Definition at line 536 of file ncmesh.cpp.

bool mfem::NCMesh::NodeSetY2 ( int  node,
int *  n 
)
inlineprotected

Definition at line 539 of file ncmesh.cpp.

bool mfem::NCMesh::NodeSetZ1 ( int  node,
int *  n 
)
inlineprotected

Definition at line 542 of file ncmesh.cpp.

bool mfem::NCMesh::NodeSetZ2 ( int  node,
int *  n 
)
inlineprotected

Definition at line 545 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 1603 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 3390 of file ncmesh.cpp.

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

Definition at line 3364 of file ncmesh.cpp.

int mfem::NCMesh::PrintMemoryDetail ( ) const

Definition at line 3566 of file ncmesh.cpp.

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

Definition at line 3591 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 3292 of file ncmesh.cpp.

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

Definition at line 231 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 1075 of file ncmesh.cpp.

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

Definition at line 666 of file ncmesh.cpp.

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

Definition at line 305 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 1721 of file ncmesh.cpp.

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

Definition at line 1353 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 3341 of file ncmesh.cpp.

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

Definition at line 96 of file ncmesh.hpp.

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

Definition at line 1874 of file ncmesh.cpp.

void mfem::NCMesh::TraverseFace ( int  vn0,
int  vn1,
int  vn2,
int  vn3,
const PointMatrix pm,
int  level 
)
protected

Definition at line 1754 of file ncmesh.cpp.

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

Definition at line 2848 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 3514 of file ncmesh.cpp.

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

Definition at line 262 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 188 of file ncmesh.cpp.

void mfem::NCMesh::UpdateElementToVertexTable ( )
inlineprotected

Definition at line 583 of file ncmesh.hpp.

void mfem::NCMesh::UpdateLeafElements ( )
protected

Definition at line 1464 of file ncmesh.cpp.

void mfem::NCMesh::UpdateVertices ( )
protectedvirtual

update Vertex::index and vertex_nodeId

Reimplemented in mfem::ParNCMesh.

Definition at line 1372 of file ncmesh.cpp.

Friends And Related Function Documentation

friend struct CompareRanks
friend

Definition at line 737 of file ncmesh.hpp.

friend class Mesh
friend

Definition at line 309 of file ncmesh.hpp.

friend class ParNCMesh
friend

Definition at line 736 of file ncmesh.hpp.

Member Data Documentation

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

subset of all faces, set by BuildFaceList

Definition at line 434 of file ncmesh.hpp.

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

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

Definition at line 675 of file ncmesh.hpp.

Table mfem::NCMesh::derefinements
protected

possible derefinements, see GetDerefinementTable

Definition at line 455 of file ncmesh.hpp.

int mfem::NCMesh::Dim
protected

Definition at line 323 of file ncmesh.hpp.

NCList mfem::NCMesh::edge_list
protected

lazy-initialized list of edges, see GetEdgeList

Definition at line 431 of file ncmesh.hpp.

Table mfem::NCMesh::element_vertex
protected

leaf-element to vertex table, see FindSetNeighbors

Definition at line 436 of file ncmesh.hpp.

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

Definition at line 399 of file ncmesh.hpp.

NCList mfem::NCMesh::face_list
protected

lazy-initialized list of faces, see GetFaceList

Definition at line 430 of file ncmesh.hpp.

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

Definition at line 397 of file ncmesh.hpp.

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

Definition at line 400 of file ncmesh.hpp.

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

Definition at line 727 of file ncmesh.hpp.

NCMesh::GeomInfo & mfem::NCMesh::gi_hex = NCMesh::GI[Geometry::CUBE]
staticprotected

Definition at line 729 of file ncmesh.hpp.

NCMesh::GeomInfo & mfem::NCMesh::gi_quad = NCMesh::GI[Geometry::SQUARE]
staticprotected

Definition at line 729 of file ncmesh.hpp.

NCMesh::GeomInfo & mfem::NCMesh::gi_tri = NCMesh::GI[Geometry::TRIANGLE]
staticprotected

Definition at line 729 of file ncmesh.hpp.

bool mfem::NCMesh::Iso
protected

true if the mesh only contains isotropic refinements

Definition at line 324 of file ncmesh.hpp.

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

Definition at line 427 of file ncmesh.hpp.

int mfem::NCMesh::NEdges
protected

Definition at line 425 of file ncmesh.hpp.

int mfem::NCMesh::NFaces
protected

Definition at line 425 of file ncmesh.hpp.

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

Definition at line 396 of file ncmesh.hpp.

int mfem::NCMesh::NVertices
protected

Definition at line 424 of file ncmesh.hpp.

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

Definition at line 660 of file ncmesh.hpp.

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

Definition at line 659 of file ncmesh.hpp.

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

Definition at line 658 of file ncmesh.hpp.

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

stack of scheduled refinements (temporary)

Definition at line 453 of file ncmesh.hpp.

int mfem::NCMesh::root_count
protected

Definition at line 403 of file ncmesh.hpp.

int mfem::NCMesh::spaceDim
protected

dimensions of the elements and the vertex coordinates

Definition at line 323 of file ncmesh.hpp.

TmpVertex* mfem::NCMesh::tmp_vertex
mutableprotected

Definition at line 690 of file ncmesh.hpp.

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

Definition at line 406 of file ncmesh.hpp.

CoarseFineTransformations mfem::NCMesh::transforms
protected

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

Definition at line 672 of file ncmesh.hpp.

NCList mfem::NCMesh::vertex_list
protected

lazy-initialized list of vertices, see GetVertexList

Definition at line 432 of file ncmesh.hpp.

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

Definition at line 428 of file ncmesh.hpp.


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