DOLFIN
DOLFIN C++ interface
CollisionPredicates.h
1 // Copyright (C) 2014-2016 Anders Logg, August Johansson and Benjamin Kehlet
2 //
3 // This file is part of DOLFIN.
4 //
5 // DOLFIN is free software: you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published by
7 // the Free Software Foundation, either version 3 of the License, or
8 // (at your option) any later version.
9 //
10 // DOLFIN is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17 //
18 // First added: 2014-02-03
19 // Last changed: 2017-09-29
20 
21 #ifndef __COLLISION_PREDICATES_H
22 #define __COLLISION_PREDICATES_H
23 
24 namespace dolfin
25 {
26 
27  // Forward declarations
28  class Point;
29  class MeshEntity;
30 
33 
35  {
36  public:
37 
38  //--- High-level collision detection predicates ---
39 
51  static bool collides(const MeshEntity& entity,
52  const Point& point);
53 
65  static bool collides(const MeshEntity& entity_0,
66  const MeshEntity& entity_1);
67 
68  //--- Low-level collision detection predicates ---
69 
71  static bool collides_segment_point(const Point& p0,
72  const Point& p1,
73  const Point& point,
74  std::size_t gdim);
75 
77  static bool collides_segment_point_1d(double p0,
78  double p1,
79  double point);
80 
82  static bool collides_segment_point_2d(const Point& p0,
83  const Point& p1,
84  const Point& point);
85 
87  static bool collides_segment_point_3d(const Point& p0,
88  const Point& p1,
89  const Point& point);
90 
92  static bool collides_segment_segment(const Point& p0,
93  const Point& p1,
94  const Point& q0,
95  const Point& q1,
96  std::size_t gdim);
97 
99  static bool collides_segment_segment_1d(double p0,
100  double p1,
101  double q0,
102  double q1);
103 
105  static bool collides_segment_segment_2d(const Point& p0,
106  const Point& p1,
107  const Point& q0,
108  const Point& q1);
109 
111  static bool collides_segment_segment_3d(const Point& p0,
112  const Point& p1,
113  const Point& q0,
114  const Point& q1);
115 
117  static bool collides_triangle_point(const Point& p0,
118  const Point& p1,
119  const Point& p2,
120  const Point& point,
121  std::size_t gdim);
122 
124  static bool collides_triangle_point_2d(const Point& p0,
125  const Point& p1,
126  const Point& p2,
127  const Point& point);
128 
130  static bool collides_triangle_point_3d(const Point& p0,
131  const Point& p1,
132  const Point& p2,
133  const Point& point);
134 
136  static bool collides_triangle_segment(const Point& p0,
137  const Point& p1,
138  const Point& p2,
139  const Point& q0,
140  const Point& q1,
141  std::size_t gdim);
142 
144  static bool collides_triangle_segment_2d(const Point& p0,
145  const Point& p1,
146  const Point& p2,
147  const Point& q0,
148  const Point& q1);
149 
151  static bool collides_triangle_segment_3d(const Point& p0,
152  const Point& p1,
153  const Point& p2,
154  const Point& q0,
155  const Point& q1);
156 
158  static bool collides_triangle_triangle(const Point& p0,
159  const Point& p1,
160  const Point& p2,
161  const Point& q0,
162  const Point& q1,
163  const Point& q2,
164  std::size_t gdim);
165 
167  static bool collides_triangle_triangle_2d(const Point& p0,
168  const Point& p1,
169  const Point& p2,
170  const Point& q0,
171  const Point& q1,
172  const Point& q2);
173 
175  static bool collides_triangle_triangle_3d(const Point& p0,
176  const Point& p1,
177  const Point& p2,
178  const Point& q0,
179  const Point& q1,
180  const Point& q2);
181 
183  static bool collides_tetrahedron_point_3d(const Point& p0,
184  const Point& p1,
185  const Point& p2,
186  const Point& p3,
187  const Point& point);
188 
190  static bool collides_tetrahedron_segment_3d(const Point& p0,
191  const Point& p1,
192  const Point& p2,
193  const Point& p3,
194  const Point& q0,
195  const Point& q1);
196 
198  static bool collides_tetrahedron_triangle_3d(const Point& p0,
199  const Point& p1,
200  const Point& p2,
201  const Point& p3,
202  const Point& q0,
203  const Point& q1,
204  const Point& q2);
205 
207  static bool collides_tetrahedron_tetrahedron_3d(const Point& p0,
208  const Point& p1,
209  const Point& p2,
210  const Point& p3,
211  const Point& q0,
212  const Point& q1,
213  const Point& q2,
214  const Point& q3);
215 
216  private:
217 
218  // Implementation of collision detection predicates
219 
220  static bool _collides_segment_point_1d(double p0,
221  double p1,
222  double point);
223 
224  static bool _collides_segment_point_2d(const Point& p0,
225  const Point& p1,
226  const Point& point);
227 
228  static bool _collides_segment_point_3d(const Point& p0,
229  const Point& p1,
230  const Point& point);
231 
232  static bool _collides_segment_segment_1d(double p0,
233  double p1,
234  double q0,
235  double q1);
236 
237  static bool _collides_segment_segment_2d(const Point& p0,
238  const Point& p1,
239  const Point& q0,
240  const Point& q1);
241 
242  static bool _collides_segment_segment_3d(const Point& p0,
243  const Point& p1,
244  const Point& q0,
245  const Point& q1);
246 
247  static bool _collides_triangle_point_2d(const Point& p0,
248  const Point& p1,
249  const Point& p2,
250  const Point& point);
251 
252  static bool _collides_triangle_point_3d(const Point& p0,
253  const Point& p1,
254  const Point& p2,
255  const Point& point);
256 
257  static bool _collides_triangle_segment_2d(const Point& p0,
258  const Point& p1,
259  const Point& p2,
260  const Point& q0,
261  const Point& q1);
262 
263  static bool _collides_triangle_segment_3d(const Point& p0,
264  const Point& p1,
265  const Point& p2,
266  const Point& q0,
267  const Point& q1);
268 
269  static bool _collides_triangle_triangle_2d(const Point& p0,
270  const Point& p1,
271  const Point& p2,
272  const Point& q0,
273  const Point& q1,
274  const Point& q2);
275 
276  static bool _collides_triangle_triangle_3d(const Point& p0,
277  const Point& p1,
278  const Point& p2,
279  const Point& q0,
280  const Point& q1,
281  const Point& q2);
282 
283  static bool _collides_tetrahedron_point_3d(const Point& p0,
284  const Point& p1,
285  const Point& p2,
286  const Point& p3,
287  const Point& point);
288 
289  static bool _collides_tetrahedron_segment_3d(const Point& p0,
290  const Point& p1,
291  const Point& p2,
292  const Point& p3,
293  const Point& q0,
294  const Point& q1);
295 
296  static bool _collides_tetrahedron_triangle_3d(const Point& p0,
297  const Point& p1,
298  const Point& p2,
299  const Point& p3,
300  const Point& q0,
301  const Point& q1,
302  const Point& q2);
303 
304  static bool _collides_tetrahedron_tetrahedron_3d(const Point& p0,
305  const Point& p1,
306  const Point& p2,
307  const Point& p3,
308  const Point& q0,
309  const Point& q1,
310  const Point& q2,
311  const Point& q3);
312  };
313 
314 }
315 
316 #endif
static bool collides_triangle_point_3d(const Point &p0, const Point &p1, const Point &p2, const Point &point)
Check whether triangle p0-p1-p2 collides with point (3D version)
Definition: CollisionPredicates.cpp:363
static bool collides_triangle_triangle_3d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2 (3D version)
Definition: CollisionPredicates.cpp:403
static bool collides_segment_segment_1d(double p0, double p1, double q0, double q1)
Check whether segment p0-p1 collides with segment q0-q1 (1D version)
Definition: CollisionPredicates.cpp:328
Definition: CollisionPredicates.h:34
static bool collides_tetrahedron_tetrahedron_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2, const Point &q3)
Check whether tetrahedron p0-p1-p2-p3 collides with tetrahedron q0-q1-q2.
Definition: CollisionPredicates.cpp:447
Definition: adapt.h:29
static bool collides_segment_segment_2d(const Point &p0, const Point &p1, const Point &q0, const Point &q1)
Check whether segment p0-p1 collides with segment q0-q1 (2D version)
Definition: CollisionPredicates.cpp:336
static bool collides_segment_point(const Point &p0, const Point &p1, const Point &point, std::size_t gdim)
Check whether segment p0-p1 collides with point.
Definition: CollisionPredicates.cpp:194
Definition: Point.h:40
static bool collides(const MeshEntity &entity, const Point &point)
Definition: CollisionPredicates.cpp:35
static bool collides_triangle_segment(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, std::size_t gdim)
Check whether triangle p0-p1-p2 collides with segment q0-q1.
Definition: CollisionPredicates.cpp:257
static bool collides_tetrahedron_point_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &point)
Check whether tetrahedron p0-p1-p2-p3 collides with point.
Definition: CollisionPredicates.cpp:414
static bool collides_triangle_point_2d(const Point &p0, const Point &p1, const Point &p2, const Point &point)
Check whether triangle p0-p1-p2 collides with point (2D version)
Definition: CollisionPredicates.cpp:354
static bool collides_tetrahedron_triangle_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)
Check whether tetrahedron p0-p1-p2-p3 collides with triangle q0-q1-q2.
Definition: CollisionPredicates.cpp:435
static bool collides_segment_point_1d(double p0, double p1, double point)
Check whether segment p0-p1 collides with point (1D version)
Definition: CollisionPredicates.cpp:304
static bool collides_triangle_segment_2d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
Check whether triangle p0-p1-p2 collides with segment q0-q1 (2D version)
Definition: CollisionPredicates.cpp:372
static bool collides_segment_point_2d(const Point &p0, const Point &p1, const Point &point)
Check whether segment p0-p1 collides with point (2D version)
Definition: CollisionPredicates.cpp:312
static bool collides_segment_point_3d(const Point &p0, const Point &p1, const Point &point)
Check whether segment p0-p1 collides with point (3D version)
Definition: CollisionPredicates.cpp:320
static bool collides_triangle_triangle_2d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2 (2D version)
Definition: CollisionPredicates.cpp:392
static bool collides_triangle_triangle(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2, std::size_t gdim)
Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2.
Definition: CollisionPredicates.cpp:278
static bool collides_segment_segment_3d(const Point &p0, const Point &p1, const Point &q0, const Point &q1)
Check whether segment p0-p1 collides with segment q0-q1 (3D version)
Definition: CollisionPredicates.cpp:345
static bool collides_triangle_point(const Point &p0, const Point &p1, const Point &p2, const Point &point, std::size_t gdim)
Check whether triangle p0-p1-p2 collides with point.
Definition: CollisionPredicates.cpp:237
Definition: MeshEntity.h:42
static bool collides_triangle_segment_3d(const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
Check whether triangle p0-p1-p2 collides with segment q0-q1 (3D version)
Definition: CollisionPredicates.cpp:382
static bool collides_segment_segment(const Point &p0, const Point &p1, const Point &q0, const Point &q1, std::size_t gdim)
Check whether segment p0-p1 collides with segment q0-q1.
Definition: CollisionPredicates.cpp:215
static bool collides_tetrahedron_segment_3d(const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1)
Check whether tetrahedron p0-p1-p2-p3 collides with segment q0-q1.
Definition: CollisionPredicates.cpp:424