1 // Copyright (c) 2010-2018   INRIA Sophia Antipolis, INRIA Nancy (France).
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/Hyperbolic_triangulation_2/include/CGAL/Hyperbolic_Delaunay_triangulation_CK_traits_2.h $
7 // $Id: Hyperbolic_Delaunay_triangulation_CK_traits_2.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 // Author(s)     : Mikhail Bogdanov
11 //                 Monique Teillaud <Monique.Teillaud@inria.fr>
12 
13 #ifndef CGAL_HYPERBOLIC_DELAUNAY_TRIANGULATION_CK_TRAITS_2_H
14 #define CGAL_HYPERBOLIC_DELAUNAY_TRIANGULATION_CK_TRAITS_2_H
15 
16 #include <CGAL/license/Hyperbolic_triangulation_2.h>
17 
18 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
19 
20 #include <CGAL/Algebraic_kernel_for_circles_2_2.h>
21 #include <CGAL/Circular_kernel_2/Intersection_traits.h>
22 #include <CGAL/Circular_kernel_2.h>
23 #include <CGAL/triangulation_assertions.h>
24 #include <CGAL/determinant.h>
25 
26 #include <boost/tuple/tuple.hpp>
27 #include <boost/variant.hpp>
28 
29 #include <CGAL/internal/Hyperbolic_Delaunay_triangulation_traits_2_functions.h>
30 
31 namespace CGAL {
32 
33 namespace internal {
34 
35   template <typename Traits>
36   class Construct_hyperbolic_circumcenter_CK_2
37   {
38     typedef typename Traits::Hyperbolic_Voronoi_point_2     Hyperbolic_Voronoi_point_2;
39     typedef typename Traits::Hyperbolic_point_2             Hyperbolic_point_2;
40     typedef typename Traits::Euclidean_circle_or_line_2     Euclidean_circle_or_line_2;
41     typedef typename Traits::Euclidean_line_2               Euclidean_line_2;
42     typedef typename Traits::Circle_2                       Circle_2;
43     typedef typename Traits::FT                             FT;
44     typedef typename Traits::Circular_arc_point_2           Circular_arc_point_2;
45 
46   public:
_gt(gt)47     Construct_hyperbolic_circumcenter_CK_2(const Traits gt = Traits()): _gt(gt) {}
48 
operator()49     Hyperbolic_Voronoi_point_2 operator()(const Hyperbolic_point_2& p,
50                                           const Hyperbolic_point_2& q,
51                                           const Hyperbolic_point_2& r) const
52     {
53 
54       Construct_circle_or_line_supporting_bisector<Traits> cclsb(_gt);
55 
56       Hyperbolic_point_2 po(CGAL::ORIGIN);
57       Circle_2 l_inf(po, FT(1));
58 
59       if( _gt.compare_distance_2_object()(po, p, q) == EQUAL &&
60           _gt.compare_distance_2_object()(po, p, r) == EQUAL)
61         return po;
62 
63       // now supporting objects cannot both be Euclidean lines
64 
65       Euclidean_circle_or_line_2 bis_pq = cclsb(p, q);
66       Euclidean_circle_or_line_2 bis_qr = cclsb(q, r);
67 
68       std::pair<Circular_arc_point_2, unsigned> pair;
69       Euclidean_line_2* l;
70       Circle_2* c;
71 
72       if(Circle_2* c_pq = boost::get<Circle_2>(&bis_pq))
73       {
74         if(Circle_2* c_qr = boost::get<Circle_2>(&bis_qr))
75         {
76           typedef typename CK2_Intersection_traits<Traits, Circle_2, Circle_2>::type Intersection_result;
77           std::vector< Intersection_result > inters;
78           intersection(*c_pq, *c_qr, std::back_inserter(inters));
79 
80           CGAL_triangulation_assertion(assign(pair, inters[0]));
81           if(pair.second == 1)
82           {
83             if(_gt.has_on_bounded_side_2_object()(l_inf, pair.first))
84               return pair.first;
85 
86             CGAL_triangulation_assertion(assign(pair, inters[1]));
87             return pair.first;
88           }
89           return pair.first;
90         }
91 
92         // here bis_qr is a line
93         l = boost::get<Euclidean_line_2>(&bis_qr);
94         c = c_pq;
95       }
96       else
97       {
98         // here bis_pq is a line, and bis_qr is necessarily a circle
99         l = boost::get<Euclidean_line_2>(&bis_pq);
100         c = boost::get<Circle_2>(&bis_qr);
101       }
102 
103       typedef typename CK2_Intersection_traits<Traits, Euclidean_line_2, Circle_2>::type Intersection_result;
104       std::vector< Intersection_result > inters;
105       intersection(*l, *c, std::back_inserter(inters));
106 
107       CGAL_triangulation_assertion(assign(pair,inters[0]));
108       if(pair.second == 1)
109       {
110         if(_gt.has_on_bounded_side_2_object()(l_inf, pair.first))
111           return pair.first;
112 
113         CGAL_triangulation_assertion(assign(pair, inters[1]));
114         return pair.first;
115       }
116       return pair.first;
117     }
118 
119   private:
120     const Traits& _gt;
121   }; // end Construct_hyperbolic_circumcenter_2
122 
123 
124   template <typename Traits>
125   class Construct_hyperbolic_bisector_CK_2
126   {
127     typedef typename Traits::Hyperbolic_segment_2         Hyperbolic_segment_2;
128     typedef typename Traits::Hyperbolic_point_2           Hyperbolic_point_2;
129     typedef typename Traits::Circle_2                     Circle_2;
130     typedef typename Traits::Euclidean_line_2             Euclidean_line_2;
131     typedef typename Traits::Euclidean_circle_or_line_2   Euclidean_circle_or_line_2;
132     typedef typename Traits::Circular_arc_2               Circular_arc_2;
133     typedef typename Traits::Line_arc_2                   Line_arc_2;
134     typedef typename Traits::Circular_arc_point_2         Circular_arc_point_2;
135     typedef typename Traits::FT                           FT;
136 
137   public:
_gt(gt)138     Construct_hyperbolic_bisector_CK_2(const Traits& gt = Traits()) : _gt(gt) {}
139 
140     // constructs a hyperbolic line
operator()141     Hyperbolic_segment_2 operator()(const Hyperbolic_point_2& p,
142                                     const Hyperbolic_point_2& q) const
143     {
144 
145       Construct_circle_or_line_supporting_bisector<Traits>  cclsb(_gt);
146 
147       Hyperbolic_point_2 po(CGAL::ORIGIN);
148       Circle_2 l_inf = Circle_2(po, FT(1));
149 
150       if(_gt.compare_distance_2_object()(po, p, q) == EQUAL)
151       {
152         Euclidean_line_2 l = _gt.construct_Euclidean_bisector_2_object()(p, q);
153         if(_gt.less_y_2_object()(p, q))
154           return Line_arc_2(l, l_inf, false, l_inf, true);
155 
156         return Line_arc_2(l, l_inf, true, l_inf, false);
157       }
158 
159       Euclidean_circle_or_line_2 bis_pq = cclsb(p, q);
160       Circle_2* c = boost::get<Circle_2>(&bis_pq);
161 
162       if(_gt.less_y_2_object()(po, c->center()))
163         return Circular_arc_2(*c, l_inf, true, l_inf, false);
164       else if(_gt.less_y_2_object()(c->center(), po))
165         return Circular_arc_2(*c, l_inf, false, l_inf, true);
166 
167       // the center of the circle is on the x-axis
168       if(_gt.less_x_2_object()(po, c->center())) { return Circular_arc_2(*c, l_inf, true, l_inf, false); }
169       return Circular_arc_2(*c, l_inf, false, l_inf, true);
170     }
171 
172     // constructs the hyperbolic bisector of segment [p, q] limited by
173     // circumcenter(p, q, r) on one side
174     // and circumcenter(p, s, q) on the other side
operator()175     Hyperbolic_segment_2 operator()(const Hyperbolic_point_2& p,
176                                     const Hyperbolic_point_2& q,
177                                     const Hyperbolic_point_2& r,
178                                     const Hyperbolic_point_2& s) const
179     {
180       CGAL_triangulation_precondition((_gt.orientation_2_object()(p, q, r) == ON_POSITIVE_SIDE) &&
181                                       (_gt.orientation_2_object()(p, s, q) == ON_POSITIVE_SIDE));
182       CGAL_triangulation_precondition((_gt.side_of_oriented_circle_2_object()(p, q, r,s) == ON_NEGATIVE_SIDE) &&
183                                       (_gt.side_of_oriented_circle_2_object()(p, s, q, r) == ON_NEGATIVE_SIDE));
184 
185       Construct_circle_or_line_supporting_bisector<Traits>  cclsb(_gt);
186       Construct_hyperbolic_circumcenter_CK_2<Traits>        chc(_gt);
187       Hyperbolic_point_2 po(CGAL::ORIGIN);
188 
189       // TODO MT this is non-optimal...
190       // the bisector is already computed here
191       // and it will be recomputed below
192       Circular_arc_point_2 a = chc(p, q, r);
193       Circular_arc_point_2 b = chc(p, s, q);
194 
195       if(_gt.compare_distance_2_object()(po, p, q) == EQUAL)
196       {
197         Euclidean_line_2 l = _gt.construct_Euclidean_bisector_2_object()(p, q);
198         return Line_arc_2(l, a, b);
199       }
200 
201       Euclidean_circle_or_line_2
202           bis_pq = cclsb(p, q);
203       Circle_2* c_pq = boost::get<Circle_2>(&bis_pq);
204 
205       if(_gt.compare_distance_2_object()(po, p, q) == POSITIVE)
206       {
207         // then p is inside the supporting circle
208         return Circular_arc_2(*c_pq, b, a);
209       }
210 
211       return Circular_arc_2(*c_pq, a, b);
212     }
213 
214     // constructs the hyperbolic bisector of segment [p, q]
215     // limited by hyperbolic circumcenter(p, q, r) on one side
216     // and going to the infinite line on the other side
operator()217     Hyperbolic_segment_2 operator()(const Hyperbolic_point_2& p,
218                                     const Hyperbolic_point_2& q,
219                                     const Hyperbolic_point_2& r) const
220     {
221       CGAL_triangulation_precondition(_gt.orientation_2_object()(p, q, r) == POSITIVE);
222 
223       Construct_circle_or_line_supporting_bisector<Traits>  cclsb(_gt);
224       Construct_hyperbolic_circumcenter_CK_2<Traits>        chc(_gt);
225 
226       Hyperbolic_point_2 po(CGAL::ORIGIN);
227       Circle_2 l_inf(po, FT(1));
228 
229       // TODO MT this is non-optimal...
230       // the bisector is computed (at least) twice
231       Circular_arc_point_2 a = chc(p, q, r);
232 
233       if(_gt.compare_distance_2_object()(po, p, q) == EQUAL)
234       {
235         Euclidean_line_2 bis_pq = _gt.construct_Euclidean_bisector_2_object()(p, q);
236         typedef typename
237         CK2_Intersection_traits<Traits, Euclidean_line_2, Circle_2>::type
238             Intersection_result;
239         std::vector< Intersection_result > inters;
240         intersection(bis_pq, l_inf, std::back_inserter(inters));
241         std::pair<Circular_arc_point_2, unsigned> pair;
242 
243         CGAL_triangulation_assertion(assign(pair,inters[0]));
244         CGAL_triangulation_assertion(pair.second == 1);
245         if(_gt.less_y_2_object()(p, q))
246           return Line_arc_2(bis_pq,a,pair.first);
247 
248         CGAL_triangulation_assertion(assign(pair,inters[1]));
249         CGAL_triangulation_assertion(pair.second == 1);
250         return Line_arc_2(bis_pq,a,pair.first);
251       }
252 
253       Euclidean_circle_or_line_2 bis_pq = cclsb(p, q);
254       Circle_2* c_pq = boost::get<Circle_2>(&bis_pq);
255 
256       Hyperbolic_point_2 approx_a(to_double(a.x()),to_double(a.y()));
257 
258       typedef typename
259       CK2_Intersection_traits<Traits, Circle_2, Circle_2>::type            Intersection_result;
260       std::vector< Intersection_result > inters;
261       intersection(*c_pq, l_inf, std::back_inserter(inters));
262       std::pair<Circular_arc_point_2, unsigned> pair;
263 
264       CGAL_triangulation_assertion(assign(pair,inters[0]));
265       CGAL_triangulation_assertion(pair.second == 1);
266 
267       Hyperbolic_point_2 approx_pinf(to_double(pair.first.x()), to_double(pair.first.y()));
268       Hyperbolic_point_2 approx_c(to_double(c_pq->center().x()),
269                                   to_double(c_pq->center().y()));
270       if(_gt.orientation_2_object()(p, q,approx_pinf) == NEGATIVE)
271       {
272         if(_gt.orientation_2_object()(approx_c,approx_a,approx_pinf) == POSITIVE)
273           return Circular_arc_2(*c_pq, a, pair.first);
274 
275         return Circular_arc_2(*c_pq, pair.first, a);
276       }
277 
278       CGAL_triangulation_assertion(assign(pair,inters[1]));
279       if(_gt.orientation_2_object()(approx_c,approx_a,approx_pinf) == POSITIVE)
280         return Circular_arc_2(*c_pq, pair.first, a);
281 
282       return Circular_arc_2(*c_pq, a, pair.first);
283     }
284 
285   private:
286     const Traits& _gt;
287   }; // end Construct_hyperbolic_bisector_2
288 
289 
290 } // end namespace internal
291 
292 
293 
294 
295 
296 template<class R = CGAL::Circular_kernel_2<
297            CGAL::Exact_predicates_inexact_constructions_kernel,
298            CGAL::Algebraic_kernel_for_circles_2_2<
299              CGAL::Exact_predicates_inexact_constructions_kernel::RT >
300            >
301          >
302 class Hyperbolic_Delaunay_triangulation_CK_traits_2
303     : public R // R is supposed to be a model of CircularKernel2
304 {
305   typedef Hyperbolic_Delaunay_triangulation_CK_traits_2<R>    Self;
306   typedef R                                                   Base;
307 
308 public:
309   typedef typename R::FT                                      FT;
310 
311   typedef typename R::Point_2                                 Hyperbolic_point_2;
312   typedef typename R::Circle_2                                Circle_2;
313   typedef typename R::Line_2                                  Euclidean_line_2;
314   typedef boost::variant<Circle_2, Euclidean_line_2>          Euclidean_circle_or_line_2;
315 
316   typedef typename R::Circular_arc_2                          Circular_arc_2;
317   typedef typename R::Line_arc_2                              Line_arc_2;
318   typedef typename R::Circular_arc_point_2                    Circular_arc_point_2;
319   typedef Circular_arc_point_2                                Hyperbolic_Voronoi_point_2;
320   typedef typename R::Segment_2                               Euclidean_segment_2; // only used internally here
321   typedef boost::variant<Circular_arc_2, Line_arc_2>          Hyperbolic_segment_2;
322 
323   typedef typename R::Triangle_2                              Hyperbolic_triangle_2;
324 
325   typedef typename R::Compare_x_2                             Compare_x_2;
326   typedef typename R::Compare_y_2                             Compare_y_2;
327   typedef typename R::Compare_distance_2                      Compare_distance_2;
328   typedef typename R::Orientation_2                           Orientation_2;
329   typedef typename R::Side_of_oriented_circle_2               Side_of_oriented_circle_2;
330 
331   // the following types are only used internally in this traits class,
332   // so they need not be documented, and they don't need _object()
333   typedef typename R::Collinear_2                             Euclidean_collinear_2;
334   typedef typename R::Construct_bisector_2                    Construct_Euclidean_bisector_2;
335   typedef typename R::Construct_midpoint_2                    Construct_Euclidean_midpoint_2;
336   typedef typename R::Compute_squared_distance_2              Compute_squared_Euclidean_distance_2;
337   typedef typename R::Has_on_bounded_side_2                   Has_on_bounded_side_2;
338 
339   typedef typename R::Less_x_2                                Less_x_2;
340   typedef typename R::Less_y_2                                Less_y_2;
341 
342   // only kept for demo to please T2graphicsitems
343   typedef Euclidean_segment_2                                 Line_segment_2;
344   typedef Hyperbolic_segment_2                                Segment_2;
345 
346   typedef internal::Construct_hyperbolic_circumcenter_CK_2<Self>   Construct_hyperbolic_circumcenter_2;
347   typedef internal::Construct_hyperbolic_bisector_CK_2<Self>       Construct_hyperbolic_bisector_2;
348   typedef internal::Is_Delaunay_hyperbolic<Self>                Is_Delaunay_hyperbolic;
349   typedef internal::Side_of_oriented_hyperbolic_segment_2<Self> Side_of_oriented_hyperbolic_segment_2;
350   typedef typename internal::Construct_circle_or_line_supporting_bisector<Self> Construct_circle_or_line_supporting_bisector;
351   typedef internal::Construct_hyperbolic_segment_2<Self>        Construct_hyperbolic_segment_2;
352   typedef typename Base::Construct_segment_2                    Construct_segment_2;
353 
354 public:
Base(kernel)355   Hyperbolic_Delaunay_triangulation_CK_traits_2(const Base& kernel = Base()) : Base(kernel) {}
356 
357 
358   Construct_hyperbolic_segment_2
construct_hyperbolic_segment_2_object()359   construct_hyperbolic_segment_2_object() const
360   { return Construct_hyperbolic_segment_2(*this); }
361 
362   Construct_segment_2
construct_segment_2_object()363   construct_segment_2_object() const
364   { return this->Base::construct_segment_2_object(); }
365 
366   Construct_hyperbolic_circumcenter_2
construct_hyperbolic_circumcenter_2_object()367   construct_hyperbolic_circumcenter_2_object() const
368   { return Construct_hyperbolic_circumcenter_2(*this); }
369 
370   Construct_hyperbolic_bisector_2
construct_hyperbolic_bisector_2_object()371   construct_hyperbolic_bisector_2_object() const
372   { return Construct_hyperbolic_bisector_2(*this); }
373 
374   Is_Delaunay_hyperbolic
is_Delaunay_hyperbolic_2_object()375   is_Delaunay_hyperbolic_2_object() const
376   { return Is_Delaunay_hyperbolic(*this); }
377 
378   Side_of_oriented_hyperbolic_segment_2
side_of_oriented_hyperbolic_segment_2_object()379   side_of_oriented_hyperbolic_segment_2_object() const
380   { return Side_of_oriented_hyperbolic_segment_2(*this); }
381 
382   Construct_Euclidean_bisector_2
construct_Euclidean_bisector_2_object()383   construct_Euclidean_bisector_2_object() const
384   { return this->Base::construct_bisector_2_object(); }
385 
386   Construct_circle_or_line_supporting_bisector
construct_circle_or_line_supporting_bisector_2_object()387   construct_circle_or_line_supporting_bisector_2_object() const
388   { return Construct_circle_or_line_supporting_bisector(*this); }
389 
390   Euclidean_collinear_2
euclidean_collinear_2_object()391   euclidean_collinear_2_object() const
392   { return this->Base::collinear_2_object(); }
393 
394   Compute_squared_Euclidean_distance_2
compute_squared_Euclidean_distance_2_object()395   compute_squared_Euclidean_distance_2_object() const
396   { return this->Base::compute_squared_distance_2_object(); }
397 
398 };
399 
400 // Take out the code below to some separate file
401 
402 #ifdef CGAL_EXACT_PREDICATES_EXACT_CONSTRUCTIONS_KERNEL_H
403 template  <>
404 struct Triangulation_structural_filtering_traits< Hyperbolic_Delaunay_triangulation_CK_traits_2<Epeck> > {
405   typedef Tag_true Use_structural_filtering_tag;
406 };
407 #endif // CGAL_EXACT_PREDICATES_EXACT_CONSTRUCTIONS_KERNEL_H
408 
409 #ifdef CGAL_EXACT_PREDICATES_INEXACT_CONSTRUCTIONS_KERNEL_H
410 template <>
411 struct Triangulation_structural_filtering_traits< Hyperbolic_Delaunay_triangulation_CK_traits_2<Epick> > {
412   typedef Tag_true Use_structural_filtering_tag;
413 };
414 #endif // CGAL_EXACT_PREDICATES_INEXACT_CONSTRUCTIONS_KERNEL_H
415 
416 } //namespace CGAL
417 
418 #endif // CGAL_HYPERBOLIC_DELAUNAY_TRIANGULATION_CK_TRAITS_2_H
419