CArl
Code Arlequin / C++ implementation
intersection_tools.h
Go to the documentation of this file.
1 /*
2  * intersection_tools.h
3  *
4  * Created on: Apr 17, 2016
5  * Author: Thiago Milanetto Schlittler
6  */
7 
8 #ifndef COMMON_INTERSECTIONS_PARALLEL_INTERSECTION_TOOLS_H_
9 #define COMMON_INTERSECTIONS_PARALLEL_INTERSECTION_TOOLS_H_
10 
11 #include "carl_headers.h"
12 #include "mesh_tables.h"
13 
14 namespace carl
15 {
16 
17 // ************************
18 // Intersection_Tools class
19 // ************************
20 
26 {
27 protected:
28 
29  // Nef polyhedrons
34 
39 
40  // Exact point vectors
41  std::vector<ExactPoint_3> m_exact_points_A;
42  std::vector<ExactPoint_3> m_exact_points_B;
43  std::vector<ExactPoint_3> m_exact_points_C;
44 
46 
47  CGAL::Convex_hull_traits_3<ExactKernel> ExactHullTraits;
48 
49  // CGAL Kernel converters
52 
53  // Dummy exact tetrahedron and triangle
56 
58 
59  // Perflog and debug variables
61  libMesh::PerfLog m_perf_log;
62 
63  // Data structures
64  /*
65  * CGAL's "do_intersect" methods only works with tetrahedrons and triangles, at most.
66  * The vectors below map an TET or HEX element to its decomposition with tetrahedrons and
67  * triangles (which is trivial for the TET elements).
68  */
69  std::vector<std::vector<unsigned int> > m_TET_tetrahedrons;
70  std::vector<std::vector<unsigned int> > m_TET_triangles;
71  std::vector<std::vector<unsigned int> > m_TET_edges;
72 
73  std::vector<std::vector<unsigned int> > m_HEX_tetrahedrons;
74  std::vector<std::vector<unsigned int> > m_HEX_triangles;
75  std::vector<std::vector<unsigned int> > m_HEX_edges;
76 
77  std::vector<std::vector<unsigned int> > * m_elem_C_tetrahedrons;
78  std::vector<std::vector<unsigned int> > * m_elem_C_triangles;
79 
80  // PROTECTED methods
90  void set_element_indexes();
91 
101  std::vector<ExactPoint_3> & elem_C_points,
102  std::vector<std::vector<unsigned int> > & elem_C_tetras,
103  std::vector<std::vector<unsigned int> > & elem_C_triangles,
104  std::vector<ExactPoint_3> & elem_D_points,
105  std::vector<std::vector<unsigned int> > & elem_D_tetras,
106  std::vector<std::vector<unsigned int> > & elem_D_triangles);
107 
108 public:
109 
110  // Constructors
114  Intersection_Tools(const libMesh::Elem * elem_C, double Min_Inter_Volume = 1E-21, bool bDoPerf_log = false) :
115  m_Min_Inter_Volume { Min_Inter_Volume },
116  MASTER_bPerfLog_intersection_tools {bDoPerf_log},
117  m_perf_log { libMesh::PerfLog("Intersection tools", MASTER_bPerfLog_intersection_tools) }
118  {
119  m_exact_points_A.resize(8);
120  m_exact_points_B.resize(8);
121  m_exact_points_C.resize(8);
122 
123  m_dummyPoly.reserve(8,18,12);
124 
126 
127  this->set_element_indexes();
128 
129  m_elem_C_tetrahedrons = NULL;
130  m_elem_C_triangles = NULL;
131  };
132 
135  Intersection_Tools(double Min_Inter_Volume = 1E-21, bool bDoPerf_log = false) :
136  m_Min_Inter_Volume { Min_Inter_Volume },
137  MASTER_bPerfLog_intersection_tools {bDoPerf_log},
138  m_perf_log { libMesh::PerfLog("Intersection tools", MASTER_bPerfLog_intersection_tools) }
139  {
140  m_exact_points_A.resize(8);
141  m_exact_points_B.resize(8);
142  m_exact_points_C.resize(8);
143 
144  m_dummyPoly.reserve(8,18,12);
145 
146  m_nef_C.clear(Nef_Polyhedron::EMPTY);
147 
148  this->set_element_indexes();
149 
150  m_elem_C_tetrahedrons = NULL;
151  m_elem_C_triangles = NULL;
152  };
153 
154  // PUBLIC methods
162  const libMesh::Elem * FindFirstIntersection( const libMesh::Elem * Query_elem,
163  std::unique_ptr<libMesh::PointLocatorBase> & point_locator,
164  bool bGuaranteeQueryIsInMesh = false);
165 
169  bool FindAllIntersection( const libMesh::Elem * Query_elem,
170  std::unique_ptr<libMesh::PointLocatorBase> & point_locator,
171  std::set<unsigned int> & Intersecting_elems);
172 
177  bool libMesh_exact_do_intersect(const libMesh::Elem * elem_A,
178  const libMesh::Elem * elem_B);
179 
182  bool libMesh_exact_intersection(const libMesh::Elem * elem_A,
183  const libMesh::Elem * elem_B,
184  std::set<Point_3> & points_out,
185  bool bCreateNewNefForA = true,
186  bool bConvertPoints = true,
187  bool bTestNeeded = true);
188 
194  bool libMesh_exact_do_intersect_inside_coupling(const libMesh::Elem * elem_A,
195  const libMesh::Elem * elem_B,
196  bool bCreateNewNefForA = true);
197 
205  bool libMesh_exact_intersection_inside_coupling(const libMesh::Elem * elem_A,
206  const libMesh::Elem * elem_B,
207  std::set<libMesh::Point> & points_out,
208  bool bCreateNewNefForA = true,
209  bool bConvertPoints = true,
210  bool bTestNeeded = false);
211 
214  void libmesh_set_coupling_nef_polyhedron(const libMesh::Elem * elem_C);
215 
218  void convert_elem_to_exact_points( const libMesh::Elem * elem_input,
219  std::vector<ExactPoint_3> & points_output);
220 
223  void convert_exact_points_to_Nef(std::vector<ExactPoint_3>::const_iterator it_begin,
224  std::vector<ExactPoint_3>::const_iterator it_end,
225  Nef_Polyhedron & nef_out);
226 };
227 }
228 
229 #endif /* COMMON_INTERSECTIONS_PARALLEL_INTERSECTION_TOOLS_H_ */
CGAL::Polyhedron_3< ExactKernel > ExactPolyhedron
Definition: CGAL_typedefs.h:99
std::vector< std::vector< unsigned int > > m_HEX_edges
Vertex indices for a HEX* element edges (not used!)
Kernel_to_ExactKernel ConvertInexactToExact
Convert inexact CGAL constructs to exact constructs.
std::vector< std::vector< unsigned int > > m_TET_edges
Vertex indices for a TET* element edges (not used!)
const libMesh::Elem * FindFirstIntersection(const libMesh::Elem *Query_elem, std::unique_ptr< libMesh::PointLocatorBase > &point_locator, bool bGuaranteeQueryIsInMesh=false)
Find an element from a mesh intersecting Query_elem, using the mesh's libMesh::PointLocatorBase.
std::vector< std::vector< unsigned int > > m_HEX_tetrahedrons
Vertex indices for a HEX* element tetrahedron decomposition.
double m_Min_Inter_Volume
Minimal volume cutoff for the intersections.
ExactKernel::Triangle_3 ExactTriangle_3
Definition: CGAL_typedefs.h:95
Intersection_Tools(double Min_Inter_Volume=1E-21, bool bDoPerf_log=false)
Constructor preallocating most objects and building the decomposition mappings.
ExactKernel_to_Kernel ConvertExactToInexact
Convert exact CGAL constructs to inexact constructs.
void set_element_indexes()
Build the tetrahedron and triangle decomposition mappings.
std::vector< ExactPoint_3 > m_exact_points_B
ExactTetrahedron m_test_tetra
std::vector< ExactPoint_3 > m_exact_points_A
CGAL Exact point vectors.
Nef_Polyhedron m_nef_I_AC
Nef polyhedron used to save the intersection between m_nef_A and m_nef_C, reducing the number of inte...
void libmesh_set_coupling_nef_polyhedron(const libMesh::Elem *elem_C)
Set the Nef polyhedron for the coupling mesh element.
bool libMesh_exact_intersection(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, std::set< Point_3 > &points_out, bool bCreateNewNefForA=true, bool bConvertPoints=true, bool bTestNeeded=true)
Build the intersection between two elements.
std::vector< std::vector< unsigned int > > m_HEX_triangles
Vertex indices for a HEX* element surface triangles decomposition.
CGAL::Cartesian_converter< ExactKernel, Kernel > ExactKernel_to_Kernel
Definition: CGAL_typedefs.h:86
void convert_elem_to_exact_points(const libMesh::Elem *elem_input, std::vector< ExactPoint_3 > &points_output)
Convert an libMesh element to a list of CGAL exact points.
bool FindAllIntersection(const libMesh::Elem *Query_elem, std::unique_ptr< libMesh::PointLocatorBase > &point_locator, std::set< unsigned int > &Intersecting_elems)
Find all elements from the mesh intersecting Query_elem, using the mesh's libMesh::PointLocatorBase.
CGAL::Nef_polyhedron_3< ExactKernel, CGAL::SNC_indexed_items > Nef_Polyhedron
CGAL::Convex_hull_traits_3< ExactKernel > ExactHullTraits
CGAL ExactKernel's convex hull traits.
bool MASTER_bPerfLog_intersection_tools
Do performance log? Default: false.
Class with a series of methods to find the intersections between libMesh's elements, using CGAL internally.
ExactPolyhedron m_dummyPoly
Intersection polyhedron.
std::vector< std::vector< unsigned int > > m_TET_triangles
Vertex indices for a TET* element surface triangles.
bool elements_do_intersect(std::vector< ExactPoint_3 > &elem_C_points, std::vector< std::vector< unsigned int > > &elem_C_tetras, std::vector< std::vector< unsigned int > > &elem_C_triangles, std::vector< ExactPoint_3 > &elem_D_points, std::vector< std::vector< unsigned int > > &elem_D_tetras, std::vector< std::vector< unsigned int > > &elem_D_triangles)
Determinate if two elements intersect, using their points and decompositions.
ExactTriangle_3 m_test_triangle
Intersection_Tools(const libMesh::Elem *elem_C, double Min_Inter_Volume=1E-21, bool bDoPerf_log=false)
Constructor with a pre-defined coupling element. It preallocates most objects and build the decomposi...
CGAL::Tetrahedron_3< ExactKernel > ExactTetrahedron
CGAL::Cartesian_converter< Kernel, ExactKernel > Kernel_to_ExactKernel
Definition: CGAL_typedefs.h:85
bool libMesh_exact_do_intersect_inside_coupling(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, bool bCreateNewNefForA=true)
Determinate if two elements intersect inside the coupling region.
std::vector< std::vector< unsigned int > > * m_elem_C_tetrahedrons
Pointer for the appropriate tetrahedron vertex index list (useless!)
std::vector< ExactPoint_3 > m_exact_points_C
bool libMesh_exact_do_intersect(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B)
Determinate if two elements intersect.
Nef_Polyhedron m_nef_A
Nef polyhedrons for the elements.
bool libMesh_exact_intersection_inside_coupling(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, std::set< libMesh::Point > &points_out, bool bCreateNewNefForA=true, bool bConvertPoints=true, bool bTestNeeded=false)
Build the intersection between two elements intersect inside the coupling region. ...
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
void convert_exact_points_to_Nef(std::vector< ExactPoint_3 >::const_iterator it_begin, std::vector< ExactPoint_3 >::const_iterator it_end, Nef_Polyhedron &nef_out)
Convert an CGAL exact point set to a CGAL Nef polyhedron.
std::vector< std::vector< unsigned int > > * m_elem_C_triangles
Pointer for the appropriate triangle vertex index list (useless!)
std::vector< std::vector< unsigned int > > m_TET_tetrahedrons
Vertex indices for a TET* element tetrahedron.