MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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]

Public Member Functions

 ParMesh (const ParMesh &pmesh, bool copy_nodes=true)
 
 ParMesh (MPI_Comm comm, Mesh &mesh, int *partitioning_=NULL, int part_method=1)
 
 ParMesh (MPI_Comm comm, std::istream &input)
 Read a parallel mesh, each MPI rank from its own file/stream. More...
 
 ParMesh (ParMesh *orig_mesh, int ref_factor, int ref_type)
 Create a uniformly refined (by any factor) version of orig_mesh. More...
 
MPI_Comm GetComm () const
 
int GetNRanks () const
 
int GetMyRank () const
 
int GetNGroups () const
 
void GenerateOffsets (int N, HYPRE_Int loc_sizes[], Array< HYPRE_Int > *offsets[]) const
 
void ExchangeFaceNbrData ()
 
void ExchangeFaceNbrNodes ()
 
int GetNFaceNeighbors () const
 
int GetFaceNbrGroup (int fn) const
 
int GetFaceNbrRank (int fn) const
 
TableGetFaceToAllElementTable () const
 
FaceElementTransformationsGetSharedFaceTransformations (int sf, bool fill2=true)
 
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 RefineGroups (const DSTable &v_to_v, int *middle)
 Update the groups after tet refinement. More...
 
void Rebalance ()
 Load balance the mesh. NC meshes only. More...
 
virtual void Print (std::ostream &out=mfem::out) const
 
virtual void PrintXG (std::ostream &out=mfem::out) const
 
void PrintAsOne (std::ostream &out=mfem::out)
 
void PrintAsOneXG (std::ostream &out=mfem::out)
 Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'. 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)
 
virtual void PrintInfo (std::ostream &out=mfem::out)
 Print various parallel mesh stats. More...
 
void ParPrint (std::ostream &out) const
 Save the mesh in a parallel mesh format. 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...
 
virtual ~ParMesh ()
 
These methods require group > 0
int GroupNVertices (int group)
 
int GroupNEdges (int group)
 
int GroupNFaces (int group)
 
int GroupVertex (int group, int i)
 
void GroupEdge (int group, int i, int &edge, int &o)
 
void GroupFace (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 (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 ()
 Finalize the construction of the secondary topology (connectivity) data of a Mesh. More...
 
void Finalize (bool refine=false, bool fix_orientation=false)
 Finalize the construction of a general Mesh. More...
 
void SetAttributes ()
 
void GetGeckoElementReordering (Array< int > &ordering)
 
void ReorderElements (const Array< int > &ordering, bool reorder_vertices=true)
 
 Mesh (int nx, int ny, int nz, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0, double sz=1.0)
 
 Mesh (int nx, int ny, Element::Type type, int generate_edges=0, double sx=1.0, double sy=1.0)
 
 Mesh (int n, double sx=1.0)
 
 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...
 
 Mesh (Mesh *orig_mesh, int ref_factor, int ref_type)
 Create a uniformly refined (by any factor) version of orig_mesh. 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...
 
long GetGlobalNE () const
 Return the total (global) number of elements. 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
 
int GetFaceBaseGeometry (int i) const
 
int GetElementBaseGeometry (int i=0) const
 
int GetBdrElementBaseGeometry (int i=0) const
 
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 > &, Array< int > &) 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 > &, Array< int > &) const
 Return the indices and the orientations of all faces of element i. More...
 
void GetBdrElementFace (int i, int *, int *) 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...
 
int GetElementType (int i) const
 Returns the type of element i. More...
 
int 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)
 
bool FaceIsInterior (int FaceNo) const
 Return true if the given face is interior. More...
 
void GetFaceElements (int Face, int *Elem1, int *Elem2)
 
void GetFaceInfos (int Face, int *Inf1, int *Inf2)
 
int GetFaceGeometryType (int Face) const
 
int GetFaceElementType (int Face) const
 
int CheckElementOrientation (bool fix_it=true)
 Check 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...
 
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)
 
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 SetCurvature (int order, bool discont=false, int space_dim=-1, int ordering=1)
 
void UniformRefinement ()
 
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 triangles_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 GetElementColoring (Array< int > &colors, int el0=0)
 
void PrintWithPartitioning (int *partitioning, std::ostream &out, int elem_attr=0) const
 
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)
 
double GetElementSize (int i, const Vector &dir)
 
double GetElementVolume (int i)
 
void GetBoundingBox (Vector &min, Vector &max, int ref=2)
 
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)
 
void MesquiteSmooth (const int mesquite_option=0)
 
virtual ~Mesh ()
 Destroys Mesh. More...
 
ElementNewElement (int geom)
 
void AddVertex (const double *)
 
void AddTri (const int *vi, int attr=1)
 
void AddTriangle (const int *vi, int attr=1)
 
void AddQuad (const int *vi, int attr=1)
 
void AddTet (const int *vi, int attr=1)
 
void AddHex (const int *vi, int attr=1)
 
void AddHexAsTets (const int *vi, int attr=1)
 
void AddElement (Element *elem)
 
void AddBdrElement (Element *elem)
 
void AddBdrSegment (const int *vi, int attr=1)
 
void AddBdrTriangle (const int *vi, int attr=1)
 
void AddBdrQuad (const int *vi, int attr=1)
 
void AddBdrQuadAsTriangles (const int *vi, 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 FinalizeHexMesh (int generate_edges=0, int refine=0, bool fix_orientation=true)
 Finalize the construction of a hexahedral Mesh. More...
 
void KnotInsert (Array< KnotVector * > &kv)
 
void DegreeElevate (int rel_degree, int degree=16)
 

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...
 

Protected Member Functions

 ParMesh ()
 
 ParMesh (const ParNCMesh &pncmesh)
 Create from a nonconforming mesh. More...
 
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...
 
int GetFaceSplittings (Element *face, const DSTable &v_to_v, int *middle)
 Return a number(0-4) identifying how the given face has been split. More...
 
void GetFaceNbrElementTransformation (int i, IsoparametricTransformation *ElTr)
 
ElementTransformationGetGhostFaceTransformation (FaceElementTransformations *FETr, int face_type, int face_geom)
 
virtual void QuadUniformRefinement ()
 Refine quadrilateral mesh. More...
 
virtual void HexUniformRefinement ()
 Refine a hexahedral mesh. 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 DeleteFaceNbrData ()
 
bool WantSkipSharedMaster (const NCMesh::Master &master) const
 
- Protected Member Functions inherited from mfem::Mesh
void Init ()
 
void InitTables ()
 
void SetEmpty ()
 
void DestroyTables ()
 
void DeleteTables ()
 
void DestroyPointers ()
 
void Destroy ()
 
void DeleteLazyTables ()
 
ElementReadElementWithoutAttr (std::istream &)
 
ElementReadElement (std::istream &)
 
void ReadMFEMMesh (std::istream &input, bool mfem_v11, 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 ReadVTKMesh (std::istream &input, int &curved, int &read_gf)
 
void ReadNURBSMesh (std::istream &input, int &curved, int &read_gf)
 
void ReadInlineMesh (std::istream &input, int generate_edges=0)
 
void ReadGmshMesh (std::istream &input)
 
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 *)
 
void Bisection (int i, const DSTable &, int *)
 
void UniformRefinement (int i, const DSTable &, int *, int *, int *)
 
void AverageVertices (int *indexes, int n, int result)
 
void InitRefinementTransforms ()
 
int FindCoarseElement (int i)
 
void UpdateNodes ()
 Update the nodes of a curved mesh after refinement. More...
 
void DerefineMesh (const Array< int > &derefinements)
 Derefine elements once a list of derefinements is known. 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 GetLocalQuadToHexTransformation (IsoparametricTransformation &loc, int i)
 Used in GetFaceElementTransformations (...) More...
 
void ApplyLocalSlaveTransformation (IsoparametricTransformation &transf, const FaceInfo &fi)
 
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 InitBaseGeom ()
 
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, int generate_edges, double sx, double sy, double sz)
 
void Make2D (int nx, int ny, Element::Type type, int generate_edges, double sx, double sy)
 
void Make1D (int n, double sx=1.0)
 Creates a 1D mesh for the interval [0,sx] divided into n equal intervals. 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 Swap (Mesh &other, bool non_geometry=false)
 
void GetElementData (const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
 

Protected Attributes

MPI_Comm MyComm
 
int NRanks
 
int MyRank
 
Array< Element * > shared_edges
 
Array< Element * > shared_faces
 
Table group_svert
 Shared objects in each group. More...
 
Table group_sedge
 
Table group_sface
 
Array< int > svert_lvert
 Shared to local index mapping. More...
 
Array< int > sedge_ledge
 
Array< int > sface_lface
 
- Protected Attributes inherited from mfem::Mesh
int Dim
 
int spaceDim
 
int NumOfVertices
 
int NumOfElements
 
int NumOfBdrElements
 
int NumOfEdges
 
int NumOfFaces
 
int BaseGeom
 
int BaseBdrGeom
 
int meshgen
 
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 FaceTransformation
 
IsoparametricTransformation EdgeTransformation
 
FaceElementTransformations FaceElemTr
 
CoarseFineTransformations CoarseFineTr
 
GridFunctionNodes
 
int own_nodes
 
MemAlloc< Tetrahedron, 1024 > TetMemory
 
Operation last_operation
 

Friends

class ParPumiMesh
 

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
 
- Static Public Member Functions inherited from mfem::Mesh
static FiniteElementGetTransformationFEforElementType (int)
 
- 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 void GetElementArrayEdgeTable (const Array< Element * > &elem_array, const DSTable &v_to_v, Table &el_to_edge)
 
static void ShiftL2R (int &, int &, int &)
 
static void Rotate3 (int &, int &, int &)
 
- Static Protected Attributes inherited from mfem::Mesh
static const int vtk_quadratic_tet [10]
 
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 ( )
inlineprotected

Definition at line 38 of file pmesh.hpp.

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

Create from a nonconforming mesh.

Definition at line 695 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 29 of file pmesh.cpp.

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

Definition at line 86 of file pmesh.cpp.

mfem::ParMesh::ParMesh ( MPI_Comm  comm,
std::istream &  input 
)

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

Definition at line 706 of file pmesh.cpp.

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

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 863 of file pmesh.cpp.

mfem::ParMesh::~ParMesh ( )
virtual

Definition at line 4598 of file pmesh.cpp.

Member Function Documentation

void mfem::ParMesh::DeleteFaceNbrData ( )
protected

Definition at line 1357 of file pmesh.cpp.

void mfem::ParMesh::ExchangeFaceNbrData ( )

Definition at line 1378 of file pmesh.cpp.

void mfem::ParMesh::ExchangeFaceNbrNodes ( )

Definition at line 1781 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 4558 of file pmesh.cpp.

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

Definition at line 1262 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 4354 of file pmesh.cpp.

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

Definition at line 4370 of file pmesh.cpp.

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

Definition at line 123 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 1203 of file pmesh.cpp.

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

Definition at line 1305 of file pmesh.cpp.

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

Definition at line 162 of file pmesh.hpp.

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

Definition at line 1839 of file pmesh.cpp.

int mfem::ParMesh::GetFaceSplittings ( Element face,
const DSTable v_to_v,
int *  middle 
)
protected

Return a number(0-4) identifying how the given face has been split.

Definition at line 1225 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 1856 of file pmesh.cpp.

ElementTransformation * mfem::ParMesh::GetGhostFaceTransformation ( FaceElementTransformations FETr,
int  face_type,
int  face_geom 
)
protected

Definition at line 1914 of file pmesh.cpp.

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

Definition at line 125 of file pmesh.hpp.

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

Definition at line 161 of file pmesh.hpp.

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

Definition at line 142 of file pmesh.hpp.

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

Definition at line 124 of file pmesh.hpp.

int mfem::ParMesh::GetNSharedFaces ( ) const

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

Definition at line 2030 of file pmesh.cpp.

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

Return the local face index for the given shared face.

Definition at line 2049 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 1945 of file pmesh.cpp.

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

Definition at line 1052 of file pmesh.cpp.

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

Definition at line 1060 of file pmesh.cpp.

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

Definition at line 146 of file pmesh.hpp.

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

Definition at line 147 of file pmesh.hpp.

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

Definition at line 145 of file pmesh.hpp.

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

Definition at line 149 of file pmesh.hpp.

void mfem::ParMesh::HexUniformRefinement ( )
protectedvirtual

Refine a hexahedral mesh.

Reimplemented from mfem::Mesh.

Definition at line 3053 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 2115 of file pmesh.cpp.

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

Reimplemented from mfem::Mesh.

Definition at line 1077 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 2698 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 2646 of file pmesh.cpp.

void mfem::ParMesh::NURBSUniformRefinement ( )
protectedvirtual

Refine NURBS mesh.

Reimplemented from mfem::Mesh.

Definition at line 3212 of file pmesh.cpp.

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

Save the mesh in a parallel mesh format.

Definition at line 4494 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 3428 of file pmesh.cpp.

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

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 3556 of file pmesh.cpp.

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

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

Definition at line 3831 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 4383 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 3220 of file pmesh.cpp.

void mfem::ParMesh::QuadUniformRefinement ( )
protectedvirtual

Refine quadrilateral mesh.

Reimplemented from mfem::Mesh.

Definition at line 2966 of file pmesh.cpp.

void mfem::ParMesh::Rebalance ( )

Load balance the mesh. NC meshes only.

Definition at line 2746 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 4487 of file pmesh.cpp.

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

Update the groups after tet refinement.

Definition at line 2779 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 2071 of file pmesh.cpp.

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

Definition at line 3415 of file pmesh.cpp.

Friends And Related Function Documentation

friend class ParPumiMesh
friend

Definition at line 35 of file pmesh.hpp.

Member Data Documentation

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

Definition at line 134 of file pmesh.hpp.

Array<int> mfem::ParMesh::face_nbr_elements_offset

Definition at line 132 of file pmesh.hpp.

Array<int> mfem::ParMesh::face_nbr_group

Definition at line 131 of file pmesh.hpp.

Array<Vertex> mfem::ParMesh::face_nbr_vertices

Definition at line 135 of file pmesh.hpp.

Array<int> mfem::ParMesh::face_nbr_vertices_offset

Definition at line 133 of file pmesh.hpp.

Table mfem::ParMesh::group_sedge
protected

Definition at line 49 of file pmesh.hpp.

Table mfem::ParMesh::group_sface
protected

Definition at line 50 of file pmesh.hpp.

Table mfem::ParMesh::group_svert
protected

Shared objects in each group.

Definition at line 48 of file pmesh.hpp.

GroupTopology mfem::ParMesh::gtopo

Definition at line 127 of file pmesh.hpp.

bool mfem::ParMesh::have_face_nbr_data

Definition at line 130 of file pmesh.hpp.

MPI_Comm mfem::ParMesh::MyComm
protected

Definition at line 41 of file pmesh.hpp.

int mfem::ParMesh::MyRank
protected

Definition at line 42 of file pmesh.hpp.

int mfem::ParMesh::NRanks
protected

Definition at line 42 of file pmesh.hpp.

ParNCMesh* mfem::ParMesh::pncmesh

Definition at line 140 of file pmesh.hpp.

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

Definition at line 54 of file pmesh.hpp.

Table mfem::ParMesh::send_face_nbr_elements

Definition at line 137 of file pmesh.hpp.

Table mfem::ParMesh::send_face_nbr_vertices

Definition at line 138 of file pmesh.hpp.

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

Definition at line 55 of file pmesh.hpp.

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

Definition at line 44 of file pmesh.hpp.

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

Definition at line 45 of file pmesh.hpp.

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

Shared to local index mapping.

Definition at line 53 of file pmesh.hpp.


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