1 // Copyright (c) 2015  Università della Svizzera italiana.
2 // All rights reserved.
3 //
4 // This file is part of CGAL (www.cgal.org).
5 //
6 // $URL: https://github.com/CGAL/cgal/blob/v5.3/Segment_Delaunay_graph_Linf_2/include/CGAL/Segment_Delaunay_graph_Linf_2/Basic_predicates_C2.h $
7 // $Id: Basic_predicates_C2.h 0779373 2020-03-26T13:31:46+01:00 Sébastien Loriot
8 // SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
9 //
10 //
11 // Author(s)     : Panagiotis Cheilaris, Sandeep Kumar Dey, Evanthia Papadopoulou
12 //philaris@gmail.com, sandeep.kr.dey@gmail.com, evanthia.papadopoulou@usi.ch
13 
14 #ifndef CGAL_SEGMENT_DELAUNAY_GRAPH_LINF_2_BASIC_PREDICATES_C2_H
15 #define CGAL_SEGMENT_DELAUNAY_GRAPH_LINF_2_BASIC_PREDICATES_C2_H
16 
17 #include <CGAL/license/Segment_Delaunay_graph_Linf_2.h>
18 
19 
20 
21 #include <CGAL/Segment_Delaunay_graph_Linf_2/basic.h>
22 #include <CGAL/Segment_Delaunay_graph_2/Are_same_points_C2.h>
23 #include <CGAL/Segment_Delaunay_graph_2/Are_same_segments_C2.h>
24 #include <CGAL/Segment_Delaunay_graph_2/Compare_x_2.h>
25 #include <CGAL/Segment_Delaunay_graph_2/Compare_y_2.h>
26 
27 #include <CGAL/Polychain_2.h>
28 
29 #include <CGAL/Segment_Delaunay_graph_2/Basic_predicates_C2.h>
30 #include <CGAL/Segment_Delaunay_graph_Linf_2/Bisector_Linf.h>
31 
32 
33 namespace CGAL {
34 
35 namespace SegmentDelaunayGraphLinf_2 {
36 
37 template<class K>
38 struct Basic_predicates_C2
39  : public CGAL::SegmentDelaunayGraph_2::Basic_predicates_C2<K>
40 {
41 public:
42   //-------------------------------------------------------------------
43   // TYPES
44   //-------------------------------------------------------------------
45 
46   typedef CGAL::SegmentDelaunayGraph_2::Basic_predicates_C2<K> Base;
47 
48   typedef typename K::RT                  RT;
49   typedef typename K::FT                  FT;
50   typedef typename K::Point_2             Point_2;
51   typedef typename K::Segment_2           Segment_2;
52   typedef typename K::Site_2              Site_2;
53   typedef typename K::Oriented_side       Oriented_side;
54   typedef typename K::Bounded_side        Bounded_side;
55   typedef typename K::Comparison_result   Comparison_result;
56   typedef typename K::Sign                Sign;
57   typedef typename K::Orientation         Orientation;
58   typedef typename K::Compute_scalar_product_2
59                         Compute_scalar_product_2;
60   typedef typename K::Boolean             Boolean;
61   typedef typename K::Direction_2         Direction_2;
62   typedef typename K::Vector_2            Vector_2;
63   typedef typename K::Compare_x_2         Compare_x_2;
64   typedef typename K::Compare_y_2         Compare_y_2;
65 
66   typedef typename CGAL::Polychainline_2<K> Polychainline_2;
67   typedef SegmentDelaunayGraph_2::Are_same_points_C2<K>  Are_same_points_2;
68   typedef SegmentDelaunayGraph_2::Are_same_segments_C2<K>
69           Are_same_segments_2;
70 
71   typedef SegmentDelaunayGraph_2::Compare_x_2<K> Compare_x_2_Sites_Type;
72   typedef SegmentDelaunayGraph_2::Compare_y_2<K> Compare_y_2_Sites_Type;
73 
74   typedef unsigned int Bearing;
75 
76 private:
77   typedef typename K::Intersections_tag ITag;
78   typedef Bisector_Linf<K>              Bisector_Linf_Type;
79 
80   typedef typename Algebraic_structure_traits<RT>::Algebraic_category RT_Category;
81   typedef typename Algebraic_structure_traits<FT>::Algebraic_category FT_Category;
82 public:
83 
84   typedef typename Base::Line_2               Line_2;
85   typedef typename Base::Homogeneous_point_2  Homogeneous_point_2;
86 
87 public:    //    compute_supporting_line(q.supporting_segment(), a1, b1, c1);
88     //    compute_supporting_line(r.supporting_segment(), a2, b2, c2);
89 
90   //-------------------------------------------------------------------
91   // BASIC CONSTRUCTIONS
92   //-------------------------------------------------------------------
93 
94   using Base::compute_supporting_line;
95 
96   // compute horizontal line l that goes through p,
97   // and leaves q on the oriented side s
98   // s: has to be either +1 or -1 (not 0)
99   // q: should not be on line l
100   static
101   Line_2
compute_horizontal_side_lineBasic_predicates_C2102   compute_horizontal_side_line(
103       const Point_2& p, const Point_2& q, Oriented_side s)
104   {
105     CGAL_precondition(s != ON_ORIENTED_BOUNDARY);
106 
107     RT b, c;
108 
109     b = RT(1);
110     c = - p.y();
111 
112     // direction is (1, 0)
113 
114     Compare_y_2 cmpy;
115     if (((cmpy(q, p) == LARGER) &&
116          (s == ON_NEGATIVE_SIDE)   ) ||
117         ((cmpy(q, p) == SMALLER) &&
118          (s == ON_POSITIVE_SIDE)   )   ) {
119       b = -b;
120       c = -c;
121     }
122     return Line_2(RT(0), b, c);
123   }
124 
125   // compute vertical line l that goes through p,
126   // and leaves q on the oriented side s
127   // s: has to be either +1 or -1 (not 0)
128   // q: should not be on line l
129   static
130   Line_2
compute_vertical_side_lineBasic_predicates_C2131   compute_vertical_side_line(
132       const Point_2& p, const Point_2& q, Oriented_side s)
133   {
134     CGAL_precondition(s != ON_ORIENTED_BOUNDARY);
135 
136     RT a, c;
137 
138     a = RT(1);
139     c = - p.x();
140 
141     // direction is (0, -1)
142 
143     Compare_x_2 cmpx;
144     if (((cmpx(q, p) == LARGER) &&
145          (s == ON_NEGATIVE_SIDE)   ) ||
146         ((cmpx(q, p) == SMALLER) &&
147          (s == ON_POSITIVE_SIDE)   )   ) {
148       a = -a;
149       c = -c;
150     }
151     return Line_2(a, RT(0), c);
152   }
153 
154 
155   //using Base::compute_projection;
156 
157   static
158   Homogeneous_point_2
compute_linf_projection_homBasic_predicates_C2159   compute_linf_projection_hom(const Line_2& l, const Point_2& p)
160   {
161     //CGAL_precondition( ! l.is_degenerate() );
162     CGAL_precondition(
163         (CGAL::sign(l.a()) != ZERO) || (CGAL::sign(l.b()) != ZERO) );
164 
165     Sign signa = CGAL::sign(l.a());
166     Sign signb = CGAL::sign(l.b());
167 
168     RT hx, hy, hw;
169 
170     if (signa == ZERO) { // l is horizontal
171       // l.a() == 0  =>  y = -c/b
172       hx = p.x() * l.b();
173       hy = - l.c();
174       hw = l.b();
175     } else if (signb == ZERO) { // l is vertical
176       // l.b() == 0  =>  x = -c/a
177       hx = - l.c();
178       hy = p.y() * l.a();
179       hw = l.a();
180     } else {
181       // here both l.a() and l.b() are non-zero
182       if ( signa == signb ) {
183         hx = l.b() * ( p.x() - p.y() ) - l.c();
184         hy = l.a() * ( p.y() - p.x() ) - l.c();
185         hw = l.a() + l.b();
186       } else { // signa != signb
187         hx = -l.b() * ( p.x() + p.y() ) - l.c();
188         hy = l.a() * ( p.x() + p.y() ) + l.c();
189         hw = l.a() - l.b();
190       }
191     }
192 
193     return Homogeneous_point_2(hx, hy, hw);
194   }
195 
196   static
197   Point_2
compute_linf_projection_nonhomBasic_predicates_C2198   compute_linf_projection_nonhom(const Line_2& l, const Point_2& p)
199   {
200     //CGAL_precondition( ! l.is_degenerate() );
201     CGAL_precondition(
202         (CGAL::sign(l.a()) != ZERO) || (CGAL::sign(l.b()) != ZERO) );
203 
204     Sign signa = CGAL::sign(l.a());
205     Sign signb = CGAL::sign(l.b());
206 
207     RT hx, hy, hw;
208 
209     if (signa == ZERO) { // l is horizontal
210       // l.a() == 0  =>  y = -c/b
211       hx = p.x() * l.b();
212       hy = - l.c();
213       hw = l.b();
214     } else if (signb == ZERO) { // l is vertical
215       // l.b() == 0  =>  x = -c/a
216       hx = - l.c();
217       hy = p.y() * l.a();
218       hw = l.a();
219     } else {
220       // here both l.a() and l.b() are non-zero
221       if ( signa == signb ) {
222         hx = l.b() * ( p.x() - p.y() ) - l.c();
223         hy = l.a() * ( p.y() - p.x() ) - l.c();
224         hw = l.a() + l.b();
225       } else { // signa != signb
226         hx = -l.b() * ( p.x() + p.y() ) - l.c();
227         hy = l.a() * ( p.x() + p.y() ) + l.c();
228         hw = l.a() - l.b();
229       }
230     }
231 
232     return Point_2(hx, hy, hw);
233   }
234 
235   static
236   Homogeneous_point_2
compute_horizontal_projection_homBasic_predicates_C2237   compute_horizontal_projection_hom(const Line_2& l, const Point_2& p)
238   {
239     //CGAL_precondition( ! l.is_horizontal() );
240 
241     CGAL_precondition(
242         (CGAL::sign(l.a()) != ZERO) );
243 
244     RT hx, hy, hw;
245 
246     hx = - l.b() * p.y() - l.c();
247     hy = p.y() * l.a();
248     hw = l.a();
249 
250     return Homogeneous_point_2(hx, hy, hw);
251   }
252 
253   static
254   Point_2
compute_horizontal_projectionBasic_predicates_C2255   compute_horizontal_projection(const Line_2& l, const Point_2& p)
256   {
257     //CGAL_precondition( ! l.is_horizontal() );
258     CGAL_precondition(
259         (CGAL::sign(l.a()) != ZERO) );
260 
261     RT hx, hy, hw;
262 
263     hx = - l.b() * p.y() - l.c();
264     hy = p.y() * l.a();
265     hw = l.a();
266 
267     return Point_2(hx, hy, hw);
268   }
269 
270   static
271   Homogeneous_point_2
compute_vertical_projection_homBasic_predicates_C2272   compute_vertical_projection_hom(const Line_2& l, const Point_2& p)
273   {
274     //CGAL_precondition( ! l.is_horizontal() );
275     CGAL_precondition(
276         (CGAL::sign(l.b()) != ZERO) );
277 
278     RT hx, hy, hw;
279 
280     hx = p.x() * l.b();
281     hy = - l.a() * p.x() - l.c();
282     hw = l.b();
283 
284     return Homogeneous_point_2(hx, hy, hw);
285   }
286 
287   static
288   Point_2
compute_vertical_projectionBasic_predicates_C2289   compute_vertical_projection(const Line_2& l, const Point_2& p)
290   {
291     //CGAL_precondition( ! l.is_horizontal() );
292     CGAL_precondition(
293         (CGAL::sign(l.b()) != ZERO) );
294 
295     RT hx, hy, hw;
296 
297     hx = p.x() * l.b();
298     hy = - l.a() * p.x() - l.c();
299     hw = l.b();
300 
301     return Point_2(hx, hy, hw);
302   }
303 
304 
305   //using Base::projection_on_line;
306 
307   using Base::midpoint;
308 
309   using Base::compute_perpendicular;
310 
311   // compute_cw_perpendicular is opposite of compute_perpendicular
312   static
compute_cw_perpendicularBasic_predicates_C2313   Line_2 compute_cw_perpendicular(const Line_2& l, const Point_2& p)
314   {
315     RT a, b, c;
316     a = l.b();
317     b = -l.a();
318     c = -l.b() * p.x() + l.a() * p.y();
319     return Line_2(a, b, c);
320   }
321 
322   static
compute_linf_perpendicularBasic_predicates_C2323   Line_2 compute_linf_perpendicular(const Line_2& l, const Point_2& p)
324   {
325     RT a, b, c;
326     a = RT( - (int) CGAL::sign(l.b()) );
327     b = RT( (int) CGAL::sign(l.a()) );
328     c = - a * p.x() - b * p.y();
329     return Line_2(a, b, c);
330   }
331 
332   using Base::opposite_line;
333 
334   // philaris: similar to compute_supporting_line
335   static
compute_line_from_toBasic_predicates_C2336   Line_2 compute_line_from_to(const Point_2& p, const Point_2&q)
337   {
338     RT a, b, c;
339     a = p.y() - q.y();
340     b = q.x() - p.x();
341 
342     CGAL_assertion((CGAL::sign(a) != ZERO) ||
343                    (CGAL::sign(b) != ZERO))   ;
344 
345     c = p.x() * q.y() - q.x() * p.y();
346 
347     return Line_2(a, b, c);
348   }
349 
350   // compute line from a point and a direction
351   inline static
compute_line_dirBasic_predicates_C2352   Line_2 compute_line_dir(
353       const Homogeneous_point_2& p, const Direction_2& d)
354   {
355     return Line_2( -d.dy()*p.hw(), d.dx()*p.hw(),
356                    -(-d.dy()*p.hx() +d.dx()*p.hy()) );
357   }
358 
359   // compute bisector of two parallel lines
360   inline static
parallel_bisBasic_predicates_C2361   Line_2 parallel_bis(const Line_2& lp, const Line_2& lq)
362   {
363     RT bisa, bisb, bisc;
364     if ( CGAL::sign ( lq.a() ) != ZERO ) {
365       bisa = RT(2) * lp.a() * lq.a();
366       bisb = RT(2) * lp.a() * lq.b();
367       bisc = lp.a() * lq.c() + lp.c() * lq.a();
368     } else {
369       bisa = RT(2) * lp.a() * lq.b();
370       bisb = RT(2) * lp.b() * lq.b();
371       bisc = lp.c() * lq.b() + lp.b() * lq.c();
372     }
373     return Line_2(bisa, bisb, bisc);
374   }
375 
376   /* use point p for y coordinate of line */
377   static
compute_horizontal_line_from_toBasic_predicates_C2378   Line_2 compute_horizontal_line_from_to(const Point_2& p, const Point_2&q)
379   {
380     RT b, c;
381     Compare_x_2 cmpx;
382     Comparison_result cmpxqp = cmpx(q,p);
383     CGAL_assertion(cmpxqp != EQUAL);
384     b = (cmpxqp == SMALLER) ? RT(-1) : RT(1);
385     c = (cmpxqp == SMALLER) ? p.y() : RT(-p.y());
386     return Line_2(RT(0), b, c);
387   }
388 
389   /* use point p for x coordinate of line */
390   static
compute_vertical_line_from_toBasic_predicates_C2391   Line_2 compute_vertical_line_from_to(const Point_2& p, const Point_2&q)
392   {
393     RT a, c;
394     Compare_y_2 cmpy;
395     Comparison_result cmpypq = cmpy(p,q);
396     CGAL_assertion(cmpypq != EQUAL);
397     a = (cmpypq == SMALLER) ? RT(-1) : RT(1);
398     //a = RT(CGAL::sign(p.y() - q.y()));
399     c = (cmpypq == SMALLER) ? p.x() : RT(-p.x());
400     return Line_2(a, RT(0), c);
401   }
402 
403 
404   // neg slope 45 degree line passing through p
405   inline
406   static
compute_neg_45_line_atBasic_predicates_C2407   Line_2 compute_neg_45_line_at(const Point_2 & p)
408   {
409     return Line_2(p.hw() , p.hw(), -p.hx()-p.hy());
410   }
411 
412   // pos slope 45 degree line passing through p
413   inline
414   static
compute_pos_45_line_atBasic_predicates_C2415   Line_2 compute_pos_45_line_at(const Point_2 & p)
416   {
417     return Line_2(RT(1),RT(-1),p.y()-p.x());
418   }
419 
420   // horizontal line passing through p
421   inline
422   static
compute_hor_line_atBasic_predicates_C2423   Line_2 compute_hor_line_at(const Point_2 & p)
424   {
425     return Line_2(RT(0), p.hw(), - p.hy());
426   }
427 
428   // vertical line passing through p
429   inline
430   static
compute_ver_line_atBasic_predicates_C2431   Line_2 compute_ver_line_at(const Point_2 & p)
432   {
433     return Line_2(p.hw(), RT(0), - p.hx());
434   }
435 
436   static
compute_linf_distanceBasic_predicates_C2437   RT compute_linf_distance(const Point_2& p, const Point_2& q)
438   {
439     return (CGAL::max)(
440               CGAL::abs(p.x() - q.x()),
441               CGAL::abs(p.y() - q.y()));
442   }
443 
444   static
445   std::pair<RT,RT>
compute_linf_distanceBasic_predicates_C2446   compute_linf_distance(const Point_2& p, const Line_2& l)
447   {
448     const RT nomin = CGAL::abs(l.a() * p.x() + l.b() * p.y() + l.c());
449     const RT denom = CGAL::abs(
450           l.a() +
451           ( CGAL::sign(l.a()) == CGAL::sign(l.b())? l.b() : RT(-l.b()) ) );
452     return std::pair<RT,RT>(nomin, denom);
453   }
454 
455   static
compute_intersection_of_linesBasic_predicates_C2456   void compute_intersection_of_lines(
457       const Line_2& l1, const Line_2& l2,
458       RT& hx, RT& hy, RT& hw)
459   {
460     hx = l1.b() * l2.c() - l1.c() * l2.b();
461     hy = l1.c() * l2.a() - l1.a() * l2.c();
462     hw = l1.a() * l2.b() - l1.b() * l2.a();
463   }
464 
465   inline
466   static
coord_atBasic_predicates_C2467   RT coord_at(const Line_2 &l, const RT & val, const bool return_y_coord)
468   {
469     return (return_y_coord)?
470       (l.a() * val + l.c()) / (- l.b()) :
471       (l.b() * val + l.c()) / (- l.a()) ;
472   }
473 
474   inline
475   static
touch_same_sideBasic_predicates_C2476   bool touch_same_side(
477       const Site_2 & p, const Site_2 & q, const Line_2 & l,
478       const bool samexpq, const bool pos_slope)
479   {
480     const RT common_coord = (samexpq) ? p.point().x() : p.point().y();
481     const RT otherp = (samexpq) ? p.point().y() : p.point().x();
482     const RT otherq = (samexpq) ? q.point().y() : q.point().x();
483     const RT lineval = coord_at(l, common_coord, samexpq);
484     return (CGAL::sign(lineval - otherp) == CGAL::sign(otherp - otherq)) ?
485       samexpq == pos_slope :
486       samexpq != pos_slope ;
487   }
488 
489   inline
490   static
is_orth_dist_smaller_than_pt_distBasic_predicates_C2491   bool is_orth_dist_smaller_than_pt_dist(
492       const RT & closest_coord, const Line_2 & l,
493       const Site_2 & p, const Site_2 & q,
494       const bool samexpq)
495   {
496     const RT lineval = coord_at(l, closest_coord, ! samexpq);
497     return CGAL::abs(lineval - ((samexpq) ? p.point().x() : p.point().y()))
498            <
499            CGAL::abs((samexpq)? p.point().y() - q.point().y() :
500                                 p.point().x() - q.point().x()  ) ;
501   }
502 
503 
504 public:
505   //-------------------------------------------------------------------
506   // BASIC PREDICATES
507   //-------------------------------------------------------------------
508   static
509   Comparison_result
compare_linf_distances_to_lineBasic_predicates_C2510   compare_linf_distances_to_line(const Line_2& l, const Point_2& p,
511                                     const Point_2& q)
512   {
513     Homogeneous_point_2 hp = compute_linf_projection_hom(l, p);
514     Homogeneous_point_2 hq = compute_linf_projection_hom(l, q);
515 
516     RT dlp = (CGAL::max)(CGAL::abs(p.x() - hp.x()),
517                          CGAL::abs(p.y() - hp.y()));
518 
519     RT dlq = (CGAL::max)(CGAL::abs(q.x() - hq.x()),
520                          CGAL::abs(q.y() - hq.y()));
521 
522     Comparison_result crude = CGAL::compare(dlp, dlq);
523 
524     if (crude != EQUAL) {
525       return crude;
526     } else {
527       CGAL_SDG_DEBUG(std::cout
528           << "compare_linf_distances_to_line refining"
529           << std::endl;);
530       return crude;
531     }
532   }
533 
534   static
535   Comparison_result
compare_linf_distances_to_linesBasic_predicates_C2536   compare_linf_distances_to_lines(const Point_2& p,
537                                      const Line_2& l1,
538                                      const Line_2& l2)
539   {
540     Homogeneous_point_2 hl1 = compute_linf_projection_hom(l1, p);
541     Homogeneous_point_2 hl2 = compute_linf_projection_hom(l2, p);
542 
543     RT dl1p = (CGAL::max)(CGAL::abs(hl1.x() - p.x()),
544                           CGAL::abs(hl1.y() - p.y()));
545 
546     RT dl2p = (CGAL::max)(CGAL::abs(hl2.x() - p.x()),
547                           CGAL::abs(hl2.y() - p.y()));
548 
549     Comparison_result crude = CGAL::compare(dl1p, dl2p);
550 
551     if (crude != EQUAL) {
552       return crude;
553     } else {
554       CGAL_SDG_DEBUG(std::cout << "compare_linf_distances_to_lines refining"
555                      << std::endl;);
556       return crude;
557     }
558   }
559 
560   using Base::oriented_side_of_line;
561 
562   using Base::is_on_positive_halfspace;
563 
564   /* compares the Linf distances pq and pr;
565      if the Linf distances pq are the same, then try to break ties
566      by comparing the minima of the dx and dy distance components
567    */
568   static
569   Comparison_result
compare_distance_to_point_linfBasic_predicates_C2570   compare_distance_to_point_linf(
571       const Point_2& p, const Point_2& q, const Point_2& r)
572   {
573     const RT pqdx = CGAL::abs(p.x()-q.x());
574     const RT pqdy = CGAL::abs(p.y()-q.y());
575     const bool pqdxlarger = CGAL::compare(pqdx, pqdy) == LARGER;
576     const RT & pqmax = pqdxlarger ? pqdx : pqdy;
577     const RT & pqmin = pqdxlarger ? pqdy : pqdx;
578 
579     const RT prdx = CGAL::abs(p.x()-r.x());
580     const RT prdy = CGAL::abs(p.y()-r.y());
581     const bool prdxlarger = CGAL::compare(prdx, prdy) == LARGER;
582     const RT & prmax = prdxlarger ? prdx : prdy;
583     const RT & prmin = prdxlarger ? prdy : prdx;
584 
585     const Comparison_result resmax = CGAL::compare(pqmax, prmax);
586 
587     if (resmax == EQUAL) {
588       CGAL_SDG_DEBUG(std::cout <<
589           "debug cmpdistlinf break ties with min" << std::endl;);
590       return CGAL::compare(pqmin, prmin);
591     } else {
592       return resmax;
593     }
594   }
595 
596   static
597   Bounded_side
bounded_side_of_bboxBasic_predicates_C2598   bounded_side_of_bbox(
599       const Point_2& p, const Point_2& q, const Point_2& r)
600   {
601     // precondition: p, q, r should be monotone points.
602     // Predicate bounded_side_of_bbox (P_bbox) returns:
603     //  0 if p = q,
604     //  0 if r = p or r = q (ON_BOUNDARY),
605     // -1 if r is strictly outside the bounding box of p,q
606     //    (ON_UNBOUNDED_SIDE),
607     // +1 if r is strictly inside the bounding box of p,q
608     //    (ON_BOUNDED_SIDE).
609     // If p and q are on the same vertical or horizontal
610     // line but are not the same point, then the bounding
611     // box of p and q degenerates to the line segment pq.
612 
613     CGAL_SDG_DEBUG(std::cout << "debug bounded_side_of_bbox (p q r) = ("
614                    << p << ") (" << q << ") (" << r << ")" << std::endl; );
615 
616     if ((CGAL::compare(p.x(), q.x()) == EQUAL) &&
617         (CGAL::compare(p.y(), q.y()) == EQUAL)    ) {
618       return ON_BOUNDARY;
619     }
620 
621     Comparison_result cmpxpr, cmpxrq, cmpypr, cmpyrq;
622 
623     cmpxpr = CGAL::compare(p.x(), r.x());
624     cmpxrq = CGAL::compare(r.x(), q.x());
625     cmpypr = CGAL::compare(p.y(), r.y());
626     cmpyrq = CGAL::compare(r.y(), q.y());
627 
628     Comparison_result comp =
629       CGAL::compare(cmpxpr*cmpxrq + cmpypr*cmpyrq, 0);
630 
631     CGAL_SDG_DEBUG(std::cout << "debug bounded_side_of_bbox returns ";);
632 
633     switch(comp) {
634       case SMALLER:
635         CGAL_SDG_DEBUG(std::cout << "ON_UNBOUNDED_SIDE" << std::endl;);
636         return ON_UNBOUNDED_SIDE;
637       case EQUAL:
638         CGAL_SDG_DEBUG(std::cout << "ON_BOUNDARY" << std::endl;);
639         return ON_BOUNDARY;
640       case LARGER:
641         CGAL_SDG_DEBUG(std::cout << "ON_BOUNDED_SIDE" << std::endl;);
642         return ON_BOUNDED_SIDE;
643       default:
644         CGAL_SDG_DEBUG(std::cout << "error: should never reach here";);
645         CGAL_assertion( false );
646         return ON_BOUNDARY;
647     }
648   }
649 
650 
651   // returns true if and only if
652   // the interior of s has non-empty intersection
653   // with the positive halfplane of oriented line l
654   static
655   Boolean
intersects_segment_positive_halfplaneBasic_predicates_C2656   intersects_segment_positive_halfplane(
657       const Site_2 & s,
658       const Line_2 & l)
659   {
660     Segment_2 seg = s.segment();
661 
662     Point_2 ssrc = seg.source();
663     Point_2 strg = seg.target();
664 
665     CGAL_SDG_DEBUG(std::cout
666       << "debug: intersects_segment_positive_halfplane "
667       << "s=" << s
668       << " l=" << l.a() << " " << l.b() << " " << l.c()
669                    << std::endl;);
670 
671     Oriented_side oslsrc = oriented_side_of_line(l, ssrc);
672     Oriented_side osltrg = oriented_side_of_line(l, strg);
673 
674       CGAL_SDG_DEBUG(std::cout
675           << "debug: intersects_segment_positive_halfplane "
676           << "oslsrc=" << oslsrc << " osltrg=" << osltrg
677                      << std::endl;);
678 
679     if ((oslsrc == ON_POSITIVE_SIDE) ||
680         (osltrg == ON_POSITIVE_SIDE)   )
681     {
682       return true;
683     } else {
684       return false;
685     }
686   }
687 
688 
689   // returns true if and only if
690   // the interior of s has non-empty intersection
691   // with the negative halfplane of oriented line l
692   static
693   Boolean
intersects_segment_negative_halfplaneBasic_predicates_C2694   intersects_segment_negative_halfplane(
695       const Site_2 & s,
696       const Line_2 & l)
697   {
698     Segment_2 seg = s.segment();
699 
700     Point_2 ssrc = seg.source();
701     Point_2 strg = seg.target();
702 
703     CGAL_SDG_DEBUG(std::cout
704         << "debug: intersects_segment_negative_halfplane "
705         << "s=" << s
706         << " l=" << l.a() << " " << l.b() << " " << l.c()
707         << std::endl;);
708 
709     Oriented_side oslsrc = oriented_side_of_line(l, ssrc);
710     Oriented_side osltrg = oriented_side_of_line(l, strg);
711 
712     CGAL_SDG_DEBUG(std::cout
713         << "debug: intersects_segment_negative_halfplane "
714         << "oslsrc=" << oslsrc << " osltrg=" << osltrg
715         << std::endl;);
716 
717     if ((oslsrc == ON_NEGATIVE_SIDE) ||
718         (osltrg == ON_NEGATIVE_SIDE)   )
719     {
720       return true;
721     } else {
722       return false;
723     }
724   }
725 
726   // returns true if and only if
727   // the interior of s has non-empty intersection
728   // with the interior of the following infinite box:
729   // the only finite corner of the infinite box is corner
730   // and if you traverse the infinite box ccw, then
731   // you meet points in that order: q, corner, p
732   static
733   Boolean
intersects_segment_interior_inf_boxBasic_predicates_C2734   intersects_segment_interior_inf_box(const Site_2 & s,
735       const Site_2 & q, const Site_2 & p,
736       const Comparison_result & cmpxpq, const Comparison_result & cmpypq)
737   {
738     CGAL_assertion(cmpxpq != EQUAL);
739     CGAL_assertion(cmpypq != EQUAL);
740     CGAL_assertion(s.is_segment());
741     const Segment_2 seg = s.segment();
742 
743     const Point_2 ssrc = seg.source();
744     const Point_2 strg = seg.target();
745 
746     const Point_2 qq = q.point();
747     const Point_2 pp = p.point();
748 
749     const bool eqcmp = cmpxpq == cmpypq;
750 
751     Are_same_points_2 same_points;
752     Compare_x_2 cmpx;
753     Compare_y_2 cmpy;
754 
755     bool is_ssrc_positive;
756     if (same_points(q, s.source_site()) ||
757         same_points(p, s.source_site())   ) {
758       is_ssrc_positive = false;
759     } else {
760       const bool conflp = eqcmp ?
761         (cmpx(pp, ssrc) == cmpxpq) : (cmpy(pp, ssrc) == cmpypq) ;
762       const bool conflq = eqcmp ?
763         (cmpy(ssrc, qq) == cmpypq) : (cmpx(ssrc, qq) == cmpxpq) ;
764       is_ssrc_positive = (conflp && conflq);
765     }
766     if (is_ssrc_positive) {
767       CGAL_SDG_DEBUG(std::cout << "debug is_segment_inside_inf_box "
768                      << "src endpoint inside" << std::endl;);
769       return true;
770     }
771 
772     bool is_strg_positive;
773     if (same_points(q, s.target_site()) ||
774         same_points(p, s.target_site())   ) {
775       is_strg_positive = false;
776     } else {
777       const bool conflp = eqcmp ?
778         (cmpx(pp, strg) == cmpxpq) : (cmpy(pp, strg) == cmpypq) ;
779       const bool conflq = eqcmp ?
780         (cmpy(strg, qq) == cmpypq) : (cmpx(strg, qq) == cmpxpq) ;
781       is_strg_positive = (conflp && conflq);
782     }
783 
784     if (is_strg_positive) {
785       CGAL_SDG_DEBUG(std::cout << "debug is_segment_inside_inf_box "
786                      << "trg endpoint inside" << std::endl;);
787       return true;
788     } else {
789       // here you have to check if the interior is inside
790 
791       CGAL_SDG_DEBUG(std::cout << "debug is_segment_inside_inf_box "
792                      << "try for interior to be inside" << std::endl;);
793 
794       const Point_2 corner = eqcmp ?
795         Point_2( pp.x(), qq.y() ) :
796         Point_2( qq.x(), pp.y() ) ;
797 
798       // in fact, here you can intersect the segment
799       // with the ray starting from corner and going to the
800       // direction of the center of the infinite box
801 
802       const RT one(1);
803 
804       const Point_2 displaced ( corner.x() + (-(int)cmpypq)*one ,
805                                 corner.y() + (int)cmpxpq * one   );
806 
807       const Line_2 l = compute_line_from_to(corner, displaced);
808 
809       const Line_2 lseg = compute_supporting_line(s.supporting_site());
810 
811       RT hx, hy, hw;
812       compute_intersection_of_lines(l, lseg, hx, hy, hw);
813 
814       if (CGAL::sign(hw) == ZERO) {
815         return false;
816       } else {
817         const Point_2 ip ( hx, hy, hw );
818         const Line_2 lqc = compute_line_from_to(qq, corner);
819         const Line_2 lcp = compute_line_from_to(corner, pp);
820         const Oriented_side os_lqc_ip = oriented_side_of_line(lqc, ip);
821         const Oriented_side os_lcp_ip = oriented_side_of_line(lcp, ip);
822 
823         const Comparison_result cmpxsrcip = cmpx(ssrc, ip);
824         const Comparison_result cmpysrcip = cmpy(ssrc, ip);
825         const Comparison_result cmpxiptrg = cmpx(ip, strg);
826         const Comparison_result cmpyiptrg = cmpy(ip, strg);
827 
828         // philaris: to check
829         Boolean is_ip_inside_segment =
830           (CGAL::sign(cmpxsrcip * cmpxiptrg +
831                       cmpysrcip * cmpyiptrg   )) == POSITIVE;
832 
833         if ((os_lqc_ip == ON_POSITIVE_SIDE) &&
834             (os_lcp_ip == ON_POSITIVE_SIDE) &&
835             is_ip_inside_segment ) {
836           return true;
837         } else {
838           return false;
839         }
840       }
841     }
842   } // end of intersects_segment_interior_inf_box
843 
844   // returns true if and only if
845   // the interior of s has non-empty intersection
846   // with the interior of the following infinite box:
847   // the only finite corner of the infinite box is corner
848   // and if you traverse the infinite box ccw, then
849   // you meet points in that order: q, corner, p
850   static
851   Boolean
intersects_segment_interior_inf_boxBasic_predicates_C2852   intersects_segment_interior_inf_box(const Site_2 & s,
853       const Site_2 & q, const Point_2 & corner,
854       const Site_2 & p)
855   {
856     CGAL_assertion(s.is_segment());
857     Segment_2 seg = s.segment();
858 
859     Point_2 ssrc = seg.source();
860     Point_2 strg = seg.target();
861 
862     Point_2 qq = q.point();
863     Point_2 pp = p.point();
864 
865     Line_2 lqc = compute_line_from_to(qq, corner);
866     Line_2 lcp = compute_line_from_to(corner, pp);
867 
868     Are_same_points_2 same_points;
869 
870     bool is_ssrc_positive;
871     if (same_points(q, s.source_site()) ||
872         same_points(p, s.source_site())   ) {
873       is_ssrc_positive = false;
874     } else {
875       Oriented_side os_lqc_ssrc = oriented_side_of_line(lqc, ssrc);
876       Oriented_side os_lcp_ssrc = oriented_side_of_line(lcp, ssrc);
877       is_ssrc_positive =
878         ((os_lqc_ssrc == ON_POSITIVE_SIDE) &&
879          (os_lcp_ssrc == ON_POSITIVE_SIDE)    ) ;
880     }
881 
882     bool is_strg_positive;
883     if (same_points(q, s.target_site()) ||
884         same_points(p, s.target_site())   ) {
885       is_strg_positive = false;
886     } else {
887       Oriented_side os_lqc_strg = oriented_side_of_line(lqc, strg);
888       Oriented_side os_lcp_strg = oriented_side_of_line(lcp, strg);
889       is_strg_positive =
890         ((os_lqc_strg == ON_POSITIVE_SIDE) &&
891          (os_lcp_strg == ON_POSITIVE_SIDE)    ) ;
892     }
893 
894     CGAL_SDG_DEBUG(std::cout << "debug qcp= (" << q << ") (" << corner
895         << ") (" << p << ")"
896         << " isssrcpos=" << is_ssrc_positive
897         << " isstrgpos=" << is_strg_positive
898         << std::endl;);
899 
900     if (is_ssrc_positive || is_strg_positive) {
901       CGAL_SDG_DEBUG(std::cout << "debug is_segment_inside_inf_box "
902                      << "endpoint inside" << std::endl;);
903       return true;
904     } else {
905       // here you have to check if the interior is inside
906 
907       CGAL_SDG_DEBUG(std::cout << "debug is_segment_inside_inf_box "
908                      << "try for interior to be inside" << std::endl;);
909 
910       // in fact, here you can intersect the segment
911       // with the ray starting from corner and going to the
912       // direction of the center of the infinite box
913 
914       Compare_x_2 cmpx;
915       Compare_y_2 cmpy;
916 
917       Comparison_result cmpxpq = cmpx(pp,qq);
918       Comparison_result cmpypq = cmpy(pp,qq);
919 
920       RT one(1);
921 
922       Point_2 displaced ( corner.x() + (-(int)cmpypq)*one ,
923                           corner.y() + (int)cmpxpq * one   );
924 
925       Line_2 l = compute_line_from_to(corner, displaced);
926 
927       Line_2 lseg = compute_supporting_line(s.supporting_site());
928 
929       RT hx, hy, hw;
930 
931       compute_intersection_of_lines(l, lseg, hx, hy, hw);
932 
933       if (CGAL::sign(hw) == ZERO) {
934         return false;
935       } else {
936         Point_2 ip ( hx/hw, hy/hw);
937         Oriented_side os_lqc_ip = oriented_side_of_line(lqc, ip);
938         Oriented_side os_lcp_ip = oriented_side_of_line(lcp, ip);
939 
940         Compare_x_2 cmpx;
941         Compare_y_2 cmpy;
942 
943         Comparison_result cmpxsrcip = cmpx(ssrc, ip);
944         Comparison_result cmpysrcip = cmpy(ssrc, ip);
945         Comparison_result cmpxiptrg = cmpx(ip, strg);
946         Comparison_result cmpyiptrg = cmpy(ip, strg);
947 
948         // philaris: to check
949         Boolean is_ip_inside_segment =
950           (CGAL::sign(cmpxsrcip * cmpxiptrg +
951                       cmpysrcip * cmpyiptrg   )) == POSITIVE;
952 
953         if ((os_lqc_ip == ON_POSITIVE_SIDE) &&
954             (os_lcp_ip == ON_POSITIVE_SIDE) &&
955             is_ip_inside_segment ) {
956           return true;
957         } else {
958           return false;
959         }
960       }
961     }
962   } // end of intersects_segment_interior_inf_box
963 
964 
965   // returns true if and only if
966   // the interior of s has non-empty intersection
967   // with the interior of the following infinite box:
968   // this infinite box is a 90 degree wedge defined
969   // by the intersection of the halfplanes
970   // with supporting lines lhor and lver, where
971   // the halfplanes are both on the positive or negative
972   // sides of the supporting lines
973   static
974   Boolean
intersects_segment_side_of_wedgeBasic_predicates_C2975   intersects_segment_side_of_wedge(const Site_2 & s,
976       const Line_2 & lhor, const Line_2 & lver,
977       Oriented_side orside)
978   {
979     CGAL_assertion(s.is_segment());
980     Segment_2 seg = s.segment();
981 
982     CGAL_SDG_DEBUG(std::cout << "debug sofw s=" << s
983                    << " orside=" << orside << std::endl;);
984 
985     Point_2 ssrc = seg.source();
986     Point_2 strg = seg.target();
987 
988     Oriented_side os_lhor_ssrc = oriented_side_of_line(lhor, ssrc);
989     Oriented_side os_lver_ssrc = oriented_side_of_line(lver, ssrc);
990 
991     Oriented_side os_lhor_strg = oriented_side_of_line(lhor, strg);
992     Oriented_side os_lver_strg = oriented_side_of_line(lver, strg);
993 
994     if (((os_lhor_ssrc == orside) &&
995          (os_lver_ssrc == orside)) ||
996         ((os_lhor_strg == orside) &&
997          (os_lver_strg == orside))   ) {
998           CGAL_SDG_DEBUG(std::cout
999               << "debug intersects_segment_side_of_wedge "
1000               << "endpoint inside" << std::endl;);
1001       return true;
1002     } else {
1003       // here we have to check if the interior is inside
1004 
1005       CGAL_SDG_DEBUG(std::cout
1006           << "debug intersects_segment_side_of_wedge "
1007           << "try for interior to be inside" << std::endl;);
1008 
1009       // in fact, here you can intersect the segment
1010       // with the ray starting from corner and going to the
1011       // direction of the center of the infinite box
1012 
1013       // corner has homogenuous coordinates cx, cy, cw
1014       RT cx, cy, cw;
1015       compute_intersection_of_lines(lhor, lver, cx, cy, cw);
1016 
1017       CGAL_assertion( CGAL::sign(cw) != ZERO );
1018 
1019       Point_2 corner ( cx, cy, cw );
1020 
1021       CGAL_SDG_DEBUG(std::cout << "debug corner=" << corner << std::endl;);
1022 
1023       RT one(1);
1024 
1025       Point_2 displaced (
1026           corner.x() + ( (+orside)*CGAL::sign(lver.a()) ) * one ,
1027           corner.y() +   (+orside)*CGAL::sign(lhor.b())   * one   );
1028 
1029       CGAL_SDG_DEBUG(std::cout
1030           << "debug displaced=" << displaced << std::endl;);
1031 
1032       Line_2 l = compute_line_from_to(corner, displaced);
1033 
1034       Line_2 lseg = compute_supporting_line(s.supporting_site());
1035 
1036       RT hx, hy, hw;
1037 
1038       CGAL_SDG_DEBUG(std::cout
1039           << "debug: intersects_segment_side_of_wedge "
1040           << " l=" << l.a() << " " << l.b() << " " << l.c()
1041           << " lseg=" << lseg.a() << " " << lseg.b() << " " << lseg.c()
1042           << std::endl;);
1043 
1044       compute_intersection_of_lines(l, lseg, hx, hy, hw);
1045 
1046       if (CGAL::sign(hw) == ZERO) {
1047         CGAL_SDG_DEBUG(std::cout
1048             << "debug l and lseg are parallel" << std::endl;);
1049         return false;
1050       } else {
1051         Point_2 ip ( hx, hy, hw );
1052         CGAL_SDG_DEBUG(std::cout << "debug ip=" << ip << std::endl;);
1053         CGAL_SDG_DEBUG(std::cout << "debug ip_hom="
1054             << hx << ' ' << hy << ' ' << hw << std::endl;);
1055         Oriented_side os_lhor_ip = oriented_side_of_line(lhor, ip);
1056         Oriented_side os_lver_ip = oriented_side_of_line(lver, ip);
1057 
1058         Compare_x_2 cmpx;
1059         Compare_y_2 cmpy;
1060 
1061         Comparison_result cmpxsrcip = cmpx(ssrc, ip);
1062         Comparison_result cmpysrcip = cmpy(ssrc, ip);
1063         Comparison_result cmpxiptrg = cmpx(ip, strg);
1064         Comparison_result cmpyiptrg = cmpy(ip, strg);
1065 
1066         // philaris: to check
1067         Boolean is_ip_inside_segment =
1068           (CGAL::sign(cmpxsrcip * cmpxiptrg +
1069                       cmpysrcip * cmpyiptrg   )) == POSITIVE;
1070 
1071         if ((os_lhor_ip == orside) &&
1072             (os_lver_ip == orside) &&
1073             is_ip_inside_segment      ) {
1074           return true;
1075         } else {
1076           return false;
1077         }
1078       }
1079     }
1080   } // end of intersects_segment_side_of_wedge
1081 
1082   // returns true if and only if
1083   // the interior of s has non-empty intersection
1084   // with the interior of the following infinite box:
1085   // this infinite box is a 90 degree wedge defined
1086   // by the intersection of the positive halfplanes
1087   // with supporting lines lhor and lver
1088   static
1089   Boolean
intersects_segment_positive_of_wedgeBasic_predicates_C21090   intersects_segment_positive_of_wedge(const Site_2 & s,
1091       const Line_2 & lhor, const Line_2 & lver)
1092   {
1093     return intersects_segment_side_of_wedge(
1094         s, lhor, lver, ON_POSITIVE_SIDE);
1095   }
1096 
1097   // returns true if and only if
1098   // the interior of s has non-empty intersection
1099   // with the interior of the following infinite box:
1100   // this infinite box is a 90 degree wedge defined
1101   // by the intersection of the positive halfplanes
1102   // with supporting lines lhor and lver
1103   static
1104   Boolean
intersects_segment_negative_of_wedgeBasic_predicates_C21105   intersects_segment_negative_of_wedge(const Site_2 & s,
1106       const Line_2 & lhor, const Line_2 & lver)
1107   {
1108     return intersects_segment_side_of_wedge(
1109         s, lhor, lver, ON_NEGATIVE_SIDE);
1110   }
1111 
1112   // returns true if and only if
1113   // the interior of t has non-empty intersection
1114   // with the interior of the following infinite box:
1115   // the only finite corner of the infinite box is p
1116   // and this infinite box is a 90 degree wedge defined
1117   // by the intersection of the halfplanes
1118   // with supporting lines lhor and lver through p, where
1119   // the halfplanes are both on the positive or negative
1120   // sides of the supporting lines. The positive or negative
1121   // side depends on the relative position of p with respect to s
1122 
1123   static
1124   Boolean
intersects_segment_interior_inf_wedge_spBasic_predicates_C21125   intersects_segment_interior_inf_wedge_sp(const Site_2 & s,
1126                                            const Site_2 & p,
1127                                            const Site_2 & t)
1128   {
1129     CGAL_assertion( t.is_segment() );
1130     CGAL_assertion( s.is_segment() );
1131     CGAL_assertion(! is_site_h_or_v(s));
1132 
1133     Segment_2 seg = s.segment();
1134 
1135     Point_2 ssrc = seg.source();
1136     Point_2 strg = seg.target();
1137 
1138     Point_2 pp = p.point();
1139 
1140     Sign dxs = CGAL::sign(strg.x() - ssrc.x());
1141     Sign dys = CGAL::sign(strg.y() - ssrc.y());
1142 
1143     Line_2 lseg = compute_supporting_line(s.supporting_site());
1144     Oriented_side os_lseg_p = oriented_side_of_line(lseg, pp);
1145 
1146     CGAL_assertion( os_lseg_p != ON_ORIENTED_BOUNDARY );
1147 
1148     //sandeep: lhor, lver remains same for left turn and right turn
1149 
1150     if (dxs == NEGATIVE && dys == NEGATIVE) {
1151 
1152       Line_2 lhor = Line_2(0,1,-pp.y());
1153       Line_2 lver = Line_2(-1,0,pp.x());
1154 
1155       return intersects_segment_side_of_wedge(t,
1156                                        lhor, lver,
1157                                        os_lseg_p);
1158     } else if (dxs == POSITIVE && dys == NEGATIVE) {
1159 
1160       Line_2 lhor = Line_2(0,-1,pp.y());
1161       Line_2 lver = Line_2(-1,0,pp.x());
1162 
1163       return intersects_segment_side_of_wedge(t,
1164                                        lhor, lver,
1165                                        os_lseg_p);
1166     } else if (dxs == POSITIVE && dys == POSITIVE) {
1167 
1168       Line_2 lhor = Line_2(0,-1,pp.y());
1169       Line_2 lver = Line_2(1,0,-pp.x());
1170 
1171       return intersects_segment_side_of_wedge(t,
1172                                        lhor, lver,
1173                                        os_lseg_p);
1174     } else {//dxs == NEGATIVE and dys == POSITIVE
1175 
1176       Line_2 lhor = Line_2(0,1,-pp.y());
1177       Line_2 lver = Line_2(1,0,-pp.x());
1178 
1179       return intersects_segment_side_of_wedge(t,
1180                                        lhor, lver,
1181                                        os_lseg_p);
1182     }
1183 
1184   } // end of intersects_segment_interior_inf_wedge_sp
1185 
1186 
1187   // returns true if and only if
1188   // the interior of s has non-empty intersection
1189   // with the interior of the bounding box of q, p
1190   // precondition: the bounding box should be non-trivial,
1191   // i.e., it should not be a segment
1192   static
1193   Boolean
intersects_segment_interior_bboxBasic_predicates_C21194   intersects_segment_interior_bbox(const Site_2 & s,
1195       const Site_2 & q,
1196       const Site_2 & p)
1197   {
1198     CGAL_precondition(s.is_segment());
1199     CGAL_precondition(p.is_point());
1200     CGAL_precondition(q.is_point());
1201 
1202     Point_2 pp = p.point();
1203     Point_2 qq = q.point();
1204 
1205     CGAL_assertion_code( Compare_x_2 cmpx; )
1206     CGAL_assertion_code( Compare_y_2 cmpy; )
1207     CGAL_assertion(cmpx(pp,qq) != EQUAL);
1208     CGAL_assertion(cmpy(pp,qq) != EQUAL);
1209 
1210     Point_2 corner1 ( pp.x(), qq.y());
1211     Point_2 corner2 ( qq.x(), pp.y());
1212 
1213     if (CGAL::orientation( qq, corner1, pp ) == LEFT_TURN) {
1214       return intersects_segment_interior_inf_box(s, q, corner1, p)
1215           && intersects_segment_interior_inf_box(s, p, corner2, q);
1216     } else {
1217       return intersects_segment_interior_inf_box(s, q, corner2, p)
1218           && intersects_segment_interior_inf_box(s, p, corner1, q);
1219     }
1220   } // end of intersects_segment_interior_bbox
1221 
1222   // returns true if and only if
1223   // the non-horizontal/non-vertical segment at site s
1224   // has positive slope
1225   static
1226   Boolean
has_positive_slopeBasic_predicates_C21227   has_positive_slope(const Site_2 & s)
1228   {
1229     CGAL_precondition(s.is_segment());
1230     CGAL_precondition(! is_site_h_or_v(s));
1231     Compare_x_2 cmpx;
1232     Compare_y_2 cmpy;
1233     Point_2 src = s.supporting_site().source();
1234     Point_2 trg = s.supporting_site().target();
1235     return cmpx(src, trg) == cmpy(src, trg);
1236   }
1237 
1238   inline
1239   static
1240   bool
has_positive_slopeBasic_predicates_C21241   has_positive_slope(const Line_2 & l) {
1242     return (CGAL::sign(l.a()) + CGAL::sign(l.b()) == 0);
1243   }
1244 
1245   inline
1246   static
1247   Boolean
have_same_slopeBasic_predicates_C21248   have_same_slope(const Site_2 & s, const Site_2 & t)
1249   {
1250     CGAL_precondition(s.is_segment());
1251     CGAL_precondition(t.is_segment());
1252     Compare_x_2 cmpx;
1253     Compare_y_2 cmpy;
1254     Point_2 ssrc = s.supporting_site().source();
1255     Point_2 strg = s.supporting_site().target();
1256     Comparison_result scmpx = cmpx(ssrc, strg);
1257     Comparison_result scmpy = cmpy(ssrc, strg);
1258     Point_2 tsrc = t.supporting_site().source();
1259     Point_2 ttrg = t.supporting_site().target();
1260     Comparison_result tcmpx = cmpx(tsrc, ttrg);
1261     Comparison_result tcmpy = cmpy(tsrc, ttrg);
1262     CGAL_SDG_DEBUG(std::cout << "debug have_same_slope"
1263         << " scmpx=" << scmpx << " scmpy=" << scmpy
1264         << " tcmpx=" << tcmpx << " tcmpy=" << tcmpy
1265         << std::endl;);
1266     if (   ((scmpx == EQUAL) && (tcmpx == EQUAL)) // vertical
1267         || ((scmpy == EQUAL) && (tcmpy == EQUAL)) // horizontal
1268         || ((scmpx == scmpy) && (tcmpx == tcmpy)) // positive
1269         || ((scmpx != EQUAL) && (scmpy != EQUAL) &&
1270             (tcmpx != EQUAL) && (tcmpy != EQUAL) &&
1271             (scmpx != scmpy) && (tcmpx != tcmpy)) // negative
1272        ) {
1273       return true;
1274     } else {
1275       return false;
1276     }
1277   }
1278 
1279   inline
1280   static
1281   Boolean
is_site_h_or_vBasic_predicates_C21282   is_site_h_or_v(const Site_2 & s)
1283   {
1284     return is_site_horizontal(s) || is_site_vertical(s);
1285   }
1286 
1287   inline
1288   static
1289   Boolean
is_site_horizontalBasic_predicates_C21290   is_site_horizontal(const Site_2 & s)
1291   {
1292     CGAL_assertion(s.is_segment());
1293     return s.supporting_site().segment().is_horizontal();
1294   }
1295 
1296   inline
1297   static
1298   Boolean
is_site_verticalBasic_predicates_C21299   is_site_vertical(const Site_2 & s)
1300   {
1301     CGAL_assertion(s.is_segment());
1302     return s.supporting_site().segment().is_vertical();
1303   }
1304 
1305   inline
1306   static
1307   Boolean
is_line_h_or_vBasic_predicates_C21308   is_line_h_or_v(const Line_2 & l)
1309   {
1310     return (CGAL::sign(l.a()) == ZERO) || (CGAL::sign(l.b()) == ZERO);
1311   }
1312 
1313   inline
1314   static
1315   Boolean
test_starBasic_predicates_C21316   test_star(const Site_2 & p, const Site_2 & u,
1317             const Site_2 & v, const Site_2 & t) {
1318     CGAL_precondition(p.is_point());
1319     CGAL_precondition(u.is_segment());
1320     CGAL_precondition(v.is_segment());
1321     CGAL_precondition(t.is_segment());
1322 
1323     Are_same_points_2 same_points;
1324 
1325     Point_2 pu =
1326       (same_points(p, u.source_site()) ? u.target_site() : u.source_site())
1327        .point();
1328     Point_2 pv =
1329       (same_points(p, v.source_site()) ? v.target_site() : v.source_site())
1330        .point();
1331     Point_2 pt =
1332       (same_points(p, t.source_site()) ? t.target_site() : t.source_site())
1333        .point();
1334 
1335     Orientation oupt = CGAL::orientation(pu, p.point(), pt);
1336     Orientation otpv = CGAL::orientation(pt, p.point(), pv);
1337 
1338     return (oupt == LEFT_TURN) && (otpv == LEFT_TURN);
1339   }
1340 
1341   static
1342   Boolean
are_in_same_open_halfspace_ofBasic_predicates_C21343   are_in_same_open_halfspace_of(
1344       const Site_2 & p, const Site_2 & q, const Site_2 & r)
1345   {
1346     CGAL_precondition(p.is_point() && q.is_point() && r.is_segment());
1347     Line_2 lseg = compute_supporting_line(r.supporting_site());
1348     Oriented_side os_lseg_p = oriented_side_of_line(lseg, p.point());
1349     if ( os_lseg_p == ON_ORIENTED_BOUNDARY ) {
1350       return false;
1351     }
1352     Oriented_side os_lseg_q = oriented_side_of_line(lseg, q.point());
1353     return os_lseg_p == os_lseg_q;
1354   }
1355 
1356   inline
1357   static
horseg_y_coordBasic_predicates_C21358   RT horseg_y_coord(const Site_2 & s) {
1359     CGAL_assertion(s.is_segment());
1360     CGAL_assertion(is_site_horizontal(s));
1361     return s.supporting_site().source_site().point().y();
1362   }
1363 
1364   inline
1365   static
verseg_x_coordBasic_predicates_C21366   RT verseg_x_coord(const Site_2 & s) {
1367     CGAL_assertion(s.is_segment());
1368     CGAL_assertion(is_site_vertical(s));
1369     return s.supporting_site().source_site().point().x();
1370   }
1371 
1372   inline
1373   static
hvseg_coordBasic_predicates_C21374   RT hvseg_coord(const Site_2 & s, const bool is_hor) {
1375     CGAL_assertion(s.is_segment());
1376     CGAL_assertion(is_site_horizontal(s) == is_hor);
1377     CGAL_assertion(is_site_vertical(s) == (! is_hor));
1378     return is_hor ? s.supporting_site().source_site().point().y() :
1379                     s.supporting_site().source_site().point().x() ;
1380   }
1381 
1382   inline static
check_if_exactBasic_predicates_C21383   bool check_if_exact(const Site_2& , unsigned int ,
1384                       const Tag_false&)
1385   {
1386     return true;
1387   }
1388 
1389   inline static
check_if_exactBasic_predicates_C21390   bool check_if_exact(const Site_2& s, unsigned int i,
1391                       const Tag_true&)
1392   {
1393     return s.is_input(i);
1394   }
1395 
1396   // determines if the segment s is on the positive halfspace as
1397   // defined by the supporting line of the segment supp; the line l
1398   // is supposed to be the supporting line of the segment supp and we
1399   // pass it so that we do not have to recompute it
1400   static bool
is_on_positive_halfspaceBasic_predicates_C21401   is_on_positive_halfspace(const Site_2& supp,
1402                            const Site_2& s, const Line_2& l)
1403   {
1404     CGAL_precondition( supp.is_segment() && s.is_segment() );
1405     Are_same_points_2 same_points;
1406     Are_same_segments_2 same_segments;
1407 
1408     if ( same_segments(supp.supporting_site(),
1409                        s.supporting_site()) ) {
1410       return false;
1411     }
1412 
1413     if ( same_points(supp.source_site(), s.source_site()) ||
1414          same_points(supp.target_site(), s.source_site()) ) {
1415       return oriented_side_of_line(l, s.target()) == ON_POSITIVE_SIDE;
1416     }
1417 
1418     if ( same_points(supp.source_site(), s.target_site()) ||
1419          same_points(supp.target_site(), s.target_site()) ) {
1420       return oriented_side_of_line(l, s.source()) == ON_POSITIVE_SIDE;
1421     }
1422 
1423     ITag itag;
1424 
1425     if ( !check_if_exact(s, 0, itag) &&
1426          same_segments(supp.supporting_site(),
1427                        s.crossing_site(0)) ) {
1428       return oriented_side_of_line(l, s.target()) == ON_POSITIVE_SIDE;
1429     }
1430 
1431     if ( !check_if_exact(s, 1, itag) &&
1432          same_segments(supp.supporting_site(),
1433                        s.crossing_site(1)) ) {
1434       return oriented_side_of_line(l, s.source()) == ON_POSITIVE_SIDE;
1435     }
1436 
1437     return Base::is_on_positive_halfspace(l, s.segment());
1438   }
1439 
1440   // The bearing of a line l is defined as a number in {0, 1, ..., 7},
1441   // encoding each permissible ordered pair (sign(l.a), sign(l.b)).
1442   // There are 8 permissible pairs; pair (ZERO, ZERO) is not allowed.
1443   // A line with pair (NEG, POS) has bearing 0 and as this line rotates
1444   // counterclockwise we get consecutive bearings. Observe that
1445   // axis-parallel lines have odd bearings (1, 3, 5, 7).
1446 
1447   inline static Bearing
bearingBasic_predicates_C21448   bearing(const Line_2 & l) {
1449     const Sign sa = CGAL::sign(l.a());
1450     const Sign sb = CGAL::sign(l.b());
1451     if (sa == NEGATIVE) {
1452       return 1-sb;
1453     } else if (sa == ZERO) {
1454       return (sb == NEGATIVE) ? 3 : 7;
1455     } else { // sa == POSITIVE
1456       return 5+sb;
1457     }
1458   }
1459 
1460 private:
1461   inline static bool
have_same_bearingBasic_predicates_C21462   have_same_bearing(const Line_2 & l1, const Line_2 & l2) {
1463     return (CGAL::sign(l1.a()) == CGAL::sign(l2.a())) &&
1464            (CGAL::sign(l1.b()) == CGAL::sign(l2.b()))    ;
1465   }
1466 
1467   inline static bool
have_opposite_bearingBasic_predicates_C21468   have_opposite_bearing(const Line_2 & l1, const Line_2 & l2) {
1469     return (CGAL::sign(l1.a()) == -CGAL::sign(l2.a())) &&
1470            (CGAL::sign(l1.b()) == -CGAL::sign(l2.b()))    ;
1471   }
1472 
1473   inline static bool
bearing_outsideBasic_predicates_C21474   bearing_outside(
1475       const Bearing bprev,
1476       const Bearing br,
1477       const Bearing bnext)
1478   {
1479     CGAL_assertion(bprev <= 7);
1480     CGAL_assertion(br <= 7);
1481     CGAL_assertion(bnext <= 7);
1482     CGAL_assertion(bprev != bnext);
1483     CGAL_assertion(bprev != br);
1484     CGAL_assertion(bnext != br);
1485     Bearing i(bprev);
1486     while (true) {
1487       i = (i + 1) % 8;
1488       if (i == bnext) {
1489         return true;
1490       } else if (i == br) {
1491         return false;
1492       }
1493     }
1494   }
1495 
1496 public:
1497   // absolute bearing difference between low and high
1498   inline static unsigned int
bearing_diffBasic_predicates_C21499   bearing_diff(const Bearing low, const Bearing high)
1500   {
1501     CGAL_assertion(low <= 7);
1502     CGAL_assertion(high <= 7);
1503     return high > low ? high - low : 8 + high - low;
1504   }
1505 
1506   // Orient the segments p, q, r so that they go counterclockwise
1507   // around the Linf square. The orientations are saved as lines
1508   // l[0], l[1], l[2] in the l array.
1509   static void
orient_lines_linfBasic_predicates_C21510   orient_lines_linf(const Site_2& p, const Site_2& q, const Site_2& r,
1511       Line_2 l[])
1512   {
1513     CGAL_precondition( p.is_segment() && q.is_segment() &&
1514                        r.is_segment() );
1515 
1516     l[0] = compute_supporting_line(p.supporting_site());
1517     l[1] = compute_supporting_line(q.supporting_site());
1518     l[2] = compute_supporting_line(r.supporting_site());
1519 
1520     bool is_oriented[3] = {false, false, false};
1521 
1522     if ( is_on_positive_halfspace(p, q, l[0]) ||
1523          is_on_positive_halfspace(p, r, l[0]) ) {
1524       is_oriented[0] = true;
1525     } else {
1526       l[0] = opposite_line(l[0]);
1527       if ( is_on_positive_halfspace(p, q, l[0]) ||
1528            is_on_positive_halfspace(p, r, l[0]) ) {
1529         is_oriented[0] = true;
1530       } else {
1531         l[0] = opposite_line(l[0]);
1532       }
1533     }
1534 
1535     if ( is_on_positive_halfspace(q, p, l[1]) ||
1536          is_on_positive_halfspace(q, r, l[1]) ) {
1537       is_oriented[1] = true;
1538     } else {
1539        l[1] = opposite_line(l[1]);
1540       if ( is_on_positive_halfspace(q, p, l[1]) ||
1541            is_on_positive_halfspace(q, r, l[1]) ) {
1542         is_oriented[1] = true;
1543       } else {
1544         l[1] = opposite_line(l[1]);
1545       }
1546     }
1547 
1548     if ( is_on_positive_halfspace(r, p, l[2]) ||
1549          is_on_positive_halfspace(r, q, l[2]) ) {
1550       is_oriented[2] = true;
1551     } else {
1552       l[2] = opposite_line(l[2]);
1553       if ( is_on_positive_halfspace(r, p, l[2]) ||
1554            is_on_positive_halfspace(r, q, l[2]) ) {
1555         is_oriented[2] = true;
1556       } else {
1557         l[2] = opposite_line(l[2]);
1558       }
1559     }
1560 
1561     if ( is_oriented[0] && is_oriented[1] && is_oriented[2] ) {
1562       return;
1563     }
1564 
1565     int i_no(-1);
1566     for (int i = 0; i < 3; i++) {
1567       if ( !is_oriented[i] ) {
1568         i_no = i;
1569         CGAL_assertion( is_oriented[(i+1)%3] && is_oriented[(i+2)%3] );
1570         break;
1571       }
1572     }
1573 
1574     CGAL_assertion( i_no != -1 );
1575 
1576     CGAL_SDG_DEBUG( std::cout << "orient_lines_linf i_no: "
1577         << p << ' ' << q << ' ' << ' ' << r << " i_no=" << i_no
1578         << std::endl;);
1579 
1580     for (int j = 1; j < 3; j++) {
1581       if (have_same_bearing(l[i_no], l[(i_no+j)%3])) {
1582         l[i_no] = opposite_line(l[i_no]);
1583         is_oriented[i_no] = true;
1584       } else if (have_opposite_bearing(l[i_no], l[(i_no+j)%3])) {
1585         is_oriented[i_no] = true;
1586       }
1587       if (is_oriented[i_no]) break;
1588     }
1589 
1590     if (is_oriented[i_no]) {
1591       return;
1592     }
1593 
1594     CGAL_SDG_DEBUG( std::cout << "orient_lines_linf lonely bearing: "
1595         << p << ' ' << q << ' ' << ' ' << r << " i_no=" << i_no
1596         << std::endl;);
1597     const unsigned int iprev = (i_no+2)%3;
1598     const unsigned int inext = (i_no+1)%3;
1599     const Bearing bprev = bearing(l[iprev]);
1600     const Bearing bnext = bearing(l[inext]);
1601     CGAL_SDG_DEBUG( std::cout << "orient_lines_linf"
1602         << " bprev=" << bprev << " bnext=" << bnext
1603         << std::endl;);
1604     CGAL_assertion(bprev != bnext);
1605     CGAL_assertion_code( const Bearing diffbear = (bnext - bprev)%8 );
1606     CGAL_assertion(diffbear != 1);
1607     CGAL_assertion(diffbear != 7);
1608     const Bearing br = bearing(l[i_no]);
1609     if ( bearing_outside(bprev, br, bnext) ) {
1610       CGAL_assertion_code( const Bearing bropp = (br+4)%8 );
1611       CGAL_assertion( ! bearing_outside(bprev, bropp, bnext) );
1612       l[i_no] = opposite_line(l[i_no]);
1613       is_oriented[i_no] = true;
1614     } else {
1615       CGAL_assertion( ! bearing_outside(bprev, br, bnext) );
1616       const Bearing bropp = (br+4)%8;
1617       if ( ! bearing_outside(bprev, bropp, bnext) ) {
1618         CGAL_assertion( diffbear == 6 );
1619         CGAL_assertion( br % 2 == 1 ); // undecided is axis-parallel
1620         const Site_2 & sprev = iprev == 0 ? p : (iprev == 1? q : r);
1621         Bearing brcorrect = (bprev+5)%8;
1622         if (Base::is_on_positive_halfspace(l[inext], sprev.segment())) {
1623           brcorrect = (bprev+1)%8;
1624         } else {
1625           CGAL_assertion_code(
1626               const Site_2 & snext = inext == 0 ? p : (inext == 1? q : r) );
1627           CGAL_assertion(
1628               Base::is_on_positive_halfspace(l[iprev], snext.segment()) );
1629         }
1630         if (brcorrect == bropp) {
1631           l[i_no] = opposite_line(l[i_no]);
1632         } else {
1633           CGAL_assertion(brcorrect == br);
1634         }
1635       }
1636       is_oriented[i_no] = true;
1637     }
1638     CGAL_assertion(is_oriented[i_no]);
1639   }
1640 
1641   inline static bool
are_parallel_linesBasic_predicates_C21642   are_parallel_lines(const Line_2 & lp, const Line_2 & lq) {
1643     return lp.a() * lq.b() == lq.a() * lp.b();
1644   }
1645 
1646   inline static Direction_2
directionBasic_predicates_C21647   direction(const Line_2 & l) {
1648     return Direction_2(l.b(), -l.a());
1649   }
1650 
1651   // compute correct bisector direction of two consecutive lines,
1652   // that are already correctly oriented from sites with orient_lines_linf
1653   static Direction_2
dir_from_linesBasic_predicates_C21654   dir_from_lines(const Line_2& lp, const Line_2& lq)
1655   {
1656     Bisector_Linf_Type linf_bisect_direction;
1657     const unsigned int bdiff = bearing_diff(bearing(lp), bearing(lq));
1658     CGAL_assertion(bdiff != 0);
1659     if (bdiff < 4) {
1660       return linf_bisect_direction(direction(lq), -direction(lp));
1661     } else if (bdiff > 4) {
1662       return linf_bisect_direction(direction(lp), -direction(lq));
1663     } else {
1664       CGAL_assertion(bdiff == 4);
1665       const Sign sgn = CGAL::sign(lp.a() * lq.b() - lq.a() * lp.b());
1666       if (sgn == POSITIVE) {
1667         return linf_bisect_direction(direction(lq), -direction(lp));
1668       } else {
1669         CGAL_assertion(sgn == NEGATIVE);
1670         return linf_bisect_direction(direction(lp), -direction(lq));
1671       }
1672     }
1673   }
1674 
1675   // Compute the bisecting line between sites p and q, such that
1676   // the sites are correctly oriented according to lines lp and lq,
1677   // respectively
1678   static Line_2
bisector_linf_lineBasic_predicates_C21679   bisector_linf_line(const Site_2& p, const Site_2& q,
1680       const Line_2 & lp, const Line_2 & lq)
1681   {
1682     if (are_parallel_lines(lp, lq)) {
1683       return parallel_bis(lp, lq);
1684     } else {
1685       return bisector_linf_line_nonpar(p, q, lp, lq);
1686     }
1687   }
1688 
1689   inline static Line_2
bisector_linf_line_nonparBasic_predicates_C21690   bisector_linf_line_nonpar(const Site_2& p, const Site_2& q,
1691       const Line_2 & lp, const Line_2 & lq)
1692   {
1693     const bool is_psrc_q = is_endpoint_of(p.source_site(), q);
1694     const bool is_ptrg_q = is_endpoint_of(p.target_site(), q);
1695     const bool have_common_pq = is_psrc_q || is_ptrg_q;
1696     Homogeneous_point_2 xpq;
1697     if (have_common_pq) {
1698       xpq = is_psrc_q ? p.source() : p.target();
1699     } else {
1700       RT hx, hy, hw;
1701       compute_intersection_of_lines(lp, lq, hx, hy, hw);
1702       CGAL_SDG_DEBUG( std::cout << "debug xpq hom="
1703         << hx << ' ' << hy << ' ' << hw << std::endl; );
1704       xpq = Homogeneous_point_2(hx, hy, hw);
1705     }
1706     const Direction_2 dirbpq = dir_from_lines(lp, lq);
1707     CGAL_SDG_DEBUG( std::cout << "debug xpq hom="
1708         << xpq.hx() << ' ' << xpq.hy() << ' ' << xpq.hw()
1709         << " dirbpq=" << dirbpq << std::endl; );
1710     return compute_line_dir(xpq, dirbpq);
1711   }
1712 
1713   // check whether the point p is an endpoint of the segment s
1714   inline static
is_endpoint_ofBasic_predicates_C21715   bool is_endpoint_of(const Site_2& p, const Site_2& s)
1716   {
1717     CGAL_precondition( p.is_point() && s.is_segment() );
1718     Are_same_points_2 same_points;
1719     return ( same_points(p, s.source_site()) ||
1720              same_points(p, s.target_site())   );
1721   }
1722 
1723   // Orient the segment s and return result as a line.
1724   // Site p is a point which is an endpoint of s and
1725   // p_before_s is true if p is just before s in the
1726   // Voronoi vertex.
1727   static Line_2
orient_line_endpBasic_predicates_C21728   orient_line_endp(const Site_2& p, const Site_2& s, const bool p_before_s)
1729   {
1730     return compute_line_from_to(
1731         p_before_s ? p.point() : other_site(p, s).point(),
1732         p_before_s ? other_site(p, s).point() : p.point()  );
1733   }
1734 
1735   // Orient the segment s and return result as a line.
1736   // Site p is a point which is not an endpoint of s.
1737   // Site p must not be on the line defined by s.
1738   static Line_2
orient_line_nonendpBasic_predicates_C21739   orient_line_nonendp(const Site_2& p, const Site_2& s)
1740   {
1741     Line_2 lseg = compute_supporting_line(s.supporting_site());
1742     Oriented_side os = oriented_side_of_line(lseg, p.point());
1743     if (os != ON_POSITIVE_SIDE) {
1744       CGAL_assertion( os == ON_NEGATIVE_SIDE );
1745       lseg = opposite_line(lseg);
1746     }
1747     return lseg;
1748   }
1749 
1750   // given that point site p is an endpoint of segment seg,
1751   // return the other endpoint of seg (as a site)
1752   inline static
other_siteBasic_predicates_C21753   const Site_2 other_site(const Site_2& p, const Site_2& seg)
1754   {
1755     CGAL_precondition( p.is_point() && seg.is_segment() );
1756     Are_same_points_2 same_points;
1757     if ( same_points(p, seg.source_site()) ){
1758       return seg.target_site();
1759     } else {
1760       CGAL_assertion(same_points(p, seg.target_site()));
1761       return seg.source_site();
1762     }
1763   }
1764 
1765   // Given corner with bearing cb (0: bottom right, 2: top right,
1766   // 4: top left, 6: bottom left) and another point p, return the
1767   // (smallest) square having that corner and passing through p.
1768   static
center_from_corner_and_ptBasic_predicates_C21769   Point_2 center_from_corner_and_pt(
1770       const Point_2 & corner, const Bearing cb, const Point_2 & p)
1771   {
1772     CGAL_precondition(cb % 2 == 0);
1773     const FT absdifx = CGAL::abs(corner.x() - p.x());
1774     const FT absdify = CGAL::abs(corner.y() - p.y());
1775     const Comparison_result cmp = CGAL::compare(absdifx, absdify);
1776     if (cmp == SMALLER) {
1777       const FT ox = corner.x() + FT((cb < 3) ? -1: +1)*absdify/FT(2);
1778       const FT oy = (corner.y() + p.y())/FT(2);
1779       return Point_2(ox, oy);
1780     } else {
1781       const FT ox = (corner.x() + p.x())/FT(2);
1782       const FT oy = corner.y() + FT((cb % 6 == 0) ? +1: -1)*absdifx/FT(2);
1783       return Point_2(ox, oy);
1784     }
1785   }
1786 
1787   inline
1788   static
center_from_opposite_cornersBasic_predicates_C21789   Point_2 center_from_opposite_corners(
1790       const Point_2 & c, const Point_2 & d)
1791   {
1792     return Point_2(c.x() + d.x(), c.y() + d.y(), RT(2));
1793   }
1794 
1795   inline
1796   static
center_from_same_side_cornersBasic_predicates_C21797   Point_2 center_from_same_side_corners(
1798       const Point_2 & c, const Point_2 & d, const Bearing bside)
1799   {
1800     CGAL_precondition(bside % 2 == 1);
1801     const FT ax = (bside % 4 == 1) ?
1802       FT(RT(2)*c.x() + c.y() - d.y()) : FT(c.x() + d.x());
1803     const FT ay = (bside % 4 == 1) ?
1804       FT(c.y() + d.y()) : FT(RT(2)*c.y() + d.x() - c.y());
1805     return Point_2(ax, ay, RT(2));
1806   }
1807 
1808   inline
1809   static
points_inside_touching_sides_vBasic_predicates_C21810   bool points_inside_touching_sides_v(
1811       const Line_2 & ls, const Site_2 & pt_site,
1812       const Site_2 & other_s, const Site_2 & t, const Point_2 & v)
1813   {
1814     CGAL_precondition(pt_site.is_point());
1815     CGAL_precondition(t.is_point());
1816     CGAL_USE(other_s);
1817     CGAL_SDG_DEBUG(std::cout << "debug points_inside_touching_sides_v "
1818         << "ls: " << ls.a() << ' ' << ls.b() << ' ' <<  ls.c()
1819         << " pt_site=" << pt_site << " other_s=" << other_s
1820         << " t=" << t << " v=" << v << std::endl;);
1821     const Point_2 corner =
1822       compute_linf_projection_nonhom(ls, v);
1823     const Line_2 ltest = has_positive_slope(ls) ?
1824       compute_pos_45_line_at(v): compute_neg_45_line_at(v);
1825     CGAL_assertion(
1826         oriented_side_of_line(ltest, v) == ON_ORIENTED_BOUNDARY);
1827     const Oriented_side ost = oriented_side_of_line(ltest, t.point());
1828     const Oriented_side osx = oriented_side_of_line(ltest, corner);
1829     CGAL_SDG_DEBUG(std::cout << "debug points_inside_touching_sides_v"
1830         << " ltest: " << ltest.a() << ' ' << ltest.b() << ' ' <<  ltest.c()
1831         << " v=" << v << " ost=" << ost
1832         << " corner=" << corner << " osx=" << osx << std::endl;);
1833     if (ost == osx) {
1834       const Point_2 & p = pt_site.point();
1835       const Oriented_side osp = oriented_side_of_line(ltest, p);
1836       if (ost == osp) {
1837         // +-pi/2 slope line through corner and v
1838         const Line_2 lcv = has_positive_slope(ls) ?
1839           compute_neg_45_line_at(v): compute_pos_45_line_at(v);
1840         const Oriented_side oslt = oriented_side_of_line(lcv, t.point());
1841         const Oriented_side oslp = oriented_side_of_line(lcv, p);
1842         if (oslt != oslp) {
1843           return true;
1844         }
1845       }
1846     }
1847     return false;
1848   }
1849 
1850   inline
1851   static
points_inside_touching_sides_vBasic_predicates_C21852   bool points_inside_touching_sides_v(
1853       const Site_2 & s, const Site_2 & pt_site,
1854       const Site_2 & other_s, const Site_2 & t, const Point_2 & v)
1855   {
1856     CGAL_precondition(t.is_point());
1857     CGAL_precondition(pt_site.is_point());
1858     CGAL_precondition(s.is_segment());
1859     CGAL_SDG_DEBUG(std::cout << "debug points_inside_touching_sides_v "
1860         << "s=" << s
1861         << " pt_site=" << pt_site << " other_s=" << other_s
1862         << " t=" << t << std::endl;);
1863     CGAL_assertion(! is_site_h_or_v(s));
1864     if (other_s.is_segment()) {
1865       // shortcut: when the point pt_site is on a corner of
1866       // the Linf square, because it is the endpoint of the
1867       // other site which is a segment; return false immediately
1868       if ((! is_site_h_or_v(other_s)) &&
1869           is_endpoint_of(pt_site, other_s)) {
1870         return false;
1871       }
1872     }
1873     const Line_2 ls = compute_supporting_line(s.supporting_site());
1874     return points_inside_touching_sides_v(ls, pt_site, other_s, t, v);
1875   }
1876 
1877   // Check if point p's Voronoi area has zero area, because point p
1878   // is sandwiched between two segments with agreeing slope and
1879   // direction. As a result: vpqr and vqps coincide in edge conflicts.
1880   inline static
1881   bool
zero_voronoi_areaBasic_predicates_C21882   zero_voronoi_area(const Site_2& p, const Site_2& r, const Site_2& s)
1883   {
1884     Are_same_points_2 same_points;
1885     if (p.is_segment()) { return false; }
1886     if (r.is_point() || s.is_point()) { return false; }
1887     const bool is_p_rsrc = same_points(p, r.source_site());
1888     const bool is_p_rtrg =
1889       (! is_p_rsrc) && same_points(p, r.target_site());
1890     const bool is_p_endp_of_r = is_p_rsrc || is_p_rtrg;
1891     if (is_p_endp_of_r) {
1892       const bool is_p_ssrc = same_points(p, s.source_site());
1893       const bool is_p_strg =
1894         (! is_p_ssrc) && same_points(p, s.target_site());
1895       const bool is_p_endp_of_s = is_p_ssrc || is_p_strg;
1896       if (is_p_endp_of_s) {
1897         if (is_site_horizontal(r) && is_site_horizontal(s)) { return true; }
1898         if (is_site_vertical(r) && is_site_vertical(s)) { return true; }
1899         if ((! is_site_h_or_v(r)) && (! is_site_h_or_v(s))) {
1900           const bool pos_r = has_positive_slope(r);
1901           const bool pos_s = has_positive_slope(s);
1902           if (pos_r == pos_s) {
1903             const Line_2 l = pos_r ? compute_neg_45_line_at(p.point()) :
1904               compute_pos_45_line_at(p.point()) ;
1905             const Oriented_side osr =
1906               oriented_side_of_line(l, is_p_rsrc ? r.target() : r.source());
1907             const Oriented_side oss =
1908               oriented_side_of_line(l, is_p_ssrc ? s.target() : s.source());
1909             if (osr != oss) { return true; }
1910           }
1911         }
1912       }
1913     }
1914     return false;
1915   }
1916 
1917   // Check if point p is on the line of the axis-parallel segment s.
1918   // It returns true only for an axis-parallel segment s argument.
is_on_hv_seg_lineBasic_predicates_C21919   static inline bool is_on_hv_seg_line(const Site_2 & p, const Site_2 & s)
1920   {
1921     CGAL_precondition(p.is_point());
1922     CGAL_precondition(s.is_segment());
1923     Compare_x_2_Sites_Type scmpx;
1924     Compare_y_2_Sites_Type scmpy;
1925     const bool is_hor = is_site_horizontal(s);
1926     const bool is_ver = (! is_hor) && is_site_vertical(s);
1927     if (is_hor || is_ver) {
1928       return ( (is_hor) ?
1929                scmpy(p, s.source_site()) : scmpx(p, s.source_site()) )
1930           == EQUAL;
1931     } else {
1932       return false;
1933     }
1934   }
1935 
1936 
1937 }; // end of struct Basic_predicates_C2
1938 
1939 
1940 } //namespace SegmentDelaunayGraphLinf_2
1941 
1942 } //namespace CGAL
1943 
1944 
1945 #endif // CGAL_SEGMENT_DELAUNAY_GRAPH_LINF_2_BASIC_PREDICATES_C2_H
1946