|
BWAPI
|
00001 // Copyright (c) 2003,2004,2005,2006 INRIA Sophia-Antipolis (France) and 00002 // Notre Dame University (U.S.A.). All rights reserved. 00003 // 00004 // This file is part of CGAL (www.cgal.org); you may redistribute it under 00005 // the terms of the Q Public License version 1.0. 00006 // See the file LICENSE.QPL distributed with CGAL. 00007 // 00008 // Licensees holding a valid commercial license may use this file in 00009 // accordance with the commercial license agreement provided with the software. 00010 // 00011 // This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE 00012 // WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE. 00013 // 00014 // $URL: svn+ssh://scm.gforge.inria.fr/svn/cgal/branches/CGAL-3.5-branch/Segment_Delaunay_graph_2/include/CGAL/Segment_Delaunay_graph_2/Finite_edge_interior_conflict_C2.h $ 00015 // $Id: Finite_edge_interior_conflict_C2.h 46227 2008-10-13 11:52:58Z afabri $ 00016 // 00017 // 00018 // Author(s) : Menelaos Karavelas <mkaravel@cse.nd.edu> 00019 00020 #ifndef CGAL_SEGMENT_DELAUNAY_GRAPH_2_FINITE_EDGE_INTERIOR_CONFLICT_C2_H 00021 #define CGAL_SEGMENT_DELAUNAY_GRAPH_2_FINITE_EDGE_INTERIOR_CONFLICT_C2_H 00022 00023 #include <CGAL/Segment_Delaunay_graph_2/Basic_predicates_C2.h> 00024 #include <CGAL/Segment_Delaunay_graph_2/Voronoi_vertex_C2.h> 00025 #include <CGAL/Segment_Delaunay_graph_2/Are_same_points_C2.h> 00026 #include <CGAL/Segment_Delaunay_graph_2/Are_same_segments_C2.h> 00027 00028 #if defined(BOOST_MSVC) 00029 # pragma warning(push) 00030 # pragma warning(disable:4800) // complaint about performance where we can't do anything 00031 #endif 00032 00033 CGAL_BEGIN_NAMESPACE 00034 00035 CGAL_SEGMENT_DELAUNAY_GRAPH_2_BEGIN_NAMESPACE 00036 00037 //----------------------------------------------------------------------------- 00038 00039 template<class K, class Method_tag> 00040 class Finite_edge_interior_conflict_C2 00041 : public Basic_predicates_C2<K> 00042 { 00043 public: 00044 00045 typedef Basic_predicates_C2<K> Base; 00046 typedef Voronoi_vertex_C2<K,Method_tag> Voronoi_vertex_2; 00047 00048 typedef typename Base::Point_2 Point_2; 00049 typedef typename Base::Segment_2 Segment_2; 00050 typedef typename Base::Line_2 Line_2; 00051 typedef typename Base::Site_2 Site_2; 00052 typedef typename Base::FT FT; 00053 typedef typename Base::RT RT; 00054 00055 typedef typename Base::Comparison_result Comparison_result; 00056 typedef typename Base::Sign Sign; 00057 typedef typename Base::Orientation Orientation; 00058 typedef typename Base::Oriented_side Oriented_side; 00059 typedef typename Base::Boolean Boolean; 00060 00061 typedef typename Base::Homogeneous_point_2 Homogeneous_point_2; 00062 00063 private: 00064 typedef Are_same_points_C2<K> Are_same_points_2; 00065 typedef Are_same_segments_C2<K> Are_same_segments_2; 00066 00067 typedef typename K::Intersections_tag ITag; 00068 00069 private: 00070 Are_same_points_2 same_points; 00071 Are_same_segments_2 same_segments; 00072 00073 private: 00074 00075 //-------------------------------------------------------------------- 00076 //-------------------------------------------------------------------- 00077 //-------------------------------------------------------------------- 00078 00079 Boolean 00080 is_interior_in_conflict_both(const Site_2& p, const Site_2& q, 00081 const Site_2& r, const Site_2& s, 00082 const Site_2& t, Method_tag tag) const 00083 { 00084 Boolean in_conflict(false); 00085 00086 if ( p.is_point() && q.is_point() ) { 00087 in_conflict = is_interior_in_conflict_both_pp(p, q, r, s, t, tag); 00088 00089 } else if ( p.is_segment() && q.is_segment() ) { 00090 00091 in_conflict = is_interior_in_conflict_both_ss(p, q, r, s, t, tag); 00092 00093 } else if ( p.is_point() && q.is_segment() ) { 00094 00095 in_conflict = is_interior_in_conflict_both_ps(p, q, r, s, t, tag); 00096 00097 } else { // p is a segment and q is a point 00098 00099 in_conflict = is_interior_in_conflict_both_sp(p, q, r, s, t, tag); 00100 } 00101 00102 return in_conflict; 00103 } 00104 00105 //-------------------------------------------------------------------- 00106 00107 bool 00108 is_interior_in_conflict_both_pp(const Site_2& sp, const Site_2& sq, 00109 const Site_2& r, const Site_2& s, 00110 const Site_2& t, Method_tag ) const 00111 { 00112 CGAL_precondition( sp.is_point() && sq.is_point() ); 00113 00114 Point_2 p = sp.point(), q = sq.point(); 00115 00116 if ( t.is_point() ) { return true; } 00117 00118 Line_2 lt = compute_supporting_line(t.supporting_site()); 00119 00120 Oriented_side op, oq; 00121 00122 if ( same_points(sp, t.source_site()) || 00123 same_points(sp, t.target_site()) ) { 00124 op = ON_ORIENTED_BOUNDARY; 00125 } else { 00126 op = oriented_side_of_line(lt, p); 00127 } 00128 00129 if ( same_points(sq, t.source_site()) || 00130 same_points(sq, t.target_site()) ) { 00131 oq = ON_ORIENTED_BOUNDARY; 00132 } else { 00133 oq = oriented_side_of_line(lt, q); 00134 } 00135 00136 00137 if ((op == ON_POSITIVE_SIDE && oq == ON_NEGATIVE_SIDE) || 00138 (op == ON_NEGATIVE_SIDE && oq == ON_POSITIVE_SIDE) || 00139 (op == ON_ORIENTED_BOUNDARY || oq == ON_ORIENTED_BOUNDARY)) { 00140 return true; 00141 } 00142 00143 Comparison_result res = 00144 compare_squared_distances_to_line(lt, p, q); 00145 00146 if ( res == EQUAL ) { return true; } 00147 00148 Voronoi_vertex_2 vpqr(sp, sq, r); 00149 Voronoi_vertex_2 vqps(sq, sp, s); 00150 00151 00152 Line_2 lperp; 00153 if ( res == SMALLER ) { 00154 // p is closer to lt than q 00155 lperp = compute_perpendicular(lt, p); 00156 } else { 00157 // q is closer to lt than p 00158 lperp = compute_perpendicular(lt, q); 00159 } 00160 00161 Oriented_side opqr = vpqr.oriented_side(lperp); 00162 Oriented_side oqps = vqps.oriented_side(lperp); 00163 00164 return ( opqr == oqps ); 00165 } 00166 00167 //-------------------------------------------------------------------- 00168 00169 bool 00170 is_interior_in_conflict_both_ss(const Site_2& p, const Site_2& q, 00171 const Site_2& , const Site_2& , 00172 const Site_2& , Method_tag) const 00173 { 00174 CGAL_precondition( p.is_segment() && q.is_segment() ); 00175 return true; 00176 } 00177 00178 //-------------------------------------------------------------------- 00179 00180 Boolean 00181 is_interior_in_conflict_both_ps(const Site_2& p, const Site_2& q, 00182 const Site_2& r, const Site_2& s, 00183 const Site_2& t, Method_tag tag) const 00184 { 00185 CGAL_precondition( p.is_point() && q.is_segment() ); 00186 00187 if ( same_points(p, q.source_site()) || 00188 same_points(p, q.target_site()) ) { 00189 return false; 00190 } 00191 00192 if ( t.is_point() ) { 00193 return is_interior_in_conflict_both_ps_p(p, q, r, s, t, tag); 00194 } 00195 return is_interior_in_conflict_both_ps_s(p, q, r, s, t, tag); 00196 } 00197 00198 //-------------------------------------------------------------------- 00199 00200 Boolean 00201 is_interior_in_conflict_both_ps_p(const Site_2& p, const Site_2& q, 00202 const Site_2& r, const Site_2& s, 00203 const Site_2& t, Method_tag ) const 00204 { 00205 CGAL_precondition( t.is_point() ); 00206 00207 // Line_2 lq = compute_supporting_line(q); 00208 Line_2 lq = compute_supporting_line(q.supporting_site()); 00209 00210 Comparison_result res = 00211 compare_squared_distances_to_line(lq, p.point(), t.point()); 00212 00213 //if ( res != SMALLER ) { return true; } 00214 if (certainly( res != SMALLER ) ) { return true; } 00215 if (! is_certain( res != SMALLER ) ) { return indeterminate<Boolean>(); } 00216 00217 Voronoi_vertex_2 vpqr(p, q, r); 00218 Voronoi_vertex_2 vqps(q, p, s); 00219 00220 Line_2 lperp = compute_perpendicular(lq, p.point()); 00221 00222 Oriented_side opqr = vpqr.oriented_side(lperp); 00223 Oriented_side oqps = vqps.oriented_side(lperp); 00224 00225 return (opqr == oqps); 00226 } 00227 00228 //-------------------------------------------------------------------- 00229 00230 bool check_if_exact(const Site_2&, const Tag_false&) const 00231 { 00232 return true; 00233 } 00234 00235 bool check_if_exact(const Site_2& t1, const Tag_true&) const 00236 { 00237 return t1.is_input(); 00238 } 00239 00240 Boolean 00241 is_interior_in_conflict_both_ps_s(const Site_2& sp, const Site_2& sq, 00242 const Site_2& r, const Site_2& s, 00243 const Site_2& st, Method_tag ) const 00244 { 00245 CGAL_precondition( st.is_segment() ); 00246 Point_2 p = sp.point(); 00247 // Segment_2 q = sq.segment(), t = st.segment(); 00248 00249 Line_2 lq = compute_supporting_line(sq.supporting_site()); 00250 00251 if ( oriented_side_of_line(lq, p) == ON_NEGATIVE_SIDE ) { 00252 lq = opposite_line(lq); 00253 } 00254 00255 if ( same_points(sp, st.source_site()) || 00256 same_points(sp, st.target_site()) ) { 00257 00258 Line_2 lqperp = compute_perpendicular(lq, p); 00259 00260 Voronoi_vertex_2 vpqr(sp, sq, r); 00261 Voronoi_vertex_2 vqps(sq, sp, s); 00262 00263 Oriented_side opqr = vpqr.oriented_side(lqperp); 00264 Oriented_side oqps = vqps.oriented_side(lqperp); 00265 00266 Boolean on_different_parabola_arcs = 00267 ((opqr == ON_NEGATIVE_SIDE) & (oqps == ON_POSITIVE_SIDE)) | 00268 ((opqr == ON_POSITIVE_SIDE) & (oqps == ON_NEGATIVE_SIDE)); 00269 00270 //if ( !on_different_parabola_arcs ) { return true; } 00271 if (certainly( !on_different_parabola_arcs ) ) { return true; } 00272 if (! is_certain( !on_different_parabola_arcs ) ) { return indeterminate<Boolean>(); } 00273 00274 Site_2 t1; 00275 if ( same_points(sp, st.source_site()) ) { 00276 t1 = st.target_site(); 00277 } else { 00278 t1 = st.source_site(); 00279 } 00280 00281 Oriented_side o_t1; 00282 00283 if ( same_points(t1, sq.source_site()) || 00284 same_points(t1, sq.target_site()) ) { 00285 o_t1 = ON_ORIENTED_BOUNDARY; 00286 } else if ( !check_if_exact(t1, ITag()) && 00287 ( same_segments(t1.supporting_site(0), 00288 sq.supporting_site()) || 00289 same_segments(t1.supporting_site(1), 00290 sq.supporting_site()) ) ) { 00291 o_t1 = ON_ORIENTED_BOUNDARY; 00292 } else { 00293 o_t1 = oriented_side_of_line(lq, t1.point()); 00294 } 00295 00296 if ( o_t1 == ON_NEGATIVE_SIDE ) { 00297 return true; 00298 } 00299 00300 Comparison_result res = 00301 compare_squared_distances_to_line(lq, p, t1.point()); 00302 00303 return ( res == LARGER ); 00304 } 00305 00306 Line_2 lt = compute_supporting_line(st.supporting_site()); 00307 00308 if ( oriented_side_of_line(lt, p) == ON_NEGATIVE_SIDE ) { 00309 lt = opposite_line(lt); 00310 } 00311 00312 Comparison_result res = 00313 CGAL::compare(lt.a() * lq.b(), lt.b() * lq.a()); 00314 bool are_parallel = (res == EQUAL); 00315 00316 if ( are_parallel ) { 00317 Sign sgn = CGAL::sign(lt.a() * lq.a() + lt.b() * lq.b()); 00318 bool have_opposite_directions = (sgn == NEGATIVE); 00319 if ( have_opposite_directions ) { lq = opposite_line(lq); } 00320 00321 if ( oriented_side_of_line(lq, p) == oriented_side_of_line(lt, p) ) { 00322 return true; 00323 } 00324 00325 if ( have_opposite_directions ) { 00326 lq = opposite_line(lq); 00327 } 00328 } 00329 00330 Line_2 l = compute_perpendicular(lt, p); 00331 00332 Voronoi_vertex_2 vpqr(sp, sq, r); 00333 Voronoi_vertex_2 vqps(sq, sp, s); 00334 00335 Oriented_side o_l_pqr = vpqr.oriented_side(l); 00336 Oriented_side o_l_qps = vqps.oriented_side(l); 00337 if (certainly( (o_l_pqr == ON_POSITIVE_SIDE) & 00338 (o_l_qps == ON_NEGATIVE_SIDE) ) ) 00339 return false; 00340 if (certainly( (o_l_pqr == ON_NEGATIVE_SIDE) & 00341 (o_l_qps == ON_POSITIVE_SIDE) ) ) 00342 return true; 00343 if (! is_certain((o_l_pqr == -o_l_qps) & (o_l_pqr != ZERO))) 00344 return indeterminate<Boolean>(); 00345 00346 //>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 00347 //>>>>>>>>>> HERE I NEED TO CHECK THE BOUNDARY CASES <<<<<< 00348 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 00349 Line_2 lqperp = compute_perpendicular(lq, p); 00350 00351 Oriented_side opqr = vpqr.oriented_side(lqperp); 00352 Oriented_side oqps = vqps.oriented_side(lqperp); 00353 00354 Boolean on_different_parabola_arcs = (opqr == -oqps) & (opqr != ZERO); 00355 00356 // if ( !on_different_parabola_arcs ) { return true; } 00357 if (certainly( !on_different_parabola_arcs ) ) { return true; } 00358 if (! is_certain( !on_different_parabola_arcs ) ) { return indeterminate<Boolean>(); } 00359 00360 Homogeneous_point_2 pv = projection_on_line(lq, p); 00361 Homogeneous_point_2 hp(p); 00362 00363 pv = Base::midpoint(pv, hp); 00364 00365 Oriented_side o_l_pv = oriented_side_of_line(l, pv); 00366 00367 CGAL_assertion( o_l_pv != ON_ORIENTED_BOUNDARY ); 00368 00369 CGAL_assertion( o_l_pqr != ON_ORIENTED_BOUNDARY || 00370 o_l_qps != ON_ORIENTED_BOUNDARY ); 00371 00372 if ( o_l_pqr == ON_ORIENTED_BOUNDARY ) { 00373 return ( o_l_qps == o_l_pv ); 00374 } else { 00375 return ( o_l_pqr == o_l_pv ); 00376 } 00377 00378 } 00379 00380 //-------------------------------------------------------------------- 00381 00382 Boolean 00383 is_interior_in_conflict_both_sp(const Site_2& p, const Site_2& q, 00384 const Site_2& r, const Site_2& s, 00385 const Site_2& t, Method_tag tag) const 00386 { 00387 return is_interior_in_conflict_both_ps(q, p, s, r, t, tag); 00388 } 00389 00390 00391 //------------------------------------------------------------------------ 00392 //------------------------------------------------------------------------ 00393 //------------------------------------------------------------------------ 00394 00395 bool 00396 is_interior_in_conflict_touch(const Site_2& p, const Site_2& q, 00397 const Site_2& r, const Site_2& s, 00398 const Site_2& t, Method_tag tag) const 00399 { 00400 // checks if interior of voronoi edge is in conflict if both extrema 00401 // of the voronoi edge touch the corresponding circles. 00402 // return true if interior is in conflict; false otherwise 00403 if ( t.is_segment() ) { return false; } 00404 00405 #if 1 00406 CGAL_assertion( p.is_segment() || q.is_segment() ); 00407 00408 Voronoi_vertex_2 vpqr(p, q, r); 00409 Voronoi_vertex_2 vqps(q, p, s); 00410 00411 if ( vpqr.incircle_no_easy(s) == ZERO && 00412 vqps.incircle_no_easy(r) == ZERO ) { 00413 return false; 00414 } 00415 00416 if ( p.is_segment() && q.is_segment() ) { 00417 return true; 00418 } 00419 #else 00420 // OLD CODE: buggy if the edge is degenerate 00421 if ( (p.is_point() && q.is_point()) || 00422 (p.is_segment() && q.is_segment()) ) { 00423 return true; 00424 } 00425 #endif 00426 00427 if ( p.is_point() && q.is_segment() ) { 00428 Line_2 lq = compute_supporting_line(q.supporting_site()); 00429 00430 Comparison_result res = 00431 compare_squared_distances_to_line(lq, p.point(), t.point()); 00432 00433 return (res != SMALLER); 00434 } 00435 00436 return is_interior_in_conflict_touch(q, p, s, r, t, tag); 00437 } 00438 00439 //------------------------------------------------------------------------ 00440 //------------------------------------------------------------------------ 00441 //------------------------------------------------------------------------ 00442 00443 00444 bool 00445 is_interior_in_conflict_none(const Site_2& p, const Site_2& q, 00446 const Site_2& r, const Site_2& s, 00447 const Site_2& t, Method_tag tag) const 00448 { 00449 if ( t.is_segment() ) { return false; } 00450 00451 bool in_conflict(false); 00452 00453 if ( p.is_point() && q.is_point() ) { 00454 in_conflict = is_interior_in_conflict_none_pp(p, q, r, s, t, tag); 00455 } else if ( p.is_point() && q.is_segment() ) { 00456 in_conflict = is_interior_in_conflict_none_ps(p, q, r, s, t, tag); 00457 } else if ( p.is_segment() && q.is_point() ) { 00458 in_conflict = is_interior_in_conflict_none_sp(p, q, r, s, t, tag); 00459 } else { // both p and q are segments 00460 in_conflict = is_interior_in_conflict_none_ss(p, q, r, s, t, tag); 00461 } 00462 00463 return in_conflict; 00464 } 00465 00466 //------------------------------------------------------------------------ 00467 00468 00469 bool 00470 is_interior_in_conflict_none_pp(const Site_2& p, const Site_2& q, 00471 const Site_2& , const Site_2& , 00472 const Site_2& t, Method_tag ) const 00473 { 00474 CGAL_precondition( p.is_point() && q.is_point() && t.is_point() ); 00475 return false; 00476 } 00477 00478 //------------------------------------------------------------------------ 00479 00480 bool 00481 is_interior_in_conflict_none_ps(const Site_2& sp, const Site_2& sq, 00482 const Site_2& r, const Site_2& s, 00483 const Site_2& st, Method_tag ) const 00484 { 00485 CGAL_precondition( sp.is_point() && sq.is_segment() && st.is_point() ); 00486 00487 if ( same_points(sp, sq.source_site()) || 00488 same_points(sp, sq.target_site()) ) { 00489 return false; 00490 } 00491 00492 Line_2 lq = compute_supporting_line(sq.supporting_site()); 00493 00494 Voronoi_vertex_2 vpqr(sp, sq, r); 00495 Voronoi_vertex_2 vqps(sq, sp, s); 00496 00497 Point_2 p = sp.point(), t = st.point(); 00498 00499 Line_2 lperp = compute_perpendicular(lq, t); 00500 00501 Oriented_side op = oriented_side_of_line(lq, p); 00502 Oriented_side ot = oriented_side_of_line(lq, t); 00503 00504 bool on_same_side = 00505 ((op == ON_POSITIVE_SIDE && ot == ON_POSITIVE_SIDE) || 00506 (op == ON_NEGATIVE_SIDE && ot == ON_NEGATIVE_SIDE)); 00507 00508 Comparison_result res = 00509 compare_squared_distances_to_line(lq, t, p); 00510 00511 Oriented_side opqr = vpqr.oriented_side(lperp); 00512 Oriented_side oqps = vqps.oriented_side(lperp); 00513 00514 bool on_different_side = 00515 ((opqr == ON_POSITIVE_SIDE && oqps == ON_NEGATIVE_SIDE) || 00516 (opqr == ON_NEGATIVE_SIDE && oqps == ON_POSITIVE_SIDE)); 00517 00518 return ( on_same_side && (res == SMALLER) && on_different_side ); 00519 } 00520 00521 //------------------------------------------------------------------------ 00522 00523 bool 00524 is_interior_in_conflict_none_sp(const Site_2& p, const Site_2& q, 00525 const Site_2& r, const Site_2& s, 00526 const Site_2& t, Method_tag tag) const 00527 { 00528 return is_interior_in_conflict_none_ps(q, p, s, r, t, tag); 00529 } 00530 00531 //------------------------------------------------------------------------ 00532 00533 bool 00534 is_interior_in_conflict_none_ss(const Site_2& p, const Site_2& q, 00535 const Site_2& r, const Site_2& s, 00536 const Site_2& t, Method_tag ) const 00537 { 00538 CGAL_precondition( p.is_segment() && q.is_segment() && t.is_point() ); 00539 00540 Voronoi_vertex_2 vpqr(p, q, r); 00541 Voronoi_vertex_2 vqps(q, p, s); 00542 00543 Line_2 lp = compute_supporting_line(p.supporting_site()); 00544 Line_2 lq = compute_supporting_line(q.supporting_site()); 00545 00546 // first orient lp according to its Voronoi vertices 00547 if ( vpqr.is_degenerate_Voronoi_circle() ) { 00548 Site_2 tpqr = Site_2::construct_site_2(vpqr.degenerate_point()); 00549 00550 if ( same_points(tpqr, p.source_site()) || 00551 same_points(tpqr, p.target_site()) ) { 00552 if ( vqps.oriented_side(lp) != ON_POSITIVE_SIDE ) { 00553 lp = opposite_line(lp); 00554 } 00555 } 00556 } else { 00557 if ( vpqr.oriented_side(lp) != ON_POSITIVE_SIDE ) { 00558 lp = opposite_line(lp); 00559 } 00560 } 00561 #if 0 // OLD CODE 00562 if ( ( vpqr.is_degenerate_Voronoi_circle() && 00563 same_points(vpqr.degenerate_point(), p.source_site()) ) || 00564 ( vpqr.is_degenerate_Voronoi_circle() && 00565 same_points(vpqr.degenerate_point(), p.target_site()) ) ) { 00566 if ( vqps.oriented_side(lp) != ON_POSITIVE_SIDE ) { 00567 lp = opposite_line(lp); 00568 } 00569 } else { 00570 if ( vpqr.oriented_side(lp) != ON_POSITIVE_SIDE ) { 00571 lp = opposite_line(lp); 00572 } 00573 } 00574 #endif 00575 00576 // then orient lq according to its Voronoi vertices 00577 if ( vpqr.is_degenerate_Voronoi_circle() ) { 00578 Site_2 tpqr = Site_2::construct_site_2(vpqr.degenerate_point()); 00579 00580 if ( same_points(tpqr, q.source_site()) || 00581 same_points(tpqr, q.target_site()) ) { 00582 if ( vqps.oriented_side(lq) != ON_POSITIVE_SIDE ) { 00583 lq = opposite_line(lq); 00584 } 00585 } 00586 } else { 00587 if ( vpqr.oriented_side(lq) != ON_POSITIVE_SIDE ) { 00588 lq = opposite_line(lq); 00589 } 00590 } 00591 #if 0 // OLD CODE 00592 if ( ( vpqr.is_degenerate_Voronoi_circle() && 00593 same_points(vpqr.degenerate_point(), q.source_site()) ) || 00594 ( vpqr.is_degenerate_Voronoi_circle() && 00595 same_points(vpqr.degenerate_point(), q.target_site()) ) ) { 00596 if ( vqps.oriented_side(lq) != ON_POSITIVE_SIDE ) { 00597 lq = opposite_line(lq); 00598 } 00599 } else { 00600 if ( vpqr.oriented_side(lq) != ON_POSITIVE_SIDE ) { 00601 lq = opposite_line(lq); 00602 } 00603 } 00604 #endif 00605 00606 Point_2 tp = t.point(); 00607 00608 // check if t is on the same side as the Voronoi vertices 00609 Oriented_side ot_lp = oriented_side_of_line(lp, tp); 00610 Oriented_side ot_lq = oriented_side_of_line(lq, tp); 00611 00612 if ( ot_lp != ON_POSITIVE_SIDE || ot_lq != ON_POSITIVE_SIDE ) { 00613 return false; 00614 } 00615 00616 Line_2 lperp; 00617 00618 Comparison_result res = 00619 compare_squared_distances_to_lines(tp, lp, lq); 00620 00621 if ( res == SMALLER ) { 00622 lperp = compute_perpendicular(lp, tp); 00623 } else { 00624 lperp = compute_perpendicular(lq, tp); 00625 } 00626 00627 CGAL_precondition( ot_lp != ON_ORIENTED_BOUNDARY && 00628 ot_lq != ON_ORIENTED_BOUNDARY ); 00629 00630 // check of lperp separates the two Voronoi vertices 00631 Oriented_side opqr_perp = vpqr.oriented_side(lperp); 00632 Oriented_side oqps_perp = vqps.oriented_side(lperp); 00633 00634 bool on_different_side = 00635 (opqr_perp == ON_POSITIVE_SIDE && 00636 oqps_perp == ON_NEGATIVE_SIDE) || 00637 (opqr_perp == ON_NEGATIVE_SIDE && 00638 oqps_perp == ON_POSITIVE_SIDE); 00639 00640 return ( on_different_side ); 00641 } 00642 00643 //------------------------------------------------------------------------ 00644 //------------------------------------------------------------------------ 00645 //------------------------------------------------------------------------ 00646 00647 public: 00648 typedef Boolean result_type; 00649 typedef Site_2 argument_type; 00650 00651 Boolean operator()(const Site_2& p, const Site_2& q, const Site_2& r, 00652 const Site_2& s, const Site_2& t, Sign sgn) const 00653 { 00654 if ( sgn == POSITIVE ) { 00655 return is_interior_in_conflict_none(p, q, r, s, t, Method_tag()); 00656 } else if ( sgn == NEGATIVE ) { 00657 return is_interior_in_conflict_both(p, q, r, s, t, Method_tag()); 00658 } else { 00659 return is_interior_in_conflict_touch(p, q, r, s, t, Method_tag()); 00660 } 00661 } 00662 00663 00664 Boolean operator()(const Site_2& p, const Site_2& q, const Site_2& , 00665 const Site_2& t, Sign sgn) const 00666 { 00667 if ( t.is_point() ) { 00668 return ( sgn == NEGATIVE ); 00669 } 00670 00671 if ( sgn != NEGATIVE ) { 00672 return false; 00673 } 00674 00675 if ( p.is_segment() || q.is_segment() ) { 00676 return false; 00677 } 00678 00679 bool p_is_endpoint = 00680 same_points(p, t.source_site()) || same_points(p, t.target_site()); 00681 bool q_is_endpoint = 00682 same_points(q, t.source_site()) || same_points(q, t.target_site()); 00683 00684 return ( p_is_endpoint && q_is_endpoint ); 00685 } 00686 00687 Boolean operator()(const Site_2& p, const Site_2& q, const Site_2& t, 00688 Sign ) const 00689 { 00690 if ( p.is_segment() || q.is_segment()) { 00691 return false; 00692 } 00693 00694 // both p and q are points 00695 if ( t.is_point() ) { 00696 RT dtpx = p.point().x() - t.point().x(); 00697 RT minus_dtpy = -p.point().y() + t.point().y(); 00698 RT dtqx = q.point().x() - t.point().x(); 00699 RT dtqy = q.point().y() - t.point().y(); 00700 00701 Sign s1 = sign_of_determinant(dtpx, minus_dtpy, dtqy, dtqx); 00702 00703 CGAL_assertion( s1 != ZERO ); 00704 return ( s1 == NEGATIVE ); 00705 } 00706 00707 bool bp = 00708 same_points(p, t.source_site()) || same_points(p, t.target_site()); 00709 bool bq = 00710 same_points(q, t.source_site()) || same_points(q, t.target_site()); 00711 00712 00713 return ( bp && bq ); 00714 } 00715 00716 }; 00717 00718 00719 //----------------------------------------------------------------------------- 00720 00721 CGAL_SEGMENT_DELAUNAY_GRAPH_2_END_NAMESPACE 00722 00723 CGAL_END_NAMESPACE 00724 00725 00726 #if defined(BOOST_MSVC) 00727 # pragma warning(pop) 00728 #endif 00729 00730 00731 #endif // CGAL_SEGMENT_DELAUNAY_GRAPH_2_FINITE_EDGE_INTERIOR_CONFLICT_C2_H
1.7.6.1