15 #include "../config/config.hpp"
19 #include "../general/communication.hpp"
20 #include "../general/globals.hpp"
44 Vert3(
int v0,
int v1,
int v2) {
v[0] = v0;
v[1] = v1;
v[2] = v2; }
45 void Set(
int v0,
int v1,
int v2) {
v[0] = v0;
v[1] = v1;
v[2] = v2; }
46 void Set(
const int *w) {
v[0] = w[0];
v[1] = w[1];
v[2] = w[2]; }
53 Vert4(
int v0,
int v1,
int v2,
int v3)
54 {
v[0] = v0;
v[1] = v1;
v[2] = v2;
v[3] = v3; }
55 void Set(
int v0,
int v1,
int v2,
int v3)
56 {
v[0] = v0;
v[1] = v1;
v[2] = v2;
v[3] = v3; }
57 void Set(
const int *w)
58 {
v[0] = w[0];
v[1] = w[1];
v[2] = w[2];
v[3] = w[3]; }
154 int nc_limit = 0)
override;
157 double threshold,
int nc_limit = 0,
158 int op = 1)
override;
178 Table* &edge_element);
192 int &nstria,
int &nsquad);
199 const Mesh &mesh,
int *partitioning,
206 const Table *edge_element);
301 ParMesh(MPI_Comm comm,
Mesh &mesh,
int *partitioning_ = NULL,
302 int part_method = 1);
312 ParMesh(MPI_Comm comm, std::istream &input,
bool refine =
true);
344 void Finalize(
bool refine =
false,
bool fix_orientation =
false)
override;
402 void GroupEdge(
int group,
int i,
int &edge,
int &o)
const;
403 void GroupTriangle(
int group,
int i,
int &face,
int &o)
const;
453 void SetCurvature(
int order,
bool discont =
false,
int space_dim = -1,
454 int ordering = 1)
override;
506 int mask = 31)
override;
556 long long ReduceInt(
int value)
const override;
583 void Save(
const char *fname,
int precision=16)
const override;
585 #ifdef MFEM_USE_ADIOS2
612 void SaveAsOne(
const char *fname,
int precision=16)
const;
622 bool high_order_output=
false,
623 int compression_level=0,
624 bool bdr=
false)
override;
627 void Load(std::istream &input,
int generate_edges = 0,
628 int refine = 1,
bool fix_orientation =
true)
override;
635 double &kappa_min,
double &kappa_max);
657 #ifdef MFEM_USE_ADIOS2
664 #endif // MFEM_USE_MPI
void UniformRefinement3D() override
Refine a mixed 3D mesh uniformly.
int GetNFaceNeighbors() const
void PrintAsOneXG(std::ostream &out=mfem::out)
Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'.
int Size() const
Return the logical size of the array.
void GroupQuadrilateral(int group, int i, int &face, int &o) const
void GetGhostFaceTransformation(FaceElementTransformations *FETr, Element::Type face_type, Geometry::Type face_geom)
void GetGlobalFaceIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
void LoadSharedEntities(std::istream &input)
long glob_offset_sequence
void UniformRefineGroups2D(int old_nv)
FaceElementTransformations * GetSharedFaceTransformationsByLocalIndex(int FaceNo, bool fill2=true)
ParMesh & operator=(ParMesh &&mesh)
Move assignment operator.
void MakeRefined_(ParMesh &orig_mesh, int ref_factor, int ref_type)
Internal function used in ParMesh::MakeRefined (and related constructor)
int FindSharedVertices(const int *partition, Table *vertex_element, ListOfIntegerSets &groups)
int GroupNVertices(int group) const
ElementTransformation * GetFaceNbrElementTransformation(int i)
void PrintInfo(std::ostream &out=mfem::out) override
Print various parallel mesh stats.
int GroupNEdges(int group) const
int GroupNTriangles(int group) const
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
int BuildLocalVertices(const Mesh &global_mesh, const int *partitioning, Array< int > &vert_global_local)
Fills out partitioned Mesh::vertices.
void UniformRefinement2D() override
Refine a mixed 2D mesh uniformly.
Array< Element * > face_nbr_elements
void GetSharedVertexCommunicator(int ordering, GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
void GetSharedVertexCommunicator(GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
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)
Data type dense matrix using column-major storage.
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
int GetLocalElementNum(long long global_element_num) const
void GetCharacteristics(double &h_min, double &h_max, double &kappa_min, double &kappa_max)
void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0) override
This function is not public anymore. Use GeneralRefinement instead.
void EnsureParNodes()
If the mesh is curved, make sure 'Nodes' is ParGridFunction.
void UniformRefineGroups3D(int old_nv, int old_nedges, const DSTable &old_v_to_v, const STable3D &old_faces, Array< int > *f2qf)
int FindPoints(DenseMatrix &point_mat, Array< int > &elem_ids, Array< IntegrationPoint > &ips, bool warn=true, InverseElementTransformation *inv_trans=NULL) override
Find the ids of the elements that contain the given points, and their corresponding reference coordin...
Abstract parallel finite element space.
Array< Vert3 > shared_trias
Array< int > face_nbr_vertices_offset
Array< int > face_nbr_group
bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1) override
NC version of GeneralDerefinement.
Data arrays will be written in ASCII format.
void RebalanceImpl(const Array< int > *partition)
void MarkTetMeshForRefinement(DSTable &v_to_v) override
void BuildFaceGroup(int ngroups, const Mesh &mesh, const Array< int > &face_group, int &nstria, int &nsquad)
void Set(int v0, int v1, int v2)
void GetSharedEdgeCommunicator(GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
STable3D * GetSharedFacesTable()
bool HasBoundaryElements() const override
Checks if any rank in the mesh has boundary elements.
Subdomain representation of a topological parent in another ParMesh.
bool WantSkipSharedMaster(const NCMesh::Master &master) const
A parallel extension of the NCMesh class.
Class for PUMI parallel meshes.
int FindSharedEdges(const Mesh &mesh, const int *partition, Table *&edge_element, ListOfIntegerSets &groups)
void ExchangeFaceNbrData()
void BuildSharedVertMapping(int nvert, const Table *vert_element, const Array< int > &vert_global_local)
Vert4(int v0, int v1, int v2, int v3)
void NURBSUniformRefinement() override
Refine NURBS mesh.
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.
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
Array< int > face_nbr_elements_offset
void Finalize(bool refine=false, bool fix_orientation=false) override
Finalize the construction of a general Mesh.
Symmetric 3D Table stored as an array of rows each of which has a stack of column, floor, number nodes. The number of the node is assigned by counting the nodes from zero as they are pushed into the table. Diagonals of any kind are not allowed so the row, column and floor must all be different for each node. Only one node is stored for all 6 symmetric entries that are indexable by unique triplets of row, column, and floor.
void GetBoundingBox(Vector &p_min, Vector &p_max, int ref=2)
Mesh GetSerialMesh(int save_rank) const
Array< Element * > shared_edges
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
void BuildSharedEdgeElems(int nedges, Mesh &mesh, const Array< int > &vert_global_local, const Table *edge_element)
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31) override
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Table send_face_nbr_vertices
void GetSharedTriCommunicator(GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
Type
Constants for the classes derived from Element.
void GetSharedQuadCommunicator(GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
VTKFormat
Data array format for VTK and VTU files.
void PrintVTU(std::string pathname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr=false) override
void Set(int v0, int v1, int v2, int v3)
int GetNFbyType(FaceType type) const override
Returns the number of local faces according to the requested type, does not count master non-conformi...
void GetGlobalVertexIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Table * face_nbr_el_to_face
long long ReduceInt(int value) const override
Utility function: sum integers from all processors (Allreduce).
STable3D * GetFaceNbrElementToFaceTable(int ret_ftbl=0)
Array< Vert4 > shared_quads
int GetFaceNbrGroup(int fn) const
int GetNFaceNeighborElements() const
void PrintXG(std::ostream &out=mfem::out) const override
int GroupNQuadrilaterals(int group) const
void FindSharedFaces(const Mesh &mesh, const int *partition, Array< int > &face_group, ListOfIntegerSets &groups)
void LocalRefinement(const Array< int > &marked_el, int type=3) override
This function is not public anymore. Use GeneralRefinement instead.
void RefineGroups(const DSTable &v_to_v, int *middle)
Update the groups after triangle refinement.
Array< Vertex > face_nbr_vertices
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void GetSharedTriCommunicator(int ordering, GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
void Swap(ParMesh &other)
void PrintSharedEntities(const char *fname_prefix) const
Debugging method.
long long glob_elem_offset
void GetFaceNbrElementFaces(int i, Array< int > &fcs, Array< int > &cor) const
ParMesh()
Default constructor. Create an empty ParMesh.
void SetAttributes() override
void GetSharedEdgeCommunicator(int ordering, GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
void BuildVertexGroup(int ngroups, const Table &vert_element)
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
void DistributeAttributes(Array< int > &attr)
Ensure that bdr_attributes and attributes agree across processors.
long long GetGlobalElementNum(int local_element_num) const
Map a local element number to a global element number.
void GroupEdge(int group, int i, int &edge, int &o) const
void Swap(Mesh &other, bool non_geometry)
void BuildEdgeGroup(int ngroups, const Table &edge_element)
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.
void GetSharedQuadCommunicator(int ordering, GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
void ComputeGlobalElementOffset() const
void ExchangeFaceNbrNodes()
int BuildLocalElements(const Mesh &global_mesh, const int *partitioning, const Array< int > &vert_global_local)
Fills out partitioned Mesh::elements.
Table group_svert
Shared objects in each group.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1) override
bool DecodeFaceSplittings(HashTable< Hashed2 > &v_to_v, const int *v, const Array< unsigned > &codes, int &pos)
double GetFaceNbrElementSize(int i, int type=0)
void GroupTriangle(int group, int i, int &face, int &o) const
void GetGlobalElementIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are supported.
int NGroups() const
Return the number of groups.
void GetGlobalEdgeIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
void SetPrintShared(bool print)
int GroupVertex(int group, int i) const
IsoparametricTransformation FaceNbrTransformation
int GetEdgeSplittings(Element *edge, const DSTable &v_to_v, int *middle)
Return a number(0-1) identifying how the given edge has been split.
static ParMesh MakeSimplicial(ParMesh &orig_mesh)
void PrintAsOne(std::ostream &out=mfem::out) const
void SetNodalFESpace(FiniteElementSpace *nfes) override
Array< int > svert_lvert
Shared to local index mapping.
void SaveAsOne(const char *fname, int precision=16) const
void Print(std::ostream &out=mfem::out) const override
Table send_face_nbr_elements
Table * GetFaceToAllElementTable() const
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void Load(std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true) override
Parallel version of Mesh::Load().
void ParPrint(std::ostream &out) const
Save the mesh in a parallel mesh format.
Class for parallel meshes.
void PrintAsSerial(std::ostream &out=mfem::out) const
Abstract data type element.
int GetFaceNbrRank(int fn) const
void Save(const char *fname, int precision=16) const override
void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[], Array< HYPRE_BigInt > *offsets[]) const
Vert3(int v0, int v1, int v2)
MFEM_DEPRECATED void ReorientTetMesh() override
See the remarks for the serial version in mesh.hpp.