Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
predicates/kernel_ftC2.h
Go to the documentation of this file.
1 // ======================================================================
2 //
3 // Copyright (c) 2000 The CGAL Consortium
4 
5 // This software and related documentation is part of the Computational
6 // Geometry Algorithms Library (CGAL).
7 // This software and documentation is provided "as-is" and without warranty
8 // of any kind. In no event shall the CGAL Consortium be liable for any
9 // damage of any kind.
10 //
11 // Every use of CGAL requires a license.
12 //
13 // Academic research and teaching license
14 // - For academic research and teaching purposes, permission to use and copy
15 // the software and its documentation is hereby granted free of charge,
16 // provided that it is not a component of a commercial product, and this
17 // notice appears in all copies of the software and related documentation.
18 //
19 // Commercial licenses
20 // - A commercial license is available through Algorithmic Solutions, who also
21 // markets LEDA (http://www.algorithmic-solutions.de).
22 // - Commercial users may apply for an evaluation license by writing to
23 // Algorithmic Solutions (contact@algorithmic-solutions.com).
24 //
25 // The CGAL Consortium consists of Utrecht University (The Netherlands),
26 // ETH Zurich (Switzerland), Free University of Berlin (Germany),
27 // INRIA Sophia-Antipolis (France), Martin-Luther-University Halle-Wittenberg
28 // (Germany), Max-Planck-Institute Saarbrucken (Germany), RISC Linz (Austria),
29 // and Tel-Aviv University (Israel).
30 //
31 // ----------------------------------------------------------------------
32 //
33 // release : CGAL-2.2
34 // release_date : 2000, September 30
35 //
36 // file : include/CGAL/predicates/kernel_ftC2.h
37 // package : C2 (4.4)
38 // revision : $Revision: 1.1.1.1 $
39 // revision_date : $Date: 2001/07/05 22:17:49 $
40 // author(s) : Herve Bronnimann
41 // coordinator : INRIA Sophia-Antipolis
42 //
43 // email : contact@cgal.org
44 // www : http://www.cgal.org
45 //
46 // ======================================================================
47 
48 #ifndef CGAL_PREDICATES_KERNEL_FTC2_H
49 #define CGAL_PREDICATES_KERNEL_FTC2_H
50 
51 #include <CGAL/number_utils.h>
54 
56 
57 template < class FT >
59 bool
60 equal_lineC2(const FT &l1a, const FT &l1b, const FT &l1c,
61  const FT &l2a, const FT &l2b, const FT &l2c)
62 {
63  if (sign_of_determinant2x2(l1a, l1b, l2a, l2b) != ZERO)
64  return false; // Not parallel.
65  CGAL::Sign s1a = CGAL_NTS sign(l1a);
66  if (s1a != ZERO)
67  return s1a == CGAL_NTS sign(l2a)
68  && sign_of_determinant2x2(l1a, l1c, l2a, l2c) == ZERO;
69  return CGAL_NTS sign(l1b) == CGAL_NTS sign(l2b)
70  && sign_of_determinant2x2(l1b, l1c, l2b, l2c) == ZERO;
71 }
72 
73 template < class FT >
76 compare_xC2(const FT &px,
77  const FT &la, const FT &lb, const FT &lc,
78  const FT &ha, const FT &hb, const FT &hc)
79 {
80  // The abscissa of the intersection point is num/den.
81  FT num = det2x2_by_formula( lb, lc, hb, hc);
82  FT den = det2x2_by_formula( la, lb, ha, hb);
83  Sign s = CGAL_NTS sign(den);
85  return Comparison_result( s * CGAL_NTS compare( px * den, num) );
86 }
87 
88 template < class FT >
91 compare_xC2(const FT &la, const FT &lb, const FT &lc,
92  const FT &h1a, const FT &h1b, const FT &h1c,
93  const FT &h2a, const FT &h2b, const FT &h2c)
94 {
95  /*
96  FT num1 = det2x2_by_formula( lb, lc, h1b, h1c);
97  FT den1 = det2x2_by_formula( la, lb, h1a, h1b);
98  FT num2 = det2x2_by_formula( lb, lc, h2b, h2c);
99  FT den2 = det2x2_by_formula( la, lb, h2a, h2b);
100  Sign s = Sign (CGAL_NTS sign(den1) * CGAL_NTS sign(den2));
101  CGAL_kernel_assertion( s != ZERO );
102  return Comparison_result( s * sign_of_determinant2x2(num1,
103  num2, den1, den2));
104  */
105  FT num1 = det2x2_by_formula( la, lc, h1a, h1c);
106  FT num2 = det2x2_by_formula( la, lc, h2a, h2c);
107  FT num = det2x2_by_formula(h1a,h1c,h2a,h2c)*lb
108  + det2x2_by_formula(num1,num2,h1b,h2b);
109  FT den1 = det2x2_by_formula( la, lb, h1a, h1b);
110  FT den2 = det2x2_by_formula( la, lb, h2a, h2b);
111  return Comparison_result( CGAL_NTS sign(lb) * CGAL_NTS sign(num) *
112  CGAL_NTS sign(den1) * CGAL_NTS sign(den2));
113 }
114 
115 template < class FT >
118 compare_xC2(const FT &l1a, const FT &l1b, const FT &l1c,
119  const FT &h1a, const FT &h1b, const FT &h1c,
120  const FT &l2a, const FT &l2b, const FT &l2c,
121  const FT &h2a, const FT &h2b, const FT &h2c)
122 {
123  FT num1 = det2x2_by_formula( l1b, l1c, h1b, h1c);
124  FT den1 = det2x2_by_formula( l1a, l1b, h1a, h1b);
125  FT num2 = det2x2_by_formula( l2b, l2c, h2b, h2c);
126  FT den2 = det2x2_by_formula( l2a, l2b, h2a, h2b);
127  Sign s = Sign (CGAL_NTS sign(den1) * CGAL_NTS sign(den2));
128  CGAL_kernel_assertion( s != ZERO );
129  return Comparison_result( s * sign_of_determinant2x2(num1, num2,
130  den1, den2));
131 }
132 
133 template < class FT >
136 compare_y_at_xC2(const FT &px, const FT &py,
137  const FT &la, const FT &lb, const FT &lc)
138 {
139  Sign s = CGAL_NTS sign(lb);
140  CGAL_kernel_assertion( s != ZERO );
141  return Comparison_result (s * CGAL_NTS sign(la*px + lb*py + lc));
142 }
143 
144 template < class FT >
147 compare_y_at_xC2(const FT &px,
148  const FT &l1a, const FT &l1b, const FT &l1c,
149  const FT &l2a, const FT &l2b, const FT &l2c)
150 {
151  Sign s = Sign (CGAL_NTS sign(l1b) * CGAL_NTS sign(l2b));
152  CGAL_kernel_assertion( s != ZERO );
153  return Comparison_result ( s * sign_of_determinant2x2(l2a*px+l2c, l2b,
154  l1a*px+l1c, l1b));
155 }
156 
157 template < class FT >
160 compare_y_at_xC2(const FT &l1a, const FT &l1b, const FT &l1c,
161  const FT &l2a, const FT &l2b, const FT &l2c,
162  const FT &ha, const FT &hb, const FT &hc)
163 {
164  Sign s = Sign (sign_of_determinant2x2(l1a, l1b, l2a, l2b) *
165  CGAL_NTS sign(hb));
166  CGAL_kernel_assertion( s != ZERO );
167  return Comparison_result( s * sign_of_determinant3x3(l1a, l1b, l1c,
168  l2a, l2b, l2c,
169  ha, hb, hc));
170 }
171 
172 template < class FT >
175 compare_y_at_xC2(const FT &l1a, const FT &l1b, const FT &l1c,
176  const FT &l2a, const FT &l2b, const FT &l2c,
177  const FT &h1a, const FT &h1b, const FT &h1c,
178  const FT &h2a, const FT &h2b, const FT &h2c)
179 {
180  // The abscissa of the intersection point is num/den.
181  FT num = det2x2_by_formula( l1b, l1c, l2b, l2c);
182  FT den = det2x2_by_formula( l1a, l1b, l2a, l2b);
183  Sign s = Sign (CGAL_NTS sign(h1b) * CGAL_NTS sign(h2b) * CGAL_NTS sign(den));
184  CGAL_kernel_assertion( s != ZERO );
185  return Comparison_result ( s * sign_of_determinant2x2(h2a*num+h2c*den, h2b,
186  h1a*num+h1c*den, h1b));
187 }
188 
189 template < class FT >
191 bool
192 equal_directionC2(const FT &dx1, const FT &dy1,
193  const FT &dx2, const FT &dy2)
194 {
195  return CGAL_NTS sign(dx1) == CGAL_NTS sign(dx2)
196  && CGAL_NTS sign(dy1) == CGAL_NTS sign(dy2)
197  && sign_of_determinant2x2(dx1, dy1, dx2, dy2) == ZERO;
198 }
199 
200 template < class FT >
201 /*CGAL_NO_FILTER*/
204 compare_angle_with_x_axisC2(const FT &dx1, const FT &dy1,
205  const FT &dx2, const FT &dy2)
206 {
207  // angles are in [-pi,pi], and the angle between Ox and d1 is compared
208  // with the angle between Ox and d2
209  int quadrant_1 = (dx1 >= FT(0)) ? ((dy1 >= FT(0))?1:4)
210  : ((dy1 >= FT(0))?2:3);
211  int quadrant_2 = (dx2 >= FT(0)) ? ((dy2 >= FT(0))?1:4)
212  : ((dy2 >= FT(0))?2:3);
213  // We can't use CGAL_NTS compare(quadrant_1,quadrant_2) because in case
214  // of tie, we need additional computation
215  if (quadrant_1 > quadrant_2)
216  return LARGER;
217  else if (quadrant_1 < quadrant_2)
218  return SMALLER;
219  return Comparison_result(-sign_of_determinant2x2(dx1,dy1,dx2,dy2));
220 }
221 
222 template < class FT >
223 inline
225 compare_deltax_deltayC2(const FT &px, const FT &qx,
226  const FT &ry, const FT &sy)
227 {
228  return CGAL_NTS compare(CGAL_NTS abs(px-qx), CGAL_NTS abs(ry-sy));
229 }
230 
231 template < class FT >
232 /*CGAL_NO_FILTER*/
233 inline
235 compare_lexicographically_xyC2(const FT &px, const FT &py,
236  const FT &qx, const FT &qy)
237 {
239  return (c != EQUAL) ? c : CGAL_NTS compare(py,qy);
240 }
241 
242 template < class FT >
243 inline
245 orientationC2(const FT &px, const FT &py,
246  const FT &qx, const FT &qy,
247  const FT &rx, const FT &ry)
248 {
249  return Orientation (sign_of_determinant2x2(px-rx, py-ry,
250  qx-rx, qy-ry));
251 }
252 
253 template < class FT >
254 /*CGAL_NO_FILTER*/
256 bool
257 collinear_are_ordered_along_lineC2(const FT &px, const FT &py,
258  const FT &qx, const FT &qy,
259  const FT &rx, const FT &ry)
260 {
261  if (px < qx) return !(rx < qx);
262  if (qx < px) return !(qx < rx);
263  if (py < qy) return !(ry < qy);
264  if (qy < py) return !(qy < ry);
265  return true; // p==q
266 }
267 
268 template < class FT >
269 /*CGAL_NO_FILTER*/
271 bool
273  const FT &qx, const FT &qy,
274  const FT &rx, const FT &ry)
275 {
276  if (px < qx) return (qx < rx);
277  if (qx < px) return (rx < qx);
278  if (py < qy) return (qy < ry);
279  if (qy < py) return (ry < qy);
280  return false;
281 }
282 
283 template < class FT >
286 side_of_oriented_circleC2(const FT &px, const FT &py,
287  const FT &qx, const FT &qy,
288  const FT &rx, const FT &ry,
289  const FT &tx, const FT &ty)
290 {
291  // Oriented_side(
292  // sign_of_determinant4x4(px, py, px*px + py*py, 1,
293  // qx, qy, qx*qx + qy*qy, 1,
294  // rx, ry, rx*rx + ry*ry, 1,
295  // tx, ty, tx*tx + ty*ty, 1));
296  // We first translate so that t is the new origin.
297  FT ptx = px-tx;
298  FT pty = py-ty;
299  FT qtx = qx-tx;
300  FT qty = qy-ty;
301  FT rtx = rx-tx;
302  FT rty = ry-ty;
303 // The usual 3x3 formula can be simplified a little bit to a 2x2.
304 // sign_of_determinant3x3(ptx, pty, square(ptx) + square(pty),
305 // qtx, qty, square(qtx) + square(qty),
306 // rtx, rty, square(rtx) + square(rty)));
308  ptx*qty - pty*qtx, qtx*(qx-px) + qty*(qy-py),
309  ptx*rty - pty*rtx, rtx*(rx-px) + rty*(ry-py)));
310 }
311 
312 template < class FT >
315 side_of_bounded_circleC2(const FT &px, const FT &py,
316  const FT &qx, const FT &qy,
317  const FT &rx, const FT &ry,
318  const FT &tx, const FT &ty)
319 {
320  // Note: if the code of these is inlined, and if they are implemented
321  // in a good way, some CSE can be done by the compiler.
322  return Bounded_side( side_of_oriented_circleC2(px,py,qx,qy,rx,ry,tx,ty)
323  * orientationC2(px,py,qx,qy,rx,ry) );
324 }
325 
326 template < class FT >
327 inline
329 cmp_dist_to_pointC2(const FT &px, const FT &py,
330  const FT &qx, const FT &qy,
331  const FT &rx, const FT &ry)
332 {
333  return CGAL_NTS compare(squared_distanceC2(px,py,qx,qy),
334  squared_distanceC2(px,py,rx,ry));
335 }
336 
337 template < class FT >
338 /*CGAL_NO_FILTER*/
339 inline
340 bool
341 has_larger_dist_to_pointC2(const FT &px, const FT &py,
342  const FT &qx, const FT &qy,
343  const FT &rx, const FT &ry)
344 {
345  return cmp_dist_to_pointC2(px,py,qx,qy,rx,ry) == LARGER;
346 }
347 
348 template < class FT >
349 /*CGAL_NO_FILTER*/
350 inline
351 bool
352 has_smaller_dist_to_pointC2(const FT &px, const FT &py,
353  const FT &qx, const FT &qy,
354  const FT &rx, const FT &ry)
355 {
356  return cmp_dist_to_pointC2(px,py,qx,qy,rx,ry) == SMALLER;
357 }
358 
359 template < class FT >
360 inline
362 cmp_signed_dist_to_directionC2(const FT &la, const FT &lb,
363  const FT &px, const FT &py,
364  const FT &qx, const FT &qy)
365 {
366  return CGAL_NTS compare(scaled_distance_to_directionC2(la,lb,px,py),
367  scaled_distance_to_directionC2(la,lb,qx,qy));
368 }
369 
370 template < class FT >
371 /*CGAL_NO_FILTER*/
372 inline
373 bool
374 has_larger_signed_dist_to_directionC2(const FT &la, const FT &lb,
375  const FT &px, const FT &py,
376  const FT &qx, const FT &qy)
377 {
378  return cmp_signed_dist_to_directionC2(la,lb,px,py,qx,qy) == LARGER;
379 }
380 
381 template < class FT >
382 /*CGAL_NO_FILTER*/
383 inline
384 bool
385 has_smaller_signed_dist_to_directionC2(const FT &la, const FT &lb,
386  const FT &px, const FT &py,
387  const FT &qx, const FT &qy)
388 {
389  return cmp_signed_dist_to_directionC2(la,lb,px,py,qx,qy) == SMALLER;
390 }
391 
392 template <class FT>
393 inline
395 cmp_signed_dist_to_lineC2(const FT &px, const FT &py,
396  const FT &qx, const FT &qy,
397  const FT &rx, const FT &ry,
398  const FT &sx, const FT &sy)
399 {
400  return CGAL_NTS compare(scaled_distance_to_lineC2(px,py,qx,qy,rx,ry),
401  scaled_distance_to_lineC2(px,py,qx,qy,sx,sy));
402 }
403 
404 template <class FT>
405 /*CGAL_NO_FILTER*/
406 inline
407 bool
408 has_larger_signed_dist_to_lineC2(const FT &px, const FT &py,
409  const FT &qx, const FT &qy,
410  const FT &rx, const FT &ry,
411  const FT &sx, const FT &sy)
412 {
413  return cmp_signed_dist_to_lineC2(px,py,qx,qy,rx,ry,sx,sy) == LARGER;
414 }
415 
416 template <class FT>
417 /*CGAL_NO_FILTER*/
418 inline
419 bool
420 has_smaller_signed_dist_to_lineC2(const FT &px, const FT &py,
421  const FT &qx, const FT &qy,
422  const FT &rx, const FT &ry,
423  const FT &sx, const FT &sy)
424 {
425  return cmp_signed_dist_to_lineC2(px,py,qx,qy,rx,ry,sx,sy) == SMALLER;
426 }
427 
428 template <class FT>
429 inline
431 side_of_oriented_lineC2(const FT &a, const FT &b, const FT &c,
432  const FT &x, const FT &y)
433 {
434  return Oriented_side(CGAL_NTS sign(a*x+b*y+c));
435 }
436 
438 
439 #ifdef CGAL_ARITHMETIC_FILTER_H
440 #include <CGAL/Arithmetic_filter/predicates/kernel_ftC2.h>
441 #endif
442 
443 #endif // CGAL_PREDICATES_KERNEL_FTC2_H
CGAL_BEGIN_NAMESPACE FT det2x2_by_formula(const FT &a00, const FT &a01, const FT &a10, const FT &a11)
Definition: determinant.h:59
static SURF_BEGIN_NAMESPACE double sign(double x)
Sign sign_of_determinant3x3(const FT &a00, const FT &a01, const FT &a02, const FT &a10, const FT &a11, const FT &a12, const FT &a20, const FT &a21, const FT &a22)
CGAL_KERNEL_MEDIUM_INLINE bool collinear_are_strictly_ordered_along_lineC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry)
Sign Orientation
Definition: enum.h:64
void int int REAL REAL * y
Definition: read.cpp:74
bool has_larger_dist_to_pointC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry)
CGAL_BEGIN_NAMESPACE Sign sign_of_determinant2x2(const FT &a00, const FT &a01, const FT &a10, const FT &a11)
double s
Definition: blastest.C:80
Orientation orientationC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry)
CGAL_KERNEL_INLINE FT squared_distanceC2(const FT &px, const FT &py, const FT &qx, const FT &qy)
Comparison_result compare_deltax_deltayC2(const FT &px, const FT &qx, const FT &ry, const FT &sy)
NT & den
bool has_smaller_dist_to_pointC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry)
Definition: enum.h:96
#define CGAL_KERNEL_MEDIUM_INLINE
Definition: kernel_basic.h:55
Comparison_result cmp_signed_dist_to_directionC2(const FT &la, const FT &lb, const FT &px, const FT &py, const FT &qx, const FT &qy)
Sign
Definition: enum.h:57
CGAL_KERNEL_MEDIUM_INLINE Comparison_result compare_angle_with_x_axisC2(const FT &dx1, const FT &dy1, const FT &dx2, const FT &dy2)
Oriented_side
Definition: enum.h:78
bool has_larger_signed_dist_to_lineC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry, const FT &sx, const FT &sy)
CGAL_BEGIN_NAMESPACE CGAL_KERNEL_MEDIUM_INLINE bool equal_lineC2(const FT &l1a, const FT &l1b, const FT &l1c, const FT &l2a, const FT &l2b, const FT &l2c)
Comparison_result compare_lexicographically_xyC2(const FT &px, const FT &py, const FT &qx, const FT &qy)
CGAL_KERNEL_MEDIUM_INLINE Comparison_result compare_y_at_xC2(const FT &px, const FT &py, const FT &la, const FT &lb, const FT &lc)
Oriented_side side_of_oriented_lineC2(const FT &a, const FT &b, const FT &c, const FT &x, const FT &y)
#define CGAL_KERNEL_LARGE_INLINE
Definition: kernel_basic.h:56
CGAL_KERNEL_LARGE_INLINE Bounded_side side_of_bounded_circleC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry, const FT &tx, const FT &ty)
void int int REAL * x
Definition: read.cpp:74
Comparison_result cmp_dist_to_pointC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry)
Comparison_result cmp_signed_dist_to_lineC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry, const FT &sx, const FT &sy)
bool has_larger_signed_dist_to_directionC2(const FT &la, const FT &lb, const FT &px, const FT &py, const FT &qx, const FT &qy)
CGAL_KERNEL_INLINE FT scaled_distance_to_directionC2(const FT &la, const FT &lb, const FT &px, const FT &py)
CGAL_KERNEL_LARGE_INLINE Oriented_side side_of_oriented_circleC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry, const FT &tx, const FT &ty)
CGAL_KERNEL_MEDIUM_INLINE bool equal_directionC2(const FT &dx1, const FT &dy1, const FT &dx2, const FT &dy2)
Comparison_result
Definition: enum.h:94
Bounded_side
Definition: enum.h:86
Definition: enum.h:98
NT abs(const NT &x)
Definition: number_utils.h:130
CGAL_KERNEL_MEDIUM_INLINE Comparison_result compare_xC2(const FT &px, const FT &la, const FT &lb, const FT &lc, const FT &ha, const FT &hb, const FT &hc)
bool has_smaller_signed_dist_to_directionC2(const FT &la, const FT &lb, const FT &px, const FT &py, const FT &qx, const FT &qy)
Definition: enum.h:97
#define CGAL_BEGIN_NAMESPACE
Definition: kdtree_d.h:86
Definition: enum.h:60
CGAL_KERNEL_INLINE Comparison_result compare(const NT &n1, const NT &n2)
Definition: number_utils.h:143
#define CGAL_NTS
bool has_smaller_signed_dist_to_lineC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry, const FT &sx, const FT &sy)
#define CGAL_kernel_assertion(EX)
#define CGAL_END_NAMESPACE
Definition: kdtree_d.h:87
CGAL_KERNEL_MEDIUM_INLINE bool collinear_are_ordered_along_lineC2(const FT &px, const FT &py, const FT &qx, const FT &qy, const FT &rx, const FT &ry)
CGAL_KERNEL_INLINE FT scaled_distance_to_lineC2(const FT &la, const FT &lb, const FT &lc, const FT &px, const FT &py)