Add a static filter for do_intersect(Segment_2, Segment_2)

This commit is contained in:
Andreas Fabri 2017-10-26 15:35:26 +01:00
parent c9a4fcc663
commit 4cd96171e2
4 changed files with 148 additions and 14 deletions

View File

@ -31,7 +31,8 @@
#include <CGAL/Filtered_kernel/Cartesian_coordinate_iterator_2.h>
#include <CGAL/Filtered_kernel/Cartesian_coordinate_iterator_3.h>
#include <CGAL/Lazy.h>
#include <CGAL/internal/Static_filters/tools.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <boost/none.hpp>
#include <boost/mpl/if.hpp>
#include <boost/mpl/or.hpp>
@ -281,6 +282,10 @@ public:
#include <CGAL/Kernel/interface_macros.h>
};
template < typename EK_, typename AK_, typename E2A_, typename Kernel_ >
class Lazy_kernel_base
: public Lazy_kernel_generic_base<EK_, AK_, E2A_, Kernel_>
@ -307,6 +312,35 @@ public:
// typedef void Compute_z_3; // to detect where .z() is called
// typedef void Construct_point_3; // to detect where the ctor is called
class Do_intersect_2 {
public:
bool operator()(const Segment_2& s, const Segment_2& t) const
{
internal::Static_filters_predicates::Get_approx<Segment_2> get_approx;
typedef Exact_predicates_inexact_constructions_kernel Fit;
typedef typename Fit::Point_2 Fpoint;
typedef typename Fit::Segment_2 Fsegment;
double ssx, ssy, stx, sty, tsx, tsy, ttx, tty;
if(fit_in_double(get_approx(s).source().x(), ssx) &&
fit_in_double(get_approx(s).source().y(), ssy) &&
fit_in_double(get_approx(s).target().x(), stx) &&
fit_in_double(get_approx(s).target().y(), sty) &&
fit_in_double(get_approx(t).source().x(), tsx) &&
fit_in_double(get_approx(t).source().y(), tsy) &&
fit_in_double(get_approx(t).target().x(), ttx) &&
fit_in_double(get_approx(t).target().y(), tty)){
return do_intersect(Fsegment(Fpoint(ssx,ssy), Fpoint(stx,sty)),
Fsegment(Fpoint(tsx,tsy), Fpoint(ttx,tty)));
}
typename Lazy_kernel_generic_base<EK_, AK_, E2A_, Kernel_>::Do_intersect_2 DI;
return DI(s,t);
}
};
Do_intersect_2 do_intersect_2_object() const
{ return Do_intersect_2(); }
Assign_2
assign_2_object() const
{ return Assign_2(); }

View File

@ -0,0 +1,90 @@
// Copyright (c) 2008 ETH Zurich (Switzerland)
// Copyright (c) 2008-2009 INRIA Sophia-Antipolis (France)
// Copyright (c) 2017 GeometryFactory Sarl (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 Lesser 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$
//
//
// Author(s) : Andreas Fabri, Laurent Rineau
#ifndef CGAL_INTERNAL_STATIC_FILTERS_DO_INTERSECT_2_H
#define CGAL_INTERNAL_STATIC_FILTERS_DO_INTERSECT_2_H
#include <iostream>
namespace CGAL {
namespace internal {
namespace Static_filters_predicates {
template < typename K_base, typename SFK >
class Do_intersect_2
: public K_base::Do_intersect_2
{
typedef typename K_base::Point_2 Point_2;
typedef typename K_base::Segment_2 Segment_2;
typedef typename K_base::Do_intersect_2 Base;
typedef K_base TA1;
typedef SFK TA2;
public:
typedef typename Base::result_type result_type;
#ifndef CGAL_CFG_MATCHING_BUG_6
using Base::operator();
#else // CGAL_CFG_MATCHING_BUG_6
template <typename T1, typename T2>
result_type
operator()(const T1& t1, const T2& t2) const
{
return Base()(t1,t2);
}
#endif // CGAL_CFG_MATCHING_BUG_6
// The internal::do_intersect(..) function
// only performs orientation tests on the vertices
// of the segment
// By calling the do_intersect function with
// the statically filtered kernel we avoid
// that doubles are put into Interval_nt
// to get taken out again with fit_in_double
result_type
operator()(const Segment_2 &s, const Segment_2& t) const
{
return internal::do_intersect(s,t, SFK());
}
result_type
operator()(const Point_2 &p, const Segment_2& t) const
{
return internal::do_intersect(p,t, SFK());
}
result_type
operator()(const Segment_2& t, const Point_2 &p) const
{
return internal::do_intersect(p,t, SFK());
}
};
} // Static_filters_predicates
} // internal
} // CGAL
#endif

View File

@ -83,7 +83,7 @@ public:
// of the triangle and the segment
// By calling the do_intersect function with
// the statically filtered kernel we avoid
// that doubles are put into Inteval_nt
// that doubles are put into Interval_nt
// to get taken out again with fit_in_double
result_type
operator()(const Segment_3 &s, const Triangle_3& t) const

View File

@ -45,7 +45,8 @@
# define CGAL_NO_COMPARE_X_2_STATIC_FILTERS 1
# define CGAL_NO_IS_DEGENERATE_3_STATIC_FILTERS 1
# define CGAL_NO_ANGLE_3_STATIC_FILTERS 1
# define CGAL_NO_DO_INTERSECT_3_STATIC_FILTERS 1
# define CGAL_NO_DO_INTERSECT_STATIC_FILTERS 1
#endif // CGAL_DISABLE_STATIC_FILTERS_ADDED_2011
@ -67,9 +68,10 @@
# include <CGAL/internal/Static_filters/Angle_3.h>
#endif // NOT CGAL_NO_ANGLE_3_STATIC_FILTERS
#ifndef CGAL_NO_DO_INTERSECT_3_STATIC_FILTERS
#ifndef CGAL_NO_DO_INTERSECT_STATIC_FILTERS
# include <CGAL/internal/Static_filters/Do_intersect_3.h>
#endif // NOT NOT CGAL_NO_DO_INTERSECT_3_STATIC_FILTERS
# include <CGAL/internal/Static_filters/Do_intersect_2.h>
#endif // NOT NOT CGAL_NO_DO_INTERSECT_STATIC_FILTERS
#include <CGAL/internal/Static_filters/Compare_y_at_x_2.h>
#include <CGAL/internal/Static_filters/Side_of_oriented_circle_2.h>
@ -139,9 +141,6 @@ public:
typedef Static_filters_predicates::Angle_3<K_base> Angle_3;
#endif // NOT CGAL_NO_ANGLE_3_STATIC_FILTERS
#ifndef CGAL_NO_DO_INTERSECT_3_STATIC_FILTERS
typedef Static_filters_predicates::Do_intersect_3<K_base,Self> Do_intersect_3;
#endif // NOT CGAL_NO_DO_INTERSECT_3_STATIC_FILTERS
typedef Static_filters_predicates::Side_of_oriented_circle_2<K_base> Side_of_oriented_circle_2;
typedef Static_filters_predicates::Side_of_oriented_sphere_3<K_base> Side_of_oriented_sphere_3;
typedef Static_filters_predicates::Compare_squared_radius_3<K_base> Compare_squared_radius_3;
@ -179,6 +178,7 @@ public:
Compare_y_2
compare_y_2_object() const
{ return Compare_y_2(); }
#endif // NOT CGAL_NO_COMPARE_Y_2_STATIC_FILTERS
#ifndef CGAL_NO_IS_DEGENERATE_3_STATIC_FILTERS
@ -193,12 +193,6 @@ Compare_y_2
{ return Angle_3(); }
#endif // NOT CGAL_NO_ANGLE_3_STATIC_FILTERS
#ifndef CGAL_NO_DO_INTERSECT_3_STATIC_FILTERS
Do_intersect_3
do_intersect_3_object() const
{ return Do_intersect_3(); }
#endif // NOT CGAL_NO_DO_INTERSECT_3_STATIC_FILTERS
Side_of_oriented_circle_2
side_of_oriented_circle_2_object() const
{ return Side_of_oriented_circle_2(); }
@ -261,6 +255,11 @@ public:
typedef CartesianKernelFunctors::Compare_xy_3<Self> Compare_xy_3;
typedef CartesianKernelFunctors::Compare_xyz_3<Self> Compare_xyz_3;
#ifndef CGAL_NO_DO_INTERSECT_STATIC_FILTERS
typedef Static_filters_predicates::Do_intersect_2<K_base,Self> Do_intersect_2;
typedef Static_filters_predicates::Do_intersect_3<K_base,Self> Do_intersect_3;
#endif // NOT CGAL_NO_DO_INTERSECT_STATIC_FILTERS
Compare_xy_2
compare_xy_2_object() const
{ return Compare_xy_2(); }
@ -329,6 +328,17 @@ public:
compare_y_at_x_2_object() const
{ return Compare_y_at_x_2(); }
#ifndef CGAL_NO_DO_INTERSECT_STATIC_FILTERS
Do_intersect_3
do_intersect_3_object() const
{ return Do_intersect_3(); }
Do_intersect_2
do_intersect_2_object() const
{ return Do_intersect_2(); }
#endif // NOT CGAL_NO_DO_INTERSECT_STATIC_FILTERS
// The two following are for degenerate cases, so I'll update them later.
//
// typedef Static_filters_predicates::Coplanar_orientation_3<Point_3, Orientation_2>