1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575
|
// Copyright (c) 2003,2004,2005,2006 INRIA Sophia-Antipolis (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
// You can redistribute it and/or modify it under the terms of the GNU
// General Public License as published by the Free Software Foundation,
// either version 3 of the License, or (at your option) any later version.
//
// Licensees holding a valid commercial license may use this file in
// accordance with the commercial license agreement provided with the software.
//
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
// SPDX-License-Identifier: GPL-3.0+
//
//
// Author(s) : Menelaos Karavelas <mkaravel@iacm.forth.gr>
#ifndef CGAL_SEGMENT_DELAUNAY_GRAPH_2_VERTEX_CONFLICT_C2_H
#define CGAL_SEGMENT_DELAUNAY_GRAPH_2_VERTEX_CONFLICT_C2_H
#include <CGAL/license/Segment_Delaunay_graph_2.h>
#include <CGAL/Segment_Delaunay_graph_2/Voronoi_vertex_C2.h>
#include <CGAL/Segment_Delaunay_graph_2/Are_same_points_C2.h>
#include <CGAL/Segment_Delaunay_graph_2/Are_same_segments_C2.h>
#ifdef CGAL_SDG_CHECK_INCIRCLE_CONSISTENCY
#ifndef CGAL_SDG_USE_OLD_INCIRCLE
#include <CGAL/Segment_Delaunay_graph_2/Voronoi_vertex_sqrt_field_C2.h>
#endif // CGAL_SDG_USE_OLD_INCIRCLE
#endif // CGAL_SDG_CHECK_INCIRCLE_CONSISTENCY
namespace CGAL {
namespace SegmentDelaunayGraph_2 {
//---------------------------------------------------------------------
template<class K, class Method_tag>
class Vertex_conflict_C2
{
private:
typedef typename K::Point_2 Point_2;
typedef typename K::Segment_2 Segment_2;
typedef typename K::Site_2 Site_2;
typedef typename K::FT FT;
typedef typename K::RT RT;
typedef typename K::Orientation Orientation;
typedef typename K::Sign Sign;
typedef Voronoi_vertex_C2<K,Method_tag> Voronoi_vertex_2;
typedef Are_same_points_C2<K> Are_same_points_2;
typedef Are_same_segments_C2<K> Are_same_segments_2;
typedef typename K::Intersections_tag ITag;
private:
Are_same_points_2 same_points;
Are_same_segments_2 same_segments;
bool is_on_common_support(const Site_2& s1, const Site_2& s2,
const Point_2& p) const
{
CGAL_precondition( !s1.is_input() && !s2.is_input() );
if ( same_segments(s1.supporting_site(0),
s2.supporting_site(0)) ||
same_segments(s1.supporting_site(0),
s2.supporting_site(1)) ) {
Site_2 support = s1.supporting_site(0);
Site_2 tp = Site_2::construct_site_2(p);
return ( same_points(support.source_site(), tp) ||
same_points(support.target_site(), tp) );
} else if ( same_segments(s1.supporting_site(1),
s2.supporting_site(1)) ||
same_segments(s1.supporting_site(1),
s2.supporting_site(0)) ) {
Site_2 support = s1.supporting_site(1);
Site_2 tp = Site_2::construct_site_2(p);
return ( same_points(support.source_site(), tp) ||
same_points(support.target_site(), tp) );
}
return false;
}
bool have_common_support(const Site_2& p, const Site_2& q) const
{
CGAL_precondition( !p.is_input() && !q.is_input() );
return
same_segments(p.supporting_site(0), q.supporting_site(0)) ||
same_segments(p.supporting_site(0), q.supporting_site(1)) ||
same_segments(p.supporting_site(1), q.supporting_site(1)) ||
same_segments(p.supporting_site(1), q.supporting_site(0));
}
bool have_common_support(const Site_2& s, const Point_2& p1,
const Point_2& p2) const
{
CGAL_precondition( !s.is_input() );
Site_2 t = Site_2::construct_site_2(p1, p2);
return ( same_segments(s.supporting_site(0), t) ||
same_segments(s.supporting_site(1), t) );
}
private:
Sign incircle_ppp(const Site_2& p, const Site_2& q,
const Site_2& t, const Tag_false&) const
{
Point_2 pp = p.point(), qp = q.point(), tp = t.point();
// MK::ERROR: here I should call a kernel object, not a
// function...; actually here (and everywhere in this class)
// use the orientation predicate for sites; it does some
// geometric filtering...
Orientation o = orientation(pp, qp, tp);
if ( o != COLLINEAR ) {
return (o == LEFT_TURN) ? POSITIVE : NEGATIVE;
}
// MK::ERROR: change the following code to use the compare_x_2
// and compare_y_2 stuff...
RT dtpx = pp.x() - tp.x();
RT dtpy = pp.y() - tp.y();
RT dtqx = qp.x() - tp.x();
RT minus_dtqy = -qp.y() + tp.y();
Sign s = sign_of_determinant(dtpx, dtpy, minus_dtqy, dtqx);
CGAL_assertion( s != ZERO );
return s;
}
Sign incircle_ppp(const Site_2& p, const Site_2& q,
const Site_2& t, const Tag_true&) const
{
Orientation o = COLLINEAR; // the initialization was done in
// order a compiler warning
// do some geometric filtering...
bool p_exact = p.is_input();
bool q_exact = q.is_input();
bool t_exact = t.is_input();
bool filtered = false;
// the following if-statement does the gometric filtering...
// maybe it is not so important since this will only be
// activated if a lot of intersection points appear on the
// convex hull
if ( !p_exact || !q_exact || !t_exact ) {
if ( !p_exact && !q_exact && !t_exact ) {
if ( have_common_support(p, q) &&
have_common_support(q, t) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( !p_exact && !q_exact && t_exact ) {
if ( is_on_common_support(p, q, t.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( !p_exact && q_exact && !t_exact ) {
if ( is_on_common_support(p, t, q.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( p_exact && !q_exact && !t_exact ) {
if ( is_on_common_support(t, q, p.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( !p_exact && q_exact && t_exact ) {
if ( have_common_support(p, q.point(), t.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( p_exact && !q_exact && t_exact ) {
if ( have_common_support(q, p.point(), t.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( p_exact && q_exact && !t_exact ) {
if ( have_common_support(t, p.point(), q.point()) ) {
o = COLLINEAR;
filtered = true;
}
}
}
Point_2 pp = p.point(), qp = q.point(), tp = t.point();
if ( !filtered ) {
// MK::ERROR: here I should call a kernel object, not a
// function...; actually here (and everywhere in this class)
// use the orientation predicate for sites; it does some
// geometric filtering...
o = orientation(pp, qp, tp);
}
if ( o != COLLINEAR ) {
return (o == LEFT_TURN) ? POSITIVE : NEGATIVE;
}
// MK::ERROR: change the following code to use the compare_x_2
// and compare_y_2 stuff...
RT dtpx = pp.x() - tp.x();
RT dtpy = pp.y() - tp.y();
RT dtqx = qp.x() - tp.x();
RT minus_dtqy = -qp.y() + tp.y();
Sign s = sign_of_determinant(dtpx, dtpy, minus_dtqy, dtqx);
CGAL_assertion( s != ZERO );
return s;
}
Sign incircle_p(const Site_2& p, const Site_2& q,
const Site_2& t) const
{
CGAL_precondition( t.is_point() );
if ( p.is_point() && q.is_point() ) {
#if 1
return incircle_ppp(p, q, t, ITag());
#else
Orientation o = COLLINEAR; // the initialization was done in
// order a compiler warning
// do some geometric filtering...
bool p_exact = p.is_input();
bool q_exact = q.is_input();
bool t_exact = t.is_input();
bool filtered = false;
// the following if-statement does the gometric filtering...
// maybe it is not so important since this will only be
// activated if a lot of intersection points appear on the
// convex hull
if ( !p_exact || !q_exact || !t_exact ) {
if ( !p_exact && !q_exact && !t_exact ) {
if ( have_common_support(p, q) &&
have_common_support(q, t) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( !p_exact && !q_exact && t_exact ) {
if ( is_on_common_support(p, q, t.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( !p_exact && q_exact && !t_exact ) {
if ( is_on_common_support(p, t, q.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( p_exact && !q_exact && !t_exact ) {
if ( is_on_common_support(t, q, p.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( !p_exact && q_exact && t_exact ) {
if ( have_common_support(p, q.point(), t.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( p_exact && !q_exact && t_exact ) {
if ( have_common_support(q, p.point(), t.point()) ) {
o = COLLINEAR;
filtered = true;
}
} else if ( p_exact && q_exact && !t_exact ) {
if ( have_common_support(t, p.point(), q.point()) ) {
o = COLLINEAR;
filtered = true;
}
}
}
Point_2 pp = p.point(), qp = q.point(), tp = t.point();
if ( !filtered ) {
// MK::ERROR: here I should call a kernel object, not a
// function...; actually here (and everywhere in this class)
// use the orientation predicate for sites; it does some
// geometric filtering...
o = orientation(pp, qp, tp);
}
if ( o != COLLINEAR ) {
return (o == LEFT_TURN) ? POSITIVE : NEGATIVE;
}
// MK::ERROR: change the following code to use the compare_x_2
// and compare_y_2 stuff...
RT dtpx = pp.x() - tp.x();
RT dtpy = pp.y() - tp.y();
RT dtqx = qp.x() - tp.x();
RT minus_dtqy = -qp.y() + tp.y();
Sign s = sign_of_determinant(dtpx, dtpy, minus_dtqy, dtqx);
CGAL_assertion( s != ZERO );
return s;
#endif
}
CGAL_assertion( p.is_point() || q.is_point() );
Orientation o;
if ( p.is_point() && q.is_segment() ) {
Point_2 pq = same_points(p, q.source_site()) ? q.target() : q.source();
o = orientation(p.point(), pq, t.point());
} else { // p is a segment and q is a point
Point_2 pp = same_points(q, p.source_site()) ? p.target() : p.source();
o = orientation(pp, q.point(), t.point());
}
if ( CGAL::is_certain(o == RIGHT_TURN) )
return CGAL::get_certain( o == RIGHT_TURN ) ? NEGATIVE : POSITIVE;
return CGAL::Uncertain<CGAL::Sign>::indeterminate();
}
//-----------------------------------------------------------------------
Sign incircle_pps(const Site_2& p, const Site_2& q,
const Site_2& t) const
{
CGAL_precondition( p.is_point() && q.is_point() );
bool is_p_tsrc = same_points(p, t.source_site());
bool is_p_ttrg = same_points(p, t.target_site());
bool is_q_tsrc = same_points(q, t.source_site());
bool is_q_ttrg = same_points(q, t.target_site());
bool is_p_on_t = is_p_tsrc || is_p_ttrg;
bool is_q_on_t = is_q_tsrc || is_q_ttrg;
if ( is_p_on_t && is_q_on_t ) {
// if t is the segment joining p and q then t must be a vertex
// on the convex hull
return NEGATIVE;
} else if ( is_p_on_t ) {
// p is an endpoint of t
// in this case the p,q,oo vertex is destroyed only if the
// other endpoint of t is beyond
Point_2 pt = is_p_tsrc ? t.target() : t.source();
Orientation o = orientation(p.point(), q.point(), pt);
return (o == RIGHT_TURN) ? NEGATIVE : POSITIVE;
} else if ( is_q_on_t ) {
Point_2 pt = is_q_tsrc ? t.target() : t.source();
Orientation o = orientation(p.point(), q.point(), pt);
return (o == RIGHT_TURN) ? NEGATIVE : POSITIVE;
} else {
// maybe here I should immediately return POSITIVE;
// since we insert endpoints of segments first, p and q cannot
// be consecutive points on the convex hull if one of the
// endpoints of t is to the right of the line pq.
Point_2 pp = p.point(), qq = q.point();
Orientation o1 = orientation(pp, qq, t.source());
Orientation o2 = orientation(pp, qq, t.target());
if ( o1 == RIGHT_TURN || o2 == RIGHT_TURN ) {
return NEGATIVE;
}
return POSITIVE;
}
}
Sign incircle_sps(const Site_2& p, const Site_2& q,
const Site_2& t) const
{
CGAL_precondition( p.is_segment() && q.is_point() );
bool is_q_tsrc = same_points(q, t.source_site());
bool is_q_ttrg = same_points(q, t.target_site());
bool is_q_on_t = is_q_tsrc || is_q_ttrg;
if ( is_q_on_t ) {
Point_2 pp = same_points(q, p.source_site()) ? p.target() : p.source();
Point_2 pt = is_q_tsrc ? t.target() : t.source();
Orientation o = orientation(pp, q.point(), pt);
return (o == RIGHT_TURN) ? NEGATIVE : POSITIVE;
} else {
return POSITIVE;
}
}
Sign incircle_pss(const Site_2& p, const Site_2& q,
const Site_2& t) const
{
CGAL_precondition( p.is_point() && q.is_segment() );
bool is_p_tsrc = same_points(p, t.source_site());
bool is_p_ttrg = same_points(p, t.target_site());
bool is_p_on_t = is_p_tsrc || is_p_ttrg;
if ( is_p_on_t ) {
Point_2 pq = same_points(p, q.source_site()) ? q.target() : q.source();
Point_2 pt = is_p_tsrc ? t.target() : t.source();
Orientation o = orientation(p.point(), pq, pt);
return (o == RIGHT_TURN) ? NEGATIVE : POSITIVE;
} else {
// if p is not an endpoint of t, then either p and q should
// not be on the convex hull or t does not affect the vertex
// of p and q.
return POSITIVE;
}
}
Sign incircle_s(const Site_2& p, const Site_2& q,
const Site_2& t) const
{
CGAL_precondition( t.is_segment() );
if ( p.is_point() && q.is_point() ) {
return incircle_pps(p, q, t);
} else if ( p.is_point() && q.is_segment() ) {
return incircle_pss(p, q, t);
} else { // p is a segment and q is a point
return incircle_sps(p, q, t);
}
}
public:
typedef Site_2 argument_type;
typedef Sign result_type;
Sign operator()(const Site_2& p, const Site_2& q,
const Site_2& r, const Site_2& t) const
{
#ifdef CGAL_PROFILE
// In case CGAL profile is called then output the sites in case of
// a filter failure
if ( Algebraic_structure_traits<FT>::Is_exact::value ) {
int np = 0;
if ( p.is_point() ) ++np;
if ( q.is_point() ) ++np;
if ( r.is_point() ) ++np;
std::string suffix("-failure-log.cin");
std::string fname;
if ( np == 3 ) {
fname = "ppp";
} else if ( np == 2 ) {
fname = "pps";
} else if ( np == 1 ) {
fname = "pss";
} else {
fname = "sss";
}
fname += suffix;
std::ofstream ofs(fname.c_str(), std::ios_base::app);
ofs.precision(16);
ofs << p << std::endl;
ofs << q << std::endl;
ofs << r << std::endl;
ofs << t << std::endl;
ofs << "=======" << std::endl;
ofs.close();
}
#endif
#ifdef CGAL_SDG_CHECK_INCIRCLE_CONSISTENCY
#ifdef CGAL_SDG_USE_OLD_INCIRCLE
typedef Voronoi_vertex_sqrt_field_new_C2<K> Alt_Voronoi_vertex_2;
#else
typedef Voronoi_vertex_sqrt_field_C2<K> Alt_Voronoi_vertex_2;
#endif
Voronoi_vertex_2 v(p, q, r);
Alt_Voronoi_vertex_2 v_alt(p, q, r);
Sign s = v.incircle(t);
Sign s_alt = v_alt.incircle(t);
if ( s != s_alt ) {
std::cerr << "different results" << std::endl;
std::cerr << p << std::endl;
std::cerr << q << std::endl;
std::cerr << r << std::endl;
std::cerr << t << std::endl;
CGAL_assertion( s == s_alt );
exit(1);
}
return s;
#else
Voronoi_vertex_2 v(p, q, r);
return v.incircle(t);
#endif // CGAL_SDG_CHECK_INCIRCLE_CONSISTENCY
}
Sign operator()(const Site_2& p, const Site_2& q,
const Site_2& t) const
{
#ifdef CGAL_PROFILE
// In case CGAL profile is called then output the sites in case of
// a filter failure
if ( Algebraic_structure_traits<FT>::Is_exact::value ) {
std::ofstream ofs("failure-log.cin", std::ios_base::app);
ofs.precision(16);
ofs << p << std::endl;
ofs << q << std::endl;
ofs << t << std::endl;
ofs << "=======" << std::endl;
ofs.close();
}
#endif
CGAL_assertion( !(p.is_segment() && q.is_segment()) );
if ( p.is_point() && q.is_segment() ) {
// p must be an endpoint of q
CGAL_assertion( same_points(p, q.source_site()) ||
same_points(p, q.target_site()) );
} else if ( p.is_segment() && q.is_point() ) {
// q must be an endpoint of p
CGAL_assertion( same_points(p.source_site(), q) ||
same_points(p.target_site(), q) );
}
if ( t.is_point() ) {
// return incircle_p(p, q, t);
return incircle_p(q, p, t);
}
// MK::ERROR: do geometric filtering when orientation is called.
// return incircle_s(p, q, t);
return incircle_s(q, p, t);
}
};
//---------------------------------------------------------------------
} //namespace SegmentDelaunayGraph_2
} //namespace CGAL
#endif // CGAL_SEGMENT_DELAUNAY_GRAPH_2_VERTEX_CONFLICT_C2_H
|