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/Arrangement_zone_2.h $ 7 // $Id: Arrangement_zone_2.h e0c8048 2020-07-02T14:08:08+03:00 Efi Fogel 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 // Efi Fogel <efif@post.tau.ac.il> 13 // (based on old version by Eyal Flato) 14 15 #ifndef CGAL_ARRANGEMENT_ZONE_2_H 16 #define CGAL_ARRANGEMENT_ZONE_2_H 17 18 #include <CGAL/license/Arrangement_on_surface_2.h> 19 20 #include <CGAL/disable_warnings.h> 21 22 23 /*! \file 24 * Defintion of the Arrangement_zone_2 class. 25 */ 26 27 #include <boost/mpl/assert.hpp> 28 #include <CGAL/Arr_tags.h> 29 #include <CGAL/Arr_accessor.h> 30 #include <CGAL/Arrangement_2/Arr_traits_adaptor_2.h> 31 #include <CGAL/Arr_point_location_result.h> 32 33 #include <list> 34 #include <map> 35 #include <set> 36 37 namespace CGAL { 38 39 /*! \class 40 * A class for computing the zone of a given $x$-monotone curve in a given 41 * arrangement. 42 * The arrangement parameter corresponds to the underlying arrangement, and 43 * the zone-visitor parameter corresponds to a visitor class which is capable 44 * of receiving notifications on the arrangment features the query curve 45 * traverses. The visitor has to support the following functions: 46 * - init(), for initializing the visitor with a given arrangement. 47 * - found_subcurve(), called when a non-intersecting x-monotone curve is 48 * computed and located in the arrangement. 49 * - found_overlap(), called when an x-monotone curve overlaps an existing 50 * halfedge in the arrangement. 51 * Both the second and the third functions return pair<Halfedge_handle, bool>, 52 * where the halfedge handle corresponds to the halfedge created or modified 53 * by the visitor (if valid), and the Boolean value indicates whether we 54 * should halt the zone-computation process. 55 */ 56 template <typename Arrangement_, typename ZoneVisitor_> 57 class Arrangement_zone_2 { 58 public: 59 typedef Arrangement_ Arrangement_2; 60 typedef typename Arrangement_2::Geometry_traits_2 Geometry_traits_2; 61 typedef typename Arrangement_2::Topology_traits Topology_traits; 62 63 protected: 64 typedef Arr_traits_adaptor_2<Geometry_traits_2> Traits_adaptor_2; 65 66 typedef typename Traits_adaptor_2::Left_side_category Left_side_category; 67 typedef typename Traits_adaptor_2::Bottom_side_category Bottom_side_category; 68 typedef typename Traits_adaptor_2::Top_side_category Top_side_category; 69 typedef typename Traits_adaptor_2::Right_side_category Right_side_category; 70 71 BOOST_MPL_ASSERT 72 ((typename Arr_sane_identified_tagging<Left_side_category, 73 Bottom_side_category, 74 Top_side_category, 75 Right_side_category>::result)); 76 77 public: 78 typedef ZoneVisitor_ Visitor; 79 80 typedef typename Arrangement_2::Vertex_handle Vertex_handle; 81 typedef typename Arrangement_2::Halfedge_handle Halfedge_handle; 82 typedef typename Arrangement_2::Face_handle Face_handle; 83 84 typedef std::pair<Halfedge_handle, bool> Visitor_result; 85 86 typedef typename Geometry_traits_2::Point_2 Point_2; 87 typedef typename Geometry_traits_2::X_monotone_curve_2 X_monotone_curve_2; 88 typedef typename Geometry_traits_2::Multiplicity Multiplicity; 89 90 protected: 91 typedef typename Arr_are_all_sides_oblivious_tag<Left_side_category, 92 Bottom_side_category, 93 Top_side_category, 94 Right_side_category>::result 95 Are_all_sides_oblivious_category; 96 97 typedef typename Arrangement_2::Vertex_const_handle Vertex_const_handle; 98 typedef typename Arrangement_2::Halfedge_const_handle Halfedge_const_handle; 99 typedef typename Arrangement_2::Face_const_handle Face_const_handle; 100 101 typedef typename Arrangement_2::Ccb_halfedge_circulator 102 Ccb_halfedge_circulator; 103 104 // Types used for caching intersection points: 105 typedef std::pair<Point_2, Multiplicity> Intersection_point; 106 typedef boost::variant<Intersection_point, X_monotone_curve_2> 107 Intersection_result; 108 typedef boost::optional<Intersection_result> Optional_intersection; 109 typedef std::list<Intersection_result> Intersect_list; 110 typedef std::map<const X_monotone_curve_2*, Intersect_list> 111 Intersect_map; 112 typedef typename Intersect_map::iterator Intersect_map_iterator; 113 114 typedef std::set<const X_monotone_curve_2*> Curves_set; 115 typedef typename Curves_set::iterator Curves_set_iterator; 116 117 typedef Arr_point_location_result<Arrangement_2> Pl_result; 118 typedef typename Pl_result::Type Pl_result_type; 119 120 // Data members: 121 Arrangement_2& m_arr; // The associated arrangement. 122 const Traits_adaptor_2* m_geom_traits; // Its associated geometry traits. 123 Arr_accessor<Arrangement_2> m_arr_access; // An accessor for the arrangement. 124 125 Visitor* m_visitor; // The zone visitor. 126 127 Intersect_map m_inter_map; // Stores all computed intersections. 128 129 const Vertex_handle m_invalid_v; // An invalid vertex handle. 130 const Halfedge_handle m_invalid_he; // An invalid halfedge handle. 131 132 X_monotone_curve_2 m_cv; // The current portion of the 133 // inserted curve. 134 Pl_result_type m_obj; // The location of the left endpoint. 135 bool m_has_left_pt; // Is the left end of the curve bounded. 136 bool m_left_on_boundary; // Is the left point on the boundary. 137 Point_2 m_left_pt; // Its current left endpoint. 138 bool m_has_right_pt; // Is the right end of the curve bounded. 139 bool m_right_on_boundary; // Is the right point on the boundary. 140 Point_2 m_right_pt; // Its right endpoint (if bounded). 141 142 Vertex_handle m_left_v; // The arrangement vertex associated 143 // with the current left_pt (if any). 144 Halfedge_handle m_left_he; // If left_v is valid, left_he is the 145 // predecessor for cv around this 146 // vertex. Otherwise, if it is valid, 147 // it is the halfedge that contains 148 // the left endpoint it its interior. 149 150 Vertex_handle m_right_v; // The arrangement vertex associated 151 // with the current m_right_pt (if any). 152 Halfedge_handle m_right_he; // If m_right_v is valid, left_he is the 153 // predecessor for cv around this 154 // vertex. Otherwise, if it is valid, 155 // it is the halfedge that contains 156 // the right endpoint it its interior. 157 158 Point_2 m_intersect_p; // The next intersection point. 159 unsigned int m_ip_multiplicity; // Its multiplicity 160 // (0 in case of an overlap). 161 bool m_found_intersect; // An intersection has been found. 162 // (or an overlap). 163 X_monotone_curve_2 m_overlap_cv; // The currently discovered overlap. 164 bool m_found_overlap; // An overlap has been found. 165 bool m_found_iso_vert; // Check whether an isolated vertex 166 // induces the next intersection. 167 Vertex_handle m_intersect_v; // The vertex that intersects cv. 168 Halfedge_handle m_intersect_he; // The halfedge that intersects cv 169 // (or overlaps it). 170 171 X_monotone_curve_2 m_sub_cv1; // Auxiliary variable (for curve split). 172 X_monotone_curve_2 m_sub_cv2; // Auxiliary variable (for curve split). 173 174 public: 175 /*! Constructor. 176 * \param _arr The arrangement for which we compute the zone. 177 * \param _visitor A pointer to a zone-visitor object. 178 */ Arrangement_zone_2(Arrangement_2 & arr,Visitor * visitor)179 Arrangement_zone_2(Arrangement_2& arr, Visitor* visitor) : 180 m_arr(arr), 181 m_arr_access(arr), 182 m_visitor(visitor), 183 m_invalid_v(), 184 m_invalid_he() 185 { 186 m_geom_traits = static_cast<const Traits_adaptor_2*>(arr.geometry_traits()); 187 CGAL_assertion(visitor != nullptr); 188 189 // Initialize the visitor. 190 visitor->init(&arr); 191 } 192 193 /*! Initialize the zone-computation process with a given curve. 194 * \param _cv The query curve. 195 * \param pl A point-location object associated with the arrangement. 196 */ 197 template <typename PointLocation> init(const X_monotone_curve_2 & cv,const PointLocation & pl)198 void init(const X_monotone_curve_2& cv, const PointLocation& pl) 199 { 200 // Set the curve and check whether its left end has boundary conditions. 201 m_cv = cv; 202 203 const Arr_parameter_space bx1 = 204 m_geom_traits->parameter_space_in_x_2_object()(m_cv, ARR_MIN_END); 205 const Arr_parameter_space by1 = 206 m_geom_traits->parameter_space_in_y_2_object()(m_cv, ARR_MIN_END); 207 208 if (bx1 == ARR_INTERIOR && by1 == ARR_INTERIOR) { 209 // The curve has a finite left endpoint with no boundary conditions: 210 // locate it in the arrangement. 211 m_has_left_pt = true; 212 m_left_on_boundary = (bx1 != ARR_INTERIOR || by1 != ARR_INTERIOR); 213 m_left_pt = m_geom_traits->construct_min_vertex_2_object()(m_cv); 214 215 m_obj = pl.locate(m_left_pt); 216 } 217 else { 218 // The left end of the curve has boundary conditions: use the topology 219 // traits use the arrangement accessor to locate it. 220 // Note that if the curve-end is unbounded, m_left_pt does not exist. 221 // Note that if the curve-end is unbounded, m_left_pt does not exist. 222 m_has_left_pt = m_geom_traits->is_closed_2_object()(m_cv, ARR_MIN_END); 223 m_left_on_boundary = true; 224 if (m_has_left_pt) 225 m_left_pt = m_geom_traits->construct_min_vertex_2_object()(m_cv); 226 m_obj = m_arr_access.locate_curve_end(m_cv, ARR_MIN_END, bx1, by1); 227 } 228 229 // Check the boundary conditions of th right curve end. 230 if (m_geom_traits->is_closed_2_object()(m_cv, ARR_MAX_END)) { 231 const Arr_parameter_space bx2 = 232 m_geom_traits->parameter_space_in_x_2_object()(m_cv, ARR_MAX_END); 233 const Arr_parameter_space by2 = 234 m_geom_traits->parameter_space_in_y_2_object()(m_cv, ARR_MAX_END); 235 236 // The right endpoint is valid. 237 m_has_right_pt = true; 238 m_right_pt = m_geom_traits->construct_max_vertex_2_object()(m_cv); 239 m_right_on_boundary = (bx2 != ARR_INTERIOR) || (by2 != ARR_INTERIOR); 240 } 241 else { 242 // The right end of the curve lies at infinity. 243 m_has_right_pt = false; 244 m_right_on_boundary = true; 245 } 246 } 247 248 /*! Initialize the zone-computation process with a given curve and an object 249 * that wraps the location of the curve's left end. 250 * \param cv The query curve. 251 * \param obj An object that represents the location of the left end of the 252 * curve. 253 */ 254 void init_with_hint(const X_monotone_curve_2& cv, Pl_result_type obj); 255 256 /*! Compute the zone of the given curve and issue the apporpriate 257 * notifications for the visitor. 258 */ 259 void compute_zone(); 260 261 private: 262 /*! Check whether two curves with a common endpoint overlap. 263 * \pre p == min_point(cv1) 264 * \pre p == min_point(cv2) 265 * \todo move this function to a more accessible place so that it can be reused 266 */ do_overlap(const X_monotone_curve_2 & cv1,const X_monotone_curve_2 & cv2,const Point_2 & p)267 bool do_overlap(const X_monotone_curve_2& cv1, const X_monotone_curve_2& cv2, 268 const Point_2& p) const 269 { return do_overlap_impl(cv1, cv2, p, Are_all_sides_oblivious_category()); } 270 271 /*! Check whether two curves with a common min endpoint overlap. 272 */ do_overlap_impl(const X_monotone_curve_2 & cv1,const X_monotone_curve_2 & cv2,const Point_2 & p,Arr_all_sides_oblivious_tag)273 bool do_overlap_impl(const X_monotone_curve_2& cv1, 274 const X_monotone_curve_2& cv2, 275 const Point_2& p, Arr_all_sides_oblivious_tag) const 276 { 277 return m_geom_traits->compare_y_at_x_right_2_object()(cv1, cv2, p) == EQUAL; 278 } 279 280 /*! Check whether two curves with a common min endpoint overlap. 281 */ 282 bool do_overlap_impl(const X_monotone_curve_2& cv1, 283 const X_monotone_curve_2& cv2, 284 const Point_2& p, Arr_not_all_sides_oblivious_tag) const; 285 286 /* Check whether the given query curve is encountered when rotating the 287 * first curve in a clockwise direction around a given point until reaching 288 * the second curve. 289 * \pre p == min_point(xcv) 290 * \pre p == min_point(xcv1) 291 * \pre p == min_point(cxv2) 292 * \pre xcv_to_right == TRUE 293 * \todo move this function to a more accessible place so that it can be reused 294 */ is_between_cw(const X_monotone_curve_2 & xcv,bool xcv_to_right,const X_monotone_curve_2 & xcv1,bool xcv1_to_right,const X_monotone_curve_2 & xcv2,bool xcv2_to_right,const Point_2 & p,bool & xcv_equal_xcv1,bool & xcv_equal_xcv2)295 bool is_between_cw(const X_monotone_curve_2& xcv, bool xcv_to_right, 296 const X_monotone_curve_2& xcv1, bool xcv1_to_right, 297 const X_monotone_curve_2& xcv2, bool xcv2_to_right, 298 const Point_2& p, 299 bool& xcv_equal_xcv1, bool& xcv_equal_xcv2) const 300 { 301 return is_between_cw_impl(xcv, xcv_to_right, 302 xcv1, xcv1_to_right, 303 xcv2, xcv2_to_right, 304 p, xcv_equal_xcv1, xcv_equal_xcv2, 305 Are_all_sides_oblivious_category()); 306 } 307 308 /* Check whether the given query curve is encountered when rotating the 309 * first curve in a clockwise direction around a given point until reaching 310 * the second curve. 311 */ is_between_cw_impl(const X_monotone_curve_2 & xcv,bool xcv_to_right,const X_monotone_curve_2 & xcv1,bool xcv1_to_right,const X_monotone_curve_2 & xcv2,bool xcv2_to_right,const Point_2 & p,bool & xcv_equal_xcv1,bool & xcv_equal_xcv2,Arr_all_sides_oblivious_tag)312 bool is_between_cw_impl(const X_monotone_curve_2& xcv, bool xcv_to_right, 313 const X_monotone_curve_2& xcv1, bool xcv1_to_right, 314 const X_monotone_curve_2& xcv2, bool xcv2_to_right, 315 const Point_2& p, 316 bool& xcv_equal_xcv1, bool& xcv_equal_xcv2, 317 Arr_all_sides_oblivious_tag) const 318 { 319 return m_geom_traits->is_between_cw_2_object()(xcv, xcv_to_right, 320 xcv1, xcv1_to_right, 321 xcv2, xcv2_to_right, 322 p, 323 xcv_equal_xcv1, 324 xcv_equal_xcv2); 325 } 326 327 /* Check whether the given query curve is encountered when rotating the 328 * first curve in a clockwise direction around a given point until reaching 329 * the second curve. 330 */ 331 bool is_between_cw_impl(const X_monotone_curve_2& xcv, bool xcv_to_right, 332 const X_monotone_curve_2& xcv1, bool xcv1_to_right, 333 const X_monotone_curve_2& xcv2, bool xcv2_to_right, 334 const Point_2& p, 335 bool& xcv_equal_xcv1, bool& xcv_equal_xcv2, 336 Arr_not_all_sides_oblivious_tag) const; 337 338 /*! Find a face containing the query curve m_cv around the given vertex. 339 * In case an overlap occurs, sets m_intersect_he to be the overlapping edge. 340 * \param v The query vertex. 341 * \param he Output: The predecessor of m_cv around the vertex. 342 * \return (true) if m_cv overlaps with the curve associated with he; 343 * (false) if there is no overlap. 344 */ 345 bool _find_prev_around_vertex(Vertex_handle v, Halfedge_handle& he); 346 347 /*! Direct the halfedge for the location of the given subcurve around a split 348 * point that occurs in the interior of a given edge, when the subcurve lies 349 * to the right of the split point. 350 * In case of overlaps, it sets also m_found_overlap and m_intersect_he. 351 * \param cv_ins The curve to be inserted, whose left endpoint coincides 352 * with the edge to be split. 353 * \param cv_left_pt The left endpoint of cv_ins. 354 * \param query_he The edge that intersects cv_ins. 355 * \pre The left endpoint of cv_ins lies in the interior of the curve 356 * associated with query_he. 357 * \return The halfedge whose incident face contains cv_ins 358 * (either query_he or its twin). 359 */ 360 Halfedge_handle 361 _direct_intersecting_edge_to_right(const X_monotone_curve_2& cv_ins, 362 const Point_2& cv_left_pt, 363 Halfedge_handle query_he); 364 365 /*! Direct the halfedge for the location of the given subcurve around a split 366 * point that occurs in the interior of a given edge, when the subcurve lies 367 * to the left of the split point. 368 * \param cv_ins The curve to be inserted, whose right endpoint coincides 369 * with the edge to be split. 370 * \param query_he The edge that intersects cv_ins. 371 * \pre The right endpoint of cv_ins lies in the interior of the curve 372 * associated with query_he. 373 * \return The halfedge whose incident face contains cv_ins 374 * (either query_he or its twin). 375 */ 376 Halfedge_handle 377 _direct_intersecting_edge_to_left(const X_monotone_curve_2& cv_ins, 378 Halfedge_handle query_he); 379 380 /*! Get the next intersection of m_cv with the given halfedge. 381 * \param he A handle to the halfedge. 382 * \param skip_first_point Should we skip the first intersection point. 383 * \param intersect_on_right_boundary Output: If an intersetion point is 384 * computed, marks whether this 385 * point coincides with the right 386 * curve-end, which lies on the 387 * surface boundary. 388 * \return An object representing the next intersection: Intersection_point 389 * in case of a simple intersection point, X_monotone_curve_2 in 390 * case of an overlap, and an empty object if there is no 391 * intersection. 392 */ 393 Optional_intersection 394 _compute_next_intersection(Halfedge_handle he, 395 bool skip_first_point, 396 bool& intersect_on_right_boundary); 397 398 /*! Remove the next intersection of m_cv with the given halfedge from the map. 399 * \param he A handle to the halfedge. 400 * \pre The list of intersections with the curve of he has already been 401 * computed, and it is not empty. 402 */ 403 void _remove_next_intersection (Halfedge_handle he); 404 405 /*! Check whether the given point lies completely to the left of the given 406 * egde. 407 * \param p The point. 408 * \param he The halfedge. 409 * \pre he is not a fictitious edge. 410 * \return Whether p lies entirely to the left of the edge. 411 */ _is_to_left(const Point_2 & p,Halfedge_handle he)412 bool _is_to_left(const Point_2& p, Halfedge_handle he) const 413 { 414 return (_is_to_left_impl(p, he, Are_all_sides_oblivious_category())); 415 } 416 _is_to_left_impl(const Point_2 & p,Halfedge_handle he,Arr_all_sides_oblivious_tag)417 bool _is_to_left_impl(const Point_2& p, Halfedge_handle he, 418 Arr_all_sides_oblivious_tag) const 419 { 420 return (((he->direction() == ARR_LEFT_TO_RIGHT) && 421 m_geom_traits->compare_xy_2_object() 422 (p, he->source()->point()) == SMALLER) || 423 (he->direction() == ARR_RIGHT_TO_LEFT && 424 m_geom_traits->compare_xy_2_object() 425 (p, he->target()->point()) == SMALLER)); 426 } 427 428 bool _is_to_left_impl(const Point_2& p, Halfedge_handle he, 429 Arr_not_all_sides_oblivious_tag) const; 430 431 /*! Check whether the given point lies completely to the right of the given 432 * egde. 433 * \param p The point. 434 * \param he The halfedge. 435 * \pre he is not a fictitious edge. 436 * \return Whether p lies entirely to the right of the edge. 437 */ _is_to_right(const Point_2 & p,Halfedge_handle he)438 bool _is_to_right(const Point_2& p, Halfedge_handle he) const 439 { 440 return (_is_to_right_impl(p, he, Are_all_sides_oblivious_category())); 441 } 442 _is_to_right_impl(const Point_2 & p,Halfedge_handle he,Arr_all_sides_oblivious_tag)443 bool _is_to_right_impl(const Point_2& p, Halfedge_handle he, 444 Arr_all_sides_oblivious_tag) const 445 { 446 return (((he->direction() == ARR_LEFT_TO_RIGHT) && 447 m_geom_traits->compare_xy_2_object()(p, he->target()->point()) == 448 LARGER) || 449 ((he->direction() == ARR_RIGHT_TO_LEFT) && 450 m_geom_traits->compare_xy_2_object()(p, he->source()->point()) == 451 LARGER)); 452 } 453 454 bool _is_to_right_impl(const Point_2& p, Halfedge_handle he, 455 Arr_not_all_sides_oblivious_tag) const; 456 457 /*! Compute the (lexicographically) leftmost intersection of the query 458 * curve with a given halfedge on the boundary of a face in the arrangement. 459 */ 460 void 461 _leftmost_intersection(Ccb_halfedge_circulator he_curr, bool on_boundary, 462 bool& leftmost_on_right_boundary); 463 464 /*! Compute the (lexicographically) leftmost intersection of the query 465 * curve with the boundary of a given face in the arrangement. 466 * The function computes sets m_intersect_p, m_intersect_he (or alternatively 467 * m_overlap_cv and m_intersect_he) and set the flags m_found_intersect and 468 * m_found_overlap accordingly. 469 * \param face A handle to the face. 470 * \param on_boundary Specifies whether the left endpoint of the curve lies 471 * on the face boundary. 472 */ 473 void _leftmost_intersection_with_face_boundary(Face_handle face, 474 bool on_boundary); 475 476 /*! Compute the zone of an x-monotone curve in a given arrangement face. 477 * The left endpoint of the curve either lies in the face interior or on 478 * the boundary of the face. 479 * This function updates m_cv and its left endpoint and also sets m_left_v 480 * and m_left_he for the remaining portion of the curve. 481 * In case of overlaps, it sets also m_overlap_cv and m_intersect_he. 482 * \param face The given face. 483 * \param on_boundary Specifies whether the left endpoint of the curve lies 484 * on the face boundary. 485 * \pre If on_boundary is (true) then m_left_he must be valid; if it is 486 * (false), then both m_left_v anf m_left_he must be invalid. 487 * \return (true) if we are done with the zone-computation process; 488 * (false) if we still have a remaining portion of m_cv to continue 489 * with. 490 */ 491 bool _zone_in_face(Face_handle face, bool on_boundary); 492 493 /*! Compute the zone of an overlapping subcurve m_overlap_cv of m_cv and the 494 * curve currently associated with m_intersect_he. 495 * This function updates m_cv and its left endpoint and also sets m_left_v 496 * and m_left_he for the remaining portion of the curve. 497 * \return (true) if we are done with the zone-computation process; 498 * (false) if we still have a remaining portion of m_cv to continue 499 * with. 500 */ 501 bool _zone_in_overlap(); 502 }; 503 504 } //namespace CGAL 505 506 // The function definitions can be found under: 507 #include <CGAL/Arrangement_2/Arrangement_zone_2_impl.h> 508 509 #include <CGAL/enable_warnings.h> 510 511 #endif 512