15 #include "../config/config.hpp"
16 #include "../general/stable3d.hpp"
17 #include "../general/globals.hpp"
22 #include "../fem/eltrans.hpp"
23 #include "../fem/coefficient.hpp"
24 #include "../general/gzstream.hpp"
34 class FiniteElementSpace;
155 #ifdef MFEM_USE_MEMALLOC
196 void ReadMFEMMesh(std::istream &input,
bool mfem_v11,
int &curved);
201 void ReadVTKMesh(std::istream &input,
int &curved,
int &read_gf);
202 void ReadNURBSMesh(std::istream &input,
int &curved,
int &read_gf);
206 #ifdef MFEM_USE_NETCDF
207 void ReadCubit(
const char *filename,
int &curved,
int &read_gf);
234 int *edge1,
int *edge2,
int *middle)
240 int *edge1,
int *edge2,
int *middle)
241 {
Bisection(i, v_to_v, edge1, edge2, middle); }
280 double threshold,
int nc_limit = 0,
342 int v0,
int v1,
int v2);
345 int v0,
int v1,
int v2,
int v3);
355 inline static void ShiftL2R(
int &,
int &,
int &);
357 inline static void Rotate3(
int &,
int &,
int &);
365 void InitMesh(
int _Dim,
int _spaceDim,
int NVert,
int NElem,
int NBdrElem);
372 void Loader(std::istream &input,
int generate_edges = 0,
373 std::string parse_tag =
"");
379 std::string section_delimiter =
"")
const;
386 double sx,
double sy,
double sz);
393 double sx,
double sy);
396 void Make1D(
int n,
double sx = 1.0);
406 void Swap(
Mesh& other,
bool non_geometry =
false);
420 explicit Mesh(
const Mesh &mesh,
bool copy_nodes =
true);
434 int *element_attributes,
int num_elements,
436 int *boundary_attributes,
int num_boundary_elements,
437 int dimension,
int space_dimension= -1);
441 Mesh(
int _Dim,
int NVert,
int NElem,
int NBdrElem = 0,
int _spaceDim = -1)
447 InitMesh(_Dim, _spaceDim, NVert, NElem, NBdrElem);
459 void AddTri(
const int *vi,
int attr = 1);
461 void AddQuad(
const int *vi,
int attr = 1);
462 void AddTet(
const int *vi,
int attr = 1);
463 void AddHex(
const int *vi,
int attr = 1);
476 bool fix_orientation =
true);
479 bool fix_orientation =
true);
482 bool fix_orientation =
true);
485 bool fix_orientation =
true);
514 void Finalize(
bool refine =
false,
bool fix_orientation =
false);
518 #ifdef MFEM_USE_GECKO
536 double sx = 1.0,
double sy = 1.0,
double sz = 1.0)
538 Make3D(nx, ny, nz, type, generate_edges, sx, sy, sz);
546 double sx = 1.0,
double sy = 1.0)
548 Make2D(nx, ny, type, generate_edges, sx, sy);
552 explicit Mesh(
int n,
double sx = 1.0)
560 Mesh(
const char *filename,
int generate_edges = 0,
int refine = 1,
561 bool fix_orientation =
true);
566 Mesh(std::istream &input,
int generate_edges = 0,
int refine = 1,
567 bool fix_orientation =
true);
570 Mesh(
Mesh *mesh_array[],
int num_pieces);
583 Mesh(
Mesh *orig_mesh,
int ref_factor,
int ref_type);
591 virtual void Load(std::istream &input,
int generate_edges = 0,
592 int refine = 1,
bool fix_orientation =
true)
594 Loader(input, generate_edges);
629 virtual long ReduceInt(
int value)
const {
return value; }
667 bool zerocopy =
false);
717 faces[i]->GetVertices(vert);
824 if (
faces_info[FaceNo].Elem2No < 0) {
return NULL; }
896 void GetNode(
int i,
double *coord);
897 void SetNode(
int i,
const double *coord);
936 void SetCurvature(
int order,
bool discont =
false,
int space_dim = -1,
951 int nonconforming = -1,
int nc_limit = 0);
956 int nonconforming = -1,
int nc_limit = 0);
960 int nonconforming = -1,
int nc_limit = 0);
964 double eps = 0.0,
int nonconforming = -1);
969 int nonconforming = -1,
int nc_limit = 0);
974 int nonconforming = -1,
int nc_limit = 0);
982 int nc_limit = 0,
int op = 1);
986 int nc_limit = 0,
int op = 1);
996 void EnsureNCMesh(
bool triangles_nonconforming =
false);
1030 void PrintVTK(std::ostream &
out,
int ref,
int field_data=0);
1039 std::ostream &
out,
int elem_attr = 0)
const;
1043 int interior_faces = 0);
1078 double &kappa_min,
double &kappa_max,
1126 std::ostream &
operator<<(std::ostream &
out,
const Mesh &mesh);
1147 Mesh *
Extrude1D(Mesh *mesh,
const int ny,
const double sy,
1148 const bool closed =
false);
1166 a = c; c = b; b = t;
Abstract class for Finite Elements.
void AddHex(const int *vi, int attr=1)
void Loader(std::istream &input, int generate_edges=0, std::string parse_tag="")
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
void GetFaceEdges(int i, Array< int > &, Array< int > &) const
void PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
void GetPointMatrix(int i, DenseMatrix &pointmat) const
static const int vtk_quadratic_hex[27]
virtual void Print(std::ostream &out=mfem::out) const
int * CartesianPartitioning(int nxyz[])
const DenseMatrix * PointMatrix
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
const double * GetVertex(int i) const
Return pointer to vertex i's coordinates.
void ScaleElements(double sf)
Class for grid function - Vector with associated FE space.
Table * GetEdgeVertexTable() const
Returns the edge-to-vertex Table (3D)
void GetBdrElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of bdr element i.
void FreeElement(Element *E)
void DerefineMesh(const Array< int > &derefinements)
Derefine elements once a list of derefinements is known.
void GetFaceInfos(int Face, int *Inf1, int *Inf2)
int CheckElementOrientation(bool fix_it=true)
Check the orientation of the elements.
void SetVertices(const Vector &vert_coord)
static const int vtk_quadratic_tet[10]
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
bool FaceIsTrueInterior(int FaceNo) const
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
long GetGlobalNE() const
Return the total (global) number of elements.
void AddHexAsTets(const int *vi, int attr=1)
void GetFaceElements(int Face, int *Elem1, int *Elem2)
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Array< Element * > boundary
Geometry::Constants< Geometry::SQUARE > quad_t
int * GeneratePartitioning(int nparts, int part_method=1)
CoarseFineTransformations CoarseFineTr
void MoveVertices(const Vector &displacements)
int GetNBE() const
Returns number of boundary elements.
IsoparametricTransformation Transformation
void GetLocalPtToSegTransformation(IsoparametricTransformation &, int)
Used in GetFaceElementTransformations (...)
void GetCharacteristics(double &h_min, double &h_max, double &kappa_min, double &kappa_max, Vector *Vh=NULL, Vector *Vk=NULL)
void PrintWithPartitioning(int *partitioning, std::ostream &out, int elem_attr=0) const
int GetFaceGeometryType(int Face) const
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
void GetElementFaces(int i, Array< int > &, Array< int > &) const
Return the indices and the orientations of all faces of element i.
void GetLocalFaceTransformation(int face_type, int elem_type, IsoparametricTransformation &Transf, int inf)
Used in GetFaceElementTransformations (...)
void AddBdrElement(Element *elem)
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
void ReadNetgen2DMesh(std::istream &input, int &curved)
void DegreeElevate(int t)
int EulerNumber2D() const
Equals 1 - num_holes.
void AddTri(const int *vi, int attr=1)
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Data type dense matrix using column-major storage.
Mesh * Extrude1D(Mesh *mesh, const int ny, const double sy, const bool closed)
Extrude a 1D mesh.
void Transform(void(*f)(const Vector &, Vector &))
int GetBdrElementEdgeIndex(int i) const
int GetNE() const
Returns number of elements.
const Element * GetFace(int i) const
Element * GetElement(int i)
void GenerateBoundaryElements()
void SetMeshGen()
Determine the mesh generator bitmask meshgen, see MeshGenerator().
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
void GetVertices(Vector &vert_coord) const
void ReadNetgen3DMesh(std::istream &input)
void AddBdrSegment(const int *vi, int attr=1)
int GetFaceElementType(int Face) const
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
void RemoveInternalBoundaries()
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
virtual long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
void GetVertexToVertexTable(DSTable &) const
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
void MarkTriMeshForRefinement()
void ReadCubit(const char *filename, int &curved, int &read_gf)
void InitRefinementTransforms()
IsoparametricTransformation FaceTransformation
Array< NCFaceInfo > nc_faces_info
Element * ReadElement(std::istream &)
void KnotInsert(Array< KnotVector * > &kv)
A parallel extension of the NCMesh class.
Element * NewElement(int geom)
FaceElementTransformations FaceElemTr
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
const GridFunction * GetNodes() const
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
void AddVertex(const double *)
Mesh(int _Dim, int NVert, int NElem, int NBdrElem=0, int _spaceDim=-1)
Init constructor: begin the construction of a Mesh object.
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Geometry::Constants< Geometry::SEGMENT > seg_t
void AddElement(Element *elem)
void GetElementJacobian(int i, DenseMatrix &J)
void MesquiteSmooth(const int mesquite_option=0)
bool Nonconforming() const
Mesh(int nx, int ny, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0)
void MoveNodes(const Vector &displacements)
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
void GetElementData(int geom, Array< int > &elem_vtx, Array< int > &attr) const
NCFaceInfo(bool slave, int master, const DenseMatrix *pm)
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a triangular Mesh.
void FinalizeTetMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a tetrahedral Mesh.
void AddBdrQuadAsTriangles(const int *vi, int attr=1)
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
int FindCoarseElement(int i)
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
static FiniteElement * GetTransformationFEforElementType(int)
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
STable3D * GetElementToFaceTable(int ret_ftbl=0)
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
static void Rotate3(int &, int &, int &)
void Make2D(int nx, int ny, Element::Type type, int generate_edges, double sx, double sy)
int MeshGenerator()
Get the mesh generator/type.
virtual void PrintInfo(std::ostream &out=mfem::out)
Type
Constants for the classes derived from Element.
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
void CheckDisplacements(const Vector &displacements, double &tmax)
void ReadInlineMesh(std::istream &input, int generate_edges=0)
Geometry::Constants< Geometry::TETRAHEDRON > tet_t
void ReadLineMesh(std::istream &input)
virtual void Load(std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
void SetNodesOwner(bool nodes_owner)
Set the mesh nodes ownership flag.
void SetLayer(const int l)
int GetElementToEdgeTable(Table &, Array< int > &)
void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const
For the given boundary element, bdr_el, return its adjacent element and its info, i...
void RandomRefinement(double prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
Refine each element with given probability. Uses GeneralRefinement.
int GetElementType(int i) const
Returns the type of element i.
const Element * GetElement(int i) const
Mesh(int n, double sx=1.0)
void GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
IsoparametricTransformation Transformation2
void SetNodalFESpace(FiniteElementSpace *nfes)
Element * GetBdrElement(int i)
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
virtual void ReorientTetMesh()
static void ShiftL2R(int &, int &, int &)
virtual void LocalRefinement(const Array< int > &marked_el, int type=3)
This function is not public anymore. Use GeneralRefinement instead.
int GetBdrElementType(int i) const
Returns the type of boundary element i.
void GetGeckoElementReordering(Array< int > &ordering)
Data type tetrahedron element.
virtual void MarkTetMeshForRefinement(DSTable &v_to_v)
double GetElementSize(int i, int type=0)
int CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
int SpaceDimension() const
void GetBdrElementFace(int i, int *, int *) const
Return the index and the orientation of the face of bdr element i. (3D)
void CheckPartitioning(int *partitioning)
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
void ReadMFEMMesh(std::istream &input, bool mfem_v11, int &curved)
void AddPointFaceElement(int lf, int gf, int el)
Used in GenerateFaces()
void Make1D(int n, double sx=1.0)
Creates a 1D mesh for the interval [0,sx] divided into n equal intervals.
const Element *const * GetElementsArray() const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void FinalizeTopology()
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
void PrintVTK(std::ostream &out)
virtual ~Mesh()
Destroys Mesh.
Geometry::Constants< Geometry::CUBE > hex_t
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
void ApplyLocalSlaveTransformation(IsoparametricTransformation &transf, const FaceInfo &fi)
static void PrintElement(const Element *, std::ostream &)
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
static void GetElementArrayEdgeTable(const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
bool IsSlaveFace(const FaceInfo &fi) const
void AddQuad(const int *vi, int attr=1)
void ReadGmshMesh(std::istream &input)
virtual void PrintXG(std::ostream &out=mfem::out) const
Print the mesh to the given stream using Netgen/Truegrid format.
virtual void QuadUniformRefinement()
Refine quadrilateral mesh.
virtual void HexUniformRefinement()
Refine hexahedral mesh.
void RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
void Swap(Mesh &other, bool non_geometry=false)
Array< Element * > elements
const Table & ElementToElementTable()
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
const CoarseFineTransformations & GetRefinementTransforms()
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Operation GetLastOperation() const
Return type of last modification of the mesh.
int GetElementBaseGeometry(int i=0) const
void AddBdrTriangle(const int *vi, int attr=1)
Table * GetFaceToElementTable() const
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
int EulerNumber() const
Equals 1 + num_holes - num_loops.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
void AddBdrQuad(const int *vi, int attr=1)
const Table & ElementToFaceTable() const
Class for integration point with weight.
Element * ReadElementWithoutAttr(std::istream &)
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
void AddTriangle(const int *vi, int attr=1)
bool OwnsNodes() const
Return the mesh nodes ownership flag.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
const FiniteElementSpace * GetNodalFESpace() const
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
bool DerefineByError(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NodeExtrudeCoefficient(const int dim, const int _n, const double _s)
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf)
virtual ~NodeExtrudeCoefficient()
Array< FaceInfo > faces_info
STable3D * GetFacesTable()
void Make3D(int nx, int ny, int nz, Element::Type type, int generate_edges, double sx, double sy, double sz)
void EnsureNCMesh(bool triangles_nonconforming=false)
NCMesh * ncmesh
Optional non-conforming mesh extension.
static void PrintElementWithoutAttr(const Element *, std::ostream &)
virtual bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NC version of GeneralDerefinement.
void Clear()
Clear the contents of the Mesh.
int GetNEdges() const
Return the number of edges.
Geometry::Constants< Geometry::TRIANGLE > tri_t
double * GetVertex(int i)
Return pointer to vertex i's coordinates.
int GetNFaces() const
Return the number of faces in a 3D mesh.
virtual int FindPoints(DenseMatrix &point_mat, Array< int > &elem_ids, Array< IntegrationPoint > &ips, bool warn=true, InverseElementTransformation *inv_trans=NULL)
Find the ids of the elements that contain the given points, and their corresponding reference coordin...
void SetAttribute(int i, int attr)
Set the attribute of element i.
void ReadTrueGridMesh(std::istream &input)
int GetFaceBaseGeometry(int i) const
void AverageVertices(int *indexes, int n, int result)
Mesh(int nx, int ny, int nz, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0, double sz=1.0)
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
void SetNode(int i, const double *coord)
void GetNode(int i, double *coord)
void GetElementColoring(Array< int > &colors, int el0=0)
void ChangeVertexDataOwnership(double *vertices, int len_vertices, bool zerocopy=false)
Set the internal Vertex array to point to the given vertices array without assuming ownership of the ...
IsoparametricTransformation EdgeTransformation
void AddTet(const int *vi, int attr=1)
MemAlloc< Tetrahedron, 1024 > TetMemory
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
void SetNodes(const Vector &node_coord)
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
named_ifgzstream(const char *mesh_name)
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
void ScaleSubdomains(double sf)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
void Bisection(int i, const DSTable &, int *, int *, int *)
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Class for parallel meshes.
Abstract data type element.
int GetAttribute(int i) const
Return the attribute of element i.
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
const Table & ElementToEdgeTable() const
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
void GenerateNCFaceInfo()
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
int GetBdrElementBaseGeometry(int i=0) const
double GetElementVolume(int i)
Array< int > attributes
A list of all unique element attributes used by the Mesh.
void InitMesh(int _Dim, int _spaceDim, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
const Element * GetBdrElement(int i) const
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
void UpdateNodes()
Update the nodes of a curved mesh after refinement.
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=mfem::out)
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
void Printer(std::ostream &out=mfem::out, std::string section_delimiter="") const
Class used to extrude the nodes of a mesh.