CELIA3D  1.0
Fluid-structure interaction using cut-cells
intersections.hpp
Go to the documentation of this file.
1 //Copyright 2017 Laurent Monasse
2 
3 /*
4  This file is part of CELIA3D.
5 
6  CELIA3D is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  CELIA3D is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with CELIA3D. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
26 #ifndef INTERSECTIONS
27 #define INTERSECTIONS
28 
29 #include <iostream>
30 #include <stdio.h>
31 #include <vector>
32 #include <math.h>
33 #include <cassert>
34 
35 #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
36 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
37 #include <CGAL/intersections.h>
38 #include <CGAL/Bbox_3.h>
39 #include <CGAL/Bbox_2.h>
40 #include <CGAL/Timer.h>
41 #include <CGAL/Triangulation_3.h>
42 #include <CGAL/Triangulation_2.h>
43 #include <CGAL/Polyhedron_3.h>
44 #include <CGAL/Tetrahedron_3.h>
45 #include <CGAL/convex_hull_3.h>
46 #include <CGAL/squared_distance_3.h>
47 #include <CGAL/box_intersection_d.h>
48 #include <CGAL/centroid.h>
49 #include <CGAL/number_utils.h>
50 
51 typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;
52 typedef CGAL::Exact_predicates_inexact_constructions_kernel IK;
53 typedef CGAL::Cartesian_converter<Kernel,IK> Exact_to_Inexact;
54 typedef CGAL::Cartesian_converter<IK,Kernel> Inexact_to_Exact;
60 typedef CGAL::Triangle_3<Kernel> Triangle_3;
61 typedef CGAL::Triangle_3<IK> InexactTriangle_3;
62 typedef CGAL::Plane_3<Kernel> Plane_3;
63 typedef std::vector<Triangle_3> Triangles;
64 typedef CGAL::Tetrahedron_3<Kernel> Tetrahedron;
65 typedef CGAL::Tetrahedron_3<IK> InexactTetrahedron;
66 typedef std::vector<Point_3> Points;
68 typedef CGAL::Bbox_3 Bbox;
69 typedef CGAL::Polyhedron_3<IK> InexactPolyhedron_3;
70 typedef CGAL::Polyhedron_3<Kernel> Polyhedron_3;
71 typedef CGAL::Triangulation_3<Kernel> ExactTriangulation;
72 typedef CGAL::Triangulation_3<IK> InexactTriangulation;
73 typedef CGAL::Aff_transformation_3<Kernel> Aff_transformation_3;
74 
75 typedef ExactTriangulation::Finite_facets_iterator ExactFinite_faces_iterator;
76 typedef ExactTriangulation::Finite_cells_iterator ExactFinite_cells_iterator;
77 typedef InexactTriangulation::Finite_facets_iterator InexactFinite_faces_iterator;
78 typedef InexactTriangulation::Finite_cells_iterator InexactFinite_cells_iterator;
89 typedef Triangles::iterator Triangle3_iterator;
90 
91 
92 typedef CGAL::Triangulation_2<Kernel> Triangulation_2;
93 
94 
96 
97 typedef CGAL::Triangle_2<Kernel> Triangle_2;
98 typedef std::vector<Triangle_2> Triangles_2;
99 typedef std::vector<Point_2> Points_2;
102 typedef Triangles_2::iterator Triangle2_iterator;
103 typedef Triangles::iterator Triangle3_iterator;
104 
105 
113 void triang_cellule(const Bbox& cel, Triangles& trianglesB){
114 
115  trianglesB.clear();
116 
117  Point_3 s1B(cel.xmin(), cel.ymin(), cel.zmin());
118  Point_3 r1B(cel.xmax(), cel.ymin(), cel.zmin());
119  Point_3 t1B(cel.xmax(), cel.ymax(), cel.zmin());
120  Point_3 v1B(cel.xmin(), cel.ymax(), cel.zmin());
121 
122  Point_3 s2B(cel.xmin(), cel.ymin(), cel.zmax());
123  Point_3 r2B(cel.xmax(), cel.ymin(), cel.zmax());
124  Point_3 t2B(cel.xmax(), cel.ymax(), cel.zmax());
125  Point_3 v2B(cel.xmin(), cel.ymax(), cel.zmax());
126 
127  //face1
128  Triangle_3 Tri1B(s1B,r1B,v1B);
129  Triangle_3 Tri2B(t1B,r1B,v1B);
130  trianglesB.push_back(Tri1B);
131  trianglesB.push_back(Tri2B);
132 
133 
134  //face2
135  Triangle_3 Tri5B(s2B,r2B,v2B);
136  Triangle_3 Tri6B(t2B,r2B,v2B);
137  trianglesB.push_back(Tri5B);
138  trianglesB.push_back(Tri6B);
139 
140  //face3
141  Triangle_3 Tri9B(s2B,s1B,v2B);
142  Triangle_3 Tri10B(v1B,s1B,v2B);
143  trianglesB.push_back(Tri9B);
144  trianglesB.push_back(Tri10B);
145 
146 
147  //face4
148  Triangle_3 Tri13B(r2B,r1B,t2B);
149  Triangle_3 Tri14B(t1B,r1B,t2B);
150  trianglesB.push_back(Tri13B);
151  trianglesB.push_back(Tri14B);
152 
153 
154  //face5
155  Triangle_3 Tri17B(v2B,v1B,t2B);
156  Triangle_3 Tri18B(t1B,v1B,t2B);
157  trianglesB.push_back(Tri17B);
158  trianglesB.push_back(Tri18B);
159 
160 
161  //face6
162  Triangle_3 Tri21B(s2B,s1B,r2B);
163  Triangle_3 Tri22B(r1B,s1B,r2B);
164  trianglesB.push_back(Tri21B);
165  trianglesB.push_back(Tri22B);
166 
167 
168 }
173 std::vector<Point_3> intersection_bis(const Segment_3& seg, const Triangle_3& t)
174 {
175  std::vector<Point_3> result;
176  Point_3 P;
177  Segment_3 s;
178  const CGAL::Object& intersec = CGAL::intersection(seg,t);
179  if(CGAL::assign(P,intersec)){
180  result.push_back(P);
181  }
182  else if(CGAL::assign(s,intersec)){
183  result.push_back(s.operator[](0));
184  result.push_back(s.operator[](1));
185  }
186 
187  return result;
188 }
189 
190 
196 std::vector<Point_3> intersection_bis(const Triangle_3& t1, const Triangle_3& t2)
197 {
198  std::vector<Point_3> result;
199  //Intersections between the segments of t1 and triangle t2
200  for(int k=0;k<3;k++){
201  int kp = (k+1)%3;
202  Point_3 P;
203  Segment_3 seg;
204  Segment_3 arete(t1.operator[](k),t1.operator[](kp));
205  const CGAL::Object& intersec = CGAL::intersection(arete,t2);
206  if(CGAL::assign(P,intersec)){
207  result.push_back(P);
208  }
209  else if(CGAL::assign(seg,intersec)){
210  result.push_back(seg.operator[](0));
211  result.push_back(seg.operator[](1));
212  }
213  }
214  //Intersections between the segments of t2 and triangle t1
215  for(int k=0;k<3;k++){
216  int kp = (k+1)%3;
217  Point_3 P;
218  Segment_3 seg;
219  Segment_3 arete(t2.operator[](k),t2.operator[](kp));
220  const CGAL::Object& intersec = CGAL::intersection(arete,t1);
221  if(CGAL::assign(P,intersec)){
222  result.push_back(P);
223  }
224  else if(CGAL::assign(seg,intersec)){
225  result.push_back(seg.operator[](0));
226  result.push_back(seg.operator[](1));
227  }
228  }
229 
230  return result;
231 }
232 
233 /* ! \brief Determine whether a point is inside a tetrahedron
234  \details Simple application of CGAL function <b> Tetrahedron_3::has_on_bounded_side <\b>
235  */
236 bool inside_tetra(const Tetrahedron &tetra, const Point_3& P){
237  return tetra.has_on_bounded_side(P);
238 }
239 
242 bool coplanar(std::vector<Point_3>::iterator begin, std::vector<Point_3>::iterator end)
243 {
244  bool test=true;
245 
246  bool test_confondus=true;
247  Point_3 P1,P2;
248  P1 = *begin;
249  std::vector<Point_3>::iterator it;
250  for(it=begin;it!=end && test_confondus;it++){
251  P2 = *it;
252  test_confondus = (P1==P2);
253  }
254  if(test_confondus){
255  return true;
256  } else {
257  bool test_collinear=true;
258  Point_3 P3;
259  std::vector<Point_3>::iterator iter;
260  for(iter=it;iter!=end && test_collinear;iter++){
261  P3 = *iter;
262  test_collinear = CGAL::collinear(P1,P2,P3);
263  }
264  if(test_collinear){
265  return true;
266  } else {
267  Point_3 P4;
268  for(std::vector<Point_3>::iterator iterat=iter;iterat!=end && test;iterat++){
269  P4 = *iterat;
270  test = CGAL::coplanar(P1,P2,P3,P4);
271  }
272  }
273  }
274 
275  return test;
276 }
277 
278 
279 #endif
Polyhedron_3::Plane_iterator Plane_iterator
Definition: intersections.hpp:82
InexactPolyhedron_3::Plane_iterator InexactPlane_iterator
Definition: intersections.hpp:87
CGAL::Triangulation_3< IK > InexactTriangulation
Definition: intersections.hpp:72
InexactPolyhedron_3::Halfedge_around_vertex_circulator InexactHalfedge_around_vertex_circulator
Definition: intersections.hpp:88
Polyhedron_3::Facet Facet
Definition: intersections.hpp:79
CGAL::Triangle_3< IK > InexactTriangle_3
Definition: intersections.hpp:61
CGAL::Triangle_2< Kernel > Triangle_2
Definition: intersections.hpp:97
InexactPolyhedron_3::Vertex_iterator InexactVertex_iterator
Definition: intersections.hpp:86
Kernel::Line_3 Line_3
Definition: intersections.hpp:58
CGAL::Exact_predicates_exact_constructions_kernel Kernel
Definition: intersections.hpp:51
Triangles_2::iterator Triangle2_iterator
Definition: intersections.hpp:102
CGAL::Bbox_2 Bbox_2
Definition: intersections.hpp:101
Polyhedron_3::Halfedge_around_vertex_circulator Halfedge_around_vertex_circulator
Definition: intersections.hpp:83
InexactPolyhedron_3::Facet InexactFacet
Definition: intersections.hpp:84
CGAL::Aff_transformation_3< Kernel > Aff_transformation_3
Definition: intersections.hpp:73
Kernel::Point_2 Point_2
Definition: intersections.hpp:95
InexactTriangulation::Finite_cells_iterator InexactFinite_cells_iterator
Definition: intersections.hpp:78
ExactTriangulation::Finite_cells_iterator ExactFinite_cells_iterator
Definition: intersections.hpp:76
InexactPolyhedron_3::Facet_iterator InexactFacet_iterator
Definition: intersections.hpp:85
std::vector< Triangle_3 > Triangles
Definition: intersections.hpp:63
CGAL::Triangle_3< Kernel > Triangle_3
Definition: intersections.hpp:60
CGAL::Polyhedron_3< Kernel > Polyhedron_3
Definition: intersections.hpp:70
Triangles::iterator Triangle3_iterator
Definition: intersections.hpp:89
bool coplanar(std::vector< Point_3 >::iterator begin, std::vector< Point_3 >::iterator end)
Test whether points are coplanar.
Definition: intersections.hpp:242
Kernel::Segment_2 Segment_2
Definition: intersections.hpp:100
Kernel::Point_3 Point_3
Definition: intersections.hpp:55
Polyhedron_3::Facet_iterator Facet_iterator
Definition: intersections.hpp:80
std::vector< Point_3 > intersection_bis(const Segment_3 &seg, const Triangle_3 &t)
Edge/triangle intersection.
Definition: intersections.hpp:173
CGAL::Tetrahedron_3< IK > InexactTetrahedron
Definition: intersections.hpp:65
std::vector< Point_2 > Points_2
Definition: intersections.hpp:99
CGAL::Cartesian_converter< IK, Kernel > Inexact_to_Exact
Definition: intersections.hpp:54
Kernel::Plane_3 Plane_3
Definition: intersections.hpp:59
CGAL::Bbox_3 Bbox
Definition: intersections.hpp:68
void triang_cellule(const Bbox &cel, Triangles &trianglesB)
Triangulation of the faces of a fluid cubic cell.
Definition: intersections.hpp:113
Kernel::Vector_3 Vector_3
Definition: intersections.hpp:57
Kernel::Segment_3 Segment_3
Definition: intersections.hpp:67
CGAL::Tetrahedron_3< Kernel > Tetrahedron
Definition: intersections.hpp:64
ExactTriangulation::Finite_facets_iterator ExactFinite_faces_iterator
Definition: intersections.hpp:75
Polyhedron_3::Vertex_iterator Vertex_iterator
Definition: intersections.hpp:81
double P(double x, double y, double z, double dx, double dy, double dz)
Definition: parametres.cpp:310
CGAL::Cartesian_converter< Kernel, IK > Exact_to_Inexact
Definition: intersections.hpp:53
InexactTriangulation::Finite_facets_iterator InexactFinite_faces_iterator
Definition: intersections.hpp:77
IK::Point_3 InexactPoint_3
Definition: intersections.hpp:56
bool inside_tetra(const Tetrahedron &tetra, const Point_3 &P)
Definition: intersections.hpp:236
CGAL::Triangulation_2< Kernel > Triangulation_2
Definition: intersections.hpp:92
CGAL::Exact_predicates_inexact_constructions_kernel IK
Definition: intersections.hpp:52
std::vector< Triangle_2 > Triangles_2
Definition: intersections.hpp:98
CGAL::Polyhedron_3< IK > InexactPolyhedron_3
Definition: intersections.hpp:69
CGAL::Triangulation_3< Kernel > ExactTriangulation
Definition: intersections.hpp:71
std::vector< Point_3 > Points
Definition: intersections.hpp:66