|
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/Voronoi_vertex_ring_C2.h $ 00015 // $Id: Voronoi_vertex_ring_C2.h 42809 2008-04-09 13:09:17Z spion $ 00016 // 00017 // 00018 // Author(s) : Menelaos Karavelas <mkaravel@cse.nd.edu> 00019 00020 00021 00022 00023 #ifndef CGAL_SEGMENT_DELAUNAY_GRAPH_2_VORONOI_VERTEX_RING_C2_H 00024 #define CGAL_SEGMENT_DELAUNAY_GRAPH_2_VORONOI_VERTEX_RING_C2_H 00025 00026 00027 #include <CGAL/Segment_Delaunay_graph_2/Basic_predicates_C2.h> 00028 #include <CGAL/Segment_Delaunay_graph_2/Are_same_points_C2.h> 00029 #include <CGAL/Segment_Delaunay_graph_2/Are_same_segments_C2.h> 00030 00031 00032 CGAL_BEGIN_NAMESPACE 00033 00034 CGAL_SEGMENT_DELAUNAY_GRAPH_2_BEGIN_NAMESPACE 00035 00036 00037 template<class K> 00038 class Voronoi_vertex_ring_C2 00039 : public Basic_predicates_C2<K> 00040 { 00041 public: 00042 typedef Basic_predicates_C2<K> Base; 00043 00044 typedef enum {PPP = 0, PPS, PSS, SSS} vertex_t; 00045 struct PPP_Type {}; 00046 struct PPS_Type {}; 00047 struct PSS_Type {}; 00048 struct SSS_Type {}; 00049 00050 typedef typename Base::Point_2 Point_2; 00051 typedef typename Base::Segment_2 Segment_2; 00052 typedef typename Base::Line_2 Line_2; 00053 typedef typename Base::Site_2 Site_2; 00054 typedef typename Base::FT FT; 00055 typedef typename Base::RT RT; 00056 00057 typedef typename Base::Sqrt_1 Sqrt_1; 00058 typedef typename Base::Sqrt_2 Sqrt_2; 00059 typedef typename Base::Sqrt_3 Sqrt_3; 00060 00061 typedef typename Base::Homogeneous_point_2 Homogeneous_point_2; 00062 00063 typedef typename Base::Orientation Orientation; 00064 typedef typename Base::Comparison_result Comparison_result; 00065 typedef typename Base::Oriented_side Oriented_side; 00066 typedef typename Base::Sign Sign; 00067 00068 private: 00069 typedef Are_same_points_C2<K> Are_same_points_2; 00070 typedef Are_same_segments_C2<K> Are_same_segments_2; 00071 00072 typedef typename K::Intersections_tag ITag; 00073 00074 Are_same_points_2 same_points; 00075 Are_same_segments_2 same_segments; 00076 00077 private: 00078 //-------------------------------------------------------------------------- 00079 00080 void 00081 compute_ppp(const Site_2& sp, const Site_2& sq, const Site_2& sr) 00082 { 00083 CGAL_precondition( sp.is_point() && sq.is_point() && 00084 sr.is_point() ); 00085 00086 Point_2 p = sp.point(), q = sq.point(), r = sr.point(); 00087 00088 v_type = PPP; 00089 00090 RT np = CGAL::square(p.x()) + CGAL::square(p.y()); 00091 RT nq = CGAL::square(q.x()) + CGAL::square(q.y()); 00092 RT nr = CGAL::square(r.x()) + CGAL::square(r.y()); 00093 00094 ux_ppp = 00095 np * (q.y() - r.y()) + nq * (r.y() - p.y()) + nr * (p.y() - q.y()); 00096 uy_ppp = 00097 -(np * (q.x() - r.x()) + nq * (r.x() - p.x()) + nr * (p.x() - q.x())); 00098 uz_ppp = RT(2) * ( (q.x() * r.y() - r.x() * q.y()) + 00099 (r.x() * p.y() - p.x() * r.y()) + 00100 (p.x() * q.y() - q.x() * p.y()) ); 00101 } 00102 00103 //-------------------------------------------------------------------------- 00104 00105 void 00106 compute_pss(const Site_2& p, const Site_2& q, const Site_2& r) 00107 { 00108 CGAL_precondition( p.is_point() && q.is_segment() && 00109 r.is_segment() ); 00110 00111 v_type = PSS; 00112 00113 bool pq = 00114 same_points(p, q.source_site()) || same_points(p, q.target_site()); 00115 bool pr = 00116 same_points(p, r.source_site()) || same_points(p, r.target_site()); 00117 00118 Point_2 pp = p.point(); 00119 00120 if ( pq && pr ) { 00121 Sqrt_1 One(RT(1), RT(0), RT(0)); 00122 00123 ux = Sqrt_3(pp.x() * One); 00124 uy = Sqrt_3(pp.y() * One); 00125 uz = Sqrt_3(One); 00126 return; 00127 } 00128 00129 00130 00131 RT a1, b1, c1, a2, b2, c2; 00132 compute_supporting_line(q.supporting_site(), a1, b1, c1); 00133 compute_supporting_line(r.supporting_site(), a2, b2, c2); 00134 00135 RT c1_ = a1 * pp.x() + b1 * pp.y() + c1; 00136 RT c2_ = a2 * pp.x() + b2 * pp.y() + c2; 00137 00138 if ( pq ) { 00139 c1_ = RT(0); 00140 } 00141 00142 if ( pr ) { 00143 c2_ = RT(0); 00144 } 00145 00146 Sign sgn_c1_ = CGAL::sign(c1_); 00147 Sign sgn_c2_ = CGAL::sign(c2_); 00148 00149 if ( sgn_c1_ == NEGATIVE ) { 00150 a1 = -a1; b1 = -b1; c1_ = -c1_; 00151 } else if ( sgn_c1_ == ZERO ) { 00152 CGAL_assertion( pq ); 00153 if ( same_points(p, q.target_site()) ) { 00154 a1 = -a1; b1 = -b1; c1_ = -c1_; 00155 } 00156 } 00157 00158 if ( sgn_c2_ == NEGATIVE ) { 00159 a2 = -a2; b2 = -b2; c2_ = -c2_; 00160 } else if ( sgn_c2_ == ZERO ) { 00161 CGAL_assertion( pr ); 00162 if ( same_points(p, r.source_site()) ) { 00163 a2 = -a2; b2 = -b2; c2_ = -c2_; 00164 } 00165 } 00166 00167 if ( pq ) { 00168 RT J = a1 * c2_; 00169 RT I = b1 * c2_; 00170 00171 RT n1 = CGAL::square(a1) + CGAL::square(b1); 00172 RT n2 = CGAL::square(a2) + CGAL::square(b2); 00173 00174 RT D1D2 = n1 * n2; 00175 00176 Sqrt_1 Zero(RT(0), RT(0), D1D2); 00177 00178 Sqrt_1 vz(-a1 * a2 - b1 * b2, RT(1), D1D2); 00179 00180 ux = Sqrt_3(J + pp.x() * vz, 00181 Zero, Zero, Zero, Zero, Zero); 00182 uy = Sqrt_3(I + pp.y() * vz, 00183 Zero, Zero, Zero, Zero, Zero); 00184 uz = Sqrt_3(vz, Zero, Zero, Zero, Zero, Zero); 00185 } else if ( pr ) { 00186 RT J = a2 * c1_; 00187 RT I = b2 * c1_; 00188 00189 RT n1 = CGAL::square(a1) + CGAL::square(b1); 00190 RT n2 = CGAL::square(a2) + CGAL::square(b2); 00191 00192 RT D1D2 = n1 * n2; 00193 00194 Sqrt_1 Zero(RT(0), RT(0), D1D2); 00195 00196 Sqrt_1 vz(-a1 * a2 - b1 * b2, RT(1), D1D2); 00197 00198 ux = Sqrt_3(J + pp.x() * vz, 00199 Zero, Zero, Zero, Zero, Zero); 00200 uy = Sqrt_3(I + pp.y() * vz, 00201 Zero, Zero, Zero, Zero, Zero); 00202 uz = Sqrt_3(vz, Zero, Zero, Zero, Zero, Zero); 00203 } else { 00204 Line_2 lq(a1, b1, c1_); 00205 Line_2 lr(a2, b2, c2_); 00206 compute_pll(pp, lq, lr); 00207 } 00208 } 00209 00210 00211 void 00212 compute_pll(const Point_2& p, const Line_2& lq, const Line_2& lr) 00213 { 00214 RT a1 = lq.a(), b1 = lq.b(), c1_ = lq.c(); 00215 RT a2 = lr.a(), b2 = lr.b(), c2_ = lr.c(); 00216 00217 CGAL_precondition( c1_ >= RT(0) ); 00218 CGAL_precondition( c2_ >= RT(0) ); 00219 00220 RT n1 = CGAL::square(a1) + CGAL::square(b1); 00221 RT n2 = CGAL::square(a2) + CGAL::square(b2); 00222 00223 RT I = b1 * c2_ + b2 * c1_; 00224 RT J = a1 * c2_ + a2 * c1_; 00225 00226 RT c1c2 = RT(2) * c1_ * c2_; 00227 RT a1a2 = a1 * a2; 00228 RT b1b2 = b1 * b2; 00229 00230 RT D1D2 = n1 * n2; 00231 00232 Sqrt_1 Zero(RT(0), RT(0), D1D2); 00233 Sqrt_1 One(RT(1), RT(0), D1D2); 00234 00235 00236 Sign s1, s2; 00237 00238 // compute sigma 00239 Sign s_sigma(ZERO); 00240 s1 = CGAL::sign(b1); 00241 s2 = CGAL::sign(-b2); 00242 if ( s1 == ZERO ) { 00243 s_sigma = s2; 00244 } else if ( s2 == ZERO ) { 00245 s_sigma = s1; 00246 } else if ( s1 == s2 ) { 00247 s_sigma = s1; 00248 } else { 00249 RT e = CGAL::square(b1) * n2 - CGAL::square(b2) * n1; 00250 s_sigma = s1 * CGAL::sign(e); 00251 } 00252 00253 Sqrt_1 sigma = Zero; 00254 if ( s_sigma == POSITIVE ) { 00255 sigma = One; 00256 } else if ( s_sigma == NEGATIVE ) { 00257 sigma = -One; 00258 } 00259 00260 // compute rho 00261 Sign s_rho(ZERO); 00262 s1 = CGAL::sign(a1); 00263 s2 = CGAL::sign(-a2); 00264 if ( s1 == ZERO ) { 00265 s_rho = s2; 00266 } else if ( s2 == ZERO ) { 00267 s_rho = s1; 00268 } else if ( s1 == s2 ) { 00269 s_rho = s1; 00270 } else { 00271 RT e = CGAL::square(a1) * n2 - CGAL::square(a2) * n1; 00272 s_rho = s1 * CGAL::sign(e); 00273 } 00274 00275 00276 Sqrt_1 rho = Zero; 00277 if ( s_rho == POSITIVE ) { 00278 rho = One; 00279 } else if ( s_rho == NEGATIVE ) { 00280 rho = -One; 00281 } 00282 00283 00284 Sqrt_1 vz(-a1a2 - b1b2, RT(1), D1D2); 00285 00286 RT A = a1a2 - b1b2; 00287 Sqrt_1 u1( c1c2 * A, c1c2, D1D2); 00288 Sqrt_1 u2(-c1c2 * A, c1c2, D1D2); 00289 00290 00291 ux = Sqrt_3(J + p.x() * vz, sigma, Zero, Zero, u1, u2); 00292 uy = Sqrt_3(I + p.y() * vz, Zero, -rho, Zero, u1, u2); 00293 uz = Sqrt_3(vz, Zero, Zero, Zero, u1, u2); 00294 } 00295 00296 00297 //-------------------------------------------------------------------------- 00298 00299 00300 void 00301 compute_pps(const Site_2& p, const Site_2& q, const Site_2& r) 00302 { 00303 CGAL_precondition( p.is_point() && q.is_point() && 00304 r.is_segment() ); 00305 00306 v_type = PPS; 00307 00308 RT a, b, c; 00309 compute_supporting_line(r.supporting_site(), a, b, c); 00310 00311 Point_2 pp = p.point(), qq = q.point(); 00312 00313 RT c_ = a * pp.x() + b * pp.y() + c; 00314 RT cq_ = a * qq.x() + b * qq.y() + c; 00315 00316 00317 if ( same_points(p, r.source_site()) || 00318 same_points(p, r.target_site()) ) { 00319 c_ = RT(0); 00320 } 00321 if ( same_points(q, r.source_site()) || 00322 same_points(q, r.target_site()) ) { 00323 cq_ = RT(0); 00324 } 00325 00326 Sign s = CGAL::sign(c_); 00327 00328 if ( s == NEGATIVE ) { 00329 a = -a; b = -b; c = -c; c_ = -c_; cq_ = -cq_; 00330 } else if ( s == ZERO ) { 00331 Sign s1 = CGAL::sign(cq_); 00332 00333 CGAL_assertion( s1 != ZERO ); 00334 if ( s1 == NEGATIVE ) { 00335 a = -a; b = -b; c = -c; c_ = -c_; cq_ = -cq_; 00336 } 00337 } 00338 00339 RT nl = CGAL::square(a) + CGAL::square(b); 00340 00341 RT x_ = qq.x() - pp.x(); 00342 RT y_ = qq.y() - pp.y(); 00343 RT n_ = CGAL::square(x_) + CGAL::square(y_); 00344 00345 00346 Comparison_result res = CGAL::compare( c_, cq_ ); 00347 00348 if ( res == EQUAL ) { 00349 RT e1 = CGAL::square(c_); 00350 RT J = nl * (a * n_ + RT(4) * c_ * x_) - RT(4) * a * e1; 00351 RT I = nl * (b * n_ + RT(4) * c_ * y_) - RT(4) * b * e1; 00352 RT X = RT(8) * nl * c_; 00353 00354 ux_pps = Sqrt_1(J + pp.x() * X); 00355 uy_pps = Sqrt_1(I + pp.y() * X); 00356 uz_pps = Sqrt_1(X); 00357 return; 00358 } 00359 00360 00361 RT e1 = a * x_ + b * y_; 00362 RT e2 = b * x_ - a * y_; 00363 RT e3 = n_ * e1; 00364 RT e4 = RT(2) * c_ * e2; 00365 00366 RT X = RT(2) * CGAL::square(e1); 00367 RT I = b * e3 + x_ * e4; 00368 RT J = a * e3 - y_ * e4; 00369 RT S = n_ * nl * c_ * cq_; 00370 00371 ux_pps = Sqrt_1(J + pp.x() * X, RT(-2) * y_, S); 00372 uy_pps = Sqrt_1(I + pp.y() * X, RT( 2) * x_, S); 00373 uz_pps = Sqrt_1(X, RT(0), S); 00374 } 00375 00376 //-------------------------------------------------------------------------- 00377 00378 bool check_if_exact(const Site_2& , unsigned int , 00379 const Tag_false&) const 00380 { 00381 return true; 00382 } 00383 00384 bool check_if_exact(const Site_2& s, unsigned int i, 00385 const Tag_true&) const 00386 { 00387 return s.is_input(i); 00388 } 00389 00390 // determines of the segment s is on the positive halfspace as 00391 // defined by the supporting line of the segment supp; the line l 00392 // is supposed to be the supporting line of the segment supp and we 00393 // pass it so that we do not have to recompute it 00394 bool 00395 is_on_positive_halfspace(const Site_2& supp, 00396 const Site_2& s, const Line_2& l) const 00397 { 00398 CGAL_precondition( supp.is_segment() && s.is_segment() ); 00399 00400 if ( same_segments(supp.supporting_site(), 00401 s.supporting_site()) ) { 00402 return false; 00403 } 00404 00405 if ( same_points(supp.source_site(), s.source_site()) || 00406 same_points(supp.target_site(), s.source_site()) ) { 00407 return oriented_side_of_line(l, s.target()) == ON_POSITIVE_SIDE; 00408 } 00409 00410 if ( same_points(supp.source_site(), s.target_site()) || 00411 same_points(supp.target_site(), s.target_site()) ) { 00412 return oriented_side_of_line(l, s.source()) == ON_POSITIVE_SIDE; 00413 } 00414 00415 ITag itag; 00416 00417 if ( !check_if_exact(s, 0, itag) && 00418 same_segments(supp.supporting_site(), 00419 s.crossing_site(0)) ) { 00420 return oriented_side_of_line(l, s.target()) == ON_POSITIVE_SIDE; 00421 } 00422 00423 if ( !check_if_exact(s, 1, itag) && 00424 same_segments(supp.supporting_site(), 00425 s.crossing_site(1)) ) { 00426 return oriented_side_of_line(l, s.source()) == ON_POSITIVE_SIDE; 00427 } 00428 00429 return Base::is_on_positive_halfspace(l, s.segment()); 00430 } 00431 00432 //-------------------------------------------------------------------------- 00433 00434 void 00435 orient_lines(const Site_2& p, const Site_2& q, 00436 const Site_2& r, RT a[], RT b[], RT c[]) const 00437 { 00438 CGAL_precondition( p.is_segment() && q.is_segment() && 00439 r.is_segment() ); 00440 00441 Line_2 l[3]; 00442 l[0] = compute_supporting_line(p.supporting_site()); 00443 l[1] = compute_supporting_line(q.supporting_site()); 00444 l[2] = compute_supporting_line(r.supporting_site()); 00445 00446 bool is_oriented[3] = {false, false, false}; 00447 00448 if ( is_on_positive_halfspace(p, q, l[0]) || 00449 is_on_positive_halfspace(p, r, l[0]) ) { 00450 is_oriented[0] = true; 00451 } else { 00452 l[0] = opposite_line(l[0]); 00453 if ( is_on_positive_halfspace(p, q, l[0]) || 00454 is_on_positive_halfspace(p, r, l[0]) ) { 00455 is_oriented[0] = true; 00456 } else { 00457 l[0] = opposite_line(l[0]); 00458 } 00459 } 00460 00461 if ( is_on_positive_halfspace(q, p, l[1]) || 00462 is_on_positive_halfspace(q, r, l[1]) ) { 00463 is_oriented[1] = true; 00464 } else { 00465 l[1] = opposite_line(l[1]); 00466 if ( is_on_positive_halfspace(q, p, l[1]) || 00467 is_on_positive_halfspace(q, r, l[1]) ) { 00468 is_oriented[1] = true; 00469 } else { 00470 l[1] = opposite_line(l[1]); 00471 } 00472 } 00473 00474 if ( is_on_positive_halfspace(r, p, l[2]) || 00475 is_on_positive_halfspace(r, q, l[2]) ) { 00476 is_oriented[2] = true; 00477 } else { 00478 l[2] = opposite_line(l[2]); 00479 if ( is_on_positive_halfspace(r, p, l[2]) || 00480 is_on_positive_halfspace(r, q, l[2]) ) { 00481 is_oriented[2] = true; 00482 } else { 00483 l[2] = opposite_line(l[2]); 00484 } 00485 } 00486 00487 if ( is_oriented[0] && is_oriented[1] && is_oriented[2] ) { 00488 for (int i = 0; i < 3; i++) { 00489 a[i] = l[i].a(); 00490 b[i] = l[i].b(); 00491 c[i] = l[i].c(); 00492 } 00493 return; 00494 } 00495 00496 int i_no(-1); 00497 for (int i = 0; i < 3; i++) { 00498 if ( !is_oriented[i] ) { 00499 i_no = i; 00500 CGAL_assertion( is_oriented[(i+1)%3] && is_oriented[(i+2)%3] ); 00501 break; 00502 } 00503 } 00504 00505 CGAL_assertion( i_no != -1 ); 00506 00507 RT d[3]; 00508 for (int i = 0; i < 3; i++) { 00509 d[i] = CGAL::square(l[i].a()) + CGAL::square(l[i].b()); 00510 } 00511 00512 RT z[3]; 00513 for (int i = 0; i < 3; i++) { 00514 z[i] = l[(i+1)%3].a() * l[(i+2)%3].b() 00515 - l[(i+2)%3].a() * l[(i+1)%3].b(); 00516 } 00517 00518 00519 Sqrt_1 Zero(RT(0), RT(0), d[0]); 00520 Sqrt_1 sqrt_D0(RT(0), RT(1), d[0]); 00521 00522 Sqrt_1 D1 = d[1] + Zero; 00523 Sqrt_1 D2 = d[2] + Zero; 00524 00525 Sqrt_3 vz(z[0] * sqrt_D0, z[1] + Zero, z[2] + Zero, Zero, D1, D2); 00526 00527 Sign s_minus_vz = CGAL::sign(vz); 00528 00529 00530 CGAL_assertion( s_minus_vz != ZERO ); 00531 00532 if ( s_minus_vz == NEGATIVE ) { 00533 l[i_no] = opposite_line(l[i_no]); 00534 00535 for (int i = 0; i < 3; i++) { 00536 a[i] = l[i].a(); 00537 b[i] = l[i].b(); 00538 c[i] = l[i].c(); 00539 } 00540 return; 00541 } 00542 00543 // now we have to check if the other orientation of l[i_no] 00544 // corresponds to a CCW triangle as well. 00545 z[(i_no+1)%3] = -z[(i_no+1)%3]; 00546 z[(i_no+2)%3] = -z[(i_no+2)%3]; 00547 00548 vz = Sqrt_3(z[0] * sqrt_D0, z[1] + Zero, z[2] + Zero, Zero, D1, D2); 00549 00550 Sign s_minus_vz_2 = CGAL::sign(vz); 00551 00552 CGAL_assertion( s_minus_vz_2 != ZERO ); 00553 00554 if ( s_minus_vz_2 == NEGATIVE ) { 00555 // the other orientation does not correspond to a CCW triangle. 00556 for (int i = 0; i < 3; i++) { 00557 a[i] = l[i].a(); 00558 b[i] = l[i].b(); 00559 c[i] = l[i].c(); 00560 } 00561 return; 00562 } 00563 00564 // now compute the Voronoi vertex; 00565 for (int i = 0; i < 3; i++) { 00566 a[i] = l[i].a(); 00567 b[i] = l[i].b(); 00568 c[i] = l[i].c(); 00569 } 00570 00571 RT x[3], y[3], w[3]; 00572 for (int i = 0; i < 3; i++) { 00573 x[i] = c[(i+1)%3] * b[(i+2)%3] - c[(i+2)%3] * b[(i+1)%3]; 00574 y[i] = -(c[(i+1)%3] * a[(i+2)%3] - c[(i+2)%3] * a[(i+1)%3]); 00575 w[i] = -(a[(i+1)%3] * b[(i+2)%3] - a[(i+2)%3] * b[(i+1)%3]); 00576 } 00577 00578 Sqrt_3 vx, vy, vw; 00579 00580 vx = Sqrt_3(x[0] * sqrt_D0, x[1] + Zero, x[2] + Zero, Zero, D1, D2); 00581 vy = Sqrt_3(y[0] * sqrt_D0, y[1] + Zero, y[2] + Zero, Zero, D1, D2); 00582 vw = Sqrt_3(w[0] * sqrt_D0, w[1] + Zero, w[2] + Zero, Zero, D1, D2); 00583 00584 Sqrt_1 a1(a[(i_no+1)%3], RT(0), d[0]); 00585 Sqrt_1 b1(b[(i_no+1)%3], RT(0), d[0]); 00586 Sqrt_1 c1(c[(i_no+1)%3], RT(0), d[0]); 00587 00588 Sqrt_3 dist = a1 * vx + b1 * vy + c1 * vw; 00589 00590 Sign s_vw = CGAL::sign(vw); 00591 00592 Sign sgn_dist = s_vw * CGAL::sign(dist); 00593 00594 CGAL_assertion( sgn_dist != ZERO ); 00595 00596 if ( sgn_dist == NEGATIVE ) { 00597 a[i_no] = -a[i_no]; 00598 b[i_no] = -b[i_no]; 00599 c[i_no] = -c[i_no]; 00600 } 00601 } 00602 00603 00604 void 00605 compute_sss(const Site_2& p, const Site_2& q, const Site_2& r) 00606 { 00607 CGAL_precondition( p.is_segment() && q.is_segment() && 00608 r.is_segment() ); 00609 00610 v_type = SSS; 00611 00612 RT a[3], b[3], c[3]; 00613 RT cx[3], cy[3], cz[3], D[3]; 00614 00615 orient_lines(p, q, r, a, b, c); 00616 00617 for (int i = 0; i < 3; i++) { 00618 cx[i] = c[(i+1)%3] * b[(i+2)%3] - c[(i+2)%3] * b[(i+1)%3]; 00619 cy[i] = -(c[(i+1)%3] * a[(i+2)%3] - c[(i+2)%3] * a[(i+1)%3]); 00620 cz[i] = -(a[(i+1)%3] * b[(i+2)%3] - a[(i+2)%3] * b[(i+1)%3]); 00621 D[i] = CGAL::square(a[i]) + CGAL::square(b[i]); 00622 00623 Sqrt_1 Zero(RT(0), RT(0), D[0]); 00624 Sqrt_1 sqrt_D0(RT(0), RT(1), D[0]); 00625 Sqrt_1 D1 = Zero + D[1]; 00626 Sqrt_1 D2 = Zero + D[2]; 00627 00628 ux = Sqrt_3(cx[0] * sqrt_D0, cx[1] + Zero, cx[2] + Zero, 00629 Zero, D1, D2); 00630 uy = Sqrt_3(cy[0] * sqrt_D0, cy[1] + Zero, cy[2] + Zero, 00631 Zero, D1, D2); 00632 uz = Sqrt_3(cz[0] * sqrt_D0, cz[1] + Zero, cz[2] + Zero, 00633 Zero, D1, D2); 00634 } 00635 } 00636 00637 //-------------------------------------------------------------------------- 00638 00639 void 00640 compute_vertex(const Site_2& s1, const Site_2& s2, const Site_2& s3) 00641 { 00642 if ( s1.is_point() && s2.is_point() && s3.is_point() ) { 00643 compute_ppp(s1, s2, s3); 00644 00645 } else if ( s1.is_segment() && s2.is_point() && s3.is_point() ) { 00646 compute_vertex(s2, s3, s1); 00647 pps_idx = 1; 00648 00649 } else if ( s1.is_point() && s2.is_segment() && s3.is_point() ) { 00650 compute_vertex(s3, s1, s2); 00651 pps_idx = 2; 00652 00653 } else if ( s1.is_point() && s2.is_point() && s3.is_segment() ) { 00654 compute_pps(s1, s2, s3); 00655 pps_idx = 0; 00656 00657 } else if ( s1.is_point() && s2.is_segment() && s3.is_segment() ) { 00658 compute_pss(s1, s2, s3); 00659 } else if ( s1.is_segment() && s2.is_point() && s3.is_segment() ) { 00660 compute_vertex(s2, s3, s1); 00661 } else if ( s1.is_segment() && s2.is_segment() && s3.is_point() ) { 00662 compute_vertex(s3, s1, s2); 00663 } else { 00664 compute_sss(s1, s2, s3); 00665 } 00666 00667 } 00668 00669 //-------------------------------------------------------------------------- 00670 00671 bool is_endpoint_of(const Site_2& p, const Site_2& s) const 00672 { 00673 CGAL_precondition( p.is_point() && s.is_segment() ); 00674 return ( same_points(p, s.source_site()) || 00675 same_points(p, s.target_site()) ); 00676 } 00677 00678 00679 //-------------------------------------------------------------------------- 00680 //-------------------------------------------------------------------------- 00681 // the orientation test 00682 //-------------------------------------------------------------------------- 00683 //-------------------------------------------------------------------------- 00684 00685 Orientation 00686 orientation(const Line_2& l, PPP_Type) const 00687 { 00688 Sign s_uz = CGAL::sign(uz_ppp); 00689 Sign s_l = 00690 CGAL::sign(l.a() * ux_ppp + l.b() * uy_ppp + l.c() * uz_ppp); 00691 00692 return s_uz * s_l; 00693 } 00694 00695 00696 //-------------------------------------------------------------------------- 00697 00698 Orientation 00699 orientation(const Line_2& l, PPS_Type) const 00700 { 00701 Sign s_uz = CGAL::sign(uz_pps); 00702 Sign s_l = CGAL::sign(l.a() * ux_pps + l.b() * uy_pps + l.c() * uz_pps); 00703 00704 return s_uz * s_l; 00705 } 00706 00707 00708 //-------------------------------------------------------------------------- 00709 00710 // the cases PSS and SSS are identical 00711 template<class Type> 00712 Orientation 00713 orientation(const Line_2& l, Type) const 00714 { 00715 Sqrt_1 Zero(RT(0), RT(0), ux.a().c()); 00716 00717 Sqrt_1 a = l.a() + Zero; 00718 Sqrt_1 b = l.b() + Zero; 00719 Sqrt_1 c = l.c() + Zero; 00720 00721 Sign s_uz = CGAL::sign(uz); 00722 Sign s_l = CGAL::sign(a * ux + b * uy + c * uz); 00723 00724 return s_uz * s_l; 00725 } 00726 00727 00728 //-------------------------------------------------------------------------- 00729 //-------------------------------------------------------------------------- 00730 // the incircle test 00731 //-------------------------------------------------------------------------- 00732 //-------------------------------------------------------------------------- 00733 00734 //-------------------------------------------------------------------------- 00735 // the incircle test when the fourth site is a point 00736 //-------------------------------------------------------------------------- 00737 00738 //-------------------------------------------------------------------------- 00739 00740 Sign check_easy_degeneracies(const Site_2& t, PPS_Type, 00741 bool& use_result) const 00742 { 00743 CGAL_precondition( t.is_point() ); 00744 00745 use_result = false; 00746 if ( ( p_.is_point() && same_points(p_, t) ) || 00747 ( q_.is_point() && same_points(q_, t) ) || 00748 ( r_.is_point() && same_points(r_, t) ) ) { 00749 use_result = true; 00750 return ZERO; 00751 } 00752 00753 if ( ( p_.is_segment() && is_endpoint_of(t, p_) ) || 00754 ( q_.is_segment() && is_endpoint_of(t, q_) ) || 00755 ( r_.is_segment() && is_endpoint_of(t, r_) ) ) { 00756 use_result = true; 00757 return POSITIVE; 00758 } 00759 00760 return ZERO; 00761 } 00762 00763 inline 00764 Sign check_easy_degeneracies(const Site_2& t, PSS_Type, 00765 bool& use_result) const 00766 { 00767 CGAL_precondition( t.is_point() ); 00768 return check_easy_degeneracies(t, PPS_Type(), use_result); 00769 } 00770 00771 inline 00772 Sign check_easy_degeneracies(const Site_2& t, SSS_Type, 00773 bool& use_result) const 00774 { 00775 CGAL_precondition( t.is_point() ); 00776 use_result = false; 00777 // ADD THE CASES WHERE t IS AN ENDPOINT OF ONE OF THE SEGMENTS 00778 return ZERO; 00779 } 00780 00781 //-------------------------------------------------------------------------- 00782 00783 template<class Type> 00784 inline 00785 Sign incircle_p(const Site_2& st, Type type) const 00786 { 00787 CGAL_precondition( st.is_point() ); 00788 00789 bool use_result(false); 00790 Sign s = check_easy_degeneracies(st, type, use_result); 00791 if ( use_result ) { return s; } 00792 00793 return incircle_p_no_easy(st, type); 00794 } 00795 00796 //-------------------------------------------------------------------------- 00797 00798 Sign incircle_p(const Site_2& st, PPP_Type) const 00799 { 00800 CGAL_precondition( st.is_point() ); 00801 00802 Point_2 t = st.point(); 00803 00804 return - side_of_oriented_circle(p_.point(), q_.point(), r_.point(), t); 00805 } 00806 00807 //-------------------------------------------------------------------------- 00808 00809 Sign incircle_p_no_easy(const Site_2& st, PPS_Type ) const 00810 { 00811 CGAL_precondition( st.is_point() ); 00812 00813 Point_2 t = st.point(); 00814 00815 Point_2 pref = p_ref().point(); 00816 00817 Sqrt_1 vx = ux_pps - pref.x() * uz_pps; 00818 Sqrt_1 vy = uy_pps - pref.y() * uz_pps; 00819 00820 Sqrt_1 Rs = CGAL::square(vx) + CGAL::square(vy); 00821 00822 Sqrt_1 Rs1 = CGAL::square(ux_pps - t.x() * uz_pps) 00823 + CGAL::square(uy_pps - t.y() * uz_pps); 00824 00825 return CGAL::sign(Rs1 - Rs); 00826 } 00827 00828 //-------------------------------------------------------------------------- 00829 00830 Sign incircle_p_no_easy(const Site_2& st, PSS_Type ) const 00831 { 00832 CGAL_precondition( st.is_point() ); 00833 Point_2 t = st.point(); 00834 00835 Sqrt_1 Zero(RT(0), RT(0), ux.a().c()); 00836 00837 Point_2 pref = p_ref().point(); 00838 00839 Sqrt_1 xref = pref.x() + Zero; 00840 Sqrt_1 yref = pref.y() + Zero; 00841 00842 Sqrt_3 vx = ux - xref * uz; 00843 Sqrt_3 vy = uy - yref * uz; 00844 00845 Sqrt_3 Rs = CGAL::square(vx) + CGAL::square(vy); 00846 00847 Sqrt_1 tx = t.x() + Zero; 00848 Sqrt_1 ty = t.y() + Zero; 00849 00850 Sqrt_3 Rs1 = 00851 CGAL::square(ux - tx * uz) + CGAL::square(uy - ty * uz); 00852 00853 Sign s_Q = CGAL::sign(Rs1 - Rs); 00854 00855 return s_Q; 00856 } 00857 00858 //-------------------------------------------------------------------------- 00859 00860 Sign incircle_p_no_easy(const Site_2& st, SSS_Type ) const 00861 { 00862 CGAL_precondition( st.is_point() ); 00863 00864 Point_2 t = st.point(); 00865 00866 Sqrt_1 Zero(RT(0), RT(0), ux.a().c()); 00867 00868 RT a1, b1, c1; 00869 compute_supporting_line(p_.supporting_site(), a1, b1, c1); 00870 00871 RT ns = CGAL::square(a1) + CGAL::square(b1); 00872 00873 Sqrt_1 a = a1 + Zero; 00874 Sqrt_1 b = b1 + Zero; 00875 Sqrt_1 c = c1 + Zero; 00876 00877 Sqrt_1 Ns = ns + Zero; 00878 Sqrt_3 Ls = CGAL::square(a * ux + b * uy + c * uz); 00879 00880 Sqrt_1 tx = t.x() + Zero; 00881 Sqrt_1 ty = t.y() + Zero; 00882 00883 Sqrt_3 R1s = CGAL::square(ux - tx * uz) 00884 + CGAL::square(uy - ty * uz); 00885 00886 return CGAL::sign(R1s * Ns - Ls); 00887 } 00888 00889 00890 00891 //-------------------------------------------------------------------------- 00892 //-------------------------------------------------------------------------- 00893 00894 Sign incircle_p(const Site_2& t) const 00895 { 00896 if ( is_degenerate_Voronoi_circle() ) { 00897 return POSITIVE; 00898 } 00899 00900 Sign s(ZERO); 00901 switch ( v_type ) { 00902 case PPP: 00903 s = incircle_p(t, PPP_Type()); 00904 break; 00905 case PPS: 00906 s = incircle_p(t, PPS_Type()); 00907 break; 00908 case PSS: 00909 s = incircle_p(t, PSS_Type()); 00910 break; 00911 case SSS: 00912 s = incircle_p(t, SSS_Type()); 00913 break; 00914 } 00915 00916 return s; 00917 } 00918 00919 Sign incircle_p_no_easy(const Site_2& t) const 00920 { 00921 Sign s(ZERO); 00922 switch ( v_type ) { 00923 case PPP: 00924 s = incircle_p(t, PPP_Type()); 00925 break; 00926 case PPS: 00927 s = incircle_p_no_easy(t, PPS_Type()); 00928 break; 00929 case PSS: 00930 s = incircle_p_no_easy(t, PSS_Type()); 00931 break; 00932 case SSS: 00933 s = incircle_p_no_easy(t, SSS_Type()); 00934 break; 00935 } 00936 00937 return s; 00938 } 00939 00940 //-------------------------------------------------------------------------- 00941 00942 //-------------------------------------------------------------------------- 00943 // the incircle test when the fourth site is a segment 00944 //-------------------------------------------------------------------------- 00945 00946 //-------------------------------------------------------------------------- 00947 00948 Oriented_side 00949 oriented_side(const Line_2& l, const Point_2& p, PPP_Type) const 00950 { 00951 Sign s_uz = CGAL::sign(uz_ppp); 00952 00953 RT px = uz_ppp * p.x() - ux_ppp; 00954 RT py = uz_ppp * p.y() - uy_ppp; 00955 00956 Sign s1 = CGAL::sign(l.b() * px - l.a() * py); 00957 00958 return s_uz * s1; 00959 } 00960 00961 Oriented_side 00962 oriented_side(const Line_2& l, const Point_2& p, PPS_Type) const 00963 { 00964 Sqrt_1 dx = ux_pps - uz_pps * p.x(); 00965 Sqrt_1 dy = uy_pps - uz_pps * p.y(); 00966 00967 return CGAL::sign(uz_pps) * CGAL::sign(dy * l.a() - dx * l.b()); 00968 } 00969 00970 // the cases PSS and SSS are identical 00971 template<class Type> 00972 Oriented_side 00973 oriented_side(const Line_2& l, const Point_2& p, Type) const 00974 { 00975 Sqrt_1 Zero(RT(0), RT(0), ux.a().c()); 00976 Sqrt_1 px = p.x() + Zero; 00977 Sqrt_1 py = p.y() + Zero; 00978 00979 Sqrt_3 dx = ux - px * uz; 00980 Sqrt_3 dy = uy - py * uz; 00981 00982 Sqrt_1 a = l.a() + Zero; 00983 Sqrt_1 b = l.b() + Zero; 00984 00985 return CGAL::sign(uz) * CGAL::sign(a * dy - b * dx); 00986 } 00987 00988 //-------------------------------------------------------------------------- 00989 //-------------------------------------------------------------------------- 00990 00991 00992 Sign incircle(const Line_2& l, PPP_Type) const 00993 { 00994 Point_2 pref = p_ref().point(); 00995 00996 RT a1 = CGAL::square(l.a()) + CGAL::square(l.b()); 00997 RT a2 = CGAL::square(ux_ppp - pref.x() * uz_ppp) + 00998 CGAL::square(uy_ppp - pref.y() * uz_ppp); 00999 01000 RT a3 = 01001 CGAL::square(l.a() * ux_ppp + l.b() * uy_ppp + l.c() * uz_ppp); 01002 01003 Comparison_result cr = CGAL::compare(a3, a1 * a2); 01004 01005 if ( cr == LARGER ) { return POSITIVE; } 01006 if ( cr == SMALLER ) { return NEGATIVE; } 01007 return ZERO; 01008 } 01009 01010 //-------------------------------------------------------------------------- 01011 01012 Sign incircle(const Line_2& l, PPS_Type) const 01013 { 01014 Point_2 pref = p_ref().point(); 01015 01016 Sqrt_1 vx = ux_pps - pref.x() * uz_pps; 01017 Sqrt_1 vy = uy_pps - pref.y() * uz_pps; 01018 01019 Sqrt_1 Rs = CGAL::square(vx) + CGAL::square(vy); 01020 01021 RT Ns = CGAL::square(l.a()) + CGAL::square(l.b()); 01022 01023 Sqrt_1 Ls = 01024 CGAL::square(l.a() * ux_pps + l.b() * uy_pps + l.c() * uz_pps); 01025 01026 return CGAL::sign(Ls - Rs * Ns); 01027 } 01028 01029 01030 //-------------------------------------------------------------------------- 01031 01032 Sign incircle(const Line_2& l, PSS_Type) const 01033 { 01034 Sqrt_1 Zero(RT(0), RT(0), ux.a().c()); 01035 01036 Point_2 pref = p_ref().point(); 01037 01038 Sqrt_3 vx = ux - (pref.x() + Zero) * uz; 01039 Sqrt_3 vy = uy - (pref.y() + Zero) * uz; 01040 01041 Sqrt_3 Rs = CGAL::square(vx) + CGAL::square(vy); 01042 01043 01044 Sqrt_1 a = l.a() + Zero; 01045 Sqrt_1 b = l.b() + Zero; 01046 Sqrt_1 c = l.c() + Zero; 01047 01048 Sqrt_1 Ns = CGAL::square(a) + CGAL::square(b); 01049 01050 Sqrt_3 Ls = CGAL::square(a * ux + b * uy + c * uz); 01051 01052 return CGAL::sign(Ls - Rs * Ns); 01053 } 01054 01055 01056 //-------------------------------------------------------------------------- 01057 01058 Sign incircle(const Line_2& l, SSS_Type) const 01059 { 01060 Sqrt_1 Zero(RT(0), RT(0), ux.a().c()); 01061 01062 RT a1, b1, c1; 01063 compute_supporting_line(p_.supporting_site(), a1, b1, c1); 01064 01065 Sqrt_1 a = a1 + Zero; 01066 Sqrt_1 b = b1 + Zero; 01067 Sqrt_1 c = c1 + Zero; 01068 01069 Sqrt_1 Ns = CGAL::square(a) + CGAL::square(b); 01070 01071 Sqrt_1 la = l.a() + Zero; 01072 Sqrt_1 lb = l.b() + Zero; 01073 Sqrt_1 lc = l.c() + Zero; 01074 01075 Sqrt_1 Ns1 = CGAL::square(la) + CGAL::square(lb); 01076 01077 Sqrt_3 Ls = CGAL::square(a * ux + b * uy + c * uz); 01078 01079 Sqrt_3 Ls1 = 01080 CGAL::square(la * ux + lb * uy + lc * uz); 01081 01082 return CGAL::sign(Ls1 * Ns - Ls * Ns1); 01083 } 01084 01085 //-------------------------------------------------------------------------- 01086 //-------------------------------------------------------------------------- 01087 01088 template<class Type> 01089 Sign incircle_s(const Site_2& t, Type type) const 01090 { 01091 CGAL_precondition( t.is_segment() ); 01092 01093 if ( v_type == PPP || v_type == PPS ) { 01094 if ( p_.is_point() && q_.is_point() && 01095 is_endpoint_of(p_, t) && is_endpoint_of(q_, t) ) { 01096 return NEGATIVE; 01097 } 01098 01099 if ( p_.is_point() && r_.is_point() && 01100 is_endpoint_of(p_, t) && is_endpoint_of(r_, t) ){ 01101 return NEGATIVE; 01102 } 01103 01104 if ( q_.is_point() && r_.is_point() && 01105 is_endpoint_of(q_, t) && is_endpoint_of(r_, t) ){ 01106 return NEGATIVE; 01107 } 01108 } 01109 01110 if ( v_type == PSS ) { 01111 if ( p_.is_segment() && 01112 same_segments(p_.supporting_site(), 01113 t.supporting_site()) ) { 01114 return POSITIVE; 01115 } 01116 if ( q_.is_segment() && 01117 same_segments(q_.supporting_site(), 01118 t.supporting_site()) ) { 01119 return POSITIVE; 01120 } 01121 if ( r_.is_segment() && 01122 same_segments(r_.supporting_site(), 01123 t.supporting_site()) ) { 01124 return POSITIVE; 01125 } 01126 } 01127 01128 return incircle_s_no_easy(t, type); 01129 } 01130 01131 template<class Type> 01132 Sign incircle_s_no_easy(const Site_2& t, Type type) const 01133 { 01134 Sign d1, d2; 01135 if ( ( p_.is_point() && same_points(p_, t.source_site()) ) || 01136 ( q_.is_point() && same_points(q_, t.source_site()) ) || 01137 ( r_.is_point() && same_points(r_, t.source_site()) ) ) { 01138 d1 = ZERO; 01139 } else { 01140 d1 = incircle_p(t.source_site()); 01141 } 01142 if ( d1 == NEGATIVE ) { return NEGATIVE; } 01143 01144 if ( ( p_.is_point() && same_points(p_, t.target_site()) ) || 01145 ( q_.is_point() && same_points(q_, t.target_site()) ) || 01146 ( r_.is_point() && same_points(r_, t.target_site()) ) ) { 01147 d2 = ZERO; 01148 } else { 01149 d2 = incircle_p(t.target_site()); 01150 } 01151 if ( d2 == NEGATIVE ) { return NEGATIVE; } 01152 01153 Line_2 l = compute_supporting_line(t.supporting_site()); 01154 Sign sl = incircle(l, type); 01155 01156 if ( sl == POSITIVE ) { return sl; } 01157 01158 if ( sl == ZERO && (d1 == ZERO || d2 == ZERO) ) { return ZERO; } 01159 01160 Oriented_side os1 = oriented_side(l, t.source(), type); 01161 Oriented_side os2 = oriented_side(l, t.target(), type); 01162 01163 if ( sl == ZERO ) { 01164 if (os1 == ON_ORIENTED_BOUNDARY || os2 == ON_ORIENTED_BOUNDARY) { 01165 return ZERO; 01166 } 01167 return ( os1 == os2 ) ? POSITIVE : ZERO; 01168 } 01169 01170 return (os1 == os2) ? POSITIVE : NEGATIVE; 01171 } 01172 01173 //-------------------------------------------------------------------------- 01174 01175 01176 01177 Sign incircle_s(const Site_2& t) const 01178 { 01179 CGAL_precondition( t.is_segment() ); 01180 01181 if ( is_degenerate_Voronoi_circle() ) { 01182 // case 1: the new segment is not adjacent to the center of the 01183 // degenerate Voronoi circle 01184 if ( !same_points( p_ref(), t.source_site() ) && 01185 !same_points( p_ref(), t.target_site() ) ) { 01186 return POSITIVE; 01187 } 01188 01189 CGAL_assertion( v_type == PSS ); 01190 01191 if ( p_.is_segment() && 01192 same_segments(p_.supporting_site(), 01193 t.supporting_site()) ) { 01194 return ZERO; 01195 } 01196 01197 if ( q_.is_segment() && 01198 same_segments(q_.supporting_site(), 01199 t.supporting_site()) ) { 01200 return ZERO; 01201 } 01202 01203 if ( r_.is_segment() && 01204 same_segments(r_.supporting_site(), 01205 t.supporting_site()) ) { 01206 return ZERO; 01207 } 01208 01209 Site_2 pr; 01210 Site_2 sp, sq; 01211 if ( p_.is_point() ) { 01212 CGAL_assertion( q_.is_segment() && r_.is_segment() ); 01213 pr = p_; 01214 sp = q_; 01215 sq = r_; 01216 } else if ( q_.is_point() ) { 01217 CGAL_assertion( r_.is_segment() && p_.is_segment() ); 01218 pr = q_; 01219 sp = r_; 01220 sq = p_; 01221 } else { 01222 CGAL_assertion( p_.is_segment() && q_.is_segment() ); 01223 pr = r_; 01224 sp = p_; 01225 sq = q_; 01226 } 01227 01228 Point_2 pq = sq.source(), pp = sp.source(), pt = t.source(); 01229 01230 if ( same_points(sp.source_site(), pr) ) { pp = sp.target(); } 01231 if ( same_points(sq.source_site(), pr) ) { pq = sq.target(); } 01232 if ( same_points( t.source_site(), pr) ) { pt = t.target(); } 01233 01234 Point_2 pr_ = pr.point(); 01235 01236 if ( CGAL::orientation(pr_, pp, pt) == LEFT_TURN && 01237 CGAL::orientation(pr_, pq, pt) == RIGHT_TURN ) { 01238 return NEGATIVE; 01239 } 01240 return ZERO; 01241 } // if ( is_degenerate_Voronoi_circle() ) 01242 01243 Sign s(ZERO); 01244 switch ( v_type ) { 01245 case PPP: 01246 s = incircle_s(t, PPP_Type()); 01247 break; 01248 case PPS: 01249 s = incircle_s(t, PPS_Type()); 01250 break; 01251 case PSS: 01252 s = incircle_s(t, PSS_Type()); 01253 break; 01254 case SSS: 01255 s = incircle_s(t, SSS_Type()); 01256 break; 01257 } 01258 01259 return s; 01260 } 01261 01262 Sign incircle_s_no_easy(const Site_2& t) const 01263 { 01264 Sign s(ZERO); 01265 switch ( v_type ) { 01266 case PPP: 01267 s = incircle_s_no_easy(t, PPP_Type()); 01268 break; 01269 case PPS: 01270 s = incircle_s_no_easy(t, PPS_Type()); 01271 break; 01272 case PSS: 01273 s = incircle_s_no_easy(t, PSS_Type()); 01274 break; 01275 case SSS: 01276 s = incircle_s_no_easy(t, SSS_Type()); 01277 break; 01278 } 01279 01280 return s; 01281 } 01282 01283 //-------------------------------------------------------------------------- 01284 // subpredicates for the incircle test 01285 //-------------------------------------------------------------------------- 01286 01287 01288 public: 01289 bool is_degenerate_Voronoi_circle() const 01290 { 01291 if ( v_type != PSS ) { return false; } 01292 01293 if ( p_.is_point() ) { 01294 return ( is_endpoint_of(p_, q_) && is_endpoint_of(p_, r_) ); 01295 } else if ( q_.is_point() ) { 01296 return ( is_endpoint_of(q_, p_) && is_endpoint_of(q_, r_) ); 01297 } else { 01298 CGAL_assertion( r_.is_point() ); 01299 return ( is_endpoint_of(r_, p_) && is_endpoint_of(r_, q_) ); 01300 } 01301 } 01302 01303 01304 //-------------------------------------------------------------------------- 01305 01306 private: 01307 01308 //-------------------------------------------------------------------------- 01309 // the reference point (valid if v_type != SSS) 01310 //-------------------------------------------------------------------------- 01311 01312 Site_2 p_ref() const 01313 { 01314 CGAL_precondition ( v_type != SSS ); 01315 01316 if ( v_type == PPS ) { 01317 if ( pps_idx == 0 ) { return p_; } 01318 if ( pps_idx == 1 ) { return q_; } 01319 return r_; 01320 } 01321 01322 if ( p_.is_point() ) { 01323 return p_; 01324 } else if ( q_.is_point() ) { 01325 return q_; 01326 } else { 01327 CGAL_assertion( r_.is_point() ); 01328 return r_; 01329 } 01330 } 01331 01332 01333 public: 01334 //-------------------------------------------------------------------------- 01335 //-------------------------------------------------------------------------- 01336 // access methods 01337 //-------------------------------------------------------------------------- 01338 //-------------------------------------------------------------------------- 01339 01340 inline FT x(Integral_domain_without_division_tag) const { 01341 return CGAL::to_double(hx()) / CGAL::to_double(hw()); 01342 } 01343 inline FT y(Integral_domain_without_division_tag) const { 01344 return CGAL::to_double(hy()) / CGAL::to_double(hw()); 01345 } 01346 01347 inline FT x(Field_tag) const { return hx() / hw(); } 01348 inline FT y(Field_tag) const { return hy() / hw(); } 01349 01350 inline FT x() const { 01351 typedef Algebraic_structure_traits<FT> AST; 01352 return x(typename AST::Algebraic_category()); 01353 } 01354 01355 inline FT y() const { 01356 typedef Algebraic_structure_traits<FT> AST; 01357 return y(typename AST::Algebraic_category()); 01358 } 01359 01360 FT hx() const { 01361 if ( v_type == PPP ) { return ux_ppp; } 01362 if ( v_type == PPS ) { return to_ft(ux_pps); } 01363 return to_ft(ux); 01364 } 01365 01366 FT hy() const { 01367 if ( v_type == PPP ) { return uy_ppp; } 01368 if ( v_type == PPS ) { return to_ft(uy_pps); } 01369 return to_ft(uy); 01370 } 01371 01372 FT hw() const { 01373 if ( v_type == PPP ) { return uz_ppp; } 01374 if ( v_type == PPS ) { return to_ft(uz_pps); } 01375 return to_ft(uz); 01376 } 01377 01378 FT squared_radius() const { 01379 switch (v_type) { 01380 case PPP: case PPS: case PSS: 01381 { 01382 Point_2 pref = p_ref().point(); 01383 FT dx2 = CGAL::square(x() - pref.x()); 01384 FT dy2 = CGAL::square(y() - pref.y()); 01385 return dx2 + dy2; 01386 } 01387 break; 01388 case SSS: 01389 { 01390 Line_2 l = compute_supporting_line(p_.supporting_site()); 01391 Homogeneous_point_2 q = compute_projection(l, point()); 01392 01393 FT dx2 = CGAL::square(x() - q.x()); 01394 FT dy2 = CGAL::square(y() - q.y()); 01395 return dx2 + dy2; 01396 } 01397 break; 01398 default: 01399 return FT(0); 01400 } 01401 } 01402 01403 Point_2 point() const { 01404 if ( is_degenerate_Voronoi_circle() ) { 01405 return degenerate_point(); 01406 } 01407 01408 return Point_2(x(), y()); 01409 } 01410 01411 01412 Point_2 degenerate_point() const 01413 { 01414 CGAL_precondition( is_degenerate_Voronoi_circle() ); 01415 return p_ref().point(); 01416 } 01417 01418 01419 typename K::Circle_2 circle() const 01420 { 01421 typedef typename K::Circle_2 K_Circle_2; 01422 return K_Circle_2(point(), squared_radius()); 01423 } 01424 01425 vertex_t type() const { return v_type; } 01426 01427 public: 01428 Voronoi_vertex_ring_C2(const Site_2& p, 01429 const Site_2& q, 01430 const Site_2& r) 01431 : p_(p), q_(q), r_(r) 01432 { 01433 compute_vertex(p, q, r); 01434 } 01435 01436 //-------------------------------------------------------------------------- 01437 01438 Sign incircle(const Site_2& t) const 01439 { 01440 Sign s; 01441 01442 if ( t.is_point() ) { 01443 s = incircle_p(t); 01444 } else { 01445 s = incircle_s(t); 01446 } 01447 01448 return s; 01449 } 01450 01451 Sign incircle_no_easy(const Site_2& t) const 01452 { 01453 Sign s; 01454 01455 if ( t.is_point() ) { 01456 s = incircle_p_no_easy(t); 01457 } else { 01458 s = incircle_s_no_easy(t); 01459 } 01460 01461 return s; 01462 } 01463 01464 //-------------------------------------------------------------------------- 01465 01466 01467 Orientation orientation(const Line_2& l) const 01468 { 01469 Orientation o(COLLINEAR); 01470 switch ( v_type ) { 01471 case PPP: 01472 o = orientation(l, PPP_Type()); 01473 break; 01474 case PPS: 01475 o = orientation(l, PPS_Type()); 01476 break; 01477 case PSS: 01478 o = orientation(l, PSS_Type()); 01479 break; 01480 case SSS: 01481 o = orientation(l, SSS_Type()); 01482 break; 01483 } 01484 01485 return o; 01486 } 01487 01488 Oriented_side oriented_side(const Line_2& l) const 01489 { 01490 Orientation o = orientation(l); 01491 01492 if ( o == COLLINEAR ) { return ON_ORIENTED_BOUNDARY; } 01493 return ( o == LEFT_TURN ) ? ON_POSITIVE_SIDE : ON_NEGATIVE_SIDE; 01494 } 01495 01496 //-------------------------------------------------------------------------- 01497 01498 private: 01499 const Site_2& p_, q_, r_; 01500 01501 vertex_t v_type; 01502 01503 // index that indicates the refence point for the case PPS 01504 short pps_idx; 01505 01506 // the case ppp 01507 RT ux_ppp, uy_ppp, uz_ppp; 01508 01509 // the case pps 01510 Sqrt_1 ux_pps, uy_pps, uz_pps; 01511 01512 // the case pss and sss 01513 Sqrt_3 ux, uy, uz; 01514 }; 01515 01516 01517 CGAL_SEGMENT_DELAUNAY_GRAPH_2_END_NAMESPACE 01518 01519 CGAL_END_NAMESPACE 01520 01521 01522 #endif // CGAL_SEGMENT_DELAUNAY_GRAPH_2_VORONOI_VERTEX_RING_C2_H
1.7.6.1