1 // Copyright (c) 2006,2007,2009,2010,2011 Tel-Aviv University (Israel).
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/Arrangement_on_surface_2/include/CGAL/Arr_linear_traits_2.h $
7 // $Id: Arr_linear_traits_2.h f55ef7d 2020-10-09T18:36:17+02:00 Mael Rouxel-Labbé
8 // SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
9 //
10 //
11 // Author(s): Ron Wein          <wein@post.tau.ac.il>
12 //            Waqar Khan        <wkhan@mpi-inf.mpg.de>
13 //            Efi fogel         <efifogel@gmail.com>
14 
15 #ifndef CGAL_ARR_LINEAR_TRAITS_2_H
16 #define CGAL_ARR_LINEAR_TRAITS_2_H
17 
18 #include <CGAL/license/Arrangement_on_surface_2.h>
19 
20 #include <CGAL/disable_warnings.h>
21 
22 /*! \file
23  * The traits-class for handling linear objects (lines, rays and segments)
24  * in the arrangement package.
25  */
26 
27 #include <fstream>
28 
29 #include <boost/variant.hpp>
30 
31 #include <CGAL/tags.h>
32 #include <CGAL/intersections.h>
33 #include <CGAL/Arr_tags.h>
34 #include <CGAL/Arr_enums.h>
35 #include <CGAL/Arr_geometry_traits/Segment_assertions.h>
36 
37 namespace CGAL {
38 
39 template <typename Kernel_> class Arr_linear_object_2;
40 
41 /*! \class
42  * A traits class for maintaining an arrangement of linear objects (lines,
43  * rays and segments), aoviding cascading of computations as much as possible.
44  */
45 template <typename Kernel_>
46 class Arr_linear_traits_2 : public Kernel_ {
47   friend class Arr_linear_object_2<Kernel_>;
48 
49 public:
50   typedef Kernel_                         Kernel;
51   typedef typename Kernel::FT             FT;
52 
53   typedef typename Algebraic_structure_traits<FT>::Is_exact
54                                           Has_exact_division;
55 
56   // Category tags:
57   typedef Tag_true                        Has_left_category;
58   typedef Tag_true                        Has_merge_category;
59   typedef Tag_false                       Has_do_intersect_category;
60 
61   typedef Arr_open_side_tag               Left_side_category;
62   typedef Arr_open_side_tag               Bottom_side_category;
63   typedef Arr_open_side_tag               Top_side_category;
64   typedef Arr_open_side_tag               Right_side_category;
65 
66   typedef typename Kernel::Line_2         Line_2;
67   typedef typename Kernel::Ray_2          Ray_2;
68   typedef typename Kernel::Segment_2      Segment_2;
69 
70   typedef CGAL::Segment_assertions<Arr_linear_traits_2<Kernel> >
71                                           Segment_assertions;
72 
73   /*!
74    * \class Representation of a linear with cached data.
75    */
76   class _Linear_object_cached_2 {
77   public:
78     typedef typename Kernel::Line_2                Line_2;
79     typedef typename Kernel::Ray_2                 Ray_2;
80     typedef typename Kernel::Segment_2             Segment_2;
81     typedef typename Kernel::Point_2               Point_2;
82 
83   protected:
84     Line_2    l;                // The supporting line.
85     Point_2   ps;               // The source point (if exists).
86     Point_2   pt;               // The target point (if exists).
87     bool      has_source;       // Is the source point valid
88                                 // (false for a line).
89     bool      has_target;       // Is the target point valid
90                                 // (false for a line and for a ray).
91     bool      is_right;         // Is the object directed to the right
92                                 // (for segments and rays).
93     bool      is_vert;          // Is this a vertical object.
94     bool      is_horiz;         // Is this a horizontal object.
95     bool      has_pos_slope;    // Does the supporting line has a positive
96                                 // slope (if all three flags is_vert, is_horiz
97                                 // and has_pos_slope are false, then the line
98                                 // has a negative slope).
99     bool      is_degen;         // Is the object degenerate (a single point).
100 
101   public:
102     /*! Default constructor.
103      */
_Linear_object_cached_2()104     _Linear_object_cached_2() :
105       has_source(true),
106       has_target(true),
107       is_vert(false),
108       is_horiz(false),
109       has_pos_slope(false),
110       is_degen(true)
111     {}
112 
113     /*! Constructor for segment from two points.
114      * \param p1 source point.
115      * \param p2 target point.
116      * \pre The two points must not be equal.
117      */
_Linear_object_cached_2(const Point_2 & source,const Point_2 & target)118     _Linear_object_cached_2(const Point_2& source, const Point_2& target) :
119       ps(source),
120       pt(target),
121       has_source(true),
122       has_target(true)
123     {
124       Kernel kernel;
125 
126       Comparison_result res = kernel.compare_xy_2_object()(source, target);
127       is_degen = (res == EQUAL);
128       is_right = (res == SMALLER);
129 
130       CGAL_precondition_msg(! is_degen,
131                             "Cannot construct a degenerate segment.");
132 
133       l = kernel.construct_line_2_object()(source, target);
134       is_vert = kernel.is_vertical_2_object()(l);
135       is_horiz = kernel.is_horizontal_2_object()(l);
136       has_pos_slope = _has_positive_slope();
137     }
138 
139     /*! Constructor from a segment.
140      * \param seg The segment.
141      * \pre The segment is not degenerate.
142      */
_Linear_object_cached_2(const Segment_2 & seg)143     _Linear_object_cached_2(const Segment_2& seg) :
144       has_source(true),
145       has_target(true)
146     {
147       Kernel kernel;
148 
149       CGAL_assertion_msg(! kernel.is_degenerate_2_object()(seg),
150                          "Cannot construct a degenerate segment.");
151 
152       auto construct_vertex = kernel.construct_vertex_2_object();
153       ps = construct_vertex(seg, 0);
154       pt = construct_vertex(seg, 1);
155 
156       Comparison_result res = kernel.compare_xy_2_object()(ps, pt);
157       CGAL_assertion(res != EQUAL);
158       is_degen = false;
159       is_right = (res == SMALLER);
160 
161       l = kernel.construct_line_2_object()(seg);
162       is_vert = kernel.is_vertical_2_object()(seg);
163       is_horiz = kernel.is_horizontal_2_object()(seg);
164       has_pos_slope = _has_positive_slope();
165     }
166 
167     /*! Constructor from a ray.
168      * \param ray The ray.
169      * \pre The ray is not degenerate.
170      */
_Linear_object_cached_2(const Ray_2 & ray)171     _Linear_object_cached_2(const Ray_2& ray) :
172       has_source(true),
173       has_target(false)
174     {
175       Kernel kernel;
176 
177       CGAL_assertion_msg(! kernel.is_degenerate_2_object()(ray),
178                          "Cannot construct a degenerate ray.");
179 
180       auto construct_vertex = kernel.construct_point_on_2_object();
181       ps = construct_vertex(ray, 0);         // The source point.
182       pt = construct_vertex(ray, 1);         // Some point on the ray.
183 
184       Comparison_result  res = kernel.compare_xy_2_object()(ps, pt);
185       CGAL_assertion(res != EQUAL);
186       is_degen = false;
187       is_right = (res == SMALLER);
188 
189       l = kernel.construct_line_2_object()(ray);
190       is_vert = kernel.is_vertical_2_object()(ray);
191       is_horiz = kernel.is_horizontal_2_object()(ray);
192       has_pos_slope = _has_positive_slope();
193     }
194 
195     /*! Constructor from a line.
196      * \param ln The line.
197      * \pre The line is not degenerate.
198      */
_Linear_object_cached_2(const Line_2 & ln)199     _Linear_object_cached_2(const Line_2& ln) :
200       l(ln),
201       has_source(false),
202       has_target(false)
203     {
204       Kernel kernel;
205 
206       CGAL_assertion_msg(! kernel.is_degenerate_2_object()(ln),
207                          "Cannot construct a degenerate line.");
208 
209       auto construct_vertex = kernel.construct_point_on_2_object();
210       ps = construct_vertex(ln, 0);         // Some point on the line.
211       pt = construct_vertex(ln, 1);         // Some point further on the line.
212 
213       Comparison_result res = kernel.compare_xy_2_object()(ps, pt);
214       CGAL_assertion(res != EQUAL);
215       is_degen = false;
216       is_right = (res == SMALLER);
217 
218       is_vert = kernel.is_vertical_2_object()(ln);
219       is_horiz = kernel.is_horizontal_2_object()(ln);
220       has_pos_slope = _has_positive_slope();
221     }
222 
223     /*! Check whether the x-coordinate of the left point is infinite.
224      * \return ARR_LEFT_BOUNDARY if the left point is near the boundary;
225      *         ARR_INTERIOR if the x-coordinate is finite.
226      */
left_infinite_in_x()227     Arr_parameter_space left_infinite_in_x() const
228     {
229       if (is_vert || is_degen) return (ARR_INTERIOR);
230 
231       return (is_right) ?
232         (has_source ? ARR_INTERIOR : ARR_LEFT_BOUNDARY) :
233         (has_target ? ARR_INTERIOR : ARR_LEFT_BOUNDARY);
234     }
235 
236     /*! Check whether the y-coordinate of the left point is infinite.
237      * \return ARR_BOTTOM_BOUNDARY if the left point is at y = -oo;
238      *         ARR_INTERIOR if the y-coordinate is finite.
239      *         ARR_TOP_BOUNDARY if the left point is at y = +oo;
240      */
left_infinite_in_y()241     Arr_parameter_space left_infinite_in_y() const
242     {
243       if (is_horiz || is_degen) return ARR_INTERIOR;
244 
245       if (is_vert) {
246         return (is_right) ?
247           (has_source ? ARR_INTERIOR : ARR_BOTTOM_BOUNDARY) :
248           (has_target ? ARR_INTERIOR : ARR_BOTTOM_BOUNDARY);
249       }
250 
251       if ((is_right && has_source) || (! is_right && has_target))
252         return ARR_INTERIOR;
253 
254       return (has_pos_slope ? ARR_BOTTOM_BOUNDARY : ARR_TOP_BOUNDARY);
255     }
256 
257     /*! Check whether the left point is finite.
258      */
has_left()259     bool has_left() const { return (is_right ? has_source : has_target); }
260 
261     /*! Obtain the (lexicographically) left endpoint.
262      * \pre The left point is finite.
263      */
left()264     const Point_2& left() const
265     {
266       CGAL_precondition(has_left());
267       return (is_right ? ps : pt);
268     }
269 
270     /*! Set the (lexicographically) left endpoint.
271      * \param p The point to set.
272      * \pre p lies on the supporting line to the left of the right endpoint.
273      */
274     void set_left(const Point_2& p,
275                   bool CGAL_assertion_code(check_validity) = true)
276     {
277       CGAL_precondition(! is_degen);
278 
279       CGAL_precondition_code(Kernel kernel);
280       CGAL_precondition
281         (Segment_assertions::_assert_is_point_on(p, l, Has_exact_division()) &&
282          (! check_validity || ! has_right() ||
283           kernel.compare_xy_2_object()(p, right()) == SMALLER));
284 
285       if (is_right) {
286         ps = p;
287         has_source = true;
288       }
289       else {
290         pt = p;
291         has_target = true;
292       }
293     }
294 
295     /*! Set the (lexicographically) left endpoint as infinite.
296      */
set_left()297     void set_left()
298     {
299       CGAL_precondition(! is_degen);
300 
301       if (is_right) has_source = false;
302       else has_target = false;
303     }
304 
305     /*! Check whether the x-coordinate of the right point is infinite.
306      * \return ARR_RIGHT_BOUNDARY if the right point is near the boundary;
307      *         ARR_INTERIOR if the x-coordinate is finite.
308      */
right_infinite_in_x()309     Arr_parameter_space right_infinite_in_x() const
310     {
311       if (is_vert || is_degen) return ARR_INTERIOR;
312 
313       return (is_right) ?
314         (has_target ? ARR_INTERIOR : ARR_RIGHT_BOUNDARY) :
315         (has_source ? ARR_INTERIOR : ARR_RIGHT_BOUNDARY);
316     }
317 
318     /*! Check whether the y-coordinate of the right point is infinite.
319      * \return ARR_BOTTOM_BOUNDARY if the right point is at y = -oo;
320      *         ARR_INTERIOR if the y-coordinate is finite.
321      *         ARR_TOP_BOUNDARY if the right point is at y = +oo;
322      */
right_infinite_in_y()323     Arr_parameter_space right_infinite_in_y() const
324     {
325       if (is_horiz || is_degen) return ARR_INTERIOR;
326 
327       if (is_vert) {
328         return (is_right) ?
329           (has_target ? ARR_INTERIOR : ARR_TOP_BOUNDARY) :
330           (has_source ? ARR_INTERIOR : ARR_TOP_BOUNDARY);
331       }
332 
333       if ((is_right && has_target) || (! is_right && has_source))
334           return ARR_INTERIOR;
335 
336       return (has_pos_slope ? ARR_TOP_BOUNDARY : ARR_BOTTOM_BOUNDARY);
337     }
338 
339     /*! Check whether the right point is finite.
340      */
has_right()341     bool has_right() const { return (is_right ? has_target : has_source); }
342 
343     /*! Obtain the (lexicographically) right endpoint.
344      * \pre The right endpoint is finite.
345      */
right()346     const Point_2& right() const
347     {
348       CGAL_precondition(has_right());
349       return (is_right ? pt : ps);
350     }
351 
352     /*! Set the (lexicographically) right endpoint.
353      * \param p The point to set.
354      * \pre p lies on the supporting line to the right of the left endpoint.
355      */
356     void set_right(const Point_2& p,
357                    bool CGAL_assertion_code(check_validity) = true)
358     {
359       CGAL_precondition(! is_degen);
360       CGAL_precondition_code(Kernel kernel);
361       CGAL_precondition
362         (Segment_assertions::_assert_is_point_on(p, l, Has_exact_division()) &&
363          (! check_validity || ! has_left() ||
364           kernel.compare_xy_2_object()(p, left()) == LARGER));
365 
366       if (is_right) {
367         pt = p;
368         has_target = true;
369       }
370       else {
371         ps = p;
372         has_source = true;
373       }
374     }
375 
376     /*! Set the (lexicographically) right endpoint as infinite.
377      */
set_right()378     void set_right()
379     {
380       CGAL_precondition (! is_degen);
381 
382       if (is_right) has_target = false;
383       else has_source = false;
384     }
385 
386     /*! Obtain the supporting line.
387      */
supp_line()388     const Line_2& supp_line() const
389     {
390       CGAL_precondition(! is_degen);
391       return (l);
392     }
393 
394     /*! Check whether the curve is vertical.
395      */
is_vertical()396     bool is_vertical() const
397     {
398       CGAL_precondition(! is_degen);
399       return (is_vert);
400     }
401 
402     /*! Check whether the curve is degenerate.
403      */
is_degenerate()404     bool is_degenerate() const { return (is_degen); }
405 
406     /*! Check whether the curve is directed lexicographic from left to right
407      */
is_directed_right()408     bool is_directed_right() const { return (is_right); }
409 
410     /*! Check whether the given point is in the x-range of the object.
411      * \param p The query point.
412      * \return (true) is in the x-range of the segment; (false) if it is not.
413      */
is_in_x_range(const Point_2 & p)414     bool is_in_x_range(const Point_2& p) const
415     {
416       Kernel kernel;
417       typename Kernel_::Compare_x_2 compare_x = kernel.compare_x_2_object();
418       Comparison_result res1;
419 
420       if (left_infinite_in_x() == ARR_INTERIOR) {
421         // Compare with some point on the curve.
422         if (left_infinite_in_y() != ARR_INTERIOR) res1 = compare_x(p, ps);
423         else res1 = compare_x(p, left());
424       }
425       else {
426         // p is obviously to the right.
427         res1 = LARGER;
428       }
429 
430       if (res1 == SMALLER) return false;
431       else if (res1 == EQUAL) return true;
432 
433       Comparison_result res2;
434 
435       if (right_infinite_in_x() == ARR_INTERIOR) {
436         // Compare with some point on the curve.
437         if (right_infinite_in_y() != ARR_INTERIOR) res2 = compare_x(p, ps);
438         else res2 = compare_x(p, right());
439       }
440       else {
441         // p is obviously to the right.
442         res2 = SMALLER;
443       }
444 
445       return (res2 != LARGER);
446     }
447 
448     /*! Check whether the given point is in the y-range of the object.
449      * \param p The query point.
450      * \pre The object is vertical.
451      * \return (true) is in the y-range of the segment; (false) if it is not.
452      */
is_in_y_range(const Point_2 & p)453     bool is_in_y_range(const Point_2& p) const
454     {
455       CGAL_precondition(is_vertical());
456 
457       Kernel kernel;
458       typename Kernel_::Compare_y_2 compare_y = kernel.compare_y_2_object();
459       Arr_parameter_space inf = left_infinite_in_y();
460       Comparison_result res1;
461 
462       CGAL_assertion(inf != ARR_TOP_BOUNDARY);
463       if (inf == ARR_INTERIOR) res1 = compare_y (p, left());
464       else res1 = LARGER;           // p is obviously above.
465 
466       if (res1 == SMALLER) return false;
467       else if (res1 == EQUAL) return true;
468 
469       Comparison_result res2;
470 
471       inf = right_infinite_in_y();
472       CGAL_assertion(inf != ARR_BOTTOM_BOUNDARY);
473       if (inf == ARR_INTERIOR) res2 = compare_y(p, right());
474       else res2 = SMALLER;          // p is obviously below.
475 
476       return (res2 != LARGER);
477     }
478 
479   private:
480     /*! Determine if the supporting line has a positive slope.
481      */
_has_positive_slope()482     bool _has_positive_slope() const
483     {
484       if (is_vert) return true;
485       if (is_horiz) return false;
486 
487       // Construct a horizontal line and compare its slope the that of l.
488       Kernel kernel;
489       Line_2 l_horiz =
490         kernel.construct_line_2_object()(Point_2(0, 0), Point_2(1, 0));
491       return (kernel.compare_slope_2_object()(l, l_horiz) == LARGER);
492     }
493   };
494 
495 public:
496   // Traits objects
497   typedef typename Kernel::Point_2              Point_2;
498   typedef Arr_linear_object_2<Kernel>           X_monotone_curve_2;
499   typedef Arr_linear_object_2<Kernel>           Curve_2;
500   typedef unsigned int                          Multiplicity;
501 
502 public:
503   /*! Default constructor.
504    */
Arr_linear_traits_2()505   Arr_linear_traits_2() {}
506 
507   /// \name Basic functor definitions.
508   //@{
509 
510   /*! A functor that compares the x-coordinates of two points */
511   class Compare_x_2 {
512   protected:
513     typedef Arr_linear_traits_2<Kernel> Traits;
514 
515     /*! The traits (in case it has state) */
516     const Traits& m_traits;
517 
518     /*! Constructor
519      * \param traits the traits (in case it has state)
520      * The constructor is declared private to allow only the functor
521      * obtaining function, which is a member of the nesting class,
522      * constructing it.
523      */
Compare_x_2(const Traits & traits)524     Compare_x_2(const Traits& traits) : m_traits(traits) {}
525 
526     //! Allow its functor obtaining function calling the private constructor.
527     friend class Arr_linear_traits_2<Kernel>;
528 
529   public:
530     /*! Compare the x-coordinates of two points.
531      * \param p1 The first point.
532      * \param p2 The second point.
533      * \return LARGER if x(p1) > x(p2);
534      *         SMALLER if x(p1) < x(p2);
535      *         EQUAL if x(p1) = x(p2).
536      */
operator()537     Comparison_result operator()(const Point_2& p1, const Point_2& p2) const
538     {
539       const Kernel& kernel = m_traits;
540       return (kernel.compare_x_2_object()(p1, p2));
541     }
542   };
543 
544   /*! Obtain a Compare_x_2 functor. */
compare_x_2_object()545   Compare_x_2 compare_x_2_object() const { return Compare_x_2(*this); }
546 
547   /*! A functor that compares the he endpoints of an $x$-monotone curve. */
548   class Compare_endpoints_xy_2{
549   public:
550     /*! Compare the endpoints of an $x$-monotone curve lexicographically.
551      * (assuming the curve has a designated source and target points).
552      * \param cv The curve.
553      * \return SMALLER if the curve is directed right;
554      *         LARGER if the curve is directed left.
555      */
operator()556      Comparison_result operator()(const X_monotone_curve_2& xcv) const
557     { return (xcv.is_directed_right()) ? (SMALLER) : (LARGER); }
558   };
559 
compare_endpoints_xy_2_object()560   Compare_endpoints_xy_2 compare_endpoints_xy_2_object() const
561   { return Compare_endpoints_xy_2(); }
562 
563   class Trim_2 {
564   protected:
565     typedef Arr_linear_traits_2<Kernel> Traits;
566 
567     /*! The traits (in case it has state) */
568     const Traits& m_traits;
569 
570     /*! Constructor
571      * \param traits the traits (in case it has state)
572      * The constructor is declared private to allow only the functor
573      * obtaining function, which is a member of the nesting class,
574      * constructing it.
575      */
Trim_2(const Traits & traits)576     Trim_2(const Traits& traits) : m_traits(traits) {}
577 
578     //! Allow its functor obtaining function calling the private constructor.
579     friend class Arr_linear_traits_2<Kernel>;
580 
581   public:
operator()582     X_monotone_curve_2 operator()(const X_monotone_curve_2 xcv,
583                                   const Point_2 src,
584                                   const Point_2 tgt)
585     {
586       /*
587        * "Line_segment, line, and ray" will become line segments
588        * when trimmed.
589        */
590       Equal_2 equal = Equal_2();
591       Compare_y_at_x_2 compare_y_at_x = m_traits.compare_y_at_x_2_object();
592 
593       //preconditions
594       //check if source and taget are distinct points and they lie on the line.
595       CGAL_precondition(!equal(src, tgt));
596       CGAL_precondition(compare_y_at_x(src, xcv) == EQUAL);
597       CGAL_precondition(compare_y_at_x(tgt, xcv) == EQUAL);
598 
599       //create trimmed line_segment
600       X_monotone_curve_2 trimmed_segment;
601 
602       if (xcv.is_directed_right() && (tgt.x() < src.x()))
603         trimmed_segment = Segment_2(tgt, src);
604       else if (! xcv.is_directed_right() && (tgt.x() > src.x()))
605         trimmed_segment = Segment_2(tgt, src);
606       else trimmed_segment = Segment_2(src, tgt);
607 
608       return trimmed_segment;
609     }
610   };
611 
trim_2_object()612   Trim_2 trim_2_object() const { return Trim_2(*this); }
613 
614   class Construct_opposite_2{
615   protected:
616     typedef Arr_linear_traits_2<Kernel> Traits;
617 
618     /*! The traits (in case it has state) */
619     const Traits& m_traits;
620 
621     /*! Constructor
622      * \param traits the traits (in case it has state)
623      * The constructor is declared private to allow only the functor
624      * obtaining function, which is a member of the nesting class,
625      * constructing it.
626      */
Construct_opposite_2(const Traits & traits)627     Construct_opposite_2(const Traits& traits) : m_traits(traits) {}
628 
629     //! Allow its functor obtaining function calling the private constructor.
630     friend class Arr_linear_traits_2<Kernel>;
631 
632   public:
operator()633     X_monotone_curve_2 operator()(const X_monotone_curve_2& xcv) const
634     {
635       CGAL_precondition(! xcv.is_degenerate());
636 
637       X_monotone_curve_2 opp_xcv;
638 
639       if (xcv.is_segment()) opp_xcv = Segment_2(xcv.target(), xcv.source());
640       if (xcv.is_line()) opp_xcv = Line_2(xcv.get_pt(), xcv.get_ps());
641       if (xcv.is_ray()) {
642         Point_2 opp_tgt = Point_2( -(xcv.get_pt().x()), -(xcv.get_pt().y()));
643         opp_xcv = Ray_2( xcv.source(),  opp_tgt);
644       }
645 
646       return opp_xcv;
647     }
648   };
649 
650   /*! Get a Construct_opposite_2 functor object. */
construct_opposite_2_object()651   Construct_opposite_2 construct_opposite_2_object() const
652   { return Construct_opposite_2(*this); }
653 
654   /*! A functor that compares the x-coordinates of two points */
655   class Compare_xy_2 {
656   public:
657     /*! Compare two points lexigoraphically: by x, then by y.
658      * \param p1 The first point.
659      * \param p2 The second point.
660      * \return LARGER if x(p1) > x(p2), or if x(p1) = x(p2) and y(p1) > y(p2);
661      *         SMALLER if x(p1) < x(p2), or if x(p1) = x(p2) and y(p1) < y(p2);
662      *         EQUAL if the two points are equal.
663      */
operator()664     Comparison_result operator()(const Point_2& p1, const Point_2& p2) const
665     {
666       Kernel kernel;
667       return (kernel.compare_xy_2_object()(p1, p2));
668     }
669   };
670 
671   /*! Obtain a Compare_xy_2 functor object. */
compare_xy_2_object()672   Compare_xy_2 compare_xy_2_object() const { return Compare_xy_2(); }
673 
674   /*! A functor that obtains the left endpoint of a segment or a ray. */
675   class Construct_min_vertex_2 {
676   public:
677     /*! Obtain the left endpoint of the x-monotone curve (segment).
678      * \param cv The curve.
679      * \pre The left end of cv is a valid (bounded) point.
680      * \return The left endpoint.
681      */
operator()682     const Point_2& operator()(const X_monotone_curve_2& cv) const
683     {
684       CGAL_precondition(! cv.is_degenerate());
685       CGAL_precondition(cv.has_left());
686 
687       return (cv.left());
688     }
689   };
690 
691   /*! Obtain a Construct_min_vertex_2 functor object. */
construct_min_vertex_2_object()692   Construct_min_vertex_2 construct_min_vertex_2_object() const
693   { return Construct_min_vertex_2(); }
694 
695   /*! A functor that obtains the right endpoint of a segment or a ray. */
696   class Construct_max_vertex_2 {
697   public:
698     /*! Obtain the right endpoint of the x-monotone curve (segment).
699      * \param cv The curve.
700      * \pre The right end of cv is a valid (bounded) point.
701      * \return The right endpoint.
702      */
operator()703     const Point_2& operator()(const X_monotone_curve_2& cv) const
704     {
705       CGAL_precondition(! cv.is_degenerate());
706       CGAL_precondition(cv.has_right());
707 
708       return (cv.right());
709     }
710   };
711 
712   /*! Obtain a Construct_max_vertex_2 functor object. */
construct_max_vertex_2_object()713   Construct_max_vertex_2 construct_max_vertex_2_object() const
714   { return Construct_max_vertex_2(); }
715 
716   /*! A functor that checks whether a given linear curve is vertical. */
717   class Is_vertical_2 {
718   public:
719     /*! Check whether the given x-monotone curve is a vertical segment.
720      * \param cv The curve.
721      * \return (true) if the curve is a vertical segment; (false) otherwise.
722      */
operator()723     bool operator()(const X_monotone_curve_2& cv) const
724     {
725       CGAL_precondition(! cv.is_degenerate());
726       return (cv.is_vertical());
727     }
728   };
729 
730   /*! Obtain an Is_vertical_2 functor object. */
is_vertical_2_object()731   Is_vertical_2 is_vertical_2_object() const { return Is_vertical_2(); }
732 
733   /*! A functor that compares the y-coordinates of a point and a line at
734    * the point x-coordinate
735    */
736   class Compare_y_at_x_2 {
737   protected:
738     typedef Arr_linear_traits_2<Kernel> Traits;
739 
740     /*! The traits (in case it has state) */
741     const Traits& m_traits;
742 
743     /*! Constructor
744      * \param traits the traits (in case it has state)
745      * The constructor is declared private to allow only the functor
746      * obtaining function, which is a member of the nesting class,
747      * constructing it.
748      */
Compare_y_at_x_2(const Traits & traits)749     Compare_y_at_x_2(const Traits& traits) : m_traits(traits) {}
750 
751     //! Allow its functor obtaining function calling the private constructor.
752     friend class Arr_linear_traits_2<Kernel>;
753 
754   public:
755     /*! Obtain the location of the given point with respect to the input curve.
756      * \param cv The curve.
757      * \param p The point.
758      * \pre p is in the x-range of cv.
759      * \return SMALLER if y(p) < cv(x(p)), i.e. the point is below the curve;
760      *         LARGER if y(p) > cv(x(p)), i.e. the point is above the curve;
761      *         EQUAL if p lies on the curve.
762      */
operator()763     Comparison_result operator()(const Point_2& p,
764                                  const X_monotone_curve_2& cv) const
765     {
766       CGAL_precondition(! cv.is_degenerate());
767       CGAL_precondition(cv.is_in_x_range(p));
768 
769       const Kernel& kernel = m_traits;
770       if (! cv.is_vertical())
771         // Compare p with the segment's supporting line.
772         return (kernel.compare_y_at_x_2_object()(p, cv.supp_line()));
773 
774       // Compare with the vertical segment's end-points.
775       typename Kernel::Compare_y_2 compare_y = kernel.compare_y_2_object();
776       const Comparison_result res1 =
777         cv.has_left() ? compare_y(p, cv.left()) : LARGER;
778       const Comparison_result res2 =
779         cv.has_right() ? compare_y(p, cv.right()) : SMALLER;
780 
781       return (res1 == res2) ? res1 : EQUAL;
782     }
783   };
784 
785   /*! Obtain a Compare_y_at_x_2 functor object. */
compare_y_at_x_2_object()786   Compare_y_at_x_2 compare_y_at_x_2_object() const
787   { return Compare_y_at_x_2(*this); }
788 
789   /*! A functor that compares compares the y-coordinates of two linear
790    * curves immediately to the left of their intersection point.
791    */
792   class Compare_y_at_x_left_2 {
793   public:
794     /*! Compare the y value of two x-monotone curves immediately to the left
795      * of their intersection point.
796      * \param cv1 The first curve.
797      * \param cv2 The second curve.
798      * \param p The intersection point.
799      * \pre The point p lies on both curves, and both of them must be also be
800      *      defined (lexicographically) to its left.
801      * \return The relative position of cv1 with respect to cv2 immdiately to
802      *         the left of p: SMALLER, LARGER or EQUAL.
803      */
operator()804     Comparison_result operator()(const X_monotone_curve_2& cv1,
805                                  const X_monotone_curve_2& cv2,
806                                  const Point_2& CGAL_precondition_code(p)) const
807     {
808       CGAL_precondition(! cv1.is_degenerate());
809       CGAL_precondition(! cv2.is_degenerate());
810 
811       Kernel                        kernel;
812 
813       // Make sure that p lies on both curves, and that both are defined to its
814       // left (so their left endpoint is lexicographically smaller than p).
815       CGAL_precondition_code(auto compare_xy = kernel.compare_xy_2_object());
816 
817       CGAL_precondition
818         (Segment_assertions::_assert_is_point_on(p, cv1,
819                                                  Has_exact_division()) &&
820          Segment_assertions::_assert_is_point_on(p, cv2, Has_exact_division()));
821 
822       CGAL_precondition((! cv1.has_left() ||
823                          compare_xy(cv1.left(), p) == SMALLER) &&
824                         (! cv2.has_left() ||
825                          compare_xy(cv2.left(), p) == SMALLER));
826 
827       // Compare the slopes of the two segments to determine thir relative
828       // position immediately to the left of q.
829       // Notice we use the supporting lines in order to compare the slopes,
830       // and that we swap the order of the curves in order to obtain the
831       // correct result to the left of p.
832       return (kernel.compare_slope_2_object()(cv2.supp_line(), cv1.supp_line()));
833     }
834   };
835 
836   /*! Obtain a Compare_y_at_x_left_2 functor object. */
compare_y_at_x_left_2_object()837   Compare_y_at_x_left_2 compare_y_at_x_left_2_object() const
838   { return Compare_y_at_x_left_2(); }
839 
840   /*! A functor that compares compares the y-coordinates of two linear
841    * curves immediately to the right of their intersection point.
842    */
843   class Compare_y_at_x_right_2 {
844   public:
845     /*! Compare the y value of two x-monotone curves immediately to the right
846      * of their intersection point.
847      * \param cv1 The first curve.
848      * \param cv2 The second curve.
849      * \param p The intersection point.
850      * \pre The point p lies on both curves, and both of them must be also be
851      *      defined (lexicographically) to its right.
852      * \return The relative position of cv1 with respect to cv2 immdiately to
853      *         the right of p: SMALLER, LARGER or EQUAL.
854      */
operator()855     Comparison_result operator()(const X_monotone_curve_2& cv1,
856                                  const X_monotone_curve_2& cv2,
857                                  const Point_2& CGAL_precondition_code(p)) const
858     {
859       CGAL_precondition(! cv1.is_degenerate());
860       CGAL_precondition(! cv2.is_degenerate());
861 
862       Kernel kernel;
863 
864       // Make sure that p lies on both curves, and that both are defined to its
865       // right (so their right endpoint is lexicographically larger than p).
866       CGAL_precondition_code(auto compare_xy = kernel.compare_xy_2_object());
867 
868       CGAL_precondition
869         (Segment_assertions::_assert_is_point_on(p, cv1,
870                                                  Has_exact_division()) &&
871          Segment_assertions::_assert_is_point_on(p, cv2, Has_exact_division()));
872 
873       CGAL_precondition((! cv1.has_right() ||
874                          compare_xy(cv1.right(), p) == LARGER) &&
875                         (! cv2.has_right() ||
876                          compare_xy(cv2.right(), p) == LARGER));
877 
878       // Compare the slopes of the two segments to determine thir relative
879       // position immediately to the left of q.
880       // Notice we use the supporting lines in order to compare the slopes.
881       return (kernel.compare_slope_2_object()(cv1.supp_line(),
882                                               cv2.supp_line()));
883     }
884   };
885 
886   /*! Obtain a Compare_y_at_x_right_2 functor object. */
compare_y_at_x_right_2_object()887   Compare_y_at_x_right_2 compare_y_at_x_right_2_object() const
888   { return Compare_y_at_x_right_2(); }
889 
890   /*! A functor that checks whether two points and two linear curves are
891    * identical.
892    */
893   class Equal_2 {
894   public:
895     /*! Check whether the two x-monotone curves are the same (have the same
896      * graph).
897      * \param cv1 The first curve.
898      * \param cv2 The second curve.
899      * \return (true) if the two curves are the same; (false) otherwise.
900      */
operator()901     bool operator()(const X_monotone_curve_2& cv1,
902                     const X_monotone_curve_2& cv2) const
903     {
904       CGAL_precondition(! cv1.is_degenerate());
905       CGAL_precondition(! cv2.is_degenerate());
906 
907       Kernel kernel;
908       typename Kernel::Equal_2 equal = kernel.equal_2_object();
909 
910       // Check that the two supporting lines are the same.
911       if (! equal(cv1.supp_line(), cv2.supp_line()) &&
912           ! equal(cv1.supp_line(),
913                   kernel.construct_opposite_line_2_object()(cv2.supp_line())))
914       {
915         return false;
916       }
917 
918       // Check that either the two left endpoints are at infinity, or they
919       // are bounded and equal.
920       if ((cv1.has_left() != cv2.has_left()) ||
921           (cv1.has_left() && ! equal(cv1.left(), cv2.left())))
922       {
923         return false;
924       }
925 
926       // Check that either the two right endpoints are at infinity, or they
927       // are bounded and equal.
928       return ((cv1.has_right() == cv2.has_right()) &&
929               (! cv1.has_right() || equal (cv1.right(), cv2.right())));
930     }
931 
932     /*! Check whether the two points are the same.
933      * \param p1 The first point.
934      * \param p2 The second point.
935      * \return (true) if the two point are the same; (false) otherwise.
936      */
operator()937     bool operator()(const Point_2& p1, const Point_2& p2) const
938     {
939       Kernel kernel;
940       return (kernel.equal_2_object()(p1, p2));
941     }
942   };
943 
944   /*! Obtain an Equal_2 functor object. */
equal_2_object()945   Equal_2 equal_2_object() const { return Equal_2(); }
946   //@}
947 
948   /// \name Functor definitions to handle boundaries
949   //@{
950 
951   /*! A function object that obtains the parameter space of a geometric
952    * entity along the x-axis
953    */
954   class Parameter_space_in_x_2 {
955   public:
956     /*! Obtains the parameter space at the end of a line along the x-axis.
957      * \param xcv the line
958      * \param ce the line end indicator:
959      *     ARR_MIN_END - the minimal end of xc or
960      *     ARR_MAX_END - the maximal end of xc
961      * \return the parameter space at the ce end of the line xcv.
962      *   ARR_LEFT_BOUNDARY  - the line approaches the identification arc from
963      *                        the right at the line left end.
964      *   ARR_INTERIOR       - the line does not approache the identification arc.
965      *   ARR_RIGHT_BOUNDARY - the line approaches the identification arc from
966      *                        the left at the line right end.
967      */
operator()968     Arr_parameter_space operator()(const X_monotone_curve_2 & xcv,
969                                    Arr_curve_end ce) const
970     {
971       CGAL_precondition(! xcv.is_degenerate());
972       return (ce == ARR_MIN_END) ?
973         xcv.left_infinite_in_x() : xcv.right_infinite_in_x();
974     }
975 
976     /*! Obtains the parameter space at a point along the x-axis.
977      * \param p the point.
978      * \return the parameter space at p.
979      */
operator()980     Arr_parameter_space operator()(const Point_2 ) const
981     { return ARR_INTERIOR; }
982   };
983 
984   /*! Obtain a Parameter_space_in_x_2 function object */
parameter_space_in_x_2_object()985   Parameter_space_in_x_2 parameter_space_in_x_2_object() const
986   { return Parameter_space_in_x_2(); }
987 
988   /*! A function object that obtains the parameter space of a geometric
989    * entity along the y-axis
990    */
991   class Parameter_space_in_y_2 {
992   public:
993     /*! Obtains the parameter space at the end of a line along the y-axis .
994      * Note that if the line end coincides with a pole, then unless the line
995      * coincides with the identification arc, the line end is considered to
996      * be approaching the boundary, but not on the boundary.
997      * If the line coincides with the identification arc, it is assumed to
998      * be smaller than any other object.
999      * \param xcv the line
1000      * \param ce the line end indicator:
1001      *     ARR_MIN_END - the minimal end of xc or
1002      *     ARR_MAX_END - the maximal end of xc
1003      * \return the parameter space at the ce end of the line xcv.
1004      *   ARR_BOTTOM_BOUNDARY  - the line approaches the south pole at the line
1005      *                          left end.
1006      *   ARR_INTERIOR         - the line does not approache a contraction point.
1007      *   ARR_TOP_BOUNDARY     - the line approaches the north pole at the line
1008      *                          right end.
1009      */
operator()1010     Arr_parameter_space operator()(const X_monotone_curve_2 & xcv,
1011                                    Arr_curve_end ce) const
1012     {
1013       CGAL_precondition (! xcv.is_degenerate());
1014 
1015       return (ce == ARR_MIN_END) ?
1016         xcv.left_infinite_in_y() : xcv.right_infinite_in_y();
1017     }
1018 
1019     /*! Obtains the parameter space at a point along the y-axis.
1020      * \param p the point.
1021      * \return the parameter space at p.
1022      */
operator()1023     Arr_parameter_space operator()(const Point_2 ) const
1024     { return ARR_INTERIOR; }
1025   };
1026 
1027   /*! Obtain a Parameter_space_in_y_2 function object */
parameter_space_in_y_2_object()1028   Parameter_space_in_y_2 parameter_space_in_y_2_object() const
1029   { return Parameter_space_in_y_2(); }
1030 
1031   /*! A function object that compares the x-limits of arc ends on the
1032    * boundary of the parameter space
1033    */
1034   class Compare_x_at_limit_2 {
1035   protected:
1036     typedef Arr_linear_traits_2<Kernel> Traits;
1037 
1038     /*! The traits (in case it has state) */
1039     const Traits& m_traits;
1040 
1041     /*! Constructor
1042      * \param traits the traits (in case it has state)
1043      * The constructor is declared private to allow only the functor
1044      * obtaining function, which is a member of the nesting class,
1045      * constructing it.
1046      */
Compare_x_at_limit_2(const Traits & traits)1047     Compare_x_at_limit_2(const Traits& traits) : m_traits(traits) {}
1048 
1049     //! Allow its functor obtaining function calling the private constructor.
1050     friend class Arr_linear_traits_2<Kernel>;
1051 
1052   public:
1053     /*! Compare the x-limit of a vertical line at a point with the x-limit of
1054      * a line end on the boundary at y = +/- oo.
1055      * \param p the point direction.
1056      * \param xcv the line, the endpoint of which is compared.
1057      * \param ce the line-end indicator -
1058      *            ARR_MIN_END - the minimal end of xc or
1059      *            ARR_MAX_END - the maximal end of xc.
1060      * \return the comparison result:
1061      *         SMALLER - x(p) < x(xc, ce);
1062      *         EQUAL   - x(p) = x(xc, ce);
1063      *         LARGER  - x(p) > x(xc, ce).
1064      * \pre p lies in the interior of the parameter space.
1065      * \pre the ce end of the line xcv lies on a boundary, implying
1066      *      that xcv1 is vertical.
1067      */
operator()1068     Comparison_result operator()(const Point_2 & p,
1069                                  const X_monotone_curve_2 & xcv,
1070                                  Arr_curve_end ) const
1071     {
1072       CGAL_precondition(! xcv.is_degenerate());
1073       CGAL_precondition(xcv.is_vertical());
1074 
1075       const Kernel& kernel = m_traits;
1076       return (kernel.compare_x_at_y_2_object()(p, xcv.supp_line()));
1077     }
1078 
1079     /*! Compare the x-limits of 2 arcs ends on the boundary of the
1080      * parameter space at y = +/- oo.
1081      * \param xcv1 the first arc.
1082      * \param ce1 the first arc end indicator -
1083      *            ARR_MIN_END - the minimal end of xcv1 or
1084      *            ARR_MAX_END - the maximal end of xcv1.
1085      * \param xcv2 the second arc.
1086      * \param ce2 the second arc end indicator -
1087      *            ARR_MIN_END - the minimal end of xcv2 or
1088      *            ARR_MAX_END - the maximal end of xcv2.
1089      * \return the second comparison result:
1090      *         SMALLER - x(xcv1, ce1) < x(xcv2, ce2);
1091      *         EQUAL   - x(xcv1, ce1) = x(xcv2, ce2);
1092      *         LARGER  - x(xcv1, ce1) > x(xcv2, ce2).
1093      * \pre the ce1 end of the line xcv1 lies on a boundary, implying
1094      *      that xcv1 is vertical.
1095      * \pre the ce2 end of the line xcv2 lies on a boundary, implying
1096      *      that xcv2 is vertical.
1097      */
operator()1098     Comparison_result operator()(const X_monotone_curve_2 & xcv1,
1099                                  Arr_curve_end /* ce1 */,
1100                                  const X_monotone_curve_2 & xcv2,
1101                                  Arr_curve_end /*! ce2 */) const
1102     {
1103       CGAL_precondition(! xcv1.is_degenerate());
1104       CGAL_precondition(! xcv2.is_degenerate());
1105       CGAL_precondition(xcv1.is_vertical());
1106       CGAL_precondition(xcv2.is_vertical());
1107 
1108       const Kernel& kernel = m_traits;
1109       const Point_2 p = kernel.construct_point_2_object()(ORIGIN);
1110       return (kernel.compare_x_at_y_2_object()(p, xcv1.supp_line(),
1111                                                xcv2.supp_line()));
1112     }
1113   };
1114 
1115   /*! Obtain a Compare_x_at_limit_2 function object */
compare_x_at_limit_2_object()1116   Compare_x_at_limit_2 compare_x_at_limit_2_object() const
1117   { return Compare_x_at_limit_2(*this); }
1118 
1119   /*! A function object that compares the x-coordinates of arc ends near the
1120    * boundary of the parameter space
1121    */
1122   class Compare_x_near_limit_2 {
1123   public:
1124     /*! Compare the x-coordinates of 2 arcs ends near the boundary of the
1125      * parameter space at y = +/- oo.
1126      * \param xcv1 the first arc.
1127      * \param ce1 the first arc end indicator -
1128      *            ARR_MIN_END - the minimal end of xcv1 or
1129      *            ARR_MAX_END - the maximal end of xcv1.
1130      * \param xcv2 the second arc.
1131      * \param ce2 the second arc end indicator -
1132      *            ARR_MIN_END - the minimal end of xcv2 or
1133      *            ARR_MAX_END - the maximal end of xcv2.
1134      * \return the second comparison result:
1135      *         SMALLER - x(xcv1, ce1) < x(xcv2, ce2);
1136      *         EQUAL   - x(xcv1, ce1) = x(xcv2, ce2);
1137      *         LARGER  - x(xcv1, ce1) > x(xcv2, ce2).
1138      * \pre the ce end of the line xcv1 lies on a boundary, implying
1139      *      that xcv1 is vertical.
1140      * \pre the ce end of the line xcv2 lies on a boundary, implying
1141      *      that xcv2 is vertical.
1142      * \pre the $x$-coordinates of xcv1 and xcv2 at their ce ends are
1143      *      equal, implying that the curves overlap!
1144      */
1145     Comparison_result
operator()1146     operator()(const X_monotone_curve_2& CGAL_precondition_code(xcv1),
1147                const X_monotone_curve_2& CGAL_precondition_code(xcv2),
1148                Arr_curve_end /*! ce2 */) const
1149     {
1150       CGAL_precondition(! xcv1.is_degenerate());
1151       CGAL_precondition(! xcv2.is_degenerate());
1152       CGAL_precondition(xcv1.is_vertical());
1153       CGAL_precondition(xcv2.is_vertical());
1154       return EQUAL;
1155     }
1156   };
1157 
1158   /*! Obtain a Compare_x_near_limit_2 function object */
compare_x_near_limit_2_object()1159   Compare_x_near_limit_2 compare_x_near_limit_2_object() const
1160   { return Compare_x_near_limit_2(); }
1161 
1162 
1163   /*! A function object that compares the y-limits of arc ends on the
1164    * boundary of the parameter space.
1165    */
1166   class Compare_y_near_boundary_2 {
1167   protected:
1168     typedef Arr_linear_traits_2<Kernel> Traits;
1169 
1170     /*! The traits (in case it has state) */
1171     const Traits& m_traits;
1172 
1173     /*! Constructor
1174      * \param traits the traits (in case it has state)
1175      * The constructor is declared private to allow only the functor
1176      * obtaining function, which is a member of the nesting class,
1177      * constructing it.
1178      */
Compare_y_near_boundary_2(const Traits & traits)1179     Compare_y_near_boundary_2(const Traits& traits) : m_traits(traits) {}
1180 
1181     //! Allow its functor obtaining function calling the private constructor.
1182     friend class Arr_linear_traits_2<Kernel>;
1183 
1184   public:
1185     /*! Compare the y-limits of 2 lines at their ends on the boundary
1186      * of the parameter space at x = +/- oo.
1187      * \param xcv1 the first arc.
1188      * \param xcv2 the second arc.
1189      * \param ce the line end indicator.
1190      * \return the second comparison result.
1191      * \pre the ce ends of the lines xcv1 and xcv2 lie either on the left
1192      * boundary or on the right boundary of the parameter space.
1193      */
operator()1194     Comparison_result operator()(const X_monotone_curve_2 & xcv1,
1195                                  const X_monotone_curve_2 & xcv2,
1196                                  Arr_curve_end ce) const
1197     {
1198       // Make sure both curves are defined at x = -oo (or at x = +oo).
1199       CGAL_precondition(! xcv1.is_degenerate());
1200       CGAL_precondition(! xcv2.is_degenerate());
1201       CGAL_precondition((ce == ARR_MIN_END &&
1202                          xcv1.left_infinite_in_x() == ARR_LEFT_BOUNDARY &&
1203                          xcv2.left_infinite_in_x() == ARR_LEFT_BOUNDARY) ||
1204                         (ce == ARR_MAX_END &&
1205                          xcv1.right_infinite_in_x() == ARR_RIGHT_BOUNDARY &&
1206                          xcv2.right_infinite_in_x() == ARR_RIGHT_BOUNDARY));
1207 
1208       // Compare the slopes of the two supporting lines.
1209       const Kernel& kernel = m_traits;
1210       const Comparison_result res_slopes =
1211         kernel.compare_slope_2_object()(xcv1.supp_line(), xcv2.supp_line());
1212 
1213       if (res_slopes == EQUAL) {
1214         // In case the two supporting line are parallel, compare their
1215         // relative position at x = 0, which is the same as their position
1216         // at infinity.
1217         const Point_2 p = kernel.construct_point_2_object()(ORIGIN);
1218         return (kernel.compare_y_at_x_2_object()(p, xcv1.supp_line(),
1219                                                  xcv2.supp_line()));
1220       }
1221 
1222       // Flip the slope result if we compare at x = -oo:
1223       return (ce == ARR_MIN_END) ? CGAL::opposite(res_slopes) : res_slopes;
1224     }
1225   };
1226 
1227   /*! Obtain a Compare_y_limit_on_boundary_2 function object */
compare_y_near_boundary_2_object()1228   Compare_y_near_boundary_2 compare_y_near_boundary_2_object() const
1229   { return Compare_y_near_boundary_2(*this); }
1230 
1231   //@}
1232 
1233   /// \name Functor definitions for supporting intersections.
1234   //@{
1235 
1236   class Make_x_monotone_2 {
1237   public:
1238     /*! Cut the given curve into x-monotone subcurves and insert them into the
1239      * given output iterator. As segments are always x_monotone, only one
1240      * object will be contained in the iterator.
1241      * \param cv The curve.
1242      * \param oi an output iterator for the result. Its dereference type is a
1243      *           variant that wraps a \c Point_2 or an \c X_monotone_curve_2
1244      *           objects.
1245      * \return The past-the-end iterator.
1246      */
1247     template <typename OutputIterator>
operator()1248     OutputIterator operator()(const Curve_2& cv, OutputIterator oi) const
1249     {
1250       // Wrap the segment with a variant.
1251       typedef boost::variant<Point_2, X_monotone_curve_2>
1252         Make_x_monotone_result;
1253       *oi++ = Make_x_monotone_result(cv);
1254       return oi;
1255     }
1256   };
1257 
1258   /*! Obtain a Make_x_monotone_2 functor object. */
make_x_monotone_2_object()1259   Make_x_monotone_2 make_x_monotone_2_object() const
1260   { return Make_x_monotone_2(); }
1261 
1262   class Split_2 {
1263   public:
1264     /*! Split a given x-monotone curve at a given point into two sub-curves.
1265      * \param cv The curve to split
1266      * \param p The split point.
1267      * \param c1 Output: The left resulting subcurve (p is its right endpoint).
1268      * \param c2 Output: The right resulting subcurve (p is its left endpoint).
1269      * \pre p lies on cv but is not one of its end-points.
1270      */
operator()1271     void operator()(const X_monotone_curve_2& cv, const Point_2& p,
1272                     X_monotone_curve_2& c1, X_monotone_curve_2& c2) const
1273     {
1274       CGAL_precondition (! cv.is_degenerate());
1275 
1276       // Make sure that p lies on the interior of the curve.
1277       CGAL_precondition_code (
1278         Kernel kernel;
1279         typename Kernel::Compare_xy_2 compare_xy = kernel.compare_xy_2_object();
1280       );
1281 
1282       CGAL_precondition
1283         (Segment_assertions::_assert_is_point_on(p, cv, Has_exact_division()) &&
1284          (! cv.has_left() || compare_xy(cv.left(), p) == SMALLER) &&
1285          (! cv.has_right() || compare_xy(cv.right(), p) == LARGER));
1286 
1287       // Perform the split.
1288       c1 = cv;
1289       c1.set_right(p);
1290 
1291       c2 = cv;
1292       c2.set_left(p);
1293     }
1294   };
1295 
1296   /*! Obtain a Split_2 functor object. */
split_2_object()1297   Split_2 split_2_object() const { return Split_2(); }
1298 
1299   class Intersect_2 {
1300   protected:
1301     typedef Arr_linear_traits_2<Kernel>        Traits;
1302 
1303     /*! The traits (in case it has state) */
1304     const Traits& m_traits;
1305 
1306     /*! Constructor
1307      * \param traits the traits (in case it has state)
1308      */
Intersect_2(const Traits & traits)1309     Intersect_2(const Traits& traits) : m_traits(traits) {}
1310 
1311     friend class Arr_linear_traits_2<Kernel>;
1312 
1313   public:
1314     /*! Find the intersections of the two given curves and insert them into the
1315      * given output iterator. As two segments may itersect only once, only a
1316      * single intersection will be contained in the iterator.
1317      * \param cv1 The first curve.
1318      * \param cv2 The second curve.
1319      * \param oi The output iterator.
1320      * \return The past-the-end iterator.
1321      */
1322     template <typename OutputIterator>
operator()1323     OutputIterator operator()(const X_monotone_curve_2& cv1,
1324                               const X_monotone_curve_2& cv2,
1325                               OutputIterator oi) const
1326     {
1327       typedef std::pair<Point_2, Multiplicity>          Intersection_point;
1328       typedef boost::variant<Intersection_point, X_monotone_curve_2>
1329                                                         Intersection_result;
1330 
1331       CGAL_precondition(! cv1.is_degenerate());
1332       CGAL_precondition(! cv2.is_degenerate());
1333 
1334       // Intersect the two supporting lines.
1335       const Kernel& kernel = m_traits;
1336       auto res = kernel.intersect_2_object()(cv1.supp_line(), cv2.supp_line());
1337 
1338       // The supporting line are parallel lines and do not intersect:
1339       if (! res) return oi;
1340 
1341       // Check whether we have a single intersection point.
1342       const Point_2* ip = boost::get<Point_2>(&*res);
1343       if (ip != nullptr) {
1344         // Check whether the intersection point ip lies on both segments.
1345         const bool ip_on_cv1 = cv1.is_vertical() ?
1346           cv1.is_in_y_range(*ip) : cv1.is_in_x_range(*ip);
1347 
1348         if (ip_on_cv1) {
1349           const bool ip_on_cv2 = cv2.is_vertical() ?
1350             cv2.is_in_y_range(*ip) : cv2.is_in_x_range(*ip);
1351 
1352           if (ip_on_cv2) {
1353             // Create a pair representing the point with its multiplicity,
1354             // which is always 1 for line segments.
1355             Intersection_point ip_mult(*ip, 1);
1356             *oi++ = Intersection_result(ip_mult);
1357           }
1358         }
1359         return oi;
1360       }
1361 
1362       // In this case, the two supporting lines overlap.
1363       // We start with the entire cv1 curve as the overlapping subcurve,
1364       // then clip it to form the true overlapping curve.
1365       auto compare_xy = kernel.compare_xy_2_object();
1366       X_monotone_curve_2 ovlp = cv1;
1367 
1368       if (cv2.has_left()) {
1369         // If the left endpoint of cv2 is to the right of cv1's left endpoint,
1370         // clip the overlapping subcurve.
1371         if (! cv1.has_left()) {
1372           ovlp.set_left (cv2.left(), false);
1373         }
1374         else {
1375           if (compare_xy(cv1.left(), cv2.left()) == SMALLER)
1376             ovlp.set_left(cv2.left(), false);
1377         }
1378       }
1379 
1380       if (cv2.has_right()) {
1381         // If the right endpoint of cv2 is to the left of cv1's right endpoint,
1382         // clip the overlapping subcurve.
1383         if (! cv1.has_right()) {
1384           ovlp.set_right(cv2.right(), false);
1385         }
1386         else {
1387           if (compare_xy(cv1.right(), cv2.right()) == LARGER)
1388             ovlp.set_right(cv2.right(), false);
1389         }
1390       }
1391 
1392       // Examine the resulting subcurve.
1393       Comparison_result cmp_res = SMALLER;
1394 
1395       if (ovlp.has_left() && ovlp.has_right())
1396         cmp_res = compare_xy(ovlp.left(), ovlp.right());
1397 
1398       if (cmp_res == SMALLER) {
1399         // We have discovered a true overlapping subcurve:
1400         *oi++ = Intersection_result(ovlp);
1401       }
1402       else if (cmp_res == EQUAL) {
1403         // The two objects have the same supporting line, but they just share
1404         // a common endpoint. Thus we have an intersection point, but we leave
1405         // the multiplicity of this point undefined.
1406         Intersection_point ip_mult(ovlp.left(), 0);
1407         *oi++ = Intersection_result(ip_mult);
1408       }
1409 
1410       return oi;
1411     }
1412   };
1413 
1414   /*! Obtain an Intersect_2 functor object. */
intersect_2_object()1415   Intersect_2 intersect_2_object () const { return Intersect_2(*this); }
1416 
1417   class Are_mergeable_2 {
1418   public:
1419     /*! Check whether it is possible to merge two given x-monotone curves.
1420      * \param cv1 The first curve.
1421      * \param cv2 The second curve.
1422      * \return (true) if the two curves are mergeable - if they are supported
1423      *         by the same line and share a common endpoint; (false) otherwise.
1424      */
operator()1425     bool operator()(const X_monotone_curve_2& cv1,
1426                      const X_monotone_curve_2& cv2) const
1427     {
1428       CGAL_precondition(! cv1.is_degenerate());
1429       CGAL_precondition(! cv2.is_degenerate());
1430 
1431       Kernel kernel;
1432       typename Kernel::Equal_2 equal = kernel.equal_2_object();
1433 
1434       // Check whether the two curves have the same supporting line.
1435       if (! equal(cv1.supp_line(), cv2.supp_line()) &&
1436           ! equal(cv1.supp_line(),
1437                   kernel.construct_opposite_line_2_object()(cv2.supp_line())))
1438         return false;
1439 
1440       // Check whether the left endpoint of one curve is the right endpoint of the
1441       // other.
1442       return ((cv1.has_right() && cv2.has_left() &&
1443                equal(cv1.right(), cv2.left())) ||
1444               (cv2.has_right() && cv1.has_left() &&
1445                equal(cv2.right(), cv1.left())));
1446     }
1447   };
1448 
1449   /*! Obtain an Are_mergeable_2 functor object. */
are_mergeable_2_object()1450   Are_mergeable_2 are_mergeable_2_object () const { return Are_mergeable_2(); }
1451 
1452   /*! \class Merge_2
1453    * A functor that merges two x-monotone arcs into one.
1454    */
1455   class Merge_2 {
1456   protected:
1457     typedef Arr_linear_traits_2<Kernel> Traits;
1458 
1459     /*! The traits (in case it has state) */
1460     const Traits& m_traits;
1461 
1462     /*! Constructor
1463      * \param traits the traits (in case it has state)
1464      */
Merge_2(const Traits & traits)1465     Merge_2(const Traits& traits) : m_traits(traits) {}
1466 
1467     friend class Arr_linear_traits_2<Kernel>;
1468 
1469   public:
1470     /*! Merge two given x-monotone curves into a single curve (segment).
1471      * \param cv1 The first curve.
1472      * \param cv2 The second curve.
1473      * \param c Output: The merged curve.
1474      * \pre The two curves are mergeable.
1475      */
operator()1476     void operator()(const X_monotone_curve_2& cv1,
1477                     const X_monotone_curve_2& cv2,
1478                     X_monotone_curve_2& c) const
1479     {
1480       CGAL_precondition(m_traits.are_mergeable_2_object()(cv2, cv1));
1481 
1482       CGAL_precondition(!cv1.is_degenerate());
1483       CGAL_precondition(!cv2.is_degenerate());
1484 
1485       Equal_2 equal = m_traits.equal_2_object();
1486 
1487       // Check which curve extends to the right of the other.
1488       if (cv1.has_right() && cv2.has_left() &&
1489           equal(cv1.right(), cv2.left()))
1490       {
1491         // cv2 extends cv1 to the right.
1492         c = cv1;
1493 
1494         if (cv2.has_right()) c.set_right(cv2.right());
1495         else c.set_right();      // Unbounded endpoint.
1496       }
1497       else {
1498         CGAL_precondition(cv2.has_right() && cv1.has_left() &&
1499                           equal(cv2.right(), cv1.left()));
1500 
1501         // cv1 extends cv2 to the right.
1502         c = cv2;
1503 
1504         if (cv1.has_right()) c.set_right(cv1.right());
1505         else c.set_right();      // Unbounded endpoint.
1506       }
1507     }
1508   };
1509 
1510   /*! Obtain a Merge_2 functor object. */
merge_2_object()1511   Merge_2 merge_2_object() const { return Merge_2(*this); }
1512   //@}
1513 
1514   /// \name Functor definitions for the landmarks point-location strategy.
1515   //@{
1516   typedef double                          Approximate_number_type;
1517 
1518   class Approximate_2 {
1519   public:
1520     /*! Obtain an approximation of a point coordinate.
1521      * \param p The exact point.
1522      * \param i The coordinate index (either 0 or 1).
1523      * \pre i is either 0 or 1.
1524      * \return An approximation of p's x-coordinate (if i == 0), or an
1525      *         approximation of p's y-coordinate (if i == 1).
1526      */
operator()1527     Approximate_number_type operator()(const Point_2& p, int i) const
1528     {
1529       CGAL_precondition((i == 0) || (i == 1));
1530       return (i == 0) ? CGAL::to_double(p.x()) : CGAL::to_double(p.y());
1531     }
1532   };
1533 
1534   /*! Obtain an Approximate_2 functor object. */
approximate_2_object()1535   Approximate_2 approximate_2_object() const { return Approximate_2(); }
1536 
1537   class Construct_x_monotone_curve_2 {
1538   public:
1539     /*! Obtain an x-monotone curve connecting the two given endpoints.
1540      * \param p The first point.
1541      * \param q The second point.
1542      * \pre p and q must not be the same.
1543      * \return A segment connecting p and q.
1544      */
operator()1545     X_monotone_curve_2 operator()(const Point_2& p, const Point_2& q) const
1546     {
1547       Kernel kernel;
1548       Segment_2 seg = kernel.construct_segment_2_object()(p, q);
1549 
1550       return (X_monotone_curve_2(seg));
1551     }
1552   };
1553 
1554   /*! Obtain a Construct_x_monotone_curve_2 functor object. */
construct_x_monotone_curve_2_object()1555   Construct_x_monotone_curve_2 construct_x_monotone_curve_2_object() const
1556   { return Construct_x_monotone_curve_2(); }
1557   //@}
1558 };
1559 
1560 /*!
1561  * \class A representation of a segment, as used by the Arr_segment_traits_2
1562  * traits-class.
1563  */
1564 template <typename Kernel_>
1565 class Arr_linear_object_2 :
1566     public Arr_linear_traits_2<Kernel_>::_Linear_object_cached_2
1567 {
1568   typedef typename Arr_linear_traits_2<Kernel_>::_Linear_object_cached_2
1569                                                                    Base;
1570 
1571 public:
1572   typedef Kernel_                                           Kernel;
1573 
1574   typedef typename Kernel::Point_2                          Point_2;
1575   typedef typename Kernel::Segment_2                        Segment_2;
1576   typedef typename Kernel::Ray_2                            Ray_2;
1577   typedef typename Kernel::Line_2                           Line_2;
1578 
1579 public:
1580   /*! Default constructor.
1581    */
Arr_linear_object_2()1582   Arr_linear_object_2() : Base() {}
1583 
1584   /*! Constructor from two points.
1585    * \param s The source point.
1586    * \param t The target point.
1587    * \pre The two points must not be the same.
1588    */
Arr_linear_object_2(const Point_2 & s,const Point_2 & t)1589   Arr_linear_object_2(const Point_2& s, const Point_2& t) : Base(s, t) {}
1590 
1591   /*! Constructor from a segment.
1592    * \param seg The segment.
1593    * \pre The segment is not degenerate.
1594    */
Arr_linear_object_2(const Segment_2 & seg)1595   Arr_linear_object_2(const Segment_2& seg) : Base(seg) {}
1596 
1597   /*! Constructor from a ray.
1598    * \param ray The segment.
1599    * \pre The ray is not degenerate.
1600    */
Arr_linear_object_2(const Ray_2 & ray)1601   Arr_linear_object_2(const Ray_2& ray) : Base(ray) {}
1602 
1603   /*! Constructor from a line.
1604    * \param line The line.
1605    * \pre The line is not degenerate.
1606    */
Arr_linear_object_2(const Line_2 & line)1607   Arr_linear_object_2(const Line_2& line) : Base(line) {}
1608 
1609   /*! Check whether the object is actually a segment.
1610    */
is_segment()1611   bool is_segment() const
1612   { return (! this->is_degen && this->has_source && this->has_target); }
1613 
1614   /*! Cast to a segment.
1615    * \pre The linear object is really a segment.
1616    */
segment()1617   Segment_2 segment() const
1618   {
1619     CGAL_precondition(is_segment());
1620 
1621     Kernel kernel;
1622     Segment_2 seg = kernel.construct_segment_2_object()(this->ps, this->pt);
1623     return seg;
1624   }
1625 
1626   /*! Check whether the object is actually a ray.
1627    */
is_ray()1628   bool is_ray() const
1629   { return (! this->is_degen && (this->has_source != this->has_target)); }
1630 
1631   /*! Cast to a ray.
1632    * \pre The linear object is really a ray.
1633    */
ray()1634   Ray_2 ray() const
1635   {
1636     CGAL_precondition(is_ray());
1637 
1638     Kernel kernel;
1639     Ray_2 ray = (this->has_source) ?
1640       kernel.construct_ray_2_object()(this->ps, this->l) :
1641       kernel.construct_ray_2_object()
1642         (this->pt, kernel.construct_opposite_line_2_object()(this->l));
1643     return ray;
1644   }
1645 
1646   /*! Check whether the object is actually a line.
1647    */
is_line()1648   bool is_line() const
1649   { return (! this->is_degen && ! this->has_source && ! this->has_target); }
1650 
1651   /*! Cast to a line.
1652    * \pre The linear object is really a line.
1653    */
line()1654   Line_2 line() const
1655   {
1656     CGAL_precondition(is_line());
1657     return (this->l);
1658   }
1659 
1660   /*! Get the supporting line.
1661    * \pre The object is not a point.
1662    */
supporting_line()1663   const Line_2& supporting_line() const
1664   {
1665     CGAL_precondition(! this->is_degen);
1666     return (this->l);
1667   }
1668 
1669   /*!
1670    * Get the source point.
1671    * \pre The object is a point, a segment or a ray.
1672    */
source()1673   const Point_2& source() const
1674   {
1675     CGAL_precondition(! is_line());
1676 
1677     if (this->is_degen) return (this->ps);      // For a point.
1678     if (this->has_source) return (this->ps);    // For a segment or a ray.
1679     else return (this->pt);                     // For a "flipped" ray.
1680   }
1681 
1682   /*! Get the target point.
1683    * \pre The object is a point or a segment.
1684    */
target()1685   const Point_2& target() const
1686   {
1687     CGAL_precondition(! is_line() && ! is_ray());
1688     return (this->pt);
1689   }
1690 
1691   /*! Create a bounding box for the linear object.
1692    */
bbox()1693   Bbox_2 bbox() const
1694   {
1695     CGAL_precondition(this->is_segment());
1696     Kernel kernel;
1697     Segment_2 seg = kernel.construct_segment_2_object()(this->ps, this->pt);
1698     return (kernel.construct_bbox_2_object()(seg));
1699   }
1700 
1701   // Introducing casting operators instead from a curve to
1702   // Kernel::Segment_2, Kernel::Ray_2, and Kernel::Line_2 creates an
1703   // umbiguity. The compiler will barf on the last one, because there are
1704   // 2 constructors of Kernel::Line_2: one from Kernel::Segment_2 and one
1705   // from Kernel::Ray_2. Together with the cast to Kernel::Line_2, the
1706   // compiler will have 3 equivalent options to choose from.
1707 };
1708 
1709 /*!
1710  * Exporter for the segment class used by the traits-class.
1711  */
1712 template <typename Kernel, typename OutputStream>
1713 OutputStream& operator<<(OutputStream& os,
1714                          const Arr_linear_object_2<Kernel>& lobj)
1715 {
1716   // Print a letter identifying the object type, then the object itself.
1717   if (lobj.is_segment()) os << " S " << lobj.segment();
1718   else if (lobj.is_ray()) os << " R " << lobj.ray();
1719   else os << " L " << lobj.line();
1720   return os;
1721 }
1722 
1723 /*! Importer for the segment class used by the traits-class.
1724  */
1725 template <typename Kernel, typename InputStream>
1726 InputStream& operator>>(InputStream& is, Arr_linear_object_2<Kernel>& lobj)
1727 {
1728   // Read the object type.
1729   char c;
1730 
1731   do {
1732     is >> c;
1733   } while ((c != 'S' && c != 's') &&
1734            (c != 'R' && c != 'r') &&
1735            (c != 'L' && c != 'l'));
1736 
1737   // Read the object accordingly.
1738   if (c == 'S' || c == 's') {
1739     typename Kernel::Segment_2  seg;
1740     is >> seg;
1741     lobj = seg;
1742   }
1743   else if (c == 'R' || c == 'r') {
1744     typename Kernel::Ray_2 ray;
1745     is >> ray;
1746     lobj = ray;
1747   }
1748   else {
1749     typename Kernel::Line_2 line;
1750     is >> line;
1751     lobj = line;
1752   }
1753 
1754   return is;
1755 }
1756 
1757 } //namespace CGAL
1758 
1759 #include <CGAL/enable_warnings.h>
1760 
1761 #endif
1762