15 #include "../config/config.hpp" 16 #include "../general/hash.hpp" 17 #include "../general/globals.hpp" 18 #include "../general/sort_pairs.hpp" 19 #include "../linalg/densemat.hpp" 22 #include "../fem/geom.hpp" 37 enum :
char {
X = 1,
Y = 2,
Z = 4,
XY = 3,
XZ = 5,
YZ = 6,
XYZ = 7 };
83 bool want_ghosts =
false)
const;
94 void Swap(CoarseFineTransformations &
a, CoarseFineTransformations &
b);
133 NCMesh(std::istream &input,
int version,
int &curved,
int &is_nc);
329 bool oriented =
true)
const;
340 int vert_index[4],
int edge_index[4],
341 int edge_orientation[4])
const;
381 void Print(std::ostream &
out)
const;
645 int n4,
int n5,
int n6,
int n7,
int attr,
646 int fattr0,
int fattr1,
int fattr2,
647 int fattr3,
int fattr4,
int fattr5);
649 int NewWedge(
int n0,
int n1,
int n2,
650 int n3,
int n4,
int n5,
int attr,
651 int fattr0,
int fattr1,
652 int fattr2,
int fattr3,
int fattr4);
655 int fattr0,
int fattr1,
int fattr2,
int fattr3);
657 int NewPyramid(
int n0,
int n1,
int n2,
int n3,
int n4,
int attr,
658 int fattr0,
int fattr1,
int fattr2,
int fattr3,
662 int eattr0,
int eattr1,
int eattr2,
int eattr3);
665 int attr,
int eattr0,
int eattr1,
int eattr2);
667 int NewSegment(
int n0,
int n1,
int attr,
int vattr1,
int vattr2);
674 bool TriFaceSplit(
int v1,
int v2,
int v3,
int mid[3] = NULL)
const;
685 int mid12,
int mid34,
int level = 0);
688 int en1,
int en2,
int en3,
int en4,
int midf);
724 int elem,
const PointMatrix &pm,
725 PointMatrix &reordered)
const;
728 const PointMatrix& pm,
int level, Face* eface[4],
731 const PointMatrix& pm,
int level,
735 void TraverseEdge(
int vn0,
int vn1,
double t0,
double t1,
int flags,
820 for (
int i = 0; i <
dim; i++)
830 for (
int i = 0; i <
dim; i++)
959 int& h_level,
int& v_level)
const;
NCList face_list
lazy-initialized list of faces, see GetFaceList
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
The PointMatrix stores the coordinates of the slave face using the master face coordinate as referenc...
void DebugDump(std::ostream &out) const
bool TriFaceSplit(int v1, int v2, int v3, int mid[3]=NULL) const
int PrintBoundary(std::ostream *out) const
void LoadLegacyFormat(std::istream &input, int &curved, int &is_nc)
Load the deprecated MFEM mesh v1.1 format for backward compatibility.
static void GridSfcOrdering3D(int width, int height, int depth, Array< int > &coords)
int faces[MaxElemFaces][4]
virtual int GetNGhostElements() const
bool IsLegacyLoaded() const
I/O: Return true if the mesh was loaded from the legacy v1.1 format.
void CheckAnisoPrism(int vn1, int vn2, int vn3, int vn4, const Refinement *refs, int nref)
void Print(std::ostream &out) const
I/O: Print the mesh in "MFEM NC mesh v1.0" format.
CoarseFineTransformations transforms
storage for data returned by Get[De]RefinementTransforms()
void DebugLeafOrder(std::ostream &out) const
static void GridSfcOrdering2D(int width, int height, Array< int > &coords)
int NewQuadrilateral(int n0, int n1, int n2, int n3, int attr, int eattr0, int eattr1, int eattr2, int eattr3)
static PointMatrix pm_seg_identity
void QuadFaceSplitLevel(int vn1, int vn2, int vn3, int vn4, int &h_level, int &v_level) const
int elem[2]
up to 2 elements sharing the face
char ref_type
refinement XYZ bit mask (7 = full isotropic)
bool ZeroRootStates() const
Return true if all root_states are zero.
int GetEdgeMaster(int v1, int v2) const
void SetDerefMatrixCodes(int parent, Array< int > &fine_coarse)
const CoarseFineTransformations & GetDerefinementTransforms()
int child[MaxElemChildren]
2-10 children (if ref_type != 0)
char flag
generic flag/marker, can be used by algorithms
void FindVertexCousins(int elem, int local, Array< int > &cousins) const
void FindEdgeElements(int vn1, int vn2, int vn3, int vn4, Array< MeshId > &prisms) const
Refinement(int index, int type=Refinement::XYZ)
void OrientedPointMatrix(const Slave &slave, DenseMatrix &oriented_matrix) const
Return the point matrix oriented according to the master and slave edges.
PointMatrix(const Point &p0, const Point &p1, const Point &p2)
This holds in one place the constants about the geometries we support.
int GetVertexRootCoord(int elem, RefCoord coord[3]) const
Array< Triple< int, int, int > > reparents
scheduled node reparents (tmp)
bool HavePyramids() const
Return true if the mesh contains pyramid elements.
Point(const Point &p0, const Point &p1)
void CollectTriFaceVertices(int v0, int v1, int v2, Array< int > &indices)
char tet_type
tetrahedron split type, currently always 0
int NewSegment(int n0, int n1, int attr, int vattr1, int vattr2)
virtual void LimitNCLevel(int max_nc_level)
virtual void Trim()
Save memory by releasing all non-essential and cached data.
int index
element number in the Mesh, -1 if refined
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
void LegacyToNewVertexOrdering(Array< int > &order) const
I/O: Return a map from old (v1.1) vertex indices to new vertex indices.
int Dimension() const
Return the dimension of the NCMesh.
void ForceRefinement(int vn1, int vn2, int vn3, int vn4)
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4, const Point &p5, const Point &p6, const Point &p7)
Lists all edges/faces in the nonconforming mesh.
long MemoryUsage() const
Return total number of bytes allocated.
bool Legacy
true if the mesh was loaded from the legacy v1.1 format
unsigned matrix
index into NCList::point_matrices[geom]
void ClearTransforms()
Free all internal data created by the above three functions.
int ReorderFacePointMat(int v0, int v1, int v2, int v3, int elem, const PointMatrix &pm, PointMatrix &reordered) const
Table derefinements
possible derefinements, see GetDerefinementTable
Data type dense matrix using column-major storage.
NCList vertex_list
lazy-initialized list of vertices, see GetVertexList
void InitDerefTransforms()
static PointMatrix pm_prism_identity
char geom
Geometry::Type of the element (char for storage only)
void InitRootState(int root_count)
int NewTetrahedron(int n0, int n1, int n2, int n3, int attr, int fattr0, int fattr1, int fattr2, int fattr3)
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4)
void CheckAnisoFace(int vn1, int vn2, int vn3, int vn4, int mid12, int mid34, int level=0)
static PointMatrix pm_tri_identity
Geometry::Type Geom() const
bool IsGhost(const Element &el) const
Return true if the Element el is a ghost element.
void CollectDerefinements(int elem, Array< Connection > &list)
Point & operator()(int i)
void GetLimitRefinements(Array< Refinement > &refinements, int max_level)
void GetMatrix(DenseMatrix &point_matrix) const
void TraverseEdge(int vn0, int vn1, double t0, double t1, int flags, int level, MatrixMap &matrix_map)
static const int MaxElemEdges
Number of edges of an element can have.
int GetNEdges() const
Return the number of edges in the NCMesh.
int spaceDim
dimensions of the elements and the vertex coordinates
Array< int > vertex_nodeId
vertex-index to node-id map, see UpdateVertices
int attribute
boundary element attribute, -1 if internal face
void FindFaceNodes(int face, int node[4])
static const PointMatrix & GetGeomIdentity(Geometry::Type geom)
bool HavePrisms() const
Return true if the mesh contains prism elements.
void DeleteAll()
Delete the whole array.
static int find_element_edge(const Element &el, int vn0, int vn1, bool abort=true)
int index
Mesh element number.
int master
master number (in Mesh numbering)
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
void OnMeshUpdated(Mesh *mesh)
A parallel extension of the NCMesh class.
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)
void CollectEdgeVertices(int v0, int v1, Array< int > &indices)
void CollectLeafElements(int elem, int state, Array< int > &ghosts, int &counter)
static int find_local_face(int geom, int a, int b, int c)
Array< Refinement > ref_stack
stack of scheduled refinements (temporary)
Point(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
void GetMeshComponents(Mesh &mesh) const
Fill Mesh::{vertices,elements,boundary} for the current finest level.
void ReferenceElement(int elem)
int GetElementDepth(int i) const
Return the distance of leaf 'i' from the root.
void GetPointMatrix(Geometry::Type geom, const char *ref_path, DenseMatrix &matrix)
void CheckIsoFace(int vn1, int vn2, int vn3, int vn4, int en1, int en2, int en3, int en4, int midf)
Face * GetFace(Element &elem, int face_no)
void TraverseRefinements(int elem, int coarse_index, std::string &ref_path, RefPathMap &map)
int GetMidEdgeNode(int node1, int node2)
Element(Geometry::Type geom, int attr)
Geometry::Type Geom() const
int index
face number in the Mesh
const Table & GetDerefinementTable()
Table element_vertex
leaf-element to vertex table, see FindSetNeighbors
int GetNFaces() const
Return the number of (2D) faces in the NCMesh.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
const double * CalcVertexPos(int node) const
int AddElement(const Element &el)
int SpaceDimension() const
Return the space dimension of the NCMesh.
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'.
virtual void Derefine(const Array< int > &derefs)
void FindNeighbors(int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
void InitGeom(Geometry::Type geom)
Array< int > coarse_elements
state of leaf_elements before Refine(), set by MarkCoarseLevel()
void UpdateVertices()
This method assigns indices to vertices (Node::vert_index) that will be seen by the Mesh class and th...
int FindMidEdgeNode(int node1, int node2) const
Master(int index, int element, int local, int geom, int sb, int se)
Array< DenseMatrix * > point_matrices[Geometry::NumGeom]
List of unique point matrices for each slave geometry.
void RegisterFaces(int elem, int *fattr=NULL)
static PointMatrix pm_pyramid_identity
virtual void BuildEdgeList()
int GetElementSizeReduction(int i) const
const CoarseFineTransformations & GetRefinementTransforms()
void LoadCoarseElements(std::istream &input)
Load the element refinement hierarchy from a legacy mesh file.
Identifies a vertex/edge/face in both Mesh and NCMesh.
int parent
parent element, -1 if this is a root element, -2 if free'd
signed char local
local number within 'element'
virtual void ElementSharesFace(int elem, int local, int face)
int GetNumRootElements()
Return the number of root elements.
bool HaveTets() const
Return true if the mesh contains tetrahedral elements.
void LoadCoordinates(std::istream &input)
Load the "coordinates" section of the mesh file.
virtual void ElementSharesVertex(int elem, int local, int vnode)
int MyRank
used in parallel, or when loading a parallel file in serial
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
int EdgeSplitLevel(int vn1, int vn2) const
int CountTopLevelNodes() const
Return the index of the last top-level node plus one.
friend struct PointMatrixHash
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
int element
NCMesh::Element containing this vertex/edge/face.
static GeomInfo GI[Geometry::NumGeom]
int slaves_end
slave faces
void ReparentNode(int node, int new_p1, int new_p2)
static const int MaxElemNodes
Number of nodes of an element can have.
int PrintMemoryDetail() const
int NewTriangle(int n0, int n1, int n2, int attr, int eattr0, int eattr1, int eattr2)
void DeleteUnusedFaces(const Array< int > &elemFaces)
Geometry::Type GetFaceGeometry(int index) const
Return face geometry type. index is the Mesh face number.
bool Iso
true if the mesh only contains isotropic refinements
Point(double x, double y, double z)
void Swap(Array< T > &, Array< T > &)
std::map< std::string, int > RefPathMap
int QuadFaceSplitType(int v1, int v2, int v3, int v4, int mid[5]=NULL) const
int GetFaceVerticesEdges(const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
unsigned ghost
For internal use: 0 if regular fine element, 1 if parallel ghost element.
Geometry::Type GetElementGeometry(int index) const
Return element geometry type. index is the Mesh element number.
virtual void BuildFaceList()
Nonconforming edge/face within a bigger edge/face.
void DerefineElement(int elem)
Derefine the element elem, does nothing on leaf elements.
int TriFaceSplitLevel(int vn1, int vn2, int vn3) const
const NCList & GetVertexList()
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void RegisterElement(int e)
NCMesh & operator=(NCMesh &)=delete
Copy assignment not supported.
Embedding(int elem, Geometry::Type geom, int matrix=0, bool ghost=false)
Array< int > leaf_sfc_index
natural tree ordering of leaf elements
const Point & operator()(int i) const
PointMatrix(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &p4, const Point &p5)
void DeleteLast()
Delete the last entry of the array.
void UpdateLeafElements()
Update the leaf elements indices in leaf_elements.
int PrintVertexParents(std::ostream *out) const
Print the "vertex_parents" section of the mesh file.
int RetrieveNode(const Element &el, int index)
Return el.node[index] correctly, even if the element is refined.
int Geoms
bit mask of element geometries present, see InitGeomFlags()
PointMatrix(const Point &p0, const Point &p1)
void BuildElementToVertexTable()
void ForgetElement(int e)
static PointMatrix pm_tet_identity
int GetMidFaceNode(int en1, int en2, int en3, int en4)
void RefineElement(int elem, char ref_type)
static const int MaxElemChildren
Number of children of an element can have.
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)
void TraverseTetEdge(int vn0, int vn1, const Point &p0, const Point &p1, MatrixMap &matrix_map)
Array< int > boundary_faces
subset of all faces, set by BuildFaceList
Array< double > coordinates
const NCList & GetNCList(int entity)
Return vertex/edge/face list (entity = 0/1/2, respectively).
int Size() const
Returns the number of TYPE I elements.
int edges[MaxElemEdges][2]
Array< int > leaf_elements
finest elements, in Mesh ordering (+ ghosts)
static PointMatrix pm_quad_identity
int node[MaxElemNodes]
element corners (if ref_type == 0)
Array< MeshId > conforming
unsigned edge_flags
orientation flags, see OrientedPointMatrix
signed char geom
Geometry::Type (faces only) (char to save RAM)
Array< int > free_element_ids
void TraverseQuadFace(int vn0, int vn1, int vn2, int vn3, const PointMatrix &pm, int level, Face *eface[4], MatrixMap &matrix_map)
int index(int i, int j, int nx, int ny)
bool operator==(const PointMatrix &pm) const
void GetElementFacesAttributes(int i, Array< int > &faces, Array< int > &fattr) const
Return the faces and face attributes of leaf element 'i'.
int parent
Coarse Element index in the coarse mesh.
Point(double x, double y)
static PointMatrix pm_hex_identity
const MeshId & LookUp(int index, int *type=NULL) const
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
int Size() const
Return the logical size of the array.
T & Last()
Return the last element in the array.
void NeighborExpand(const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
BlockArray< Element > elements
void FindSetNeighbors(const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
Slave(int index, int element, int local, int geom)
mfem::Element * NewMeshElement(int geom) const
void InitRootElements()
Count root elements and initialize root_state.
int GetEdgeNCOrientation(const MeshId &edge_id) const
HashTable< Node > shadow
temporary storage for reparented nodes
virtual void BuildVertexList()
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
MeshId(int index, int element, int local, int geom=-1)
void UpdateElementToVertexTable()
void LoadBoundary(std::istream &input)
Load the "boundary" section of the mesh file.
int NewPyramid(int n0, int n1, int n2, int n3, int n4, int attr, int fattr0, int fattr1, int fattr2, int fattr3, int fattr4)
void PrintCoordinates(std::ostream &out) const
Print the "coordinates" section of the mesh file.
void CopyElements(int elem, const BlockArray< Element > &tmp_elements)
virtual void ElementSharesEdge(int elem, int local, int enode)
int FindNodeExt(const Element &el, int node, bool abort=true)
Extended version of find_node: works if 'el' is refined.
int GetNVertices() const
Return the number of vertices in the NCMesh.
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
static const int MaxElemFaces
Number of faces of an element can have.
void CollectQuadFaceVertices(int v0, int v1, int v2, int v3, Array< int > &indices)
static int find_node(const Element &el, int node)
void UnreferenceElement(int elem, Array< int > &elemFaces)
void LoadVertexParents(std::istream &input)
Load the vertex parent hierarchy from a mesh file.
Point points[MaxElemNodes]
Rank 3 tensor (array of matrices)
Abstract data type element.
void CollectIncidentElements(int elem, const RefCoord coord[3], Array< int > &list) const
void CountSplits(int elem, int splits[3]) const
Point & operator=(const Point &src)
int GetSingleElement() const
Return one of elem[0] or elem[1] and make sure the other is -1.
int rank
processor number (ParNCMesh), -1 if undefined/unknown
bool TraverseTriFace(int vn0, int vn1, int vn2, const PointMatrix &pm, int level, MatrixMap &matrix_map)
Defines the position of a fine element within a coarse element.
virtual void Refine(const Array< Refinement > &refinements)
Array< char > face_geom
face geometry by face index, set by OnMeshUpdated