MFEM v4.7.0
Finite element discretization library
|
#include <mesh.hpp>
Classes | |
struct | FaceInfo |
This structure stores the low level information necessary to interpret the configuration of elements on a specific face. This information can be accessed using methods like GetFaceElements(), GetFaceInfos(), FaceIsInterior(), etc. More... | |
struct | FaceInformation |
This structure is used as a human readable output format that deciphers the information contained in Mesh::FaceInfo when using the Mesh::GetFaceInformation() method. More... | |
class | GeometryList |
List of mesh geometries stored as Array<Geometry::Type>. More... | |
struct | NCFaceInfo |
Public Types | |
enum | Operation { NONE , REFINE , DEREFINE , REBALANCE } |
enum class | FaceTopology { Boundary , Conforming , Nonconforming , NA } |
enum class | ElementLocation { Local , FaceNbr , NA } |
enum class | ElementConformity { Coincident , Superset , Subset , NA } |
enum class | FaceInfoTag { Boundary , LocalConforming , LocalSlaveNonconforming , SharedConforming , SharedSlaveNonconforming , MasterNonconforming , GhostSlave , GhostMaster } |
typedef Geometry::Constants< Geometry::SEGMENT > | seg_t |
typedef Geometry::Constants< Geometry::TRIANGLE > | tri_t |
typedef Geometry::Constants< Geometry::SQUARE > | quad_t |
typedef Geometry::Constants< Geometry::TETRAHEDRON > | tet_t |
typedef Geometry::Constants< Geometry::CUBE > | hex_t |
typedef Geometry::Constants< Geometry::PRISM > | pri_t |
typedef Geometry::Constants< Geometry::PYRAMID > | pyr_t |
Public Member Functions | |
Standard Mesh constructors and related methods | |
These constructors and assignment operators accept mesh information in a variety of common forms. For more specialized constructors see Named mesh constructors. | |
Mesh () | |
Mesh (const Mesh &mesh, bool copy_nodes=true) | |
Mesh (Mesh &&mesh) | |
Move constructor, useful for using a Mesh as a function return value. | |
Mesh & | operator= (Mesh &&mesh) |
Move assignment operator. | |
Mesh & | operator= (const Mesh &mesh)=delete |
Explicitly delete the copy assignment operator. | |
Mesh (real_t *vertices, int num_vertices, int *element_indices, Geometry::Type element_type, int *element_attributes, int num_elements, int *boundary_indices, Geometry::Type boundary_type, int *boundary_attributes, int num_boundary_elements, int dimension, int space_dimension=-1) | |
Construct a Mesh from the given primary data. | |
Mesh (int Dim_, int NVert, int NElem, int NBdrElem=0, int spaceDim_=-1) | |
Init constructor: begin the construction of a Mesh object. | |
Mesh (const std::string &filename, int generate_edges=0, int refine=1, bool fix_orientation=true) | |
Mesh (std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true) | |
Mesh (Mesh *mesh_array[], int num_pieces) | |
virtual void | Load (std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true) |
void | Swap (Mesh &other, bool non_geometry) |
void | Clear () |
Clear the contents of the Mesh. | |
virtual | ~Mesh () |
Destroys Mesh. | |
Methods for piecewise Mesh construction. | |
These methods are intended to be used with the init constructor. | |
Element * | NewElement (int geom) |
int | AddVertex (real_t x, real_t y=0.0, real_t z=0.0) |
int | AddVertex (const real_t *coords) |
int | AddVertex (const Vector &coords) |
void | AddVertexParents (int i, int p1, int p2) |
Mark vertex i as nonconforming, with parent vertices p1 and p2. | |
int | AddVertexAtMeanCenter (const int *vi, const int nverts, int dim=3) |
int | AddSegment (int v1, int v2, int attr=1) |
Adds a segment to the mesh given by 2 vertices v1 and v2. | |
int | AddSegment (const int *vi, int attr=1) |
Adds a segment to the mesh given by 2 vertices vi. | |
int | AddTriangle (int v1, int v2, int v3, int attr=1) |
Adds a triangle to the mesh given by 3 vertices v1 through v3. | |
int | AddTriangle (const int *vi, int attr=1) |
Adds a triangle to the mesh given by 3 vertices vi. | |
int | AddTri (const int *vi, int attr=1) |
Adds a triangle to the mesh given by 3 vertices vi. | |
int | AddQuad (int v1, int v2, int v3, int v4, int attr=1) |
Adds a quadrilateral to the mesh given by 4 vertices v1 through v4. | |
int | AddQuad (const int *vi, int attr=1) |
Adds a quadrilateral to the mesh given by 4 vertices vi. | |
int | AddTet (int v1, int v2, int v3, int v4, int attr=1) |
Adds a tetrahedron to the mesh given by 4 vertices v1 through v4. | |
int | AddTet (const int *vi, int attr=1) |
Adds a tetrahedron to the mesh given by 4 vertices vi. | |
int | AddWedge (int v1, int v2, int v3, int v4, int v5, int v6, int attr=1) |
Adds a wedge to the mesh given by 6 vertices v1 through v6. | |
int | AddWedge (const int *vi, int attr=1) |
Adds a wedge to the mesh given by 6 vertices vi. | |
int | AddPyramid (int v1, int v2, int v3, int v4, int v5, int attr=1) |
Adds a pyramid to the mesh given by 5 vertices v1 through v5. | |
int | AddPyramid (const int *vi, int attr=1) |
Adds a pyramid to the mesh given by 5 vertices vi. | |
int | AddHex (int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1) |
Adds a hexahedron to the mesh given by 8 vertices v1 through v8. | |
int | AddHex (const int *vi, int attr=1) |
Adds a hexahedron to the mesh given by 8 vertices vi. | |
void | AddHexAsTets (const int *vi, int attr=1) |
Adds 6 tetrahedrons to the mesh by splitting a hexahedron given by 8 vertices vi. | |
void | AddHexAsWedges (const int *vi, int attr=1) |
Adds 2 wedges to the mesh by splitting a hexahedron given by 8 vertices vi. | |
void | AddHexAsPyramids (const int *vi, int attr=1) |
Adds 6 pyramids to the mesh by splitting a hexahedron given by 8 vertices vi. | |
void | AddHexAs24TetsWithPoints (int *vi, std::map< std::array< int, 4 >, int > &hex_face_verts, int attr=1) |
Adds 24 tetrahedrons to the mesh by splitting a hexahedron. | |
void | AddQuadAs4TrisWithPoints (int *vi, int attr=1) |
Adds 4 triangles to the mesh by splitting a quadrilateral given by 4 vertices vi. | |
void | AddQuadAs5QuadsWithPoints (int *vi, int attr=1) |
Adds 5 quadrilaterals to the mesh by splitting a quadrilateral given by 4 vertices vi. | |
int | AddElement (Element *elem) |
int | AddBdrElement (Element *elem) |
int | AddBdrSegment (int v1, int v2, int attr=1) |
int | AddBdrSegment (const int *vi, int attr=1) |
int | AddBdrTriangle (int v1, int v2, int v3, int attr=1) |
int | AddBdrTriangle (const int *vi, int attr=1) |
int | AddBdrQuad (int v1, int v2, int v3, int v4, int attr=1) |
int | AddBdrQuad (const int *vi, int attr=1) |
void | AddBdrQuadAsTriangles (const int *vi, int attr=1) |
int | AddBdrPoint (int v, int attr=1) |
virtual void | GenerateBoundaryElements () |
void | FinalizeTriMesh (int generate_edges=0, int refine=0, bool fix_orientation=true) |
Finalize the construction of a triangular Mesh. | |
void | FinalizeQuadMesh (int generate_edges=0, int refine=0, bool fix_orientation=true) |
Finalize the construction of a quadrilateral Mesh. | |
void | FinalizeTetMesh (int generate_edges=0, int refine=0, bool fix_orientation=true) |
Finalize the construction of a tetrahedral Mesh. | |
void | FinalizeWedgeMesh (int generate_edges=0, int refine=0, bool fix_orientation=true) |
Finalize the construction of a wedge Mesh. | |
void | FinalizeHexMesh (int generate_edges=0, int refine=0, bool fix_orientation=true) |
Finalize the construction of a hexahedral Mesh. | |
void | FinalizeMesh (int refine=0, bool fix_orientation=true) |
Finalize the construction of any type of Mesh. | |
Mesh consistency methods | |
void | FinalizeTopology (bool generate_bdr=true) |
Finalize the construction of the secondary topology (connectivity) data of a Mesh. | |
virtual void | Finalize (bool refine=false, bool fix_orientation=false) |
Finalize the construction of a general Mesh. | |
virtual void | SetAttributes () |
Determine the sets of unique attribute values in domain and boundary elements. | |
int | CheckElementOrientation (bool fix_it=true) |
Check (and optionally attempt to fix) the orientation of the elements. | |
int | CheckBdrElementOrientation (bool fix_it=true) |
Check the orientation of the boundary elements. | |
virtual MFEM_DEPRECATED void | ReorientTetMesh () |
void | RemoveUnusedVertices () |
Remove unused vertices and rebuild mesh connectivity. | |
void | RemoveInternalBoundaries () |
Element ordering methods | |
real_t | GetGeckoElementOrdering (Array< int > &ordering, int iterations=4, int window=4, int period=2, int seed=0, bool verbose=false, real_t time_limit=0) |
void | GetHilbertElementOrdering (Array< int > &ordering) |
void | ReorderElements (const Array< int > &ordering, bool reorder_vertices=true) |
Deprecated mesh constructors | |
These constructors have been deprecated in favor of Named mesh constructors. | |
MFEM_DEPRECATED | Mesh (int nx, int ny, int nz, Element::Type type, bool generate_edges=false, real_t sx=1.0, real_t sy=1.0, real_t sz=1.0, bool sfc_ordering=true) |
Deprecated: see MakeCartesian3D. | |
MFEM_DEPRECATED | Mesh (int nx, int ny, Element::Type type, bool generate_edges=false, real_t sx=1.0, real_t sy=1.0, bool sfc_ordering=true) |
Deprecated: see MakeCartesian2D. | |
MFEM_DEPRECATED | Mesh (int n, real_t sx=1.0) |
Deprecated: see MakeCartesian1D. | |
MFEM_DEPRECATED | Mesh (Mesh *orig_mesh, int ref_factor, int ref_type) |
Deprecated: see MakeRefined. | |
Information about the mesh as a whole | |
int | Dimension () const |
Dimension of the reference space used within the elements. | |
int | SpaceDimension () const |
Dimension of the physical space containing the mesh. | |
int | EulerNumber () const |
Equals 1 + num_holes - num_loops. | |
int | EulerNumber2D () const |
Equals 1 - num_holes. | |
int | MeshGenerator () const |
Get the mesh generator/type. | |
virtual bool | HasBoundaryElements () const |
Checks if the mesh has boundary elements. | |
bool | HasGeometry (Geometry::Type geom) const |
Return true iff the given geom is encountered in the mesh. Geometries of dimensions lower than Dimension() are counted as well. | |
int | GetNumGeometries (int dim) const |
Return the number of geometries of the given dimension present in the mesh. | |
void | GetGeometries (int dim, Array< Geometry::Type > &el_geoms) const |
Return all element geometries of the given dimension present in the mesh. | |
void | GetBoundingBox (Vector &min, Vector &max, int ref=2) |
Returns the minimum and maximum corners of the mesh bounding box. | |
void | GetCharacteristics (real_t &h_min, real_t &h_max, real_t &kappa_min, real_t &kappa_max, Vector *Vh=NULL, Vector *Vk=NULL) |
Information concerning numbers of mesh entities | |
int | GetNV () const |
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them in the lowest-order mesh. | |
int | GetNE () const |
Returns number of elements. | |
int | GetNBE () const |
Returns number of boundary elements. | |
int | GetNEdges () const |
Return the number of edges. | |
int | GetNFaces () const |
Return the number of faces in a 3D mesh. | |
int | GetNumFaces () const |
Return the number of faces (3D), edges (2D) or vertices (1D). | |
int | GetNumFacesWithGhost () const |
Return the number of faces (3D), edges (2D) or vertices (1D) including ghost faces. | |
virtual int | GetNFbyType (FaceType type) const |
Returns the number of faces according to the requested type, does not count master nonconforming faces. | |
long long | GetGlobalNE () const |
Return the total (global) number of elements. | |
Access to individual mesh entities | |
const real_t * | GetVertex (int i) const |
Return pointer to vertex i's coordinates. | |
real_t * | GetVertex (int i) |
Return pointer to vertex i's coordinates. | |
const Element * | GetElement (int i) const |
Return pointer to the i'th element object. | |
Element * | GetElement (int i) |
Return pointer to the i'th element object. | |
const Element * | GetBdrElement (int i) const |
Return pointer to the i'th boundary element object. | |
Element * | GetBdrElement (int i) |
Return pointer to the i'th boundary element object. | |
const Element * | GetFace (int i) const |
Return pointer to the i'th face element object. | |
Access to groups of mesh entities | |
const Element *const * | GetElementsArray () const |
void | GetElementData (int geom, Array< int > &elem_vtx, Array< int > &attr) const |
void | GetBdrElementData (int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const |
Access information concerning individual mesh entites | |
int | GetAttribute (int i) const |
Return the attribute of element i. | |
void | SetAttribute (int i, int attr) |
Set the attribute of element i. | |
int | GetBdrAttribute (int i) const |
Return the attribute of boundary element i. | |
void | SetBdrAttribute (int i, int attr) |
Set the attribute of boundary element i. | |
int | GetPatchAttribute (int i) const |
Return the attribute of patch i, for a NURBS mesh. | |
void | SetPatchAttribute (int i, int attr) |
Set the attribute of patch i, for a NURBS mesh. | |
int | GetPatchBdrAttribute (int i) const |
Return the attribute of patch boundary element i, for a NURBS mesh. | |
void | SetPatchBdrAttribute (int i, int attr) |
Set the attribute of patch boundary element i, for a NURBS mesh. | |
Element::Type | GetElementType (int i) const |
Returns the type of element i. | |
Element::Type | GetBdrElementType (int i) const |
Returns the type of boundary element i. | |
MFEM_DEPRECATED Geometry::Type | GetFaceGeometryType (int Face) const |
Deprecated in favor of Mesh::GetFaceGeometry. | |
Element::Type | GetFaceElementType (int Face) const |
Geometry::Type | GetFaceGeometry (int i) const |
Return the Geometry::Type associated with face i. | |
Geometry::Type | GetElementGeometry (int i) const |
Geometry::Type | GetBdrElementGeometry (int i) const |
MFEM_DEPRECATED Geometry::Type | GetFaceBaseGeometry (int i) const |
Deprecated in favor of Mesh::GetFaceGeometry. | |
Geometry::Type | GetElementBaseGeometry (int i) const |
Geometry::Type | GetBdrElementBaseGeometry (int i) const |
bool | FaceIsInterior (int FaceNo) const |
Return true if the given face is interior. | |
real_t | GetElementSize (int i, int type=0) |
Get the size of the i-th element relative to the perfect reference element. | |
real_t | GetElementSize (int i, const Vector &dir) |
real_t | GetElementSize (ElementTransformation *T, int type=0) const |
real_t | GetElementVolume (int i) |
void | GetElementCenter (int i, Vector ¢er) |
void | GetElementJacobian (int i, DenseMatrix &J, const IntegrationPoint *ip=NULL) |
Access connectivity for individual mesh entites | |
void | GetElementVertices (int i, Array< int > &v) const |
Returns the indices of the vertices of element i. | |
void | GetBdrElementVertices (int i, Array< int > &v) const |
Returns the indices of the vertices of boundary element i. | |
void | GetElementEdges (int i, Array< int > &edges, Array< int > &cor) const |
Return the indices and the orientations of all edges of element i. | |
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 | GetFaceEdges (int i, Array< int > &edges, Array< int > &o) const |
void | GetFaceVertices (int i, Array< int > &vert) const |
Returns the indices of the vertices of face i. | |
void | GetEdgeVertices (int i, Array< int > &vert) const |
Returns the indices of the vertices of edge i. | |
void | GetElementFaces (int i, Array< int > &faces, Array< int > &ori) const |
Return the indices and the orientations of all faces of element i. | |
Array< int > | FindFaceNeighbors (const int elem) const |
Returns the sorted, unique indices of elements sharing a face with element elem, including elem. | |
void | GetBdrElementFace (int i, int *f, int *o) 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.e. 64*local_bdr_index+bdr_orientation. | |
MFEM_DEPRECATED void | GetBdrElementAdjacentElement2 (int bdr_el, int &el, int &info) const |
Deprecated. | |
int | GetBdrElementFaceIndex (int be_idx) const |
Return the local face (codimension-1) index for the given boundary element index. | |
MFEM_DEPRECATED int | GetBdrFace (int i) const |
Deprecated in favor of GetBdrElementFaceIndex(). | |
MFEM_DEPRECATED int | GetBdrElementEdgeIndex (int i) const |
Access connectivity data | |
Table * | GetVertexToElementTable () |
Table * | GetFaceToElementTable () const |
Table * | GetFaceEdgeTable () const |
Table * | GetEdgeVertexTable () const |
void | GetVertexToVertexTable (DSTable &) const |
const Table & | ElementToElementTable () |
const Table & | ElementToFaceTable () const |
const Table & | ElementToEdgeTable () const |
Array< int > | GetFaceToBdrElMap () const |
Access the coordinate transformation for individual elements | |
See also the methods related to Geometric Factors for accessing information cached at quadrature points. | |
void | GetElementTransformation (int i, IsoparametricTransformation *ElTr) const |
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and will be owned by the caller. | |
ElementTransformation * | GetElementTransformation (int i) |
Returns a pointer to the transformation defining the i-th element. | |
void | GetElementTransformation (int i, const Vector &nodes, IsoparametricTransformation *ElTr) const |
Builds the transformation defining the i-th element in ElTr assuming position of the vertices/nodes are given by nodes. ElTr must be allocated in advance and will be owned by the caller. | |
ElementTransformation * | GetBdrElementTransformation (int i) |
Returns a pointer to the transformation defining the i-th boundary element. | |
void | GetBdrElementTransformation (int i, IsoparametricTransformation *ElTr) const |
Builds the transformation defining the i-th boundary element in ElTr. ElTr must be allocated in advance and will be owned by the caller. | |
ElementTransformation * | GetFaceTransformation (int FaceNo) |
Returns a pointer to the transformation defining the given face element. | |
void | GetFaceTransformation (int i, IsoparametricTransformation *FTr) const |
Builds the transformation defining the i-th face element in FTr. FTr must be allocated in advance and will be owned by the caller. | |
void | GetLocalFaceTransformation (int face_type, int elem_type, IsoparametricTransformation &Transf, int info) const |
A helper method that constructs a transformation from the reference space of a face to the reference space of an element. | |
void | GetEdgeTransformation (int i, IsoparametricTransformation *EdTr) const |
Builds the transformation defining the i-th edge element in EdTr. EdTr must be allocated in advance and will be owned by the caller. | |
ElementTransformation * | GetEdgeTransformation (int EdgeNo) |
Returns a pointer to the transformation defining the given edge element. | |
virtual FaceElementTransformations * | GetFaceElementTransformations (int FaceNo, int mask=31) |
virtual void | GetFaceElementTransformations (int FaceNo, FaceElementTransformations &FElTr, IsoparametricTransformation &ElTr1, IsoparametricTransformation &ElTr2, int mask=31) const |
Variant of GetFaceElementTransformations using a user allocated FaceElementTransformations object. | |
FaceElementTransformations * | GetInteriorFaceTransformations (int FaceNo) |
See GetFaceElementTransformations(). | |
void | GetInteriorFaceTransformations (int FaceNo, FaceElementTransformations &FElTr, IsoparametricTransformation &ElTr1, IsoparametricTransformation &ElTr2) const |
Variant of GetInteriorFaceTransformations using a user allocated FaceElementTransformations object. | |
FaceElementTransformations * | GetBdrFaceTransformations (int BdrElemNo) |
Builds the transformation defining the given boundary face. | |
void | GetBdrFaceTransformations (int BdrElemNo, FaceElementTransformations &FElTr, IsoparametricTransformation &ElTr1, IsoparametricTransformation &ElTr2) const |
Variant of GetBdrFaceTransformations using a user allocated FaceElementTransformations object. | |
Access the coordinate transformation at quadrature points | |
See also methods related to Element-wise coordinate transformation. | |
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. | |
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. | |
void | DeleteGeometricFactors () |
Destroy all GeometricFactors stored by the Mesh. | |
More advanced entity information access methods | |
void | GetPointMatrix (int i, DenseMatrix &pointmat) const |
void | GetBdrPointMatrix (int i, DenseMatrix &pointmat) const |
FaceInformation | GetFaceInformation (int f) const |
void | GetFaceElements (int Face, int *Elem1, int *Elem2) const |
void | GetFaceInfos (int Face, int *Inf1, int *Inf2) const |
void | GetFaceInfos (int Face, int *Inf1, int *Inf2, int *NCFace) const |
Methods related to mesh partitioning | |
int * | CartesianPartitioning (int nxyz[]) |
int * | GeneratePartitioning (int nparts, int part_method=1) |
void | CheckPartitioning (int *partitioning_) |
Methods related to accessing/altering mesh coordinates | |
See also Coordinates as a GridFunction. | |
void | MoveVertices (const Vector &displacements) |
void | GetVertices (Vector &vert_coord) const |
void | SetVertices (const Vector &vert_coord) |
void | ChangeVertexDataOwnership (real_t *vertices, int len_vertices, bool zerocopy=false) |
Set the internal Vertex array to point to the given vertices array without assuming ownership of the pointer. | |
void | GetNode (int i, real_t *coord) const |
void | SetNode (int i, const real_t *coord) |
void | MoveNodes (const Vector &displacements) |
void | GetNodes (Vector &node_coord) const |
void | SetNodes (const Vector &node_coord) |
Updates the vertex/node locations. Invokes NodesUpdated(). | |
void | ScaleSubdomains (real_t sf) |
void | ScaleElements (real_t sf) |
void | Transform (void(*f)(const Vector &, Vector &)) |
void | Transform (VectorCoefficient &deformation) |
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(). | |
Methods related to nodal coordinates stored as a GridFunction | |
See also Mesh Transformations. | |
GridFunction * | GetNodes () |
Return a pointer to the internal node GridFunction (may be NULL). | |
const GridFunction * | GetNodes () const |
bool | OwnsNodes () const |
Return the mesh nodes ownership flag. | |
void | SetNodesOwner (bool nodes_owner) |
Set the mesh nodes ownership flag. | |
void | NewNodes (GridFunction &nodes, bool make_owner=false) |
Replace the internal node GridFunction with the given GridFunction. | |
void | SwapNodes (GridFunction *&nodes, int &own_nodes_) |
Swap the internal node GridFunction pointer and ownership flag members with the given ones. | |
void | GetNodes (GridFunction &nodes) const |
Return the mesh nodes/vertices projected on the given GridFunction. | |
virtual void | SetNodalFESpace (FiniteElementSpace *nfes) |
void | SetNodalGridFunction (GridFunction *nodes, bool make_owner=false) |
const FiniteElementSpace * | GetNodalFESpace () const |
void | EnsureNodes () |
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element grid function (even if it is a low-order mesh with straight edges). | |
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. | |
Methods related to mesh refinement | |
void | UniformRefinement (int ref_algo=0) |
Refine all mesh elements. | |
virtual void | NURBSUniformRefinement (int rf=2, real_t tol=1.0e-12) |
Refine NURBS mesh, with an optional refinement factor, generally anisotropic. | |
virtual void | NURBSUniformRefinement (const Array< int > &rf, real_t tol=1.e-12) |
void | NURBSCoarsening (int cf=2, real_t tol=1.0e-12) |
void | GeneralRefinement (const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0) |
void | GeneralRefinement (const Array< int > &el_to_refine, int nonconforming=-1, int nc_limit=0) |
void | RandomRefinement (real_t prob, bool aniso=false, int nonconforming=-1, int nc_limit=0) |
Refine each element with given probability. Uses GeneralRefinement. | |
void | RefineAtVertex (const Vertex &vert, real_t eps=0.0, int nonconforming=-1) |
Refine elements sharing the specified vertex. Uses GeneralRefinement. | |
bool | RefineByError (const Array< real_t > &elem_error, real_t threshold, int nonconforming=-1, int nc_limit=0) |
bool | RefineByError (const Vector &elem_error, real_t threshold, int nonconforming=-1, int nc_limit=0) |
bool | DerefineByError (Array< real_t > &elem_error, real_t threshold, int nc_limit=0, int op=1) |
bool | DerefineByError (const Vector &elem_error, real_t threshold, int nc_limit=0, int op=1) |
Same as DerefineByError for an error vector. | |
void | EnsureNCMesh (bool simplices_nonconforming=false) |
bool | Conforming () const |
bool | Nonconforming () const |
const CoarseFineTransformations & | GetRefinementTransforms () const |
Operation | GetLastOperation () const |
Return type of last modification of the mesh. | |
long | GetSequence () const |
long | GetNodesSequence () const |
Return the nodes update counter. | |
NURBS mesh refinement methods | |
void | RefineNURBSFromFile (std::string ref_file) |
void | KnotInsert (Array< KnotVector * > &kv) |
For NURBS meshes, insert the new knots in kv, for each direction. | |
void | KnotInsert (Array< Vector * > &kv) |
For NURBS meshes, insert the knots in kv, for each direction. | |
void | KnotRemove (Array< Vector * > &kv) |
For NURBS meshes, remove the knots in kv, for each direction. | |
void | DegreeElevate (int rel_degree, int degree=16) |
Miscellaneous or undocumented methods | |
std::vector< int > | CreatePeriodicVertexMapping (const std::vector< Vector > &translations, real_t tol=1e-8) const |
Creates a mapping v2v from the vertex indices of the mesh such that coincident vertices under the given translations are identified. | |
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 coordinates. | |
void | GetGeometricParametersFromJacobian (const DenseMatrix &J, real_t &volume, Vector &aspr, Vector &skew, Vector &ori) const |
Computes geometric parameters associated with a Jacobian matrix in 2D/3D. These parameters are (1) Area/Volume, (2) Aspect-ratio (1 in 2D, and 2 non-dimensional and 2 dimensional parameters in 3D. Dimensional parameters are used for target construction in TMOP), (3) skewness (1 in 2D and 3 in 3D), and finally (4) orientation (1 in 2D and 3 in 3D). | |
virtual long long | ReduceInt (int value) const |
Utility function: sum integers from all processors (Allreduce). | |
void | GetElementColoring (Array< int > &colors, int el0=0) |
void | MesquiteSmooth (const int mesquite_option=0) |
void | CheckDisplacements (const Vector &displacements, real_t &tmax) |
Static Public Member Functions | |
static FiniteElement * | GetTransformationFEforElementType (Element::Type) |
Return FiniteElement for reference element of the specified type. | |
static IntegrationPoint | TransformBdrElementToFace (Geometry::Type geom, int o, const IntegrationPoint &ip) |
For the vertex (1D), edge (2D), or face (3D) of a boundary element with the orientation o, return the transformation of the boundary element integration point @ ip to the face element. In 2D, the the orientation is 0 or 1 as returned by GetBdrElementFace, not +/-1. Supports both internal and external boundaries. | |
static int | DecodeFaceInfoOrientation (int info) |
Given a "face info int", return the face orientation. | |
static int | DecodeFaceInfoLocalIndex (int info) |
Given a "face info int", return the local face index. | |
static int | EncodeFaceInfo (int local_face_index, int orientation) |
Given local_face_index and orientation, return the corresponding encoded "face info int". | |
Named mesh constructors. | |
Each of these constructors uses the move constructor, and can be used as the right-hand side of an assignment when creating new meshes. For more general mesh constructors see Standard mesh constructors. | |
static Mesh | LoadFromFile (const std::string &filename, int generate_edges=0, int refine=1, bool fix_orientation=true) |
static Mesh | MakeCartesian1D (int n, real_t sx=1.0) |
Creates 1D mesh, divided into n equal intervals. | |
static Mesh | MakeCartesian2D (int nx, int ny, Element::Type type, bool generate_edges=false, real_t sx=1.0, real_t sy=1.0, bool sfc_ordering=true) |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if type = TRIANGLE. | |
static Mesh | MakeCartesian3D (int nx, int ny, int nz, Element::Type type, real_t sx=1.0, real_t sy=1.0, real_t sz=1.0, bool sfc_ordering=true) |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedra if type = HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if type = TETRAHEDRON. | |
static Mesh | MakeCartesian3DWith24TetsPerHex (int nx, int ny, int nz, real_t sx=1.0, real_t sy=1.0, real_t sz=1.0) |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz*24 tetrahedrons. | |
static Mesh | MakeCartesian2DWith4TrisPerQuad (int nx, int ny, real_t sx=1.0, real_t sy=1.0) |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny*4 triangles. | |
static Mesh | MakeCartesian2DWith5QuadsPerQuad (int nx, int ny, real_t sx=1.0, real_t sy=1.0) |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny*5 quadrilaterals. | |
static Mesh | MakeRefined (Mesh &orig_mesh, int ref_factor, int ref_type) |
Create a refined (by any factor) version of orig_mesh. | |
static Mesh | MakeRefined (Mesh &orig_mesh, const Array< int > &ref_factors, int ref_type) |
refined ref_factors[i] times in each dimension. | |
static Mesh | MakeSimplicial (const Mesh &orig_mesh) |
static Mesh | MakePeriodic (const Mesh &orig_mesh, const std::vector< int > &v2v) |
Create a periodic mesh by identifying vertices of orig_mesh. | |
Public Attributes | |
Array< int > | attributes |
A list of all unique element attributes used by the Mesh. | |
Array< int > | bdr_attributes |
A list of all unique boundary attributes used by the Mesh. | |
AttributeSets | attribute_sets |
Named sets of element attributes. | |
AttributeSets | bdr_attribute_sets |
Named sets of boundary element attributes. | |
NURBSExtension * | NURBSext |
Optional NURBS mesh extension. | |
NCMesh * | ncmesh |
Optional nonconforming mesh extension. | |
Array< GeometricFactors * > | geom_factors |
Optional geometric factors. | |
Array< FaceGeometricFactors * > | face_geom_factors |
Static Public Attributes | |
static bool | remove_unused_vertices = true |
Protected Member Functions | |
void | Init () |
void | InitTables () |
void | SetEmpty () |
void | DestroyTables () |
void | DeleteTables () |
void | DestroyPointers () |
void | Destroy () |
void | ResetLazyData () |
Element * | ReadElementWithoutAttr (std::istream &input) |
Element * | ReadElement (std::istream &input) |
void | ReadMFEMMesh (std::istream &input, int version, int &curved) |
void | ReadLineMesh (std::istream &input) |
void | ReadNetgen2DMesh (std::istream &input, int &curved) |
void | ReadNetgen3DMesh (std::istream &input) |
void | ReadTrueGridMesh (std::istream &input) |
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) |
void | ReadVTKMesh (std::istream &input, int &curved, int &read_gf, bool &finalize_topo) |
void | ReadXML_VTKMesh (std::istream &input, int &curved, int &read_gf, bool &finalize_topo, const std::string &xml_prefix="") |
void | ReadNURBSMesh (std::istream &input, int &curved, int &read_gf, bool spacing=false) |
void | ReadInlineMesh (std::istream &input, bool generate_edges=false) |
void | ReadGmshMesh (std::istream &input, int &curved, int &read_gf) |
void | ReadCubit (const std::string &filename, int &curved, int &read_gf) |
Load a mesh from a Genesis file. | |
void | SetMeshGen () |
Determine the mesh generator bitmask meshgen, see MeshGenerator(). | |
real_t | GetLength (int i, int j) const |
Return the length of the segment from node i to node j. | |
void | MarkForRefinement () |
void | MarkTriMeshForRefinement () |
void | GetEdgeOrdering (const DSTable &v_to_v, Array< int > &order) |
virtual void | MarkTetMeshForRefinement (const DSTable &v_to_v) |
void | PrepareNodeReorder (DSTable **old_v_to_v, Table **old_elem_vert) |
void | DoNodeReorder (DSTable *old_v_to_v, Table *old_elem_vert) |
STable3D * | GetFacesTable () |
STable3D * | GetElementToFaceTable (int ret_ftbl=0) |
void | RedRefinement (int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle) |
void | GreenRefinement (int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle) |
void | Bisection (int i, const DSTable &, int *, int *, int *) |
Bisect a triangle: element with index i is bisected. | |
void | Bisection (int i, HashTable< Hashed2 > &) |
Bisect a tetrahedron: element with index i is bisected. | |
void | BdrBisection (int i, const HashTable< Hashed2 > &) |
Bisect a boundary triangle: boundary element with index i is bisected. | |
void | UniformRefinement (int i, const DSTable &, int *, int *, int *) |
void | AverageVertices (const int *indexes, int n, int result) |
Averages the vertices with given indexes and saves the result in vertices[result]. | |
void | InitRefinementTransforms () |
int | FindCoarseElement (int i) |
void | UpdateNodes () |
Update the nodes of a curved mesh after the topological part of a Mesh::Operation, such as refinement, has been performed. | |
void | SetVerticesFromNodes (const GridFunction *nodes) |
Helper to set vertex coordinates given a high-order curvature function. | |
void | UniformRefinement2D_base (bool update_nodes=true) |
virtual void | UniformRefinement2D () |
Refine a mixed 2D mesh uniformly. | |
void | UniformRefinement3D_base (Array< int > *f2qf=NULL, DSTable *v_to_v_p=NULL, bool update_nodes=true) |
virtual void | UniformRefinement3D () |
Refine a mixed 3D mesh uniformly. | |
virtual void | LocalRefinement (const Array< int > &marked_el, int type=3) |
This function is not public anymore. Use GeneralRefinement instead. | |
virtual void | NonconformingRefinement (const Array< Refinement > &refinements, int nc_limit=0) |
This function is not public anymore. Use GeneralRefinement instead. | |
virtual bool | NonconformingDerefinement (Array< real_t > &elem_error, real_t threshold, int nc_limit=0, int op=1) |
NC version of GeneralDerefinement. | |
real_t | AggregateError (const Array< real_t > &elem_error, const int *fine, int nfine, int op) |
Derefinement helper. | |
void | LoadPatchTopo (std::istream &input, Array< int > &edge_to_knot) |
Read NURBS patch/macro-element mesh. | |
void | UpdateNURBS () |
void | PrintTopo (std::ostream &os, const Array< int > &e_to_k, const int version, const std::string &comment="") const |
Write the beginning of a NURBS mesh to os, specifying the NURBS patch topology. Optional file comments can be provided in comments. | |
void | GetLocalPtToSegTransformation (IsoparametricTransformation &, int i) const |
Used in GetFaceElementTransformations (...) | |
void | GetLocalSegToTriTransformation (IsoparametricTransformation &loc, int i) const |
void | GetLocalSegToQuadTransformation (IsoparametricTransformation &loc, int i) const |
void | GetLocalTriToTetTransformation (IsoparametricTransformation &loc, int i) const |
void | GetLocalTriToWdgTransformation (IsoparametricTransformation &loc, int i) const |
void | GetLocalTriToPyrTransformation (IsoparametricTransformation &loc, int i) const |
void | GetLocalQuadToHexTransformation (IsoparametricTransformation &loc, int i) const |
void | GetLocalQuadToWdgTransformation (IsoparametricTransformation &loc, int i) const |
void | GetLocalQuadToPyrTransformation (IsoparametricTransformation &loc, int i) const |
void | ApplyLocalSlaveTransformation (FaceElementTransformations &FT, const FaceInfo &fi, bool is_ghost) const |
bool | IsSlaveFace (const FaceInfo &fi) const |
int | GetElementToEdgeTable (Table &) |
void | AddPointFaceElement (int lf, int gf, int el) |
Used in GenerateFaces() | |
void | AddSegmentFaceElement (int lf, int gf, int el, int v0, int v1) |
void | AddTriangleFaceElement (int lf, int gf, int el, int v0, int v1, int v2) |
void | AddQuadFaceElement (int lf, int gf, int el, int v0, int v1, int v2, int v3) |
bool | FaceIsTrueInterior (int FaceNo) const |
void | FreeElement (Element *E) |
void | GenerateFaces () |
void | GenerateNCFaceInfo () |
void | InitMesh (int Dim_, int spaceDim_, int NVert, int NElem, int NBdrElem) |
Begin construction of a mesh. | |
void | FinalizeCheck () |
void | Loader (std::istream &input, int generate_edges=0, std::string parse_tag="") |
void | Printer (std::ostream &os=mfem::out, std::string section_delimiter="", const std::string &comments="") const |
void | Make3D (int nx, int ny, int nz, Element::Type type, real_t sx, real_t sy, real_t sz, bool sfc_ordering) |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedra if type = HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if type = TETRAHEDRON. | |
void | Make3D24TetsFromHex (int nx, int ny, int nz, real_t sx, real_t sy, real_t sz) |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz*24 tetrahedrons. | |
void | Make2D4TrisFromQuad (int nx, int ny, real_t sx, real_t sy) |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny*4 triangles. | |
void | Make2D5QuadsFromQuad (int nx, int ny, real_t sx, real_t sy) |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny*5 quadrilaterals. | |
void | Make2D (int nx, int ny, Element::Type type, real_t sx, real_t sy, bool generate_edges, bool sfc_ordering) |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if type = TRIANGLE. | |
void | Make1D (int n, real_t sx=1.0) |
void | MakeRefined_ (Mesh &orig_mesh, const Array< int > &ref_factors, int ref_type) |
Internal function used in Mesh::MakeRefined. | |
void | InitFromNCMesh (const NCMesh &ncmesh) |
Initialize vertices/elements/boundary/tables from a nonconforming mesh. | |
Mesh (const NCMesh &ncmesh) | |
Create from a nonconforming mesh. | |
void | GetElementData (const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const |
void | MakeSimplicial_ (const Mesh &orig_mesh, int *vglobal) |
Static Protected Member Functions | |
static void | PrintElementWithoutAttr (const Element *el, std::ostream &os) |
static void | PrintElement (const Element *el, std::ostream &os) |
static int | GetTriOrientation (const int *base, const int *test) |
Returns the orientation of "test" relative to "base". | |
static int | InvertTriOrientation (int ori) |
static int | ComposeTriOrientations (int ori_a_b, int ori_b_c) |
static int | GetQuadOrientation (const int *base, const int *test) |
Returns the orientation of "test" relative to "base". | |
static int | InvertQuadOrientation (int ori) |
static int | ComposeQuadOrientations (int ori_a_b, int ori_b_c) |
static int | GetTetOrientation (const int *base, const int *test) |
Returns the orientation of "test" relative to "base". | |
static void | GetElementArrayEdgeTable (const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge) |
Static Protected Attributes | |
static const int | vtk_quadratic_tet [10] |
static const int | vtk_quadratic_pyramid [13] |
static const int | vtk_quadratic_wedge [18] |
static const int | vtk_quadratic_hex [27] |
Friends | |
class | NCMesh |
class | NURBSExtension |
class | ParMesh |
class | ParNCMesh |
class | adios2stream |
class | Tetrahedron |
Print/Save/Export methods | |
virtual void | PrintXG (std::ostream &os=mfem::out) const |
Print the mesh to the given stream using Netgen/Truegrid format. | |
virtual void | Print (std::ostream &os=mfem::out, const std::string &comments="") const |
virtual void | Save (const std::string &fname, int precision=16) const |
virtual void | Print (adios2stream &os) const |
Print the mesh to the given stream using the adios2 bp format. | |
void | PrintVTK (std::ostream &os) |
void | PrintVTK (std::ostream &os, int ref, int field_data=0) |
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) |
virtual void | PrintVTU (std::string fname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr=false) |
void | PrintBdrVTU (std::string fname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0) |
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 of subdomain i has boundary attribute i+1. | |
void | PrintElementsWithPartitioning (int *partitioning, std::ostream &os, int interior_faces=0) |
void | PrintSurfaces (const Table &Aface_face, std::ostream &os) const |
Print set of disjoint surfaces: | |
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. | |
virtual void | PrintInfo (std::ostream &os=mfem::out) |
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the parallel decomposition is also printed. | |
void | DebugDump (std::ostream &os) const |
Output an NCMesh-compatible debug dump. | |
static void | PrintElementsByGeometry (int dim, const Array< int > &num_elems_by_geom, std::ostream &os) |
Auxiliary method used by PrintCharacteristics(). | |
|
strong |
This enumerated type describes the topological relation of an element to a face:
Enumerator | |
---|---|
Coincident | |
Superset | |
Subset | |
NA |
|
strong |
This enumerated type describes the location of the two elements sharing a face, Local meaning that the element is local to the MPI rank, FaceNbr meaning that the element is distributed on a different MPI rank, this typically means that methods with FaceNbr should be used to access the relevant information, e.g., ParFiniteElementSpace::GetFaceNbrElementVDofs.
Enumerator | |
---|---|
Local | |
FaceNbr | |
NA |
|
strong |
This enumerated type describes the corresponding FaceInfo internal representation (encoded cases), c.f. FaceInfo's documentation: Classification of a local (non-ghost) face based on its FaceInfo:
Enumerator | |
---|---|
Boundary | |
LocalConforming | |
LocalSlaveNonconforming | |
SharedConforming | |
SharedSlaveNonconforming | |
MasterNonconforming | |
GhostSlave | |
GhostMaster |
|
strong |
This enumerated type describes the three main face topologies:
Enumerator | |
---|---|
Boundary | |
Conforming | |
Nonconforming | |
NA |
|
explicitprotected |
|
explicit |
mfem::Mesh::Mesh | ( | Mesh && | mesh | ) |
mfem::Mesh::Mesh | ( | real_t * | vertices, |
int | num_vertices, | ||
int * | element_indices, | ||
Geometry::Type | element_type, | ||
int * | element_attributes, | ||
int | num_elements, | ||
int * | boundary_indices, | ||
Geometry::Type | boundary_type, | ||
int * | boundary_attributes, | ||
int | num_boundary_elements, | ||
int | dimension, | ||
int | space_dimension = -1 ) |
Construct a Mesh from the given primary data.
The array vertices is used as external data, i.e. the Mesh does not copy the data and will not delete the pointer.
The data from the other arrays is copied into the internal Mesh data structures.
This method calls the method FinalizeTopology(). The method Finalize() may be called after this constructor and after optionally setting the Mesh nodes.
|
inline |
Init constructor: begin the construction of a Mesh object.
Construct a shell of a mesh object allocating space to store pointers to the vertices, elements, and boundary elements. The vertices and elements themselves can later be added using methods from the Mesh construction group.
|
explicit |
Creates mesh by reading a file in MFEM, Netgen, or VTK format. If generate_edges = 0 (default) edges are not generated, if 1 edges are generated. See also Mesh::LoadFromFile. See Mesh::Finalize for the meaning of refine.
|
explicit |
mfem::Mesh::Mesh | ( | Mesh * | mesh_array[], |
int | num_pieces ) |
|
inline |
|
inline |
|
inlineexplicit |
mfem::Mesh::Mesh | ( | Mesh * | orig_mesh, |
int | ref_factor, | ||
int | ref_type ) |
int mfem::Mesh::AddBdrElement | ( | Element * | elem | ) |
The parameter elem should be allocated using the NewElement() method
int mfem::Mesh::AddBdrQuad | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddBdrQuad | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | v4, | ||
int | attr = 1 ) |
void mfem::Mesh::AddBdrQuadAsTriangles | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddBdrSegment | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddBdrSegment | ( | int | v1, |
int | v2, | ||
int | attr = 1 ) |
int mfem::Mesh::AddBdrTriangle | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddBdrTriangle | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | attr = 1 ) |
int mfem::Mesh::AddElement | ( | Element * | elem | ) |
The parameter elem should be allocated using the NewElement() method
int mfem::Mesh::AddHex | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddHex | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | v4, | ||
int | v5, | ||
int | v6, | ||
int | v7, | ||
int | v8, | ||
int | attr = 1 ) |
void mfem::Mesh::AddHexAs24TetsWithPoints | ( | int * | vi, |
std::map< std::array< int, 4 >, int > & | hex_face_verts, | ||
int | attr = 1 ) |
Adds 24 tetrahedrons to the mesh by splitting a hexahedron.
vi are the 8 vertices of the hexahedron, hex_face_verts has the map from the 4 vertices of each face of the hexahedron to the index of the point created at the center of the face, and attr is the attribute of the new elements. See Make3D24TetsFromHex for usage.
void mfem::Mesh::AddHexAsPyramids | ( | const int * | vi, |
int | attr = 1 ) |
void mfem::Mesh::AddHexAsTets | ( | const int * | vi, |
int | attr = 1 ) |
void mfem::Mesh::AddHexAsWedges | ( | const int * | vi, |
int | attr = 1 ) |
|
protected |
Used in GenerateFaces()
int mfem::Mesh::AddPyramid | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddPyramid | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | v4, | ||
int | v5, | ||
int | attr = 1 ) |
int mfem::Mesh::AddQuad | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddQuad | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | v4, | ||
int | attr = 1 ) |
void mfem::Mesh::AddQuadAs4TrisWithPoints | ( | int * | vi, |
int | attr = 1 ) |
void mfem::Mesh::AddQuadAs5QuadsWithPoints | ( | int * | vi, |
int | attr = 1 ) |
|
protected |
int mfem::Mesh::AddSegment | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddSegment | ( | int | v1, |
int | v2, | ||
int | attr = 1 ) |
|
protected |
int mfem::Mesh::AddTet | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddTet | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | v4, | ||
int | attr = 1 ) |
|
inline |
int mfem::Mesh::AddTriangle | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddTriangle | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | attr = 1 ) |
|
protected |
int mfem::Mesh::AddVertexAtMeanCenter | ( | const int * | vi, |
const int | nverts, | ||
int | dim = 3 ) |
void mfem::Mesh::AddVertexParents | ( | int | i, |
int | p1, | ||
int | p2 ) |
int mfem::Mesh::AddWedge | ( | const int * | vi, |
int | attr = 1 ) |
int mfem::Mesh::AddWedge | ( | int | v1, |
int | v2, | ||
int | v3, | ||
int | v4, | ||
int | v5, | ||
int | v6, | ||
int | attr = 1 ) |
|
protected |
|
protected |
|
protected |
int * mfem::Mesh::CartesianPartitioning | ( | int | nxyz[] | ) |
void mfem::Mesh::ChangeVertexDataOwnership | ( | real_t * | vertices, |
int | len_vertices, | ||
bool | zerocopy = false ) |
int mfem::Mesh::CheckBdrElementOrientation | ( | bool | fix_it = true | ) |
int mfem::Mesh::CheckElementOrientation | ( | bool | fix_it = true | ) |
Check (and optionally attempt to fix) the orientation of the elements.
[in] | fix_it | If true , attempt to fix the orientations of some elements: triangles, quads, and tets. |
true
.
|
inline |
|
staticprotected |
|
staticprotected |
std::vector< int > mfem::Mesh::CreatePeriodicVertexMapping | ( | const std::vector< Vector > & | translations, |
real_t | tol = 1e-8 ) const |
Creates a mapping v2v from the vertex indices of the mesh such that coincident vertices under the given translations are identified.
Each Vector in translations should be of size sdim (the spatial dimension of the mesh). Two vertices are considered coincident if the translated coordinates of one vertex are within the given tolerance (tol, relative to the mesh diameter) of the coordinates of the other vertex.
|
protected |
Definition at line 396 of file mesh_readers.cpp.
void mfem::Mesh::DebugDump | ( | std::ostream & | os | ) | const |
|
inlinestatic |
|
inlinestatic |
void mfem::Mesh::DegreeElevate | ( | int | rel_degree, |
int | degree = 16 ) |
void mfem::Mesh::DeleteGeometricFactors | ( | ) |
Destroy all GeometricFactors stored by the Mesh.
This method can be used to force recomputation of the GeometricFactors, for example, after the mesh nodes are modified externally.
bool mfem::Mesh::DerefineByError | ( | Array< real_t > & | elem_error, |
real_t | threshold, | ||
int | nc_limit = 0, | ||
int | op = 1 ) |
Derefine the mesh based on an error measure associated with each element. A derefinement is performed if the sum of errors of its fine elements is smaller than 'threshold'. If 'nc_limit' > 0, derefinements that would increase the maximum level of hanging nodes of the mesh are skipped. Returns true if the mesh changed, false otherwise.
|
inline |
|
inlinestatic |
void mfem::Mesh::EnsureNCMesh | ( | bool | simplices_nonconforming = false | ) |
void mfem::Mesh::EnsureNodes | ( | ) |
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element grid function (even if it is a low-order mesh with straight edges).
|
inline |
|
inline |
|
inline |
Return true if the given face is interior.
|
inlineprotected |
|
virtual |
Finalize the construction of a general Mesh.
This method will:
[in] | refine | If true, prepare the Mesh for conforming refinement of triangular or tetrahedral meshes. |
[in] | fix_orientation | If true, fix the orientation of inverted mesh elements by permuting their vertices. |
Reimplemented in mfem::ParMesh.
void mfem::Mesh::FinalizeHexMesh | ( | int | generate_edges = 0, |
int | refine = 0, | ||
bool | fix_orientation = true ) |
void mfem::Mesh::FinalizeMesh | ( | int | refine = 0, |
bool | fix_orientation = true ) |
Finalize the construction of any type of Mesh.
This method calls FinalizeTopology() and Finalize().
void mfem::Mesh::FinalizeQuadMesh | ( | int | generate_edges = 0, |
int | refine = 0, | ||
bool | fix_orientation = true ) |
void mfem::Mesh::FinalizeTetMesh | ( | int | generate_edges = 0, |
int | refine = 0, | ||
bool | fix_orientation = true ) |
void mfem::Mesh::FinalizeTopology | ( | bool | generate_bdr = true | ) |
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
This method does not require any actual coordinate data (either vertex coordinates for linear meshes or node coordinates for meshes with nodes) to be available. However, the data generated by this method is generally required by the FiniteElementSpace class.
After calling this method, setting the Mesh vertices or nodes, it may be appropriate to call the method Finalize().
void mfem::Mesh::FinalizeTriMesh | ( | int | generate_edges = 0, |
int | refine = 0, | ||
bool | fix_orientation = true ) |
void mfem::Mesh::FinalizeWedgeMesh | ( | int | generate_edges = 0, |
int | refine = 0, | ||
bool | fix_orientation = true ) |
Array< int > mfem::Mesh::FindFaceNeighbors | ( | const int | elem | ) | const |
|
virtual |
Find the ids of the elements that contain the given points, and their corresponding reference coordinates.
The DenseMatrix point_mat describes the given points - one point for each column; it should have SpaceDimension() rows.
The InverseElementTransformation object, inv_trans, is used to attempt the element transformation inversion. If NULL pointer is given, the method will use a default constructed InverseElementTransformation. Note that the algorithms in the base class InverseElementTransformation can be completely overwritten by deriving custom classes that override the Transform() method.
If no element is found for the i-th point, elem_ids[i] is set to -1.
In the ParMesh implementation, the point_mat is expected to be the same on all ranks. If the i-th point is found by multiple ranks, only one of them will mark that point as found, i.e. set its elem_ids[i] to a non-negative number; the other ranks will set their elem_ids[i] to -2 to indicate that the point was found but assigned to another rank.
Reimplemented in mfem::ParMesh.
void mfem::Mesh::GeneralRefinement | ( | const Array< int > & | el_to_refine, |
int | nonconforming = -1, | ||
int | nc_limit = 0 ) |
void mfem::Mesh::GeneralRefinement | ( | const Array< Refinement > & | refinements, |
int | nonconforming = -1, | ||
int | nc_limit = 0 ) |
Refine selected mesh elements. Refinement type can be specified for each element. The function can do conforming refinement of triangles and tetrahedra and nonconforming refinement (i.e., with hanging-nodes) of triangles, quadrilaterals and hexahedra. If 'nonconforming' = -1, suitable refinement method is selected automatically (namely, conforming refinement for triangles). Use nonconforming = 0/1 to force the method. For nonconforming refinements, nc_limit optionally specifies the maximum level of hanging nodes (unlimited by default).
|
virtual |
Reimplemented in mfem::ParMesh.
int * mfem::Mesh::GeneratePartitioning | ( | int | nparts, |
int | part_method = 1 ) |
|
inline |
|
inline |
|
inline |
|
inline |
Return pointer to the i'th boundary element object.
The index i should be in the range [0, Mesh::GetNBE())
In parallel, i is the local boundary element index which is in the same range mentioned above.
void mfem::Mesh::GetBdrElementAdjacentElement | ( | int | bdr_el, |
int & | el, | ||
int & | info ) const |
For the given boundary element, bdr_el, return its adjacent element and its info, i.e. 64*local_bdr_index+bdr_orientation.
The returned bdr_orientation is that of the boundary element relative to the respective face element.
void mfem::Mesh::GetBdrElementAdjacentElement2 | ( | int | bdr_el, |
int & | el, | ||
int & | info ) const |
Deprecated.
For the given boundary element, bdr_el, return its adjacent element and its info, i.e. 64*local_bdr_index+inverse_bdr_orientation.
The returned inverse_bdr_orientation is the inverse of the orientation of the boundary element relative to the respective face element. In other words this is the orientation of the face element relative to the boundary element.
|
inline |
|
inline |
Return the vertex index of boundary element i. (1D) Return the edge index of boundary element i. (2D) Return the face index of boundary element i. (3D)
Deprecated in favor of GetBdrElementFaceIndex().
void mfem::Mesh::GetBdrElementFace | ( | int | i, |
int * | f, | ||
int * | o ) const |
Return the index and the orientation of the vertex of bdr element i. (1D) Return the index and the orientation of the edge of bdr element i. (2D) Return the index and the orientation of the face of bdr element i. (3D)
In 2D, the returned edge orientation is 0 or 1, not +/-1 as returned by GetElementEdges/GetBdrElementEdges.
|
inline |
|
inline |
ElementTransformation * mfem::Mesh::GetBdrElementTransformation | ( | int | i | ) |
Returns a pointer to the transformation defining the i-th boundary element.
void mfem::Mesh::GetBdrElementTransformation | ( | int | i, |
IsoparametricTransformation * | ElTr ) const |
Builds the transformation defining the i-th boundary element in ElTr. ElTr must be allocated in advance and will be owned by the caller.
Element::Type mfem::Mesh::GetBdrElementType | ( | int | i | ) | const |
|
inline |
|
inline |
Deprecated in favor of GetBdrElementFaceIndex().
FaceElementTransformations * mfem::Mesh::GetBdrFaceTransformations | ( | int | BdrElemNo | ) |
void mfem::Mesh::GetBdrFaceTransformations | ( | int | BdrElemNo, |
FaceElementTransformations & | FElTr, | ||
IsoparametricTransformation & | ElTr1, | ||
IsoparametricTransformation & | ElTr2 ) const |
Variant of GetBdrFaceTransformations using a user allocated FaceElementTransformations object.
void mfem::Mesh::GetBdrPointMatrix | ( | int | i, |
DenseMatrix & | pointmat ) const |
ElementTransformation * mfem::Mesh::GetEdgeTransformation | ( | int | EdgeNo | ) |
Returns a pointer to the transformation defining the given edge element.
void mfem::Mesh::GetEdgeTransformation | ( | int | i, |
IsoparametricTransformation * | EdTr ) const |
Builds the transformation defining the i-th edge element in EdTr. EdTr must be allocated in advance and will be owned by the caller.
Table * mfem::Mesh::GetEdgeVertexTable | ( | ) | const |
void mfem::Mesh::GetEdgeVertices | ( | int | i, |
Array< int > & | vert ) const |
|
inline |
|
inline |
Return pointer to the i'th element object.
The index i should be in the range [0, Mesh::GetNE())
In parallel, i is the local element index which is in the same range mentioned above.
|
inline |
void mfem::Mesh::GetElementCenter | ( | int | i, |
Vector & | center ) |
void mfem::Mesh::GetElementColoring | ( | Array< int > & | colors, |
int | el0 = 0 ) |
|
inline |
void mfem::Mesh::GetElementJacobian | ( | int | i, |
DenseMatrix & | J, | ||
const IntegrationPoint * | ip = NULL ) |
|
inline |
real_t mfem::Mesh::GetElementSize | ( | ElementTransformation * | T, |
int | type = 0 ) const |
real_t mfem::Mesh::GetElementSize | ( | int | i, |
int | type = 0 ) |
|
protected |
Return element to edge table and the indices for the boundary edges. The entries in the table are ordered according to the order of the nodes in the elements. For example, if T is the element to edge table T(i, 0) gives the index of edge in element i that connects vertex 0 to vertex 1, etc. Returns the number of the edges.
|
protected |
ElementTransformation * mfem::Mesh::GetElementTransformation | ( | int | i | ) |
void mfem::Mesh::GetElementTransformation | ( | int | i, |
const Vector & | nodes, | ||
IsoparametricTransformation * | ElTr ) const |
Builds the transformation defining the i-th element in ElTr assuming position of the vertices/nodes are given by nodes. ElTr must be allocated in advance and will be owned by the caller.
void mfem::Mesh::GetElementTransformation | ( | int | i, |
IsoparametricTransformation * | ElTr ) const |
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and will be owned by the caller.
Element::Type mfem::Mesh::GetElementType | ( | int | i | ) | const |
|
inline |
|
inline |
Return pointer to the i'th face element object.
The index i should be in the range [0, Mesh::GetNFaces())
|
inline |
Deprecated in favor of Mesh::GetFaceGeometry.
Table * mfem::Mesh::GetFaceEdgeTable | ( | ) | const |
void mfem::Mesh::GetFaceElements | ( | int | Face, |
int * | Elem1, | ||
int * | Elem2 ) const |
|
virtual |
Variant of GetFaceElementTransformations using a user allocated FaceElementTransformations object.
Reimplemented in mfem::ParMesh.
|
virtual |
Returns (a pointer to an object containing) the following data:
1) Elem1No - the index of the first element that contains this face this is the element that has the same outward unit normal vector as the face;
2) Elem2No - the index of the second element that contains this face this element has outward unit normal vector as the face multiplied with -1;
3) Elem1, Elem2 - pointers to the ElementTransformation's of the first and the second element respectively;
4) Face - pointer to the ElementTransformation of the face;
5) Loc1, Loc2 - IntegrationPointTransformation's mapping the face coordinate system to the element coordinate system (both in their reference elements). Used to transform IntegrationPoints from face to element. More formally, let: TL1, TL2 be the transformations represented by Loc1, Loc2, TE1, TE2 - the transformations represented by Elem1, Elem2, TF - the transformation represented by Face, then TF(x) = TE1(TL1(x)) = TE2(TL2(x)) for all x in the reference face.
6) FaceGeom - the base geometry for the face.
The mask specifies which fields in the structure to return: mask & 1 - Elem1, mask & 2 - Elem2 mask & 4 - Loc1, mask & 8 - Loc2, mask & 16 - Face. These mask values are defined in the ConfigMasks enum type as part of the FaceElementTransformations class in fem/eltrans.hpp.
Reimplemented in mfem::ParMesh.
Element::Type mfem::Mesh::GetFaceElementType | ( | int | Face | ) | const |
const FaceGeometricFactors * mfem::Mesh::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.
The IntegrationRule used with GetFaceGeometricFactors needs to remain valid until the internally stored FaceGeometricFactors objects are destroyed (by either calling Mesh::DeleteGeometricFactors(), Mesh::NodesUpdated(), or the Mesh destructor).
If the device MemoryType parameter d_mt is specified, then the returned object will use that type unless it was previously allocated with a different type.
The returned pointer points to an internal object that may be invalidated by mesh operations such as refinement, vertex/node movement, etc. Since not all such modifications can be tracked by the Mesh class (e.g. when using the pointer returned by GetNodes() to change the nodes) one needs to account for such changes by calling the method NodesUpdated() which, in particular, will call DeleteGeometricFactors().
Geometry::Type mfem::Mesh::GetFaceGeometry | ( | int | i | ) | const |
Return the Geometry::Type associated with face i.
|
inline |
Deprecated in favor of Mesh::GetFaceGeometry.
Mesh::FaceInformation mfem::Mesh::GetFaceInformation | ( | int | f | ) | const |
This method aims to provide face information in a deciphered format, i.e. Mesh::FaceInformation, compared to the raw encoded information returned by Mesh::GetFaceElements() and Mesh::GetFaceInfos().
void mfem::Mesh::GetFaceInfos | ( | int | Face, |
int * | Inf1, | ||
int * | Inf2 ) const |
void mfem::Mesh::GetFaceInfos | ( | int | Face, |
int * | Inf1, | ||
int * | Inf2, | ||
int * | NCFace ) const |
Table * mfem::Mesh::GetFaceToElementTable | ( | ) | const |
ElementTransformation * mfem::Mesh::GetFaceTransformation | ( | int | FaceNo | ) |
Returns a pointer to the transformation defining the given face element.
void mfem::Mesh::GetFaceTransformation | ( | int | i, |
IsoparametricTransformation * | FTr ) const |
Builds the transformation defining the i-th face element in FTr. FTr must be allocated in advance and will be owned by the caller.
|
inline |
real_t mfem::Mesh::GetGeckoElementOrdering | ( | Array< int > & | ordering, |
int | iterations = 4, | ||
int | window = 4, | ||
int | period = 2, | ||
int | seed = 0, | ||
bool | verbose = false, | ||
real_t | time_limit = 0 ) |
This is our integration with the Gecko library. The method finds an element ordering that will increase memory coherency by putting elements that are in physical proximity closer in memory. It can also be used to obtain a space-filling curve ordering for ParNCMesh partitioning.
[out] | ordering | Output element ordering. |
iterations | Total number of V cycles. The ordering may improve with more iterations. The best iteration is returned at the end. | |
window | Initial window size. This determines the number of permutations tested at each multigrid level and strongly influences the quality of the result, but the cost of increasing 'window' is exponential. | |
period | The window size is incremented every 'period' iterations. | |
seed | Seed for initial random ordering (0 = skip random reorder). | |
verbose | Print the progress of the optimization to mfem::out. | |
time_limit | Optional time limit for the optimization, in seconds. When reached, ordering from the best iteration so far is returned (0 = no limit). |
const GeometricFactors * mfem::Mesh::GetGeometricFactors | ( | const IntegrationRule & | ir, |
const int | flags, | ||
MemoryType | d_mt = MemoryType::DEFAULT ) |
Return the mesh geometric factors corresponding to the given integration rule.
The IntegrationRule used with GetGeometricFactors needs to remain valid until the internally stored GeometricFactors objects are destroyed (by calling Mesh::DeleteGeometricFactors(), Mesh::NodesUpdated(), or the Mesh destructor).
If the device MemoryType parameter d_mt is specified, then the returned object will use that type unless it was previously allocated with a different type.
The returned pointer points to an internal object that may be invalidated by mesh operations such as refinement, vertex/node movement, etc. Since not all such modifications can be tracked by the Mesh class (e.g. when using the pointer returned by GetNodes() to change the nodes) one needs to account for such changes by calling the method NodesUpdated() which, in particular, will call DeleteGeometricFactors().
void mfem::Mesh::GetGeometricParametersFromJacobian | ( | const DenseMatrix & | J, |
real_t & | volume, | ||
Vector & | aspr, | ||
Vector & | skew, | ||
Vector & | ori ) const |
Computes geometric parameters associated with a Jacobian matrix in 2D/3D. These parameters are (1) Area/Volume, (2) Aspect-ratio (1 in 2D, and 2 non-dimensional and 2 dimensional parameters in 3D. Dimensional parameters are used for target construction in TMOP), (3) skewness (1 in 2D and 3 in 3D), and finally (4) orientation (1 in 2D and 3 in 3D).
void mfem::Mesh::GetGeometries | ( | int | dim, |
Array< Geometry::Type > & | el_geoms ) const |
|
inline |
void mfem::Mesh::GetHilbertElementOrdering | ( | Array< int > & | ordering | ) |
Return an ordering of the elements that approximately follows the Hilbert curve. The method performs a spatial (Hilbert) sort on the centers of all elements and returns the resulting sequence, which can then be passed to ReorderElements. This is a cheap alternative to GetGeckoElementOrdering.
FaceElementTransformations * mfem::Mesh::GetInteriorFaceTransformations | ( | int | FaceNo | ) |
See GetFaceElementTransformations().
void mfem::Mesh::GetInteriorFaceTransformations | ( | int | FaceNo, |
FaceElementTransformations & | FElTr, | ||
IsoparametricTransformation & | ElTr1, | ||
IsoparametricTransformation & | ElTr2 ) const |
Variant of GetInteriorFaceTransformations using a user allocated FaceElementTransformations object.
|
inline |
|
protected |
void mfem::Mesh::GetLocalFaceTransformation | ( | int | face_type, |
int | elem_type, | ||
IsoparametricTransformation & | Transf, | ||
int | info ) const |
A helper method that constructs a transformation from the reference space of a face to the reference space of an element.
The local index of the face as a face in the element and its orientation are given by the input parameter info, as info = 64*loc_face_idx + loc_face_orientation.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
inline |
|
inline |
|
inline |
|
inline |
|
virtual |
Returns the number of faces according to the requested type, does not count master nonconforming faces.
If type==Boundary returns only the number of true boundary faces contrary to GetNBE() that returns all "boundary" elements which may include actual interior faces. Similarly, if type==Interior, only the true interior faces are counted excluding all master nonconforming faces.
Reimplemented in mfem::ParMesh.
const FiniteElementSpace * mfem::Mesh::GetNodalFESpace | ( | ) | const |
Return the FiniteElementSpace on which the current mesh nodes are defined or NULL if the mesh does not have nodes.
|
inline |
Return a pointer to the internal node GridFunction (may be NULL).
If the mesh is straight-sided (low-order), it may not have a GridFunction for the nodes, in which case this function returns NULL. To ensure that the nodal GridFunction exists, first call EnsureNodes().
|
inline |
void mfem::Mesh::GetNodes | ( | GridFunction & | nodes | ) | const |
Return the mesh nodes/vertices projected on the given GridFunction.
|
inline |
Return the nodes update counter.
This counter starts at zero, and is incremented every time the geometric factors must be recomputed (e.g. on calls to Mesh::Transform, Mesh::NodesUpdated, etc.)
int mfem::Mesh::GetNumFaces | ( | ) | const |
int mfem::Mesh::GetNumFacesWithGhost | ( | ) | const |
int mfem::Mesh::GetNumGeometries | ( | int | dim | ) | const |
|
inline |
int mfem::Mesh::GetPatchAttribute | ( | int | i | ) | const |
int mfem::Mesh::GetPatchBdrAttribute | ( | int | i | ) | const |
void mfem::Mesh::GetPointMatrix | ( | int | i, |
DenseMatrix & | pointmat ) const |
|
staticprotected |
const CoarseFineTransformations & mfem::Mesh::GetRefinementTransforms | ( | ) | const |
|
inline |
Return update counter. The counter starts at zero and is incremented each time refinement, derefinement, or rebalancing method is called. It is used for checking proper sequence of Space:: and GridFunction:: Update() calls.
|
staticprotected |
|
static |
Return FiniteElement for reference element of the specified type.
|
staticprotected |
|
inline |
Return pointer to vertex i's coordinates.
|
inline |
Table * mfem::Mesh::GetVertexToElementTable | ( | ) |
void mfem::Mesh::GetVertexToVertexTable | ( | DSTable & | v_to_v | ) | const |
Return vertex to vertex table. The connections stored in the table are from smaller to bigger vertex index, i.e. if i<j and (i, j) is in the table, then (j, i) is not stored.
|
inlineprotected |
|
inlinevirtual |
Checks if the mesh has boundary elements.
Reimplemented in mfem::ParMesh.
|
inline |
Return true iff the given geom is encountered in the mesh. Geometries of dimensions lower than Dimension() are counted as well.
|
protected |
|
protected |
|
staticprotected |
|
staticprotected |
|
protected |
void mfem::Mesh::KnotInsert | ( | Array< KnotVector * > & | kv | ) |
|
inlinevirtual |
This is similar to the mesh constructor with the same arguments, but here the current mesh is destroyed and another one created based on the data stream again given in MFEM, Netgen, or VTK format. If generate_edges = 0 (default) edges are not generated, if 1 edges are generated.
Reimplemented in mfem::ParMesh, and mfem::PumiMesh.
|
protected |
|
static |
|
protected |
|
protectedvirtual |
This function is not public anymore. Use GeneralRefinement instead.
Reimplemented in mfem::ParMesh.
|
protected |
|
protected |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if type = TRIANGLE.
If generate_edges = 0 (default) edges are not generated, if 1 edges are generated. The parameter sfc_ordering controls how the elements (when type = QUADRILATERAL) are ordered: true - use space-filling curve ordering, or false - use lexicographic ordering.
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny*5 quadrilaterals.
The mesh is generated by taking nx*ny quadrilaterals and splitting each quadrilateral into 5 quadrilaterals. Each quadrilateral is projected inwards and connected to the original quadrilateral.
|
protected |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedra if type = HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if type = TETRAHEDRON.
The parameter sfc_ordering controls how the elements (when type = HEXAHEDRON) are ordered: true - use space-filling curve ordering, or false - use lexicographic ordering.
|
protected |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz*24 tetrahedrons.
The mesh is generated by taking nx*ny*nz hexahedra and splitting each hexahedron into 24 tetrahedrons. Each face of the hexahedron is split into 4 triangles (face edges are connected to a face-centered point), and the triangles are connected to a hex-centered point.
|
static |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATERAL or into 2*nx*ny triangles if type = TRIANGLE.
If generate_edges = 0 (default) edges are not generated, if 1 edges are generated. The parameter sfc_ordering controls how the elements (when type = QUADRILATERAL) are ordered: true - use space-filling curve ordering, or false - use lexicographic ordering.
|
static |
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny*5 quadrilaterals.
The mesh is generated by taking nx*ny quadrilaterals and splitting each quadrilateral into 5 quadrilaterals. Each quadrilateral is projected inwards and connected to the original quadrilateral.
|
static |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedra if type = HEXAHEDRON or into 6*nx*ny*nz tetrahedrons if type = TETRAHEDRON.
The parameter sfc_ordering controls how the elements (when type = HEXAHEDRON) are ordered: true - use space-filling curve ordering, or false - use lexicographic ordering.
|
static |
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz*24 tetrahedrons.
The mesh is generated by taking nx*ny*nz hexahedra and splitting each hexahedron into 24 tetrahedrons. Each face of the hexahedron is split into 4 triangles (face edges are connected to a face-centered point), and the triangles are connected to a hex-centered point.
Create a periodic mesh by identifying vertices of orig_mesh.
Each vertex i will be mapped to vertex v2v[i], such that all vertices that are coincident under the periodic mapping get mapped to the same index. The mapping v2v can be generated from translation vectors using Mesh::CreatePeriodicVertexMapping.
|
static |
refined ref_factors[i] times in each dimension.
Create a refined mesh, where each element of the original mesh may be refined by a different factor.
[in] | orig_mesh | The starting coarse mesh. |
[in] | ref_factors | An array of integers whose size is the number of elements of orig_mesh. The ith element of orig_mesh is refined by refinement factor ref_factors[i]. |
[in] | ref_type | Specify the positions of the new vertices. The options are BasisType::ClosedUniform or BasisType::GaussLobatto. |
The refinement data which can be accessed with GetRefinementTransforms() is set to reflect the performed refinements.
Create a refined (by any factor) version of orig_mesh.
[in] | orig_mesh | The starting coarse mesh. |
[in] | ref_factor | The refinement factor, an integer > 1. |
[in] | ref_type | Specify the positions of the new vertices. The options are BasisType::ClosedUniform or BasisType::GaussLobatto. |
The refinement data which can be accessed with GetRefinementTransforms() is set to reflect the performed refinements.
|
protected |
Internal function used in Mesh::MakeRefined.
Create a mesh by splitting each element of orig_mesh into simplices. Quadrilaterals are split into two triangles, prisms are split into 3 tetrahedra, and hexahedra are split into either 5 or 6 tetrahedra depending on the configuration.
|
protected |
|
protectedvirtual |
Reimplemented in mfem::ParMesh.
|
inline |
Get the mesh generator/type.
The purpose of this is to be able to quickly tell what type of elements one has in the mesh. Examination of this bitmask along with knowledge of the mesh dimension can be used to identify which element types are present.
In parallel, the result takes into account elements on all processors.
void mfem::Mesh::MesquiteSmooth | ( | const int | mesquite_option = 0 | ) |
void mfem::Mesh::MoveNodes | ( | const Vector & | displacements | ) |
void mfem::Mesh::MoveVertices | ( | const Vector & | displacements | ) |
Element * mfem::Mesh::NewElement | ( | int | geom | ) |
void mfem::Mesh::NewNodes | ( | GridFunction & | nodes, |
bool | make_owner = false ) |
Replace the internal node GridFunction with the given GridFunction.
Invokes NodesUpdated().
|
inline |
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().
It deletes internal quantities derived from the node coordinates, such as the (Face)GeometricFactors.
|
protectedvirtual |
NC version of GeneralDerefinement.
Reimplemented in mfem::ParMesh.
|
protectedvirtual |
This function is not public anymore. Use GeneralRefinement instead.
Reimplemented in mfem::ParMesh.
void mfem::Mesh::NURBSCoarsening | ( | int | cf = 2, |
real_t | tol = 1.0e-12 ) |
Reimplemented in mfem::ParMesh.
|
virtual |
Refine NURBS mesh, with an optional refinement factor, generally anisotropic.
[in] | rf | Optional refinement factor. If scalar, the factor is used for all dimensions. If an array, factors can be specified for each dimension. The factor multiplies the number of elements in each dimension. Some factors can be 1. |
[in] | tol | NURBS geometry deviation tolerance, cf. Algorithm A5.8 of "The NURBS Book", 2nd ed, Piegl and Tiller. |
Reimplemented in mfem::ParMesh.
Explicitly delete the copy assignment operator.
|
inline |
|
virtual |
Print the mesh to the given stream using the adios2 bp format.
Reimplemented in mfem::ParMesh.
|
inlinevirtual |
Print the mesh to the given stream using the default MFEM mesh format.
Reimplemented in mfem::ParMesh.
void mfem::Mesh::PrintBdrVTU | ( | std::string | fname, |
VTKFormat | format = VTKFormat::ASCII, | ||
bool | high_order_output = false, | ||
int | compression_level = 0 ) |
void mfem::Mesh::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.
If Vh or Vk are not NULL, return the element sizes and aspect ratios for all elements in the given Vectors.
|
staticprotected |
|
static |
Auxiliary method used by PrintCharacteristics().
It is also used in the mesh-explorer
miniapp.
void mfem::Mesh::PrintElementsWithPartitioning | ( | int * | partitioning, |
std::ostream & | os, | ||
int | interior_faces = 0 ) |
|
staticprotected |
|
protected |
If NURBS mesh, write NURBS format. If NCMesh, write mfem v1.1 format. If section_delimiter is empty, write mfem v1.0 format. Otherwise, write mfem v1.2 format with the given section_delimiter at the end. If comments is non-empty, it will be printed after the first line of the file, and each line should begin with '#'.
|
inlinevirtual |
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the parallel decomposition is also printed.
Reimplemented in mfem::ParMesh.
void mfem::Mesh::PrintSurfaces | ( | const Table & | Aface_face, |
std::ostream & | os ) const |
|
protected |
Write the beginning of a NURBS mesh to os, specifying the NURBS patch topology. Optional file comments can be provided in comments.
[in] | os | Output stream to which to write. |
[in] | e_to_k | Map from edge to signed knotvector indices. |
[in] | version | NURBS mesh version number times 10 (e.g. 11 for v1.1). |
[in] | comment | Optional comment string, written after version line. |
void mfem::Mesh::PrintVTK | ( | std::ostream & | os | ) |
Print the mesh in VTK format (linear and quadratic meshes only).
void mfem::Mesh::PrintVTK | ( | std::ostream & | os, |
int | ref, | ||
int | field_data = 0 ) |
Print the mesh in VTK format. The parameter ref > 0 specifies an element subdivision number (useful for high order fields and curved meshes). If the optional field_data is set, we also add a FIELD section in the beginning of the file with additional dataset information.
void mfem::Mesh::PrintVTU | ( | std::ostream & | os, |
int | ref = 1, | ||
VTKFormat | format = VTKFormat::ASCII, | ||
bool | high_order_output = false, | ||
int | compression_level = 0, | ||
bool | bdr_elements = false ) |
Print the mesh in VTU format. The parameter ref > 0 specifies an element subdivision number (useful for high order fields and curved meshes). If bdr_elements is true, then output (only) the boundary elements, otherwise output only the non-boundary elements.
|
virtual |
Print the mesh in VTU format with file name fname.
Reimplemented in mfem::ParMesh.
void mfem::Mesh::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 of subdomain i has boundary attribute i+1.
|
virtual |
Print the mesh to the given stream using Netgen/Truegrid format.
Reimplemented in mfem::ParMesh.
void mfem::Mesh::RandomRefinement | ( | real_t | prob, |
bool | aniso = false, | ||
int | nonconforming = -1, | ||
int | nc_limit = 0 ) |
|
protected |
Load a mesh from a Genesis file.
Definition at line 3570 of file mesh_readers.cpp.
|
protected |
|
protected |
|
protected |
The following mappings convert the Gmsh node orderings for high order elements to MFEM's L2 degree of freedom ordering. To support more options examine Gmsh's ordering and read off the indices in MFEM's order. For example 2nd order Gmsh quadrilaterals use the following ordering:
3–6–2 | | | 7 8 5 | | | 0–4–1
(from https://gmsh.info/doc/texinfo/gmsh.html#Node-ordering)
Whereas MFEM uses a tensor product ordering with the horizontal axis cycling fastest so we would read off:
0 4 1 7 8 5 3 6 2
This corresponds to the quad9 mapping below.
Definition at line 1509 of file mesh_readers.cpp.
|
protected |
Definition at line 1349 of file mesh_readers.cpp.
|
protected |
Definition at line 132 of file mesh_readers.cpp.
|
protected |
Definition at line 39 of file mesh_readers.cpp.
|
protected |
Definition at line 166 of file mesh_readers.cpp.
|
protected |
Definition at line 228 of file mesh_readers.cpp.
|
protected |
Definition at line 1307 of file mesh_readers.cpp.
|
protected |
Definition at line 282 of file mesh_readers.cpp.
|
protected |
Definition at line 1176 of file mesh_readers.cpp.
|
protected |
Definition at line 1058 of file mesh_readers.cpp.
|
inlineprotected |
|
inlinevirtual |
Utility function: sum integers from all processors (Allreduce).
Reimplemented in mfem::ParMesh.
void mfem::Mesh::RefineNURBSFromFile | ( | std::string | ref_file | ) |
Refine a NURBS mesh with the knots specified in the file named ref_file. The file has the number of knot vectors on the first line. It is the same number of knot vectors specified in the NURBS mesh in the section edges. Then for each knot vector specified in the section edges (with the same ordering), a line describes (in this order): 1) an integer giving the number of knots inserted, 2) the knots inserted as a double. The advantage of this method is that it is possible to specifically refine a coarse NURBS mesh without changing the mesh file itself. Examples in miniapps/nurbs/meshes.
void mfem::Mesh::RemoveInternalBoundaries | ( | ) |
void mfem::Mesh::RemoveUnusedVertices | ( | ) |
void mfem::Mesh::ReorderElements | ( | const Array< int > & | ordering, |
bool | reorder_vertices = true ) |
|
virtual |
This method modifies a tetrahedral mesh so that Nedelec spaces of order greater than 1 can be defined on the mesh. Specifically, we 1) rotate all tets in the mesh so that the vertices {v0, v1, v2, v3} satisfy: v0 < v1 < min(v2, v3). 2) rotate all boundary triangles so that the vertices {v0, v1, v2} satisfy: v0 < min(v1, v2).
Reimplemented in mfem::ParMesh.
|
virtual |
Save the mesh to a file using Mesh::Print. The given precision will be used for ASCII output.
Reimplemented in mfem::ParMesh.
|
inline |
|
virtual |
Determine the sets of unique attribute values in domain and boundary elements.
Separately scan the domain and boundary elements to generate unique, sorted sets of the element attribute values present in the mesh and store these in the Mesh::attributes and Mesh::bdr_attributes arrays.
Reimplemented in mfem::ParMesh.
|
inline |
|
virtual |
Set the curvature of the mesh nodes using the given polynomial degree.
Creates a nodal GridFunction if one doesn't already exist.
[in] | order | Polynomial degree of the nodal FE space. If this value is <= 0 then the method will remove the nodal GridFunction and the Mesh will use the vertices array instead; the other arguments are ignored in this case. |
[in] | discont | Whether to use a discontinuous or continuous finite element space (continuous is default). |
[in] | space_dim | The space dimension (optional). |
[in] | ordering | The Ordering of the finite element space (Ordering::byVDIM is the default). |
Reimplemented in mfem::ParMesh.
|
protected |
Determine the mesh generator bitmask meshgen, see MeshGenerator().
Also, initializes mesh_geoms.
|
virtual |
Replace the internal node GridFunction with a new GridFunction defined on the given FiniteElementSpace. The new node coordinates are projected (derived) from the current nodes/vertices.
Reimplemented in mfem::ParMesh.
void mfem::Mesh::SetNodalGridFunction | ( | GridFunction * | nodes, |
bool | make_owner = false ) |
Replace the internal node GridFunction with the given GridFunction. The given GridFunction is updated with node coordinates projected (derived) from the current nodes/vertices.
void mfem::Mesh::SetNodes | ( | const Vector & | node_coord | ) |
Updates the vertex/node locations. Invokes NodesUpdated().
|
inline |
void mfem::Mesh::SetPatchAttribute | ( | int | i, |
int | attr ) |
void mfem::Mesh::SetPatchBdrAttribute | ( | int | i, |
int | attr ) |
|
protected |
|
inline |
void mfem::Mesh::Swap | ( | Mesh & | other, |
bool | non_geometry ) |
void mfem::Mesh::SwapNodes | ( | GridFunction *& | nodes, |
int & | own_nodes_ ) |
Swap the internal node GridFunction pointer and ownership flag members with the given ones.
Invokes NodesUpdated().
void mfem::Mesh::Transform | ( | VectorCoefficient & | deformation | ) |
|
static |
For the vertex (1D), edge (2D), or face (3D) of a boundary element with the orientation o, return the transformation of the boundary element integration point @ ip to the face element. In 2D, the the orientation is 0 or 1 as returned by GetBdrElementFace, not +/-1. Supports both internal and external boundaries.
|
protected |
Uniform Refinement. Element with index i is refined uniformly.
void mfem::Mesh::UniformRefinement | ( | int | ref_algo = 0 | ) |
Refine all mesh elements.
[in] | ref_algo | Refinement algorithm. Currently used only for pure tetrahedral meshes. If set to zero (default), a tet mesh will be refined using algorithm A, that produces elements with better quality compared to algorithm B used when the parameter is non-zero. |
For tetrahedral meshes, after using algorithm A, the mesh cannot be refined locally using methods like GeneralRefinement() unless it is re-finalized using Finalize() with the parameter refine set to true. Note that calling Finalize() in this way will generally invalidate any FiniteElementSpaces and GridFunctions defined on the mesh.
|
inlineprotectedvirtual |
Refine a mixed 2D mesh uniformly.
Reimplemented in mfem::ParMesh.
|
protected |
|
inlineprotectedvirtual |
Refine a mixed 3D mesh uniformly.
Reimplemented in mfem::ParMesh.
|
protected |
Update the nodes of a curved mesh after the topological part of a Mesh::Operation, such as refinement, has been performed.
If Nodes GridFunction is defined, i.e. not NULL, this method calls NodesUpdated().
|
friend |
|
friend |
|
friend |
AttributeSets mfem::Mesh::attribute_sets |
Array<int> mfem::Mesh::attributes |
AttributeSets mfem::Mesh::bdr_attribute_sets |
Array<int> mfem::Mesh::bdr_attributes |
|
protected |
|
mutableprotected |
|
protected |
Array<FaceGeometricFactors*> mfem::Mesh::face_geom_factors |
|
protected |
|
protected |
Array<GeometricFactors*> mfem::Mesh::geom_factors |
|
mutableprotected |
|
protected |
NCMesh* mfem::Mesh::ncmesh |
|
protected |
|
protected |
NURBSExtension* mfem::Mesh::NURBSext |
|
protected |
|
protected |
|
protected |
|
staticprotected |
|
staticprotected |
|
staticprotected |
|
staticprotected |