15 #include "../config/config.hpp" 16 #include "../general/stable3d.hpp" 17 #include "../general/globals.hpp" 23 #include "../fem/eltrans.hpp" 24 #include "../fem/coefficient.hpp" 25 #include "../general/zstr.hpp" 26 #ifdef MFEM_USE_ADIOS2 27 #include "../general/adios2stream.hpp" 36 class GeometricFactors;
37 class FaceGeometricFactors;
40 class FiniteElementSpace;
61 #ifdef MFEM_USE_ADIOS2 253 #ifdef MFEM_USE_MEMALLOC 308 void ReadMFEMMesh(std::istream &input,
int version,
int &curved);
317 int &curved,
int &read_gf,
bool &finalize_topo);
318 void ReadVTKMesh(std::istream &input,
int &curved,
int &read_gf,
319 bool &finalize_topo);
321 bool &finalize_topo,
const std::string &xml_prefix=
"");
322 void ReadNURBSMesh(std::istream &input,
int &curved,
int &read_gf);
323 void ReadInlineMesh(std::istream &input,
bool generate_edges =
false);
324 void ReadGmshMesh(std::istream &input,
int &curved,
int &read_gf);
326 #ifdef MFEM_USE_NETCDF 327 void ReadCubit(
const char *filename,
int &curved,
int &read_gf);
356 int *edge1,
int *edge2,
int *middle)
362 int *edge1,
int *edge2,
int *middle)
363 {
Bisection(i, v_to_v, edge1, edge2, middle); }
407 bool update_nodes =
true);
424 double threshold,
int nc_limit = 0,
428 const int *fine,
int nfine,
int op);
465 const FaceInfo &fi,
bool is_ghost);
515 int v0,
int v1,
int v2);
518 int v0,
int v1,
int v2,
int v3);
533 void InitMesh(
int Dim_,
int spaceDim_,
int NVert,
int NElem,
int NBdrElem);
538 void Loader(std::istream &input,
int generate_edges = 0,
539 std::string parse_tag =
"");
545 std::string section_delimiter =
"")
const;
553 double sx,
double sy,
double sz,
bool sfc_ordering);
562 bool generate_edges,
bool sfc_ordering);
565 void Make1D(
int n,
double sx = 1.0);
601 explicit Mesh(
const Mesh &mesh,
bool copy_nodes =
true);
624 int *element_attributes,
int num_elements,
626 int *boundary_attributes,
int num_boundary_elements,
627 int dimension,
int space_dimension = -1);
637 Mesh(
int Dim_,
int NVert,
int NElem,
int NBdrElem = 0,
int spaceDim_ = -1)
639 if (spaceDim_ == -1) { spaceDim_ = Dim_; }
640 InitMesh(Dim_, spaceDim_, NVert, NElem, NBdrElem);
646 explicit Mesh(
const std::string &filename,
int generate_edges = 0,
647 int refine = 1,
bool fix_orientation =
true);
652 explicit Mesh(std::istream &input,
int generate_edges = 0,
int refine = 1,
653 bool fix_orientation =
true);
658 Mesh(
Mesh *mesh_array[],
int num_pieces);
666 virtual void Load(std::istream &input,
int generate_edges = 0,
667 int refine = 1,
bool fix_orientation =
true)
669 Loader(input, generate_edges);
675 void Swap(
Mesh& other,
bool non_geometry);
701 int generate_edges = 0,
int refine = 1,
702 bool fix_orientation =
true);
713 int nx,
int ny,
Element::Type type,
bool generate_edges =
false,
714 double sx = 1.0,
double sy = 1.0,
bool sfc_ordering =
true);
722 double sx = 1.0,
double sy = 1.0,
double sz = 1.0,
723 bool sfc_ordering =
true);
787 int AddVertex(
double x,
double y = 0.0,
double z = 0.0);
796 int AddTriangle(
int v1,
int v2,
int v3,
int attr = 1);
800 int AddQuad(
int v1,
int v2,
int v3,
int v4,
int attr = 1);
801 int AddQuad(
const int *vi,
int attr = 1);
803 int AddTet(
int v1,
int v2,
int v3,
int v4,
int attr = 1);
804 int AddTet(
const int *vi,
int attr = 1);
806 int AddWedge(
int v1,
int v2,
int v3,
int v4,
int v5,
int v6,
int attr = 1);
807 int AddWedge(
const int *vi,
int attr = 1);
809 int AddPyramid(
int v1,
int v2,
int v3,
int v4,
int v5,
int attr = 1);
812 int AddHex(
int v1,
int v2,
int v3,
int v4,
int v5,
int v6,
int v7,
int v8,
814 int AddHex(
const int *vi,
int attr = 1);
832 int AddBdrQuad(
int v1,
int v2,
int v3,
int v4,
int attr = 1);
841 bool fix_orientation =
true);
844 bool fix_orientation =
true);
847 bool fix_orientation =
true);
850 bool fix_orientation =
true);
853 bool fix_orientation =
true);
856 void FinalizeMesh(
int refine = 0,
bool fix_orientation =
true);
888 virtual void Finalize(
bool refine =
false,
bool fix_orientation =
false);
960 int iterations = 4,
int window = 4,
961 int period = 2,
int seed = 0,
962 bool verbose =
false,
double time_limit = 0);
986 double sx = 1.0,
double sy = 1.0,
double sz = 1.0,
987 bool sfc_ordering =
true)
989 Make3D(nx, ny, nz, type, sx, sy, sz, sfc_ordering);
996 double sx = 1.0,
double sy = 1.0,
bool sfc_ordering =
true)
998 Make2D(nx, ny, type, sx, sy, generate_edges, sfc_ordering);
1004 explicit Mesh(
int n,
double sx = 1.0)
1012 Mesh(
Mesh *orig_mesh,
int ref_factor,
int ref_type);
1073 double &kappa_min,
double &kappa_max,
1230 return elements[i]->GetGeometryType();
1235 return boundary[i]->GetGeometryType();
1319 faces[i]->GetVertices(vert);
1519 if (
faces_info[FaceNo].Elem2No < 0) {
return NULL; }
1770 void GetFaceInfos (
int Face,
int *Inf1,
int *Inf2)
const;
1771 void GetFaceInfos (
int Face,
int *Inf1,
int *Inf2,
int *NCFace)
const;
1805 bool zerocopy =
false);
1810 void GetNode(
int i,
double *coord)
const;
1811 void SetNode(
int i,
const double *coord);
1896 virtual void SetCurvature(
int order,
bool discont =
false,
int space_dim = -1,
1926 int nonconforming = -1,
int nc_limit = 0);
1931 int nonconforming = -1,
int nc_limit = 0);
1935 int nonconforming = -1,
int nc_limit = 0);
1939 double eps = 0.0,
int nonconforming = -1);
1944 int nonconforming = -1,
int nc_limit = 0);
1949 int nonconforming = -1,
int nc_limit = 0);
1957 int nc_limit = 0,
int op = 1);
1961 int nc_limit = 0,
int op = 1);
1966 void EnsureNCMesh(
bool simplices_nonconforming =
false);
2015 virtual void Save(
const std::string &fname,
int precision=16)
const;
2018 #ifdef MFEM_USE_ADIOS2 2029 void PrintVTK(std::ostream &os,
int ref,
int field_data=0);
2037 bool high_order_output=
false,
2038 int compression_level=0,
2039 bool bdr_elements=
false);
2041 virtual void PrintVTU(std::string fname,
2043 bool high_order_output=
false,
2044 int compression_level=0,
2050 bool high_order_output=
false,
2051 int compression_level=0);
2058 std::ostream &os,
int elem_attr = 0)
const;
2062 int interior_faces = 0);
2113 const std::vector<Vector> &translations,
double tol = 1e-8)
const;
2160 virtual long long ReduceInt(
int value)
const {
return value; }
2176 std::ostream &
operator<<(std::ostream &
out,
const Mesh &mesh);
2305 Mesh *
Extrude1D(Mesh *mesh,
const int ny,
const double sy,
2306 const bool closed =
false);
2309 Mesh *
Extrude2D(Mesh *mesh,
const int nz,
const double sz);
2315 a = c; c =
b;
b =
t;
2319 std::ostream&
operator<<(std::ostream& os,
const Mesh::FaceInformation& info);
static Mesh MakePeriodic(const Mesh &orig_mesh, const std::vector< int > &v2v)
Create a periodic mesh by identifying vertices of orig_mesh.
Abstract class for all finite elements.
const IntegrationRule * IntRule
MFEM_DEPRECATED Mesh(int n, double sx=1.0)
Deprecated: see MakeCartesian1D.
void Loader(std::istream &input, int generate_edges=0, std::string parse_tag="")
double GetGeckoElementOrdering(Array< int > &ordering, int iterations=4, int window=4, int period=2, int seed=0, bool verbose=false, double time_limit=0)
Table * GetEdgeVertexTable() const
Geometry::Constants< Geometry::PYRAMID > pyr_t
Array< int > GetFaceToBdrElMap() const
static const int vtk_quadratic_hex[27]
int * CartesianPartitioning(int nxyz[])
Class for an integration rule - an Array of IntegrationPoint.
void SetPatchAttribute(int i, int attr)
Set the attribute of patch i, for a NURBS mesh.
const FiniteElementSpace * GetNodalFESpace() const
const DenseMatrix * PointMatrix
void ScaleElements(double sf)
const GridFunction * GetNodes() const
Class for grid function - Vector with associated FE space.
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo)
void FreeElement(Element *E)
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
int CheckElementOrientation(bool fix_it=true)
Check (and optionally attempt to fix) the orientation of the elements.
static bool remove_unused_vertices
void SetVertices(const Vector &vert_coord)
static const int vtk_quadratic_tet[10]
void Make2D(int nx, int ny, Element::Type type, double sx, double sy, bool generate_edges, bool sfc_ordering)
Vector J
Jacobians of the element transformations at all quadrature points.
Vector X
Mapped (physical) coordinates of all quadrature points.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
static Mesh MakeSimplicial(const Mesh &orig_mesh)
Base class for vector Coefficients that optionally depend on time and space.
void UniformRefinement3D_base(Array< int > *f2qf=NULL, DSTable *v_to_v_p=NULL, bool update_nodes=true)
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
int EulerNumber2D() const
Equals 1 - num_holes.
Geometry::Type GetElementBaseGeometry(int i) const
int GetBdrElementEdgeIndex(int i) const
static void PrintElementsByGeometry(int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
Auxiliary method used by PrintCharacteristics().
void AddHexAsTets(const int *vi, int attr=1)
static int ComposeQuadOrientations(int ori_a_b, int ori_b_c)
void AddHexAsWedges(const int *vi, int attr=1)
virtual void PrintInfo(std::ostream &os=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
const Table & ElementToFaceTable() const
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
Array< Element * > boundary
Geometry::Constants< Geometry::SQUARE > quad_t
int Dimension() const
Dimension of the reference space used within the elements.
int * GeneratePartitioning(int nparts, int part_method=1)
CoarseFineTransformations CoarseFineTr
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
void MoveVertices(const Vector &displacements)
bool HasGeometry(Geometry::Type geom) const
Return true iff the given geom is encountered in the mesh. Geometries of dimensions lower than Dimens...
static Mesh LoadFromFile(const std::string &filename, int generate_edges=0, int refine=1, bool fix_orientation=true)
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)
static FiniteElement * GetTransformationFEforElementType(Element::Type)
Return FiniteElement for reference element of the specified type.
FaceGeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, FaceType type, MemoryType d_mt=MemoryType::DEFAULT)
static Mesh MakeCartesian3D(int nx, int ny, int nz, Element::Type type, double sx=1.0, double sy=1.0, double sz=1.0, bool sfc_ordering=true)
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
virtual void UniformRefinement2D()
Refine a mixed 2D mesh uniformly.
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones...
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
static int InvertTriOrientation(int ori)
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
void ReadNetgen2DMesh(std::istream &input, int &curved)
void ShiftRight(int &a, int &b, int &c)
void GetVertexToVertexTable(DSTable &) const
static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type)
Create a refined (by any factor) version of orig_mesh.
void CreateVTKMesh(const Vector &points, const Array< int > &cell_data, const Array< int > &cell_offsets, const Array< int > &cell_types, const Array< int > &cell_attributes, int &curved, int &read_gf, bool &finalize_topo)
int AddTriangle(int v1, int v2, int v3, int attr=1)
Vector detJ
Determinants of the Jacobians at all quadrature points.
bool Nonconforming() const
int AddBdrPoint(int v, int attr=1)
Data type dense matrix using column-major storage.
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
Mesh * Extrude1D(Mesh *mesh, const int ny, const double sy, const bool closed)
Extrude a 1D mesh.
void Transform(void(*f)(const Vector &, Vector &))
bool IsSlaveFace(const FaceInfo &fi) const
void ApplyLocalSlaveTransformation(FaceElementTransformations &FT, const FaceInfo &fi, bool is_ghost)
void GetElementJacobian(int i, DenseMatrix &J, const IntegrationPoint *ip=NULL)
void AverageVertices(const int *indexes, int n, int result)
Averages the vertices with given indexes and saves the result in vertices[result].
void FinalizeWedgeMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a wedge Mesh.
Element * GetElement(int i)
Return pointer to the i'th element object.
virtual void Save(const std::string &fname, int precision=16) const
void GenerateBoundaryElements()
void GetElementData(const Array< Element *> &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
int GetNEdges() const
Return the number of edges.
const Element * GetElement(int i) const
Return pointer to the i'th element object.
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &os=mfem::out)
Compute and print mesh characteristics such as number of vertices, number of elements, number of boundary elements, minimal and maximal element sizes, minimal and maximal element aspect ratios, etc.
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
NodeExtrudeCoefficient(const int dim, const int n_, const double s_)
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 ReadNetgen3DMesh(std::istream &input)
Data arrays will be written in ASCII format.
const Element *const * GetElementsArray() const
void ReadInlineMesh(std::istream &input, bool generate_edges=false)
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
void RemoveInternalBoundaries()
Mesh * Extrude2D(Mesh *mesh, const int nz, const double sz)
Extrude a 2D mesh.
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
std::vector< int > CreatePeriodicVertexMapping(const std::vector< Vector > &translations, double tol=1e-8) const
Creates a mapping v2v from the vertex indices of the mesh such that coincident vertices under the giv...
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
void GetPointMatrix(int i, DenseMatrix &pointmat) const
int GetNBE() const
Returns number of boundary elements.
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)
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
void InitRefinementTransforms()
IsoparametricTransformation FaceTransformation
void UniformRefinement2D_base(bool update_nodes=true)
Array< FaceGeometricFactors * > face_geom_factors
Array< NCFaceInfo > nc_faces_info
Element * ReadElement(std::istream &)
void GetLocalTriToPyrTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
A parallel extension of the NCMesh class.
Element * NewElement(int geom)
FaceElementTransformations FaceElemTr
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
void SetVerticesFromNodes(const GridFunction *nodes)
Helper to set vertex coordinates given a high-order curvature function.
Structure for storing face geometric factors: coordinates, Jacobians, determinants of the Jacobians...
ElementTransformation * GetBdrElementTransformation(int i)
std::function< double(const Vector &)> f(double mass_coeff)
void DebugDump(std::ostream &out) const
Output an NCMesh-compatible debug dump.
int AddBdrTriangle(int v1, int v2, int v3, int attr=1)
void SetPatchBdrAttribute(int i, int attr)
Set the attribute of patch boundary element i, for a NURBS mesh.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
void MakeRefined_(Mesh &orig_mesh, const Array< int > ref_factors, int ref_type)
Internal function used in Mesh::MakeRefined.
Geometry::Constants< Geometry::SEGMENT > seg_t
int GetAttribute(int i) const
Return the attribute of element i.
void MesquiteSmooth(const int mesquite_option=0)
static int GetTetOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
void GetBdrElementAdjacentElement2(int bdr_el, int &el, int &info) const
For the given boundary element, bdr_el, return its adjacent element and its info, i...
MFEM_DEPRECATED Mesh(int nx, int ny, int nz, Element::Type type, bool generate_edges=false, double sx=1.0, double sy=1.0, double sz=1.0, bool sfc_ordering=true)
Deprecated: see MakeCartesian3D.
int AddTri(const int *vi, int attr=1)
int AddVertex(double x, double y=0.0, double z=0.0)
void GetBdrElementFace(int i, int *f, int *o) const
Return the index and the orientation of the face of bdr element i. (3D)
Vector J
Jacobians of the element transformations at all quadrature points.
long long GetGlobalNE() const
Return the total (global) number of elements.
void MoveNodes(const Vector &displacements)
void GetEdgeOrdering(DSTable &v_to_v, Array< int > &order)
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
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.
double GetElementSize(ElementTransformation *T, int type=0)
Symmetric 3D Table stored as an array of rows each of which has a stack of column, floor, number nodes. The number of the node is assigned by counting the nodes from zero as they are pushed into the table. Diagonals of any kind are not allowed so the row, column and floor must all be different for each node. Only one node is stored for all 6 symmetric entries that are indexable by unique triplets of row, column, and floor.
void GetVertices(Vector &vert_coord) const
void PrintVTK(std::ostream &os)
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 Make3D(int nx, int ny, int nz, Element::Type type, double sx, double sy, double sz, bool sfc_ordering)
void AddSegmentFaceElement(int lf, int gf, int el, int v0, int v1)
virtual void SetAttributes()
Determine the sets of unique attribute values in domain and boundary elements.
void EnsureNCMesh(bool simplices_nonconforming=false)
GeometryList(Mesh &mesh, int dim)
Construct a GeometryList of all geometries of dimension dim in mesh.
void GetNode(int i, double *coord) const
Vector detJ
Determinants of the Jacobians at all quadrature points.
void AddHexAsPyramids(const int *vi, int attr=1)
int FindCoarseElement(int i)
GeometryList(Mesh &mesh)
Construct a GeometryList of all element geometries in mesh.
int AddTet(int v1, int v2, int v3, int v4, int attr=1)
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
virtual void PrintXG(std::ostream &os=mfem::out) const
Print the mesh to the given stream using Netgen/Truegrid format.
void GetLocalSegToTriTransformation(IsoparametricTransformation &loc, int i)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void AddQuadFaceElement(int lf, int gf, int el, int v0, int v1, int v2, int v3)
void CheckPartitioning(int *partitioning_)
void KnotInsert(Array< KnotVector *> &kv)
void LoadPatchTopo(std::istream &input, Array< int > &edge_to_knot)
Read NURBS patch/macro-element mesh.
void ReadGmshMesh(std::istream &input, int &curved, int &read_gf)
STable3D * GetElementToFaceTable(int ret_ftbl=0)
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
int MeshGenerator()
Get the mesh generator/type.
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Type
Constants for the classes derived from Element.
void CheckDisplacements(const Vector &displacements, double &tmax)
This structure stores the low level information necessary to interpret the configuration of elements ...
const double * GetVertex(int i) const
Return pointer to vertex i's coordinates.
VTKFormat
Data array format for VTK and VTU files.
void GetElementData(int geom, Array< int > &elem_vtx, Array< int > &attr) const
int AddBdrSegment(int v1, int v2, int attr=1)
Geometry::Constants< Geometry::TETRAHEDRON > tet_t
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
int AddElement(Element *elem)
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 GetGeometricParametersFromJacobian(const DenseMatrix &J, double &volume, Vector &aspr, Vector &skew, Vector &ori) const
Computes geometric parameters associated with a Jacobian matrix in 2D/3D. These parameters are (1) Ar...
void SetLayer(const int l)
Table * GetFaceToElementTable() const
Array< int > FindFaceNeighbors(const int elem) const
Returns the sorted, unique indices of elements sharing a face with element elem, including elem...
int GetElementToEdgeTable(Table &, Array< int > &)
void GetHilbertElementOrdering(Array< int > &ordering)
void RandomRefinement(double prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
Refine each element with given probability. Uses GeneralRefinement.
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 GetLocalQuadToHexTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
static const int vtk_quadratic_pyramid[13]
bool OwnsNodes() const
Return the mesh nodes ownership flag.
void GetElementCenter(int i, Vector ¢er)
MFEM_DEPRECATED Geometry::Type GetFaceBaseGeometry(int i) const
Deprecated in favor of Mesh::GetFaceGeometry.
IsoparametricTransformation BdrTransformation
void GetLocalFaceTransformation(int face_type, int elem_type, IsoparametricTransformation &Transf, int info)
A helper method that constructs a transformation from the reference space of a face to the reference ...
IsoparametricTransformation Transformation2
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
Element * GetBdrElement(int i)
Return pointer to the i'th boundary element object.
Vector normal
Normals at all quadrature points.
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Geometry::Constants< Geometry::PRISM > pri_t
virtual MFEM_DEPRECATED void ReorientTetMesh()
int GetBdrFace(int BdrElemNo) const
Return the local face index for the given boundary face.
void ReadMFEMMesh(std::istream &input, int version, int &curved)
virtual void LocalRefinement(const Array< int > &marked_el, int type=3)
This function is not public anymore. Use GeneralRefinement instead.
Table * GetFaceEdgeTable() const
Data type tetrahedron element.
virtual void MarkTetMeshForRefinement(DSTable &v_to_v)
int GetPatchAttribute(int i) const
Return the attribute of patch i, for a NURBS mesh.
Array< Triple< int, int, int > > tmp_vertex_parents
MFEM_DEPRECATED Mesh(int nx, int ny, Element::Type type, bool generate_edges=false, double sx=1.0, double sy=1.0, bool sfc_ordering=true)
Deprecated: see MakeCartesian2D.
List of mesh geometries stored as Array<Geometry::Type>.
void PrintVTU(std::ostream &os, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr_elements=false)
int CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Geometry::Type GetElementGeometry(int i) const
GeometricFactors(const Mesh *mesh, const IntegrationRule &ir, int flags, MemoryType d_mt=MemoryType::DEFAULT)
void PrintElementsWithPartitioning(int *partitioning, std::ostream &out, int interior_faces=0)
int GetPatchBdrAttribute(int i) const
Return the attribute of patch boundary element i, for a NURBS mesh.
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
void GetLocalQuadToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
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.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
virtual ~Mesh()
Destroys Mesh.
void FinalizeMesh(int refine=0, bool fix_orientation=true)
Finalize the construction of any type of Mesh.
const Table & ElementToEdgeTable() const
int AddSegment(int v1, int v2, int attr=1)
bool FaceIsTrueInterior(int FaceNo) const
Geometry::Constants< Geometry::CUBE > hex_t
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Geometry::Type geom_buf[Geometry::NumGeom]
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
static void PrintElement(const Element *, std::ostream &)
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
const IntegrationRule * IntRule
MemoryType
Memory types supported by MFEM.
Vector X
Mapped (physical) coordinates of all quadrature points.
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
GridFunction * GetNodes()
Return a pointer to the internal node GridFunction (may be NULL).
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void PrepareNodeReorder(DSTable **old_v_to_v, Table **old_elem_vert)
int AddBdrQuad(int v1, int v2, int v3, int v4, int attr=1)
void PrintTopo(std::ostream &out, const Array< int > &e_to_k) const
static Mesh MakeCartesian2D(int nx, int ny, Element::Type type, bool generate_edges=false, double sx=1.0, double sy=1.0, bool sfc_ordering=true)
void InitMesh(int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem)
Begin construction of a mesh.
void SetBdrAttribute(int i, int attr)
Set the attribute of boundary element i.
void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo, const std::string &xml_prefix="")
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
void RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
static Mesh MakeCartesian1D(int n, double sx=1.0)
Array< Element * > elements
const Table & ElementToElementTable()
int SpaceDimension() const
Dimension of the physical space containing the mesh.
const CoarseFineTransformations & GetRefinementTransforms()
void RefineNURBSFromFile(std::string ref_file)
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
void GetGeometries(int dim, Array< Geometry::Type > &el_geoms) const
Return all element geometries of the given dimension present in the mesh.
Geometry::Type GetBdrElementGeometry(int i) const
int GetNE() const
Returns number of elements.
void GetLocalTriToTetTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
void GetFaceEdges(int i, Array< int > &edges, Array< int > &o) const
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void PrintBdrVTU(std::string fname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0)
int EulerNumber() const
Equals 1 + num_holes - num_loops.
Class for integration point with weight.
Element::Type GetElementType(int i) const
Returns the type of element i.
void Swap(Mesh &other, bool non_geometry)
Element * ReadElementWithoutAttr(std::istream &)
static const int vtk_quadratic_wedge[18]
virtual long long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
void GetLocalSegToQuadTransformation(IsoparametricTransformation &loc, int i)
virtual int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type, does not count master nonconforming face...
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
const Element * GetFace(int i) const
virtual 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)
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
void GetLocalQuadToPyrTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
void DegreeElevate(int rel_degree, int degree=16)
virtual bool HasBoundaryElements() const
Checks if the mesh has boundary elements.
Element::Type GetFaceElementType(int Face) const
virtual ~NodeExtrudeCoefficient()
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
void MakeSimplicial_(const Mesh &orig_mesh, int *vglobal)
Array< FaceInfo > faces_info
static int InvertQuadOrientation(int ori)
STable3D * GetFacesTable()
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
NCMesh * ncmesh
Optional nonconforming mesh extension.
virtual void UniformRefinement3D()
Refine a mixed 3D mesh uniformly.
int AddBdrElement(Element *elem)
static void PrintElementWithoutAttr(const Element *, std::ostream &)
Mesh(int Dim_, int NVert, int NElem, int NBdrElem=0, int spaceDim_=-1)
Init constructor: begin the construction of a Mesh object.
int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr=1)
FaceInformation GetFaceInformation(int f) const
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.
Geometry::Constants< Geometry::TRIANGLE > tri_t
double AggregateError(const Array< double > &elem_error, const int *fine, int nfine, int op)
Derefinement helper.
double * GetVertex(int i)
Return pointer to vertex i's coordinates.
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.
Mesh & operator=(Mesh &&mesh)
Move assignment operator.
MFEM_DEPRECATED Geometry::Type GetFaceGeometryType(int Face) const
Deprecated in favor of Mesh::GetFaceGeometry.
int GetNFaces() const
Return the number of faces in a 3D mesh.
void ReadTrueGridMesh(std::istream &input)
virtual void Print(std::ostream &os=mfem::out) const
const FaceGeometricFactors * GetFaceGeometricFactors(const IntegrationRule &ir, const int flags, FaceType type, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors for the faces corresponding to the given integration rule...
static void GetElementArrayEdgeTable(const Array< Element *> &elem_array, const DSTable &v_to_v, Table &el_to_edge)
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 PrintSurfaces(const Table &Aface_face, std::ostream &out) const
Print set of disjoint surfaces:
void GetLocalTriToWdgTransformation(IsoparametricTransformation &loc, int i)
Used in GetFaceElementTransformations (...)
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 ...
void DeleteGeometricFactors()
Destroy all GeometricFactors stored by the Mesh.
int GetNumFacesWithGhost() const
Return the number of faces (3D), edges (2D) or vertices (1D) including ghost faces.
IsoparametricTransformation EdgeTransformation
void Printer(std::ostream &out=mfem::out, std::string section_delimiter="") const
MemAlloc< Tetrahedron, 1024 > TetMemory
void DoNodeReorder(DSTable *old_v_to_v, Table *old_elem_vert)
void SetNodes(const Vector &node_coord)
Updates the vertex/node locations. Invokes NodesUpdated().
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
void SetNodalGridFunction(GridFunction *nodes, bool make_owner=false)
Operation GetLastOperation() const
Return type of last modification of the mesh.
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
void ScaleSubdomains(double sf)
int AddPyramid(int v1, int v2, int v3, int v4, int v5, int attr=1)
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
void Bisection(int i, const DSTable &, int *, int *, int *)
Bisect a triangle: element with index i is bisected.
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Class for parallel meshes.
Geometry::Type GetBdrElementBaseGeometry(int i) const
Abstract data type element.
void AddTriangleFaceElement(int lf, int gf, int el, int v0, int v1, int v2)
Array< GeometricFactors * > geom_factors
Optional geometric factors.
void BdrBisection(int i, const HashTable< Hashed2 > &)
Bisect a boundary triangle: boundary element with index i is bisected.
static int ComposeTriOrientations(int ori_a_b, int ori_b_c)
void GenerateNCFaceInfo()
double GetLength(int i, int j) const
Return the length of the segment from node i to node j.
double GetElementVolume(int i)
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Table * GetVertexToElementTable()
The returned Table should be deleted by the caller.
void EnsureNodes()
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element gr...
void PrintWithPartitioning(int *partitioning, std::ostream &os, int elem_attr=0) const
Prints the mesh with boundary elements given by the boundary of the subdomains, so that the boundary ...
void UpdateNodes()
Update the nodes of a curved mesh after the topological part of a Mesh::Operation, such as refinement, has been performed.
void GetBdrPointMatrix(int i, DenseMatrix &pointmat) const
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 GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
void NodesUpdated()
This function should be called after the mesh node coordinates have been updated externally, e.g. by modifying the internal nodal GridFunction returned by GetNodes().
Class used to extrude the nodes of a mesh.
void AddVertexParents(int i, int p1, int p2)
Mark vertex i as nonconforming, with parent vertices p1 and p2.