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

Class for parallel meshes. More...

#include <pmesh.hpp>

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

Classes

struct  Vert3
 
struct  Vert4
 

Public Member Functions

 ParMesh ()
 Default constructor. Create an empty ParMesh. More...
 
 ParMesh (MPI_Comm comm, Mesh &mesh, int *partitioning_=NULL, int part_method=1)
 Create a parallel mesh by partitioning a serial Mesh. More...
 
 ParMesh (const ParMesh &pmesh, bool copy_nodes=true)
 
 ParMesh (MPI_Comm comm, std::istream &input, bool refine=true)
 Read a parallel mesh, each MPI rank from its own file/stream. More...
 
MFEM_DEPRECATED ParMesh (ParMesh *orig_mesh, int ref_factor, int ref_type)
 Deprecated: see ParMesh::MakeRefined. More...
 
 ParMesh (ParMesh &&mesh)
 Move constructor. Used for named constructors. More...
 
ParMeshoperator= (ParMesh &&mesh)
 Move assignment operator. More...
 
ParMeshoperator= (ParMesh &mesh)=delete
 Explicitly delete the copy assignment operator. More...
 
virtual void Finalize (bool refine=false, bool fix_orientation=false)
 Finalize the construction of a general Mesh. More...
 
virtual void SetAttributes ()
 
MPI_Comm GetComm () const
 
int GetNRanks () const
 
int GetMyRank () const
 
int GetLocalElementNum (long global_element_num) const
 
long GetGlobalElementNum (int local_element_num) const
 Map a local element number to a global element number. More...
 
void GetGlobalVertexIndices (Array< HYPRE_Int > &gi) const
 AMR meshes are not supported. More...
 
void GetGlobalEdgeIndices (Array< HYPRE_Int > &gi) const
 AMR meshes are not supported. More...
 
void GetGlobalFaceIndices (Array< HYPRE_Int > &gi) const
 AMR meshes are not supported. More...
 
void GetGlobalElementIndices (Array< HYPRE_Int > &gi) const
 AMR meshes are supported. More...
 
int GetNGroups () const
 
void GenerateOffsets (int N, HYPRE_BigInt loc_sizes[], Array< HYPRE_BigInt > *offsets[]) const
 
void ExchangeFaceNbrData ()
 
void ExchangeFaceNbrNodes ()
 
virtual void SetCurvature (int order, bool discont=false, int space_dim=-1, int ordering=1)
 
int GetNFaceNeighbors () const
 
int GetNFaceNeighborElements () const
 
int GetFaceNbrGroup (int fn) const
 
int GetFaceNbrRank (int fn) const
 
TableGetFaceToAllElementTable () const
 
FaceElementTransformationsGetSharedFaceTransformations (int sf, bool fill2=true)
 
ElementTransformationGetFaceNbrElementTransformation (int i)
 
double GetFaceNbrElementSize (int i, int type=0)
 
int GetNSharedFaces () const
 Return the number of shared faces (3D), edges (2D), vertices (1D) More...
 
int GetSharedFace (int sface) const
 Return the local face index for the given shared face. More...
 
virtual void ReorientTetMesh ()
 See the remarks for the serial version in mesh.hpp. More...
 
virtual long ReduceInt (int value) const
 Utility function: sum integers from all processors (Allreduce). More...
 
void Rebalance ()
 
void Rebalance (const Array< int > &partition)
 
void ParPrint (std::ostream &out) const
 Save the mesh in a parallel mesh format. More...
 
virtual void Print (std::ostream &out=mfem::out) const
 
virtual void Save (const char *fname, int precision=16) const
 
virtual void Print (adios2stream &out) const
 
virtual void PrintXG (std::ostream &out=mfem::out) const
 
void PrintAsOne (std::ostream &out=mfem::out) const
 
void SaveAsOne (const char *fname, int precision=16) const
 
void PrintAsOneXG (std::ostream &out=mfem::out)
 Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'. More...
 
virtual void PrintVTU (std::string pathname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr=false)
 
virtual void Load (std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
 Parallel version of Mesh::Load(). More...
 
void GetBoundingBox (Vector &p_min, Vector &p_max, int ref=2)
 
void GetCharacteristics (double &h_min, double &h_max, double &kappa_min, double &kappa_max)
 
void Swap (ParMesh &other)
 
virtual void PrintInfo (std::ostream &out=mfem::out)
 Print various parallel mesh stats. 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 PrintSharedEntities (const char *fname_prefix) const
 Debugging method. More...
 
virtual ~ParMesh ()
 
These methods require group > 0
int GroupNVertices (int group)
 
int GroupNEdges (int group)
 
int GroupNTriangles (int group)
 
int GroupNQuadrilaterals (int group)
 
int GroupVertex (int group, int i)
 
void GroupEdge (int group, int i, int &edge, int &o)
 
void GroupTriangle (int group, int i, int &face, int &o)
 
void GroupQuadrilateral (int group, int i, int &face, int &o)
 
- 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= (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...
 
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...
 
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 GetNFbyType (FaceType type) const
 Returns the number of faces according to the requested type. More...
 
long GetGlobalNE () const
 Return the total (global) number of elements. More...
 
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)
 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...
 
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
 
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...
 
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...
 
FaceElementTransformationsGetFaceElementTransformations (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...
 
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
 
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
 
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...
 
void SetNodalFESpace (FiniteElementSpace *nfes)
 
void SetNodalGridFunction (GridFunction *nodes, bool make_owner=false)
 
const FiniteElementSpaceGetNodalFESpace () const
 
void EnsureNodes ()
 
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
 
void PrintVTK (std::ostream &out)
 
void PrintVTK (std::ostream &out, int ref, int field_data=0)
 
void PrintVTU (std::ostream &out, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr_elements=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 &out, 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 &out=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...
 
void MesquiteSmooth (const int mesquite_option=0)
 
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 non-conforming, 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 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)
 
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 ParMesh MakeRefined (ParMesh &orig_mesh, int ref_factor, int ref_type)
 Create a uniformly refined (by any factor) version of orig_mesh. More...
 
static ParMesh MakeSimplicial (ParMesh &orig_mesh)
 
- 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...
 

Public Attributes

GroupTopology gtopo
 
bool have_face_nbr_data
 
Array< int > face_nbr_group
 
Array< int > face_nbr_elements_offset
 
Array< int > face_nbr_vertices_offset
 
Array< Element * > face_nbr_elements
 
Array< Vertexface_nbr_vertices
 
Table send_face_nbr_elements
 
Table send_face_nbr_vertices
 
ParNCMeshpncmesh
 
- 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 non-conforming mesh extension. More...
 
Array< GeometricFactors * > geom_factors
 Optional geometric factors. More...
 
Array< FaceGeometricFactors * > face_geom_factors
 Optional face geometric factors. More...
 

Protected Member Functions

void ComputeGlobalElementOffset () const
 
 ParMesh (const ParNCMesh &pncmesh)
 Create from a nonconforming mesh. More...
 
void ReduceMeshGen ()
 
void FinalizeParTopo ()
 
virtual void MarkTetMeshForRefinement (DSTable &v_to_v)
 
int GetEdgeSplittings (Element *edge, const DSTable &v_to_v, int *middle)
 Return a number(0-1) identifying how the given edge has been split. More...
 
void GetFaceSplittings (const int *fv, const HashTable< Hashed2 > &v_to_v, Array< unsigned > &codes)
 Append codes identifying how the given face has been split to codes. More...
 
bool DecodeFaceSplittings (HashTable< Hashed2 > &v_to_v, const int *v, const Array< unsigned > &codes, int &pos)
 
void GetFaceNbrElementTransformation (int i, IsoparametricTransformation *ElTr)
 
void GetGhostFaceTransformation (FaceElementTransformations *FETr, Element::Type face_type, Geometry::Type face_geom)
 
void RefineGroups (const DSTable &v_to_v, int *middle)
 Update the groups after triangle refinement. More...
 
void RefineGroups (int old_nv, const HashTable< Hashed2 > &v_to_v)
 Update the groups after tetrahedron refinement. More...
 
void UniformRefineGroups2D (int old_nv)
 
void UniformRefineGroups3D (int old_nv, int old_nedges, const DSTable &old_v_to_v, const STable3D &old_faces, Array< int > *f2qf)
 
void ExchangeFaceNbrData (Table *gr_sface, int *s2l_face)
 
virtual void UniformRefinement2D ()
 Refine a mixed 2D mesh uniformly. More...
 
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...
 
void RebalanceImpl (const Array< int > *partition)
 
void DeleteFaceNbrData ()
 
bool WantSkipSharedMaster (const NCMesh::Master &master) const
 
int BuildLocalVertices (const Mesh &global_mesh, const int *partitioning, Array< int > &vert_global_local)
 Fills out partitioned Mesh::vertices. More...
 
int BuildLocalElements (const Mesh &global_mesh, const int *partitioning, const Array< int > &vert_global_local)
 Fills out partitioned Mesh::elements. More...
 
int BuildLocalBoundary (const Mesh &global_mesh, const int *partitioning, const Array< int > &vert_global_local, Array< bool > &activeBdrElem, Table *&edge_element)
 Fills out partitioned Mesh::boundary. More...
 
void FindSharedFaces (const Mesh &mesh, const int *partition, Array< int > &face_group, ListOfIntegerSets &groups)
 
int FindSharedEdges (const Mesh &mesh, const int *partition, Table *&edge_element, ListOfIntegerSets &groups)
 
int FindSharedVertices (const int *partition, Table *vertex_element, ListOfIntegerSets &groups)
 
void BuildFaceGroup (int ngroups, const Mesh &mesh, const Array< int > &face_group, int &nstria, int &nsquad)
 
void BuildEdgeGroup (int ngroups, const Table &edge_element)
 
void BuildVertexGroup (int ngroups, const Table &vert_element)
 
void BuildSharedFaceElems (int ntri_faces, int nquad_faces, const Mesh &mesh, int *partitioning, const STable3D *faces_tbl, const Array< int > &face_group, const Array< int > &vert_global_local)
 
void BuildSharedEdgeElems (int nedges, Mesh &mesh, const Array< int > &vert_global_local, const Table *edge_element)
 
void BuildSharedVertMapping (int nvert, const Table *vert_element, const Array< int > &vert_global_local)
 
void DistributeAttributes (Array< int > &attr)
 Ensure that bdr_attributes and attributes agree across processors. More...
 
void LoadSharedEntities (std::istream &input)
 
void EnsureParNodes ()
 If the mesh is curved, make sure 'Nodes' is ParGridFunction. More...
 
void MakeRefined_ (ParMesh &orig_mesh, int ref_factor, int ref_type)
 Internal function used in ParMesh::MakeRefined (and related constructor) More...
 
void Destroy ()
 
- 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)
 
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)
 
void UniformRefinement3D_base (Array< int > *f2qf=NULL, DSTable *v_to_v_p=NULL, bool update_nodes=true)
 
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 GetLocalQuadToHexTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void GetLocalQuadToWdgTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void ApplyLocalSlaveTransformation (FaceElementTransformations &FT, const FaceInfo &fi, bool is_ghost)
 
bool IsSlaveFace (const FaceInfo &fi) const
 
void GetVertexToVertexTable (DSTable &) 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)
 

Protected Attributes

MPI_Comm MyComm
 
int NRanks
 
int MyRank
 
Array< Element * > shared_edges
 
Array< Vert3shared_trias
 
Array< Vert4shared_quads
 
Table group_svert
 Shared objects in each group. More...
 
Table group_sedge
 
Table group_stria
 
Table group_squad
 
Array< int > svert_lvert
 Shared to local index mapping. More...
 
Array< int > sedge_ledge
 
Array< int > sface_lface
 
IsoparametricTransformation FaceNbrTransformation
 
long glob_elem_offset
 
long glob_offset_sequence
 
- 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
 

Friends

class ParNCMesh
 
class ParPumiMesh
 
class adios2stream
 

Additional Inherited Members

- Public Types inherited from mfem::Mesh
enum  Operation { NONE, REFINE, DEREFINE, REBALANCE }
 
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
 
- Static Public Attributes inherited from mfem::Mesh
static bool remove_unused_vertices = true
 
- 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)
 
- Static Protected Attributes inherited from mfem::Mesh
static const int vtk_quadratic_tet [10]
 
static const int vtk_quadratic_wedge [18]
 
static const int vtk_quadratic_hex [27]
 

Detailed Description

Class for parallel meshes.

Definition at line 32 of file pmesh.hpp.

Constructor & Destructor Documentation

mfem::ParMesh::ParMesh ( const ParNCMesh pncmesh)
protected

Create from a nonconforming mesh.

Definition at line 851 of file pmesh.cpp.

mfem::ParMesh::ParMesh ( )
inline

Default constructor. Create an empty ParMesh.

Definition at line 219 of file pmesh.hpp.

mfem::ParMesh::ParMesh ( MPI_Comm  comm,
Mesh mesh,
int *  partitioning_ = NULL,
int  part_method = 1 
)

Create a parallel mesh by partitioning a serial Mesh.

The mesh is partitioned automatically or using external partitioning data (the optional parameter 'partitioning_[i]' contains the desired MPI rank for element 'i'). Automatic partitioning uses METIS for conforming meshes and quick space-filling curve equipartitioning for nonconforming meshes (elements of nonconforming meshes should ideally be ordered as a sequence of face-neighbors).

Definition at line 106 of file pmesh.cpp.

mfem::ParMesh::ParMesh ( const ParMesh pmesh,
bool  copy_nodes = true 
)
explicit

Copy constructor. Performs a deep copy of (almost) all data, so that the source mesh can be modified (e.g. deleted, refined) without affecting the new mesh. If 'copy_nodes' is false, use a shallow (pointer) copy for the nodes, if present.

Definition at line 31 of file pmesh.cpp.

mfem::ParMesh::ParMesh ( MPI_Comm  comm,
std::istream &  input,
bool  refine = true 
)

Read a parallel mesh, each MPI rank from its own file/stream.

The refine parameter is passed to the method Mesh::Finalize().

Definition at line 920 of file pmesh.cpp.

mfem::ParMesh::ParMesh ( ParMesh orig_mesh,
int  ref_factor,
int  ref_type 
)

Deprecated: see ParMesh::MakeRefined.

Definition at line 1122 of file pmesh.cpp.

mfem::ParMesh::ParMesh ( ParMesh &&  mesh)

Move constructor. Used for named constructors.

Definition at line 95 of file pmesh.cpp.

mfem::ParMesh::~ParMesh ( )
virtual

Definition at line 5944 of file pmesh.cpp.

Member Function Documentation

void mfem::ParMesh::BuildEdgeGroup ( int  ngroups,
const Table edge_element 
)
protected

Definition at line 666 of file pmesh.cpp.

void mfem::ParMesh::BuildFaceGroup ( int  ngroups,
const Mesh mesh,
const Array< int > &  face_group,
int &  nstria,
int &  nsquad 
)
protected

Definition at line 620 of file pmesh.cpp.

int mfem::ParMesh::BuildLocalBoundary ( const Mesh global_mesh,
const int *  partitioning,
const Array< int > &  vert_global_local,
Array< bool > &  activeBdrElem,
Table *&  edge_element 
)
protected

Fills out partitioned Mesh::boundary.

Definition at line 381 of file pmesh.cpp.

int mfem::ParMesh::BuildLocalElements ( const Mesh global_mesh,
const int *  partitioning,
const Array< int > &  vert_global_local 
)
protected

Fills out partitioned Mesh::elements.

Definition at line 348 of file pmesh.cpp.

int mfem::ParMesh::BuildLocalVertices ( const Mesh global_mesh,
const int *  partitioning,
Array< int > &  vert_global_local 
)
protected

Fills out partitioned Mesh::vertices.

Definition at line 300 of file pmesh.cpp.

void mfem::ParMesh::BuildSharedEdgeElems ( int  nedges,
Mesh mesh,
const Array< int > &  vert_global_local,
const Table edge_element 
)
protected

Definition at line 795 of file pmesh.cpp.

void mfem::ParMesh::BuildSharedFaceElems ( int  ntri_faces,
int  nquad_faces,
const Mesh mesh,
int *  partitioning,
const STable3D faces_tbl,
const Array< int > &  face_group,
const Array< int > &  vert_global_local 
)
protected

Definition at line 718 of file pmesh.cpp.

void mfem::ParMesh::BuildSharedVertMapping ( int  nvert,
const Table vert_element,
const Array< int > &  vert_global_local 
)
protected

Definition at line 832 of file pmesh.cpp.

void mfem::ParMesh::BuildVertexGroup ( int  ngroups,
const Table vert_element 
)
protected

Definition at line 692 of file pmesh.cpp.

void mfem::ParMesh::ComputeGlobalElementOffset ( ) const
protected

Definition at line 865 of file pmesh.cpp.

bool mfem::ParMesh::DecodeFaceSplittings ( HashTable< Hashed2 > &  v_to_v,
const int *  v,
const Array< unsigned > &  codes,
int &  pos 
)
protected

Definition at line 1800 of file pmesh.cpp.

void mfem::ParMesh::DeleteFaceNbrData ( )
protected

Definition at line 1937 of file pmesh.cpp.

void mfem::ParMesh::Destroy ( )
protected

Definition at line 5930 of file pmesh.cpp.

void mfem::ParMesh::DistributeAttributes ( Array< int > &  attr)
protected

Ensure that bdr_attributes and attributes agree across processors.

Definition at line 1538 of file pmesh.cpp.

void mfem::ParMesh::EnsureParNodes ( )
protected

If the mesh is curved, make sure 'Nodes' is ParGridFunction.

Note that this method is not related to the public 'Mesh::EnsureNodes`.

Definition at line 1978 of file pmesh.cpp.

void mfem::ParMesh::ExchangeFaceNbrData ( Table gr_sface,
int *  s2l_face 
)
protected

Definition at line 2085 of file pmesh.cpp.

void mfem::ParMesh::ExchangeFaceNbrData ( )

Definition at line 2000 of file pmesh.cpp.

void mfem::ParMesh::ExchangeFaceNbrNodes ( )

Definition at line 2454 of file pmesh.cpp.

void mfem::ParMesh::Finalize ( bool  refine = false,
bool  fix_orientation = false 
)
virtual

Finalize the construction of a general Mesh.

This method will:

  • check and optionally fix the orientation of regular elements
  • check and fix the orientation of boundary elements
  • assume that vertices are defined, if Nodes == NULL
  • assume that Nodes are defined, if Nodes != NULL.
    Parameters
    [in]refineIf true, prepare the Mesh for conforming refinement of triangular or tetrahedral meshes.
    [in]fix_orientationIf true, fix the orientation of inverted mesh elements by permuting their vertices.
    Before calling this method, call FinalizeTopology() and ensure that the Mesh vertices or nodes are set.

Reimplemented from mfem::Mesh.

Definition at line 1510 of file pmesh.cpp.

void mfem::ParMesh::FinalizeParTopo ( )
protected

Definition at line 883 of file pmesh.cpp.

int mfem::ParMesh::FindPoints ( DenseMatrix point_mat,
Array< int > &  elem_ids,
Array< IntegrationPoint > &  ips,
bool  warn = true,
InverseElementTransformation inv_trans = NULL 
)
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.

Returns
The total number of points that were found.
Note
This method is not 100 percent reliable, i.e. it is not guaranteed to find a point, even if it lies inside a mesh element.

Reimplemented from mfem::Mesh.

Definition at line 5675 of file pmesh.cpp.

int mfem::ParMesh::FindSharedEdges ( const Mesh mesh,
const int *  partition,
Table *&  edge_element,
ListOfIntegerSets groups 
)
protected

Definition at line 529 of file pmesh.cpp.

void mfem::ParMesh::FindSharedFaces ( const Mesh mesh,
const int *  partition,
Array< int > &  face_group,
ListOfIntegerSets groups 
)
protected

Definition at line 502 of file pmesh.cpp.

int mfem::ParMesh::FindSharedVertices ( const int *  partition,
Table vertex_element,
ListOfIntegerSets groups 
)
protected

Definition at line 583 of file pmesh.cpp.

void mfem::ParMesh::GenerateOffsets ( int  N,
HYPRE_BigInt  loc_sizes[],
Array< HYPRE_BigInt > *  offsets[] 
) const

Definition at line 1836 of file pmesh.cpp.

void mfem::ParMesh::GetBoundingBox ( Vector p_min,
Vector p_max,
int  ref = 2 
)

Returns the minimum and maximum corners of the mesh bounding box. For high-order meshes, the geometry is refined first "ref" times.

Definition at line 5377 of file pmesh.cpp.

void mfem::ParMesh::GetCharacteristics ( double &  h_min,
double &  h_max,
double &  kappa_min,
double &  kappa_max 
)

Definition at line 5393 of file pmesh.cpp.

MPI_Comm mfem::ParMesh::GetComm ( ) const
inline

Definition at line 276 of file pmesh.hpp.

int mfem::ParMesh::GetEdgeSplittings ( Element edge,
const DSTable v_to_v,
int *  middle 
)
protected

Return a number(0-1) identifying how the given edge has been split.

Definition at line 1750 of file pmesh.cpp.

double mfem::ParMesh::GetFaceNbrElementSize ( int  i,
int  type = 0 
)

Get the size of the i-th face neighbor element relative to the reference element.

Definition at line 1932 of file pmesh.cpp.

void mfem::ParMesh::GetFaceNbrElementTransformation ( int  i,
IsoparametricTransformation ElTr 
)
protected

Definition at line 1879 of file pmesh.cpp.

ElementTransformation* mfem::ParMesh::GetFaceNbrElementTransformation ( int  i)
inline

Definition at line 355 of file pmesh.hpp.

int mfem::ParMesh::GetFaceNbrGroup ( int  fn) const
inline

Definition at line 341 of file pmesh.hpp.

int mfem::ParMesh::GetFaceNbrRank ( int  fn) const

Definition at line 2514 of file pmesh.cpp.

void mfem::ParMesh::GetFaceSplittings ( const int *  fv,
const HashTable< Hashed2 > &  v_to_v,
Array< unsigned > &  codes 
)
protected

Append codes identifying how the given face has been split to codes.

Definition at line 1765 of file pmesh.cpp.

Table * mfem::ParMesh::GetFaceToAllElementTable ( ) const

Similar to Mesh::GetFaceToElementTable with added face-neighbor elements with indices offset by the local number of elements.

Definition at line 2531 of file pmesh.cpp.

void mfem::ParMesh::GetGhostFaceTransformation ( FaceElementTransformations FETr,
Element::Type  face_type,
Geometry::Type  face_geom 
)
protected

Definition at line 2589 of file pmesh.cpp.

void mfem::ParMesh::GetGlobalEdgeIndices ( Array< HYPRE_Int > &  gi) const

AMR meshes are not supported.

Definition at line 5827 of file pmesh.cpp.

void mfem::ParMesh::GetGlobalElementIndices ( Array< HYPRE_Int > &  gi) const

AMR meshes are supported.

Definition at line 5878 of file pmesh.cpp.

long mfem::ParMesh::GetGlobalElementNum ( int  local_element_num) const

Map a local element number to a global element number.

Definition at line 1532 of file pmesh.cpp.

void mfem::ParMesh::GetGlobalFaceIndices ( Array< HYPRE_Int > &  gi) const

AMR meshes are not supported.

Definition at line 5850 of file pmesh.cpp.

void mfem::ParMesh::GetGlobalVertexIndices ( Array< HYPRE_Int > &  gi) const

AMR meshes are not supported.

The following functions define global indices for all local vertices, edges, faces, or elements. The global indices have no meaning or significance for ParMesh, but can be used for purposes beyond this class.

Definition at line 5811 of file pmesh.cpp.

int mfem::ParMesh::GetLocalElementNum ( long  global_element_num) const

Map a global element number to a local element number. If the global element is not on this processor, return -1.

Definition at line 1524 of file pmesh.cpp.

int mfem::ParMesh::GetMyRank ( ) const
inline

Definition at line 278 of file pmesh.hpp.

int mfem::ParMesh::GetNFaceNeighborElements ( ) const
inline

Definition at line 340 of file pmesh.hpp.

int mfem::ParMesh::GetNFaceNeighbors ( ) const
inline

Definition at line 339 of file pmesh.hpp.

int mfem::ParMesh::GetNGroups ( ) const
inline

Definition at line 315 of file pmesh.hpp.

int mfem::ParMesh::GetNRanks ( ) const
inline

Definition at line 277 of file pmesh.hpp.

int mfem::ParMesh::GetNSharedFaces ( ) const

Return the number of shared faces (3D), edges (2D), vertices (1D)

Definition at line 2731 of file pmesh.cpp.

int mfem::ParMesh::GetSharedFace ( int  sface) const

Return the local face index for the given shared face.

Definition at line 2750 of file pmesh.cpp.

FaceElementTransformations * mfem::ParMesh::GetSharedFaceTransformations ( int  sf,
bool  fill2 = true 
)

Get the FaceElementTransformations for the given shared face (edge 2D). In the returned object, 1 and 2 refer to the local and the neighbor elements, respectively.

Definition at line 2622 of file pmesh.cpp.

void mfem::ParMesh::GroupEdge ( int  group,
int  i,
int &  edge,
int &  o 
)

Definition at line 1592 of file pmesh.cpp.

int mfem::ParMesh::GroupNEdges ( int  group)
inline

Definition at line 319 of file pmesh.hpp.

int mfem::ParMesh::GroupNQuadrilaterals ( int  group)
inline

Definition at line 321 of file pmesh.hpp.

int mfem::ParMesh::GroupNTriangles ( int  group)
inline

Definition at line 320 of file pmesh.hpp.

int mfem::ParMesh::GroupNVertices ( int  group)
inline

Definition at line 318 of file pmesh.hpp.

void mfem::ParMesh::GroupQuadrilateral ( int  group,
int  i,
int &  face,
int &  o 
)

Definition at line 1611 of file pmesh.cpp.

void mfem::ParMesh::GroupTriangle ( int  group,
int  i,
int &  face,
int &  o 
)

Definition at line 1600 of file pmesh.cpp.

int mfem::ParMesh::GroupVertex ( int  group,
int  i 
)
inline

Definition at line 323 of file pmesh.hpp.

void mfem::ParMesh::Load ( std::istream &  input,
int  generate_edges = 0,
int  refine = 1,
bool  fix_orientation = true 
)
virtual

Parallel version of Mesh::Load().

Reimplemented from mfem::Mesh.

Definition at line 937 of file pmesh.cpp.

void mfem::ParMesh::LoadSharedEntities ( std::istream &  input)
protected

Definition at line 973 of file pmesh.cpp.

void mfem::ParMesh::LocalRefinement ( const Array< int > &  marked_el,
int  type = 3 
)
protectedvirtual

This function is not public anymore. Use GeneralRefinement instead.

Reimplemented from mfem::Mesh.

Definition at line 2988 of file pmesh.cpp.

ParMesh mfem::ParMesh::MakeRefined ( ParMesh orig_mesh,
int  ref_factor,
int  ref_type 
)
static

Create a uniformly refined (by any factor) version of orig_mesh.

Parameters
[in]orig_meshThe starting coarse mesh.
[in]ref_factorThe refinement factor, an integer > 1.
[in]ref_typeSpecify 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.

Note
The constructed ParMesh is linear, i.e. it does not have nodes.

Definition at line 1361 of file pmesh.cpp.

void mfem::ParMesh::MakeRefined_ ( ParMesh orig_mesh,
int  ref_factor,
int  ref_type 
)
protected

Internal function used in ParMesh::MakeRefined (and related constructor)

Definition at line 1127 of file pmesh.cpp.

ParMesh mfem::ParMesh::MakeSimplicial ( ParMesh orig_mesh)
static

Create a mesh by splitting each element of orig_mesh into simplices. See Mesh::MakeSimplicial for more details.

Definition at line 1368 of file pmesh.cpp.

void mfem::ParMesh::MarkTetMeshForRefinement ( DSTable v_to_v)
protectedvirtual

Reimplemented from mfem::Mesh.

Definition at line 1622 of file pmesh.cpp.

bool mfem::ParMesh::NonconformingDerefinement ( Array< double > &  elem_error,
double  threshold,
int  nc_limit = 0,
int  op = 1 
)
protectedvirtual

NC version of GeneralDerefinement.

Reimplemented from mfem::Mesh.

Definition at line 3544 of file pmesh.cpp.

void mfem::ParMesh::NonconformingRefinement ( const Array< Refinement > &  refinements,
int  nc_limit = 0 
)
protectedvirtual

This function is not public anymore. Use GeneralRefinement instead.

Reimplemented from mfem::Mesh.

Definition at line 3492 of file pmesh.cpp.

void mfem::ParMesh::NURBSUniformRefinement ( )
protectedvirtual

Refine NURBS mesh.

Reimplemented from mfem::Mesh.

Definition at line 4160 of file pmesh.cpp.

ParMesh & mfem::ParMesh::operator= ( ParMesh &&  mesh)

Move assignment operator.

Definition at line 100 of file pmesh.cpp.

ParMesh& mfem::ParMesh::operator= ( ParMesh mesh)
delete

Explicitly delete the copy assignment operator.

void mfem::ParMesh::ParPrint ( std::ostream &  out) const

Save the mesh in a parallel mesh format.

Definition at line 5525 of file pmesh.cpp.

void mfem::ParMesh::Print ( std::ostream &  out = mfem::out) const
virtual

Print the part of the mesh in the calling processor adding the interface as boundary (for visualization purposes) using the mfem v1.0 format.

Reimplemented from mfem::Mesh.

Definition at line 4382 of file pmesh.cpp.

void mfem::ParMesh::Print ( adios2stream out) const
virtual

Print the part of the mesh in the calling processor using adios2 bp format.

Reimplemented from mfem::Mesh.

Definition at line 4509 of file pmesh.cpp.

void mfem::ParMesh::PrintAsOne ( std::ostream &  out = mfem::out) const

Write the mesh to the stream 'out' on Process 0 in a form suitable for visualization: the mesh is written as a disjoint mesh and the shared boundary is added to the actual boundary; both the element and boundary attributes are set to the processor number.

Definition at line 4527 of file pmesh.cpp.

void mfem::ParMesh::PrintAsOneXG ( std::ostream &  out = mfem::out)

Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'.

Definition at line 4847 of file pmesh.cpp.

void mfem::ParMesh::PrintInfo ( std::ostream &  out = mfem::out)
virtual

Print various parallel mesh stats.

Reimplemented from mfem::Mesh.

Definition at line 5406 of file pmesh.cpp.

void mfem::ParMesh::PrintSharedEntities ( const char *  fname_prefix) const

Debugging method.

Definition at line 5724 of file pmesh.cpp.

void mfem::ParMesh::PrintVTU ( std::string  pathname,
VTKFormat  format = VTKFormat::ASCII,
bool  high_order_output = false,
int  compression_level = 0,
bool  bdr = false 
)
virtual

Print the mesh in parallel PVTU format. The PVTU and VTU files will be stored in the directory specified by pathname. If the directory does not exist, it will be created.

Reimplemented from mfem::Mesh.

Definition at line 5608 of file pmesh.cpp.

void mfem::ParMesh::PrintXG ( std::ostream &  out = mfem::out) const
virtual

Print the part of the mesh in the calling processor adding the interface as boundary (for visualization purposes) using Netgen/Truegrid format .

Reimplemented from mfem::Mesh.

Definition at line 4168 of file pmesh.cpp.

void mfem::ParMesh::Rebalance ( )

Load balance the mesh by equipartitioning the global space-filling sequence of elements. Works for nonconforming meshes only.

Definition at line 3602 of file pmesh.cpp.

void mfem::ParMesh::Rebalance ( const Array< int > &  partition)

Load balance a nonconforming mesh using a user-defined partition. Each local element 'i' is migrated to processor rank 'partition[i]', for 0 <= i < GetNE().

Definition at line 3607 of file pmesh.cpp.

void mfem::ParMesh::RebalanceImpl ( const Array< int > *  partition)
protected

Definition at line 3612 of file pmesh.cpp.

long mfem::ParMesh::ReduceInt ( int  value) const
virtual

Utility function: sum integers from all processors (Allreduce).

Reimplemented from mfem::Mesh.

Definition at line 5518 of file pmesh.cpp.

void mfem::ParMesh::ReduceMeshGen ( )
protected

Definition at line 877 of file pmesh.cpp.

void mfem::ParMesh::RefineGroups ( const DSTable v_to_v,
int *  middle 
)
protected

Update the groups after triangle refinement.

Definition at line 3650 of file pmesh.cpp.

void mfem::ParMesh::RefineGroups ( int  old_nv,
const HashTable< Hashed2 > &  v_to_v 
)
protected

Update the groups after tetrahedron refinement.

Definition at line 3723 of file pmesh.cpp.

void mfem::ParMesh::ReorientTetMesh ( )
virtual

See the remarks for the serial version in mesh.hpp.

Reimplemented from mfem::Mesh.

Definition at line 2798 of file pmesh.cpp.

void mfem::ParMesh::Save ( const char *  fname,
int  precision = 16 
) const
virtual

Save the ParMesh to files (one for each MPI rank). The files will be given suffixes according to the MPI rank. The mesh will be written to the files using ParMesh::Print. The given precision will be used for ASCII output.

Reimplemented from mfem::Mesh.

Definition at line 4499 of file pmesh.cpp.

void mfem::ParMesh::SaveAsOne ( const char *  fname,
int  precision = 16 
) const

Save the mesh as a single file (using ParMesh::PrintAsOne). The given precision is used for ASCII output.

Definition at line 4836 of file pmesh.cpp.

void mfem::ParMesh::SetAttributes ( )
virtual

Reimplemented from mfem::Mesh.

Definition at line 1574 of file pmesh.cpp.

void mfem::ParMesh::SetCurvature ( int  order,
bool  discont = false,
int  space_dim = -1,
int  ordering = 1 
)
virtual

Set the curvature of the mesh nodes using the given polynomial degree, 'order', and optionally: discontinuous or continuous FE space, 'discont', new space dimension, 'space_dim' (if != -1), and 'ordering'.

Reimplemented from mfem::Mesh.

Definition at line 1958 of file pmesh.cpp.

void mfem::ParMesh::Swap ( ParMesh other)

Swaps internal data with another ParMesh, including non-geometry members. See Mesh::Swap

Definition at line 5891 of file pmesh.cpp.

void mfem::ParMesh::UniformRefineGroups2D ( int  old_nv)
protected

Definition at line 3907 of file pmesh.cpp.

void mfem::ParMesh::UniformRefineGroups3D ( int  old_nv,
int  old_nedges,
const DSTable old_v_to_v,
const STable3D old_faces,
Array< int > *  f2qf 
)
protected

Definition at line 3958 of file pmesh.cpp.

void mfem::ParMesh::UniformRefinement2D ( )
protectedvirtual

Refine a mixed 2D mesh uniformly.

Reimplemented from mfem::Mesh.

Definition at line 4107 of file pmesh.cpp.

void mfem::ParMesh::UniformRefinement3D ( )
protectedvirtual

Refine a mixed 3D mesh uniformly.

Reimplemented from mfem::Mesh.

Definition at line 4131 of file pmesh.cpp.

bool mfem::ParMesh::WantSkipSharedMaster ( const NCMesh::Master master) const
protected

Definition at line 4369 of file pmesh.cpp.

Friends And Related Function Documentation

friend class adios2stream
friend

Definition at line 464 of file pmesh.hpp.

friend class ParNCMesh
friend

Definition at line 459 of file pmesh.hpp.

friend class ParPumiMesh
friend

Definition at line 461 of file pmesh.hpp.

Member Data Documentation

Array<Element *> mfem::ParMesh::face_nbr_elements

Definition at line 307 of file pmesh.hpp.

Array<int> mfem::ParMesh::face_nbr_elements_offset

Definition at line 305 of file pmesh.hpp.

Array<int> mfem::ParMesh::face_nbr_group

Definition at line 304 of file pmesh.hpp.

Array<Vertex> mfem::ParMesh::face_nbr_vertices

Definition at line 308 of file pmesh.hpp.

Array<int> mfem::ParMesh::face_nbr_vertices_offset

Definition at line 306 of file pmesh.hpp.

IsoparametricTransformation mfem::ParMesh::FaceNbrTransformation
protected

Definition at line 78 of file pmesh.hpp.

long mfem::ParMesh::glob_elem_offset
mutableprotected

Definition at line 81 of file pmesh.hpp.

long mfem::ParMesh::glob_offset_sequence
mutableprotected

Definition at line 81 of file pmesh.hpp.

Table mfem::ParMesh::group_sedge
protected

Definition at line 68 of file pmesh.hpp.

Table mfem::ParMesh::group_squad
protected

Definition at line 70 of file pmesh.hpp.

Table mfem::ParMesh::group_stria
protected

Definition at line 69 of file pmesh.hpp.

Table mfem::ParMesh::group_svert
protected

Shared objects in each group.

Definition at line 67 of file pmesh.hpp.

GroupTopology mfem::ParMesh::gtopo

Definition at line 300 of file pmesh.hpp.

bool mfem::ParMesh::have_face_nbr_data

Definition at line 303 of file pmesh.hpp.

MPI_Comm mfem::ParMesh::MyComm
protected

Definition at line 35 of file pmesh.hpp.

int mfem::ParMesh::MyRank
protected

Definition at line 36 of file pmesh.hpp.

int mfem::ParMesh::NRanks
protected

Definition at line 36 of file pmesh.hpp.

ParNCMesh* mfem::ParMesh::pncmesh

Definition at line 313 of file pmesh.hpp.

Array<int> mfem::ParMesh::sedge_ledge
protected

Definition at line 74 of file pmesh.hpp.

Table mfem::ParMesh::send_face_nbr_elements

Definition at line 310 of file pmesh.hpp.

Table mfem::ParMesh::send_face_nbr_vertices

Definition at line 311 of file pmesh.hpp.

Array<int> mfem::ParMesh::sface_lface
protected

Definition at line 76 of file pmesh.hpp.

Array<Element *> mfem::ParMesh::shared_edges
protected

Definition at line 59 of file pmesh.hpp.

Array<Vert4> mfem::ParMesh::shared_quads
protected

Definition at line 64 of file pmesh.hpp.

Array<Vert3> mfem::ParMesh::shared_trias
protected

Definition at line 63 of file pmesh.hpp.

Array<int> mfem::ParMesh::svert_lvert
protected

Shared to local index mapping.

Definition at line 73 of file pmesh.hpp.


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