MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
mfem::SubMesh Class Reference

Subdomain representation of a topological parent in another Mesh. More...

#include <submesh.hpp>

Inheritance diagram for mfem::SubMesh:
[legend]
Collaboration diagram for mfem::SubMesh:
[legend]

Public Types

enum  From { Domain, Boundary }
 Indicator from which part of the parent Mesh the SubMesh is created. More...
 
- Public Types inherited from mfem::Mesh
enum  Operation { NONE, REFINE, DEREFINE, REBALANCE }
 
enum  FaceTopology { FaceTopology::Boundary, FaceTopology::Conforming, FaceTopology::Nonconforming, FaceTopology::NA }
 
enum  ElementLocation { ElementLocation::Local, ElementLocation::FaceNbr, ElementLocation::NA }
 
enum  ElementConformity { ElementConformity::Coincident, ElementConformity::Superset, ElementConformity::Subset, ElementConformity::NA }
 
enum  FaceInfoTag {
  FaceInfoTag::Boundary, FaceInfoTag::LocalConforming, FaceInfoTag::LocalSlaveNonconforming, FaceInfoTag::SharedConforming,
  FaceInfoTag::SharedSlaveNonconforming, FaceInfoTag::MasterNonconforming, FaceInfoTag::GhostSlave, FaceInfoTag::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

 SubMesh ()=delete
 
const MeshGetParent () const
 Get the parent Mesh object. More...
 
From GetFrom () const
 Get the From indicator. More...
 
const Array< int > & GetParentElementIDMap () const
 Get the parent element id map. More...
 
const Array< int > & GetParentFaceIDMap () const
 Get the face id map. More...
 
const Array< int > & GetParentVertexIDMap () const
 Get the parent vertex id map. More...
 
 ~SubMesh ()
 
- Public Member Functions inherited from mfem::Mesh
 Mesh ()
 
 Mesh (const Mesh &mesh, bool copy_nodes=true)
 
 Mesh (Mesh &&mesh)
 Move constructor, useful for using a Mesh as a function return value. More...
 
Meshoperator= (Mesh &&mesh)
 Move assignment operstor. More...
 
Meshoperator= (const Mesh &mesh)=delete
 Explicitly delete the copy assignment operator. More...
 
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 given translations are identified. More...
 
 Mesh (double *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. More...
 
 Mesh (int Dim_, int NVert, int NElem, int NBdrElem=0, int spaceDim_=-1)
 Init constructor: begin the construction of a Mesh object. More...
 
void FinalizeTopology (bool generate_bdr=true)
 Finalize the construction of the secondary topology (connectivity) data of a Mesh. More...
 
virtual void Finalize (bool refine=false, bool fix_orientation=false)
 Finalize the construction of a general Mesh. More...
 
virtual void SetAttributes ()
 
double GetGeckoElementOrdering (Array< int > &ordering, int iterations=4, int window=4, int period=2, int seed=0, bool verbose=false, double time_limit=0)
 
void GetHilbertElementOrdering (Array< int > &ordering)
 
void ReorderElements (const Array< int > &ordering, bool reorder_vertices=true)
 
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. More...
 
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. More...
 
MFEM_DEPRECATED Mesh (int n, double sx=1.0)
 Deprecated: see MakeCartesian1D. More...
 
 Mesh (const char *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)
 Create a disjoint mesh from the given mesh array. More...
 
MFEM_DEPRECATED Mesh (Mesh *orig_mesh, int ref_factor, int ref_type)
 Deprecated: see MakeRefined. More...
 
virtual void Load (std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
 
void Clear ()
 Clear the contents of the Mesh. More...
 
int MeshGenerator ()
 Get the mesh generator/type. More...
 
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. More...
 
int GetNE () const
 Returns number of elements. More...
 
int GetNBE () const
 Returns number of boundary elements. More...
 
int GetNEdges () const
 Return the number of edges. More...
 
int GetNFaces () const
 Return the number of faces in a 3D mesh. More...
 
int GetNumFaces () const
 Return the number of faces (3D), edges (2D) or vertices (1D). More...
 
int GetNumFacesWithGhost () const
 Return the number of faces (3D), edges (2D) or vertices (1D) including ghost faces. More...
 
virtual int GetNFbyType (FaceType type) const
 Returns the number of faces according to the requested type, does not count master nonconforming faces. More...
 
virtual long long ReduceInt (int value) const
 Utility function: sum integers from all processors (Allreduce). More...
 
long long GetGlobalNE () const
 Return the total (global) number of elements. More...
 
void GetVertexToVertexTable (DSTable &) const
 
const GeometricFactorsGetGeometricFactors (const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
 Return the mesh geometric factors corresponding to the given integration rule. More...
 
const FaceGeometricFactorsGetFaceGeometricFactors (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. More...
 
void DeleteGeometricFactors ()
 Destroy all GeometricFactors stored by the Mesh. More...
 
void NodesUpdated ()
 This function should be called after the mesh node coordinates have changed, e.g. after the mesh has moved. More...
 
int EulerNumber () const
 Equals 1 + num_holes - num_loops. More...
 
int EulerNumber2D () const
 Equals 1 - num_holes. More...
 
int Dimension () const
 
int SpaceDimension () const
 
const double * GetVertex (int i) const
 Return pointer to vertex i's coordinates. More...
 
double * GetVertex (int i)
 Return pointer to vertex i's coordinates. More...
 
void GetElementData (int geom, Array< int > &elem_vtx, Array< int > &attr) const
 
virtual bool HasBoundaryElements () const
 Checks if the mesh has boundary elements. More...
 
void GetBdrElementData (int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
 
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 pointer. More...
 
const Element *const * GetElementsArray () const
 
const ElementGetElement (int i) const
 
ElementGetElement (int i)
 
const ElementGetBdrElement (int i) const
 
ElementGetBdrElement (int i)
 
const ElementGetFace (int i) const
 
Geometry::Type GetFaceGeometry (int i) const
 
Geometry::Type GetElementGeometry (int i) const
 
Geometry::Type GetBdrElementGeometry (int i) const
 
Geometry::Type GetFaceBaseGeometry (int i) const
 
Geometry::Type GetElementBaseGeometry (int i) const
 
Geometry::Type GetBdrElementBaseGeometry (int i) const
 
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. More...
 
int GetNumGeometries (int dim) const
 Return the number of geometries of the given dimension present in the mesh. More...
 
void GetGeometries (int dim, Array< Geometry::Type > &el_geoms) const
 Return all element geometries of the given dimension present in the mesh. More...
 
void GetElementVertices (int i, Array< int > &v) const
 Returns the indices of the vertices of element i. More...
 
void GetBdrElementVertices (int i, Array< int > &v) const
 Returns the indices of the vertices of boundary element i. More...
 
void GetElementEdges (int i, Array< int > &edges, Array< int > &cor) const
 Return the indices and the orientations of all edges of element i. More...
 
void GetBdrElementEdges (int i, Array< int > &edges, Array< int > &cor) const
 Return the indices and the orientations of all edges of bdr element i. More...
 
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. More...
 
void GetEdgeVertices (int i, Array< int > &vert) const
 Returns the indices of the vertices of edge i. More...
 
TableGetFaceEdgeTable () const
 Returns the face-to-edge Table (3D) More...
 
TableGetEdgeVertexTable () const
 Returns the edge-to-vertex Table (3D) More...
 
void GetElementFaces (int i, Array< int > &faces, Array< int > &ori) const
 Return the indices and the orientations of all faces of element i. More...
 
void GetBdrElementFace (int i, int *f, int *o) const
 Return the index and the orientation of the face of bdr element i. (3D) More...
 
int GetBdrElementEdgeIndex (int i) 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. More...
 
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.e. 64*local_bdr_index+inverse_bdr_orientation. More...
 
Element::Type GetElementType (int i) const
 Returns the type of element i. More...
 
Element::Type GetBdrElementType (int i) const
 Returns the type of boundary element i. More...
 
void GetPointMatrix (int i, DenseMatrix &pointmat) const
 
void GetBdrPointMatrix (int i, DenseMatrix &pointmat) const
 
void GetElementTransformation (int i, IsoparametricTransformation *ElTr)
 
ElementTransformationGetElementTransformation (int i)
 Returns the transformation defining the i-th element. More...
 
void GetElementTransformation (int i, const Vector &nodes, IsoparametricTransformation *ElTr)
 
ElementTransformationGetBdrElementTransformation (int i)
 Returns the transformation defining the i-th boundary element. More...
 
void GetBdrElementTransformation (int i, IsoparametricTransformation *ElTr)
 
void GetFaceTransformation (int i, IsoparametricTransformation *FTr)
 Returns the transformation defining the given face element in a user-defined variable. More...
 
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 space of an element. More...
 
ElementTransformationGetFaceTransformation (int FaceNo)
 Returns the transformation defining the given face element. More...
 
void GetEdgeTransformation (int i, IsoparametricTransformation *EdTr)
 
ElementTransformationGetEdgeTransformation (int EdgeNo)
 Returns the transformation defining the given face element. More...
 
virtual
FaceElementTransformations
GetFaceElementTransformations (int FaceNo, int mask=31)
 
FaceElementTransformationsGetInteriorFaceTransformations (int FaceNo)
 
FaceElementTransformationsGetBdrFaceTransformations (int BdrElemNo)
 
int GetBdrFace (int BdrElemNo) const
 Return the local face index for the given boundary face. More...
 
bool FaceIsInterior (int FaceNo) const
 Return true if the given face is interior. More...
 
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
 
Geometry::Type GetFaceGeometryType (int Face) const
 
Element::Type GetFaceElementType (int Face) const
 
Array< int > GetFaceToBdrElMap () const
 
int CheckElementOrientation (bool fix_it=true)
 Check (and optionally attempt to fix) the orientation of the elements. More...
 
int CheckBdrElementOrientation (bool fix_it=true)
 Check the orientation of the boundary elements. More...
 
int GetAttribute (int i) const
 Return the attribute of element i. More...
 
void SetAttribute (int i, int attr)
 Set the attribute of element i. More...
 
int GetBdrAttribute (int i) const
 Return the attribute of boundary element i. More...
 
void SetBdrAttribute (int i, int attr)
 Set the attribute of boundary element i. More...
 
const TableElementToElementTable ()
 
const TableElementToFaceTable () const
 
const TableElementToEdgeTable () const
 
TableGetVertexToElementTable ()
 The returned Table must be destroyed by the caller. More...
 
TableGetFaceToElementTable () const
 
virtual MFEM_DEPRECATED void ReorientTetMesh ()
 
int * CartesianPartitioning (int nxyz[])
 
int * GeneratePartitioning (int nparts, int part_method=1)
 
void CheckPartitioning (int *partitioning_)
 
void CheckDisplacements (const Vector &displacements, double &tmax)
 
void MoveVertices (const Vector &displacements)
 
void GetVertices (Vector &vert_coord) const
 
void SetVertices (const Vector &vert_coord)
 
void GetNode (int i, double *coord) const
 
void SetNode (int i, const double *coord)
 
void MoveNodes (const Vector &displacements)
 
void GetNodes (Vector &node_coord) const
 
void SetNodes (const Vector &node_coord)
 
GridFunctionGetNodes ()
 Return a pointer to the internal node GridFunction (may be NULL). More...
 
const GridFunctionGetNodes () const
 
bool OwnsNodes () const
 Return the mesh nodes ownership flag. More...
 
void SetNodesOwner (bool nodes_owner)
 Set the mesh nodes ownership flag. More...
 
void NewNodes (GridFunction &nodes, bool make_owner=false)
 Replace the internal node GridFunction with the given GridFunction. More...
 
void SwapNodes (GridFunction *&nodes, int &own_nodes_)
 
void GetNodes (GridFunction &nodes) const
 Return the mesh nodes/vertices projected on the given GridFunction. More...
 
virtual void SetNodalFESpace (FiniteElementSpace *nfes)
 
void SetNodalGridFunction (GridFunction *nodes, bool make_owner=false)
 
const FiniteElementSpaceGetNodalFESpace () const
 
void EnsureNodes ()
 
virtual void SetCurvature (int order, bool discont=false, int space_dim=-1, int ordering=1)
 
void UniformRefinement (int ref_algo=0)
 Refine all mesh elements. More...
 
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 (double prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
 Refine each element with given probability. Uses GeneralRefinement. More...
 
void RefineAtVertex (const Vertex &vert, double eps=0.0, int nonconforming=-1)
 Refine elements sharing the specified vertex. Uses GeneralRefinement. More...
 
bool RefineByError (const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
 
bool RefineByError (const Vector &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
 
bool DerefineByError (Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
 
bool DerefineByError (const Vector &elem_error, double threshold, int nc_limit=0, int op=1)
 Same as DerefineByError for an error vector. More...
 
void EnsureNCMesh (bool simplices_nonconforming=false)
 
bool Conforming () const
 
bool Nonconforming () const
 
const CoarseFineTransformationsGetRefinementTransforms ()
 
Operation GetLastOperation () const
 Return type of last modification of the mesh. More...
 
long GetSequence () const
 
virtual void PrintXG (std::ostream &os=mfem::out) const
 Print the mesh to the given stream using Netgen/Truegrid format. More...
 
virtual void Print (std::ostream &os=mfem::out) const
 
virtual void Save (const char *fname, int precision=16) const
 
virtual void Print (adios2stream &os) const
 Print the mesh to the given stream using the adios2 bp format. More...
 
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 GetElementColoring (Array< int > &colors, int el0=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. More...
 
void PrintElementsWithPartitioning (int *partitioning, std::ostream &out, int interior_faces=0)
 
void PrintSurfaces (const Table &Aface_face, std::ostream &out) const
 Print set of disjoint surfaces: More...
 
void ScaleSubdomains (double sf)
 
void ScaleElements (double sf)
 
void Transform (void(*f)(const Vector &, Vector &))
 
void Transform (VectorCoefficient &deformation)
 
void RemoveUnusedVertices ()
 Remove unused vertices and rebuild mesh connectivity. More...
 
void RemoveInternalBoundaries ()
 
double GetElementSize (int i, int type=0)
 Get the size of the i-th element relative to the perfect reference element. More...
 
double GetElementSize (int i, const Vector &dir)
 
double GetElementVolume (int i)
 
void GetElementCenter (int i, Vector &center)
 
void GetBoundingBox (Vector &min, Vector &max, int ref=2)
 Returns the minimum and maximum corners of the mesh bounding box. More...
 
void GetCharacteristics (double &h_min, double &h_max, double &kappa_min, double &kappa_max, Vector *Vh=NULL, Vector *Vk=NULL)
 
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. More...
 
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. More...
 
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. More...
 
void Swap (Mesh &other, bool non_geometry)
 
virtual ~Mesh ()
 Destroys Mesh. More...
 
void DebugDump (std::ostream &out) const
 Output an NCMesh-compatible debug dump. More...
 
ElementNewElement (int geom)
 
int AddVertex (double x, double y=0.0, double z=0.0)
 
int AddVertex (const double *coords)
 
void AddVertexParents (int i, int p1, int p2)
 Mark vertex i as nonconforming, with parent vertices p1 and p2. More...
 
int AddSegment (int v1, int v2, int attr=1)
 
int AddSegment (const int *vi, int attr=1)
 
int AddTriangle (int v1, int v2, int v3, int attr=1)
 
int AddTriangle (const int *vi, int attr=1)
 
int AddTri (const int *vi, int attr=1)
 
int AddQuad (int v1, int v2, int v3, int v4, int attr=1)
 
int AddQuad (const int *vi, int attr=1)
 
int AddTet (int v1, int v2, int v3, int v4, int attr=1)
 
int AddTet (const int *vi, int attr=1)
 
int AddWedge (int v1, int v2, int v3, int v4, int v5, int v6, int attr=1)
 
int AddWedge (const int *vi, int attr=1)
 
int AddPyramid (int v1, int v2, int v3, int v4, int v5, int attr=1)
 
int AddPyramid (const int *vi, int attr=1)
 
int AddHex (int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
 
int AddHex (const int *vi, int attr=1)
 
void AddHexAsTets (const int *vi, int attr=1)
 
void AddHexAsWedges (const int *vi, int attr=1)
 
void AddHexAsPyramids (const int *vi, int attr=1)
 
int AddElement (Element *elem)
 The parameter elem should be allocated using the NewElement() method. More...
 
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)
 
void GenerateBoundaryElements ()
 
void FinalizeTriMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 Finalize the construction of a triangular Mesh. More...
 
void FinalizeQuadMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 Finalize the construction of a quadrilateral Mesh. More...
 
void FinalizeTetMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 Finalize the construction of a tetrahedral Mesh. More...
 
void FinalizeWedgeMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 Finalize the construction of a wedge Mesh. More...
 
void FinalizeHexMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 Finalize the construction of a hexahedral Mesh. More...
 
void FinalizeMesh (int refine=0, bool fix_orientation=true)
 Finalize the construction of any type of Mesh. More...
 
void KnotInsert (Array< KnotVector * > &kv)
 
void KnotInsert (Array< Vector * > &kv)
 
void DegreeElevate (int rel_degree, int degree=16)
 

Static Public Member Functions

static SubMesh CreateFromDomain (const Mesh &parent, Array< int > domain_attributes)
 Create a domain SubMesh from its parent. More...
 
static SubMesh CreateFromBoundary (const Mesh &parent, Array< int > boundary_attributes)
 Create a surface SubMesh from its parent. More...
 
static void Transfer (const GridFunction &src, GridFunction &dst)
 Transfer the dofs of a GridFunction. More...
 
static TransferMap CreateTransferMap (const GridFunction &src, const GridFunction &dst)
 Create a Transfer Map object. More...
 
static bool IsSubMesh (const Mesh *m)
 Check if Mesh m is a SubMesh. More...
 
- Static Public Member Functions inherited from mfem::Mesh
static FiniteElementGetTransformationFEforElementType (Element::Type)
 
static void PrintElementsByGeometry (int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
 Auxiliary method used by PrintCharacteristics(). More...
 
static Mesh LoadFromFile (const char *filename, int generate_edges=0, int refine=1, bool fix_orientation=true)
 
static Mesh MakeCartesian1D (int n, double sx=1.0)
 
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)
 
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)
 
static Mesh MakeRefined (Mesh &orig_mesh, int ref_factor, int ref_type)
 Create a refined (by any factor) version of orig_mesh. More...
 
static Mesh MakeRefined (Mesh &orig_mesh, const Array< int > &ref_factors, int ref_type)
 refined ref_factors[i] times in each dimension. More...
 
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. More...
 

Static Public Attributes

static const int GENERATED_ATTRIBUTE = 900
 
- Static Public Attributes inherited from mfem::Mesh
static bool remove_unused_vertices = true
 

Additional Inherited Members

- Public Attributes inherited from mfem::Mesh
Array< int > attributes
 A list of all unique element attributes used by the Mesh. More...
 
Array< int > bdr_attributes
 A list of all unique boundary attributes used by the Mesh. More...
 
NURBSExtensionNURBSext
 Optional NURBS mesh extension. More...
 
NCMeshncmesh
 Optional nonconforming mesh extension. More...
 
Array< GeometricFactors * > geom_factors
 Optional geometric factors. More...
 
Array< FaceGeometricFactors * > face_geom_factors
 Optional face geometric factors. More...
 
- Protected Member Functions inherited from mfem::Mesh
void Init ()
 
void InitTables ()
 
void SetEmpty ()
 
void DestroyTables ()
 
void DeleteTables ()
 
void DestroyPointers ()
 
void Destroy ()
 
void ResetLazyData ()
 
ElementReadElementWithoutAttr (std::istream &)
 
ElementReadElement (std::istream &)
 
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)
 
void ReadInlineMesh (std::istream &input, bool generate_edges=false)
 
void ReadGmshMesh (std::istream &input, int &curved, int &read_gf)
 
void ReadCubit (const char *filename, int &curved, int &read_gf)
 
void SetMeshGen ()
 Determine the mesh generator bitmask meshgen, see MeshGenerator(). More...
 
double GetLength (int i, int j) const
 Return the length of the segment from node i to node j. More...
 
void GetElementJacobian (int i, DenseMatrix &J)
 
void MarkForRefinement ()
 
void MarkTriMeshForRefinement ()
 
void GetEdgeOrdering (DSTable &v_to_v, Array< int > &order)
 
virtual void MarkTetMeshForRefinement (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)
 
STable3DGetFacesTable ()
 
STable3DGetElementToFaceTable (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. More...
 
void Bisection (int i, HashTable< Hashed2 > &)
 Bisect a tetrahedron: element with index i is bisected. More...
 
void BdrBisection (int i, const HashTable< Hashed2 > &)
 Bisect a boundary triangle: boundary element with index i is bisected. More...
 
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]. More...
 
void InitRefinementTransforms ()
 
int FindCoarseElement (int i)
 
void UpdateNodes ()
 Update the nodes of a curved mesh after refinement. More...
 
void SetVerticesFromNodes (const GridFunction *nodes)
 Helper to set vertex coordinates given a high-order curvature function. More...
 
void UniformRefinement2D_base (bool update_nodes=true)
 
virtual void UniformRefinement2D ()
 Refine a mixed 2D mesh uniformly. More...
 
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. More...
 
virtual void NURBSUniformRefinement ()
 Refine NURBS mesh. More...
 
virtual void LocalRefinement (const Array< int > &marked_el, int type=3)
 This function is not public anymore. Use GeneralRefinement instead. More...
 
virtual void NonconformingRefinement (const Array< Refinement > &refinements, int nc_limit=0)
 This function is not public anymore. Use GeneralRefinement instead. More...
 
virtual bool NonconformingDerefinement (Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
 NC version of GeneralDerefinement. More...
 
double AggregateError (const Array< double > &elem_error, const int *fine, int nfine, int op)
 Derefinement helper. More...
 
void LoadPatchTopo (std::istream &input, Array< int > &edge_to_knot)
 Read NURBS patch/macro-element mesh. More...
 
void UpdateNURBS ()
 
void PrintTopo (std::ostream &out, const Array< int > &e_to_k) const
 
void GetLocalPtToSegTransformation (IsoparametricTransformation &, int)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalSegToTriTransformation (IsoparametricTransformation &loc, int i)
 
void GetLocalSegToQuadTransformation (IsoparametricTransformation &loc, int i)
 
void GetLocalTriToTetTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalTriToWdgTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalTriToPyrTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalQuadToHexTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalQuadToWdgTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalQuadToPyrTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void ApplyLocalSlaveTransformation (FaceElementTransformations &FT, const FaceInfo &fi, bool is_ghost)
 
bool IsSlaveFace (const FaceInfo &fi) const
 
int GetElementToEdgeTable (Table &, Array< int > &)
 
void AddPointFaceElement (int lf, int gf, int el)
 Used in GenerateFaces() More...
 
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. More...
 
void FinalizeCheck ()
 
void Loader (std::istream &input, int generate_edges=0, std::string parse_tag="")
 
void Printer (std::ostream &out=mfem::out, std::string section_delimiter="") const
 
void Make3D (int nx, int ny, int nz, Element::Type type, double sx, double sy, double sz, bool sfc_ordering)
 
void Make2D (int nx, int ny, Element::Type type, double sx, double sy, bool generate_edges, bool sfc_ordering)
 
void Make1D (int n, double sx=1.0)
 Creates a 1D mesh for the interval [0,sx] divided into n equal intervals. More...
 
void MakeRefined_ (Mesh &orig_mesh, const Array< int > ref_factors, int ref_type)
 Internal function used in Mesh::MakeRefined. More...
 
void InitFromNCMesh (const NCMesh &ncmesh)
 Initialize vertices/elements/boundary/tables from a nonconforming mesh. More...
 
 Mesh (const NCMesh &ncmesh)
 Create from a nonconforming mesh. More...
 
void GetElementData (const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
 
double GetElementSize (ElementTransformation *T, int type=0)
 
void MakeSimplicial_ (const Mesh &orig_mesh, int *vglobal)
 
- Static Protected Member Functions inherited from mfem::Mesh
static void PrintElementWithoutAttr (const Element *, std::ostream &)
 
static void PrintElement (const Element *, std::ostream &)
 
static int GetTriOrientation (const int *base, const int *test)
 Returns the orientation of "test" relative to "base". More...
 
static int GetQuadOrientation (const int *base, const int *test)
 Returns the orientation of "test" relative to "base". More...
 
static int GetTetOrientation (const int *base, const int *test)
 Returns the orientation of "test" relative to "base". More...
 
static void GetElementArrayEdgeTable (const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
 
- Protected Attributes inherited from mfem::Mesh
int Dim
 
int spaceDim
 
int NumOfVertices
 
int NumOfElements
 
int NumOfBdrElements
 
int NumOfEdges
 
int NumOfFaces
 
int nbInteriorFaces
 
int nbBoundaryFaces
 
int meshgen
 
int mesh_geoms
 
long sequence
 
Array< Element * > elements
 
Array< Vertexvertices
 
Array< Element * > boundary
 
Array< Element * > faces
 
Array< FaceInfofaces_info
 
Array< NCFaceInfonc_faces_info
 
Tableel_to_edge
 
Tableel_to_face
 
Tableel_to_el
 
Array< int > be_to_edge
 
Tablebel_to_edge
 
Array< int > be_to_face
 
Tableface_edge
 
Tableedge_vertex
 
IsoparametricTransformation Transformation
 
IsoparametricTransformation Transformation2
 
IsoparametricTransformation BdrTransformation
 
IsoparametricTransformation FaceTransformation
 
IsoparametricTransformation EdgeTransformation
 
FaceElementTransformations FaceElemTr
 
CoarseFineTransformations CoarseFineTr
 
GridFunctionNodes
 
int own_nodes
 
MemAlloc< Tetrahedron, 1024 > TetMemory
 
Array< Triple< int, int, int > > tmp_vertex_parents
 
Operation last_operation
 
- Static Protected Attributes inherited from mfem::Mesh
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]
 

Detailed Description

Subdomain representation of a topological parent in another Mesh.

SubMesh is a subdomain representation of a Mesh defined on its parents attributes. The current implementation creates either a domain or surface subset of the parents Mesh and reuses the parallel distribution.

The attributes are taken from the parent. That means if a volume is extracted from a volume, it has the same domain attribute as the parent. Its boundary attributes are generated (there will be one boundary attribute 1 for all of the boundaries).

If a surface is extracted from a volume, the boundary attribute from the parent is assigned to be the new domain attribute. Its boundary attributes are generated (there will be one boundary attribute 1 for all of the boundaries).

For more customized boundary attributes, the resulting SubMesh has to be postprocessed.

Definition at line 42 of file submesh.hpp.

Member Enumeration Documentation

Indicator from which part of the parent Mesh the SubMesh is created.

Enumerator
Domain 
Boundary 

Definition at line 46 of file submesh.hpp.

Constructor & Destructor Documentation

mfem::SubMesh::SubMesh ( )
delete
mfem::SubMesh::~SubMesh ( )

Definition at line 104 of file submesh.cpp.

Member Function Documentation

SubMesh mfem::SubMesh::CreateFromBoundary ( const Mesh parent,
Array< int >  boundary_attributes 
)
static

Create a surface SubMesh from its parent.

The SubMesh object expects the parent Mesh object to be valid for the entire object lifetime. The boundary_attributes have to mark exactly one connected subset of the parent Mesh.

Parameters
[in]parentParent Mesh
[in]boundary_attributesBoundary attributes to extract

Definition at line 25 of file submesh.cpp.

SubMesh mfem::SubMesh::CreateFromDomain ( const Mesh parent,
Array< int >  domain_attributes 
)
static

Create a domain SubMesh from its parent.

The SubMesh object expects the parent Mesh object to be valid for the entire object lifetime. The domain_attributes have to mark exactly one connected subset of the parent Mesh.

Parameters
[in]parentParent Mesh
[in]domain_attributesDomain attributes to extract

Definition at line 19 of file submesh.cpp.

TransferMap mfem::SubMesh::CreateTransferMap ( const GridFunction src,
const GridFunction dst 
)
static

Create a Transfer Map object.

The src GridFunction can either be defined on a Mesh or a SubMesh and is transferred appropriately.

Note
Either src or dst has to be defined on a SubMesh.

Definition at line 112 of file submesh.cpp.

From mfem::SubMesh::GetFrom ( ) const
inline

Get the From indicator.

Indicates whether the SubMesh has been created from a domain or surface.

Definition at line 98 of file submesh.hpp.

const Mesh* mfem::SubMesh::GetParent ( ) const
inline

Get the parent Mesh object.

Definition at line 87 of file submesh.hpp.

const Array<int>& mfem::SubMesh::GetParentElementIDMap ( ) const
inline

Get the parent element id map.

SubMesh element id (array index) to parent Mesh element id.

Definition at line 108 of file submesh.hpp.

const Array<int>& mfem::SubMesh::GetParentFaceIDMap ( ) const
inline

Get the face id map.

SubMesh element id (array index) to parent Mesh face id.

Definition at line 118 of file submesh.hpp.

const Array<int>& mfem::SubMesh::GetParentVertexIDMap ( ) const
inline

Get the parent vertex id map.

SubMesh vertex id (array index) to parent Mesh vertex id.

Definition at line 128 of file submesh.hpp.

static bool mfem::SubMesh::IsSubMesh ( const Mesh m)
inlinestatic

Check if Mesh m is a SubMesh.

Parameters
mThe input Mesh

Definition at line 162 of file submesh.hpp.

void mfem::SubMesh::Transfer ( const GridFunction src,
GridFunction dst 
)
static

Transfer the dofs of a GridFunction.

The src GridFunction can either be defined on a Mesh or a SubMesh and is transferred appropriately.

Note
Either src or dst has to be defined on a SubMesh.
Parameters
[in]src
[out]dst

Definition at line 106 of file submesh.cpp.

Member Data Documentation

const int mfem::SubMesh::GENERATED_ATTRIBUTE = 900
static

Definition at line 52 of file submesh.hpp.


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