8 #ifndef COMMON_INTERSECTIONS_PARALLEL_INTERSECTION_MESH_LIBMESH_H_
9 #define COMMON_INTERSECTIONS_PARALLEL_INTERSECTION_MESH_LIBMESH_H_
39 const libMesh::Parallel::Communicator&
m_comm;
142 int map_preallocation = 1E6,
long grid_n_min = static_cast<long>(1E9),
bool debugOutput =
false) :
143 m_comm { mesh.comm() },
144 m_nodes { m_comm.size() },
145 m_rank { m_comm.rank() },
146 m_global_comm { mesh_A.comm() },
147 m_global_nodes { m_global_comm.size() },
148 m_global_rank { m_global_comm.rank() },
149 m_libMesh_Mesh { mesh },
150 m_libMesh_PolyhedronMesh { libMesh::Mesh(m_comm) },
151 m_TetGenInterface { libMesh::TetGenMeshInterface(m_libMesh_PolyhedronMesh) },
152 m_nb_of_intersections { 0 },
154 m_GridN { std::vector<long >(3,-1) },
155 m_dummy_discrete_point { std::vector<long >(3,-1) },
156 m_Grid_MinPoint { libMesh::Point(0,0,0) },
157 m_Grid_MaxPoint { libMesh::Point(1,1,1) },
158 m_GridN_min { grid_n_min },
159 m_nb_of_elements { 0 },
160 m_nb_of_vertices { 0 },
161 m_nb_of_points { 0 },
162 m_bMeshFinalized {
false },
163 m_bPrintDebug { debugOutput },
164 m_MeshingMethod { MeshingMethod }
167 m_intersection_point_indexes.resize(24);
171 m_libMesh_Mesh.reserve_nodes(map_preallocation);
172 m_libMesh_Mesh.reserve_elem(10*map_preallocation);
174 m_intersection_pairs.reserve(map_preallocation);
175 m_intersection_couplings.reserve(map_preallocation);
176 m_intersection_element_range.reserve(map_preallocation);
178 m_libMesh_Mesh.allow_renumbering(
false);
180 switch(m_MeshingMethod)
183 std::cout <<
" -> Using TetGen to generate mesh" << std::endl;
186 std::cout <<
" -> Using CGAL to generate mesh" << std::endl;
192 const libMesh::ReplicatedMesh &
mesh();
208 void set_grid_constraints(
const libMesh::Mesh & mesh_A,
const libMesh::Mesh & mesh_B,
double vol_tol = -1);
215 unsigned int elem_idx_A,
216 unsigned int elem_idx_B,
217 unsigned int elem_idx_C);
220 void export_intersection_data(
const std::string & filename_base,
const std::string & mesh_format = std::string(
".e"));
double get_total_volume()
Calculate the intersection mesh's total volume.
libMesh::Point & min_point()
Returns the minimal point of the discrete grid.
unsigned int m_nb_of_points
(Unused!)
libMesh::ReplicatedMesh & m_libMesh_Mesh
Address of the intersection mesh.
void update_intersection_mesh()
Update the intersection mesh with the current polyhedron mesh data.
const libMesh::Parallel::Communicator & m_comm
Intersection mesh communicator.
unsigned int m_nb_of_intersections
libMesh::TetGenMeshInterface m_TetGenInterface
TetGen interface.
libMesh::Point & max_point()
Returns the maximum point of the discrete grid.
void convert_to_discrete(const libMesh::Point &iPoint, std::vector< long > &oPoint)
Convert a real valued point to a discrete point.
const libMesh::ReplicatedMesh & mesh()
Returns the address of the intersection mesh.
std::unordered_map< std::vector< long >, unsigned int, PointHash_3D, PointHash_3D_Equal > m_discrete_vertices
Unordered map (or hash) for the discrete points (key).
long grid_min_size()
Returns the minimum discrete grid size.
const unsigned int m_nodes
Intersection mesh number of processors (useless!).
void increase_intersection_mesh(const std::set< libMesh::Point > &input_points, unsigned int elem_idx_A, unsigned int elem_idx_B, unsigned int elem_idx_C)
Increase the intersection mesh with the current polyhedron mesh and update its data structures...
void export_intersection_data(const std::string &filename_base, const std::string &mesh_format=std::string(".e"))
Export both the intersection mesh its data structures.
std::unordered_map< unsigned int, std::pair< unsigned int, unsigned int > > m_intersection_pairs
Unordered map containing the element pair associated to each intersection (key).
IntersectionMeshingMethod m_MeshingMethod
Choice of meshing algorithm.
CGAL::Cartesian_converter< ExactKernel, Kernel > ExactKernel_to_Kernel
libMesh::Point m_Grid_MinPoint
Minimal point of the discrete grid.
void set_grid_constraints(const libMesh::Mesh &mesh_A, const libMesh::Mesh &mesh_B, double vol_tol=-1)
Set the boundaries of the discrete points grid, taking into account the geometry of meshes mesh_A and...
std::unordered_map< unsigned int, unsigned int > m_intersection_couplings
Unordered map associating each intersection (key) to a coupling element (useless!) ...
double eps()
Returns the discrete grid step.
void triangulate_intersection(const std::set< libMesh::Point > &input_points)
Triangulate an intersection defined by a set of libMesh::Points, using either CGAL or Tetgen...
ExactKernel_to_Kernel ConvertExactToInexact
Convert exact CGAL constructs to inexact constructs.
const unsigned int m_global_nodes
Full system number of processors.
unsigned int m_nb_of_vertices
Number of vertices of the intersection mesh.
CGAL::Delaunay_triangulation_3< Kernel, Tds_3 > DT_3
Mesh_Intersection(libMesh::ReplicatedMesh &mesh, const libMesh::Mesh &mesh_A, const libMesh::Mesh &mesh_B, IntersectionMeshingMethod MeshingMethod=IntersectionMeshingMethod::CGAL, int map_preallocation=1E6, long grid_n_min=static_cast< long >(1E9), bool debugOutput=false)
IntersectionMeshingMethod
std::unordered_map< unsigned int, std::pair< unsigned int, unsigned int > > m_intersection_element_range
Unordered map associating each intersection (key) to a range of elements in the intersection mesh...
bool m_bPrintDebug
Print debug information? Default: false.
Kernel_to_ExactKernel ConvertInexactToExact
Convert inexact CGAL constructs to exact constructs.
double min_vol()
Returns the minimal polyhedron volume.
double m_eps
Precision of the discrete point grid.
bool m_bMeshFinalized
Boolean controlling if the intersection mesh was finished or not.
void initialize()
Initialize the data structures.
std::vector< long > & grid_sizes()
Returns the discrete grid sizes.
std::vector< long > m_dummy_discrete_point
libMesh::Mesh m_libMesh_PolyhedronMesh
libMesh Mesh containing the tetrahedrization of an intersection polyhedron (the polyhedron mesh)...
double get_intersection_volume(std::set< libMesh::Point > &input_points)
Calculate the volume of an polyhedron defined by the point set.
CGAL::Cartesian_converter< Kernel, ExactKernel > Kernel_to_ExactKernel
Mesh_Intersection()
(Protected) default constructor
Class used to construct the intersection meshes.
void preallocate_grid(int map_preallocation)
Preallocate the discrete points grid.
std::vector< unsigned int > m_intersection_point_indexes
(Unused!)
void update_intersection_pairs(unsigned int elem_idx_A, unsigned int elem_idx_B, unsigned int inter_id)
Update the intersection pair map with a new pair associated to the intersection no. inter_id.
double m_vol_tol
Minimal volume for the intersection polyhedrons.
long m_GridN_min
Minimal number of integers over each dimension of the discrete grid.
DT_3 m_CGAL_PolyhedronMesh
Auxiliary CGAL mesh, used if CGAL is chosen for the tetrahedrization.
void prepare_for_use()
Prepare the intersection mesh for use.
const unsigned int m_global_rank
Full system processor rank.
const unsigned int m_rank
Intersection mesh processor rank (useless!).
std::vector< long > m_GridN
Dimensions of the discrete point grid.
const libMesh::Parallel::Communicator & m_global_comm
Full system communicator.
unsigned int m_nb_of_elements
Number of elements of the intersection mesh.
void update_intersection_element_range(unsigned int range_start, unsigned int range_end, unsigned int inter_id)
Update the intersection element range map with a new range associated to the intersection no...
libMesh::Point m_Grid_MaxPoint
Maximum point of the discrete grid.