52 Vert3(
int v0,
int v1,
int v2) {
v[0] = v0;
v[1] = v1;
v[2] = v2; }
53 void Set(
int v0,
int v1,
int v2) {
v[0] = v0;
v[1] = v1;
v[2] = v2; }
54 void Set(
const int *w) {
v[0] = w[0];
v[1] = w[1];
v[2] = w[2]; }
61 Vert4(
int v0,
int v1,
int v2,
int v3)
62 {
v[0] = v0;
v[1] = v1;
v[2] = v2;
v[3] = v3; }
63 void Set(
int v0,
int v1,
int v2,
int v3)
64 {
v[0] = v0;
v[1] = v1;
v[2] = v2;
v[3] = v3; }
65 void Set(
const int *w)
66 {
v[0] = w[0];
v[1] = w[1];
v[2] = w[2];
v[3] = w[3]; }
149 const std::unique_ptr<STable3D> &shared_faces,
150 int elem,
int start,
int end,
const int fverts[][N]);
159 MFEM_ASSERT(FElTr,
"Missing FaceElementTransformations object!");
199 int nc_limit = 0)
override;
202 real_t threshold,
int nc_limit = 0,
203 int op = 1)
override;
223 Table* &edge_element);
237 int &nstria,
int &nsquad);
244 const Mesh &mesh,
int *partitioning,
251 const Table *edge_element);
346 ParMesh(MPI_Comm comm,
Mesh &mesh,
int *partitioning_ = NULL,
347 int part_method = 1);
362 ParMesh(MPI_Comm comm, std::istream &input,
bool refine =
true,
363 int generate_edges = 1,
bool fix_orientation =
true);
395 void Finalize(
bool refine =
false,
bool fix_orientation =
false)
override;
453 void GroupEdge(
int group,
int i,
int &edge,
int &o)
const;
454 void GroupTriangle(
int group,
int i,
int &face,
int &o)
const;
508 void SetCurvature(
int order,
bool discont =
false,
int space_dim = -1,
509 int ordering = 1)
override;
574 int mask = 31)
const override;
594 bool fill2 =
true)
const;
614 bool fill2 =
true)
const;
650 { MFEM_ABORT(
"Generation of boundary elements works properly only on serial meshes."); }
656 long long ReduceInt(
int value)
const override;
670 void ParPrint(std::ostream &
out,
const std::string &comments =
"")
const;
682 const std::string &comments =
"")
const override;
688 void Save(
const std::string &fname,
int precision=16)
const override;
690#ifdef MFEM_USE_ADIOS2
707 const std::string &comments =
"")
const;
714 const std::string &comments =
"")
const;
722 void SaveAsOne(
const std::string &fname,
int precision=16)
const;
732 bool high_order_output=
false,
733 int compression_level=0,
734 bool bdr=
false)
override;
737 void Load(std::istream &input,
int generate_edges = 0,
738 int refine = 1,
bool fix_orientation =
true)
override;
int Size() const
Return the logical size of the array.
Data type dense matrix using column-major storage.
Abstract data type element.
Type
Constants for the classes derived from Element.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
int NGroups() const
Return the number of groups.
bool FaceIsTrueInterior(int FaceNo) const
void Swap(Mesh &other, bool non_geometry)
Abstract parallel finite element space.
Class for parallel meshes.
Mesh GetSerialMesh(int save_rank) const
void GetCharacteristics(real_t &h_min, real_t &h_max, real_t &kappa_min, real_t &kappa_max)
void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0) override
This function is not public anymore. Use GeneralRefinement instead.
int GroupNQuadrilaterals(int group) const
ElementTransformation * GetFaceNbrElementTransformation(int FaceNo)
Returns a pointer to the transformation defining the i-th face neighbor.
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 GetGhostFaceTransformation(FaceElementTransformations *FElTr, Element::Type face_type, Geometry::Type face_geom) const
Table send_face_nbr_elements
Array< int > face_nbr_vertices_offset
void BuildVertexGroup(int ngroups, const Table &vert_element)
void PrintVTU(std::string pathname, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr=false) override
void NURBSUniformRefinement(int rf=2, real_t tol=1.0e-12) override
Refine NURBS mesh, with an optional refinement factor.
void PrintXG(std::ostream &out=mfem::out) const override
Array< Element * > shared_edges
int GetEdgeSplittings(Element *edge, const DSTable &v_to_v, int *middle)
Return a number(0-1) identifying how the given edge has been split.
void RefineGroups(const DSTable &v_to_v, int *middle)
Update the groups after triangle refinement.
void ParPrint(std::ostream &out, const std::string &comments="") const
void GetSharedTriCommunicator(int ordering, GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
bool NonconformingDerefinement(Array< real_t > &elem_error, real_t threshold, int nc_limit=0, int op=1) override
NC version of GeneralDerefinement.
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
void SaveAsOne(const std::string &fname, int precision=16) const
void ExchangeFaceNbrData()
void BuildEdgeGroup(int ngroups, const Table &edge_element)
Table group_svert
Shared objects in each group.
MFEM_DEPRECATED void ReorientTetMesh() override
See the remarks for the serial version in mesh.hpp.
int GroupVertex(int group, int i) const
void UniformRefinement2D() override
Refine a mixed 2D mesh uniformly.
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)
bool WantSkipSharedMaster(const NCMesh::Master &master) const
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31) override
void BuildSharedVertMapping(int nvert, const Table *vert_element, const Array< int > &vert_global_local)
long long GetGlobalElementNum(int local_element_num) const
Map a local element number to a global element number.
Table * GetFaceToAllElementTable() const
void GetGlobalElementIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are supported.
void AddTriFaces(const Array< int > &v, const std::unique_ptr< STable3D > &faces, const std::unique_ptr< STable3D > &shared_faces, int elem, int start, int end, const int fverts[][N])
Helper function for adding triangle face neighbor element to face table entries. Have to use a templa...
void FindSharedFaces(const Mesh &mesh, const int *partition, Array< int > &face_group, ListOfIntegerSets &groups)
long long glob_elem_offset
void GetSharedEdgeCommunicator(int ordering, GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
void PrintInfo(std::ostream &out=mfem::out) override
Print various parallel mesh stats.
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.
int GetNFaceNeighborElements() const
bool HasBoundaryElements() const override
Checks if any rank in the mesh has boundary elements.
void LocalRefinement(const Array< int > &marked_el, int type=3) override
This function is not public anymore. Use GeneralRefinement instead.
void GetFaceNbrElementFaces(int i, Array< int > &faces, Array< int > &orientation) const
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1) override
Set the curvature of the mesh nodes using the given polynomial degree.
void SetPrintShared(bool print)
void GetSharedTriCommunicator(GroupCommunicator &stria_comm) const
Get the shared face triangles GroupCommunicator.
long glob_offset_sequence
int GetLocalElementNum(long long global_element_num) const
void SetAttributes() override
Determine the sets of unique attribute values in domain and boundary elements.
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...
void GetSharedQuadCommunicator(int ordering, GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
void BuildFaceNbrElementToFaceTable()
void MakeRefined_(ParMesh &orig_mesh, int ref_factor, int ref_type)
Internal function used in ParMesh::MakeRefined (and related constructor)
void GetGlobalFaceIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
Array< Vertex > face_nbr_vertices
void UniformRefineGroups2D(int old_nv)
void GetGlobalEdgeIndices(Array< HYPRE_BigInt > &gi) const
AMR meshes are not supported.
void ExchangeFaceNbrNodes()
void SetNodalFESpace(FiniteElementSpace *nfes) override
void UniformRefinement3D() override
Refine a mixed 3D mesh uniformly.
long long ReduceInt(int value) const override
Utility function: sum integers from all processors (Allreduce).
Table send_face_nbr_vertices
ParMesh()
Default constructor. Create an empty ParMesh.
int BuildLocalElements(const Mesh &global_mesh, const int *partitioning, const Array< int > &vert_global_local)
Fills out partitioned Mesh::elements.
bool DecodeFaceSplittings(HashTable< Hashed2 > &v_to_v, const int *v, const Array< unsigned > &codes, int &pos)
static ParMesh MakeSimplicial(ParMesh &orig_mesh)
Array< Vert4 > shared_quads
void GenerateBoundaryElements() override
void LoadSharedEntities(std::istream &input)
void Finalize(bool refine=false, bool fix_orientation=false) override
Finalize the construction of a general Mesh.
int GetFaceNbrGroup(int fn) const
void GetSharedEdgeCommunicator(GroupCommunicator &sedge_comm) const
Get the shared edges GroupCommunicator.
void BuildFaceGroup(int ngroups, const Mesh &mesh, const Array< int > &face_group, int &nstria, int &nsquad)
void BuildSharedEdgeElems(int nedges, Mesh &mesh, const Array< int > &vert_global_local, const Table *edge_element)
Array< int > svert_lvert
Shared to local index mapping.
Array< Element * > face_nbr_elements
FaceElementTransformations * GetSharedFaceTransformationsByLocalIndex(int FaceNo, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the face index FaceNo....
void PrintSharedEntities(const std::string &fname_prefix) const
Debugging method.
int GroupNTriangles(int group) const
void GetSharedVertexCommunicator(int ordering, GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
void EnsureParNodes()
If the mesh is curved, make sure 'Nodes' is ParGridFunction.
int GroupNEdges(int group) const
void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[], Array< HYPRE_BigInt > *offsets[]) const
void MarkTetMeshForRefinement(const DSTable &v_to_v) override
Array< int > face_nbr_group
Array< int > face_nbr_elements_offset
ParMesh & operator=(ParMesh &&mesh)
Move assignment operator.
void Load(std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true) override
Parallel version of Mesh::Load().
int GetNFaceNeighbors() const
void UniformRefineGroups3D(int old_nv, int old_nedges, const DSTable &old_v_to_v, const STable3D &old_faces, Array< int > *f2qf)
void GetGhostFaceTransformation(FaceElementTransformations &FElTr, Element::Type face_type, Geometry::Type face_geom) const
void RebalanceImpl(const Array< int > *partition)
int FindSharedEdges(const Mesh &mesh, const int *partition, Table *&edge_element, ListOfIntegerSets &groups)
void PrintAsOneXG(std::ostream &out=mfem::out)
Old mesh format (Netgen/Truegrid) version of 'PrintAsOne'.
std::unique_ptr< Table > face_nbr_el_to_face
Array< Vert3 > shared_trias
void GroupQuadrilateral(int group, int i, int &face, int &o) const
void Save(const std::string &fname, int precision=16) const override
void DistributeAttributes(Array< int > &attr)
Ensure that bdr_attributes and attributes agree across processors.
ParMesh & operator=(const ParMesh &mesh)=delete
Explicitly delete the copy assignment operator.
STable3D * GetSharedFacesTable()
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
real_t GetFaceNbrElementSize(int i, int type=0)
std::unique_ptr< Table > face_nbr_el_ori
orientations for each face (from nbr processor)
bool FaceIsTrueInterior(int FaceNo) const
int GetFaceNbrRank(int fn) const
int FindSharedVertices(const int *partition, Table *vertex_element, ListOfIntegerSets &groups)
void GroupTriangle(int group, int i, int &face, int &o) const
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the shared face index sf...
void PrintAsSerial(std::ostream &out=mfem::out, const std::string &comments="") const
void GetSharedVertexCommunicator(GroupCommunicator &svert_comm) const
Get the shared vertices GroupCommunicator.
void GetSharedQuadCommunicator(GroupCommunicator &squad_comm) const
Get the shared face quadrilaterals GroupCommunicator.
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.
int BuildLocalVertices(const Mesh &global_mesh, const int *partitioning, Array< int > &vert_global_local)
Fills out partitioned Mesh::vertices.
void ComputeGlobalElementOffset() const
void GetBoundingBox(Vector &p_min, Vector &p_max, int ref=2)
void PrintAsOne(std::ostream &out=mfem::out, const std::string &comments="") const
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
void Swap(ParMesh &other)
void GroupEdge(int group, int i, int &edge, int &o) const
int GroupNVertices(int group) const
A parallel extension of the NCMesh class.
Class for PUMI parallel meshes.
Subdomain representation of a topological parent in another ParMesh.
Symmetric 3D Table stored as an array of rows each of which has a stack of column,...
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
VTKFormat
Data array format for VTK and VTU files.
@ ASCII
Data arrays will be written in ASCII format.
void Set(int v0, int v1, int v2)
Vert3(int v0, int v1, int v2)
void Set(int v0, int v1, int v2, int v3)
Vert4(int v0, int v1, int v2, int v3)