Merge pull request #5667 from afabri/Kernel-compare_slopes-GF

Kernel: Add compare_slope with 4 points
This commit is contained in:
Sebastien Loriot 2021-08-12 10:12:26 +02:00 committed by GitHub
commit 92d90a4a11
No known key found for this signature in database
GPG Key ID: 4AEE18F83AFDEB23
7 changed files with 85 additions and 26 deletions

View File

@ -606,6 +606,7 @@ namespace CartesianKernelFunctors {
template <typename K> template <typename K>
class Compare_slope_2 class Compare_slope_2
{ {
typedef typename K::Point_2 Point_2;
typedef typename K::Line_2 Line_2; typedef typename K::Line_2 Line_2;
typedef typename K::Segment_2 Segment_2; typedef typename K::Segment_2 Segment_2;
public: public:
@ -625,6 +626,15 @@ namespace CartesianKernelFunctors {
s2.source().x(), s2.source().y(), s2.source().x(), s2.source().y(),
s2.target().x(), s2.target().y()); s2.target().x(), s2.target().y());
} }
result_type
operator()(const Point_2& s1s, const Point_2& s1t, const Point_2& s2s, const Point_2& s2t) const
{
return compare_slopesC2(s1s.x(), s1s.y(),
s1t.x(), s1t.y(),
s2s.x(), s2s.y(),
s2t.x(), s2t.y());
}
}; };
template <typename K> template <typename K>

View File

@ -810,6 +810,7 @@ namespace HomogeneousKernelFunctors {
template <typename K> template <typename K>
class Compare_slope_2 class Compare_slope_2
{ {
typedef typename K::Point_2 Point_2;
typedef typename K::Line_2 Line_2; typedef typename K::Line_2 Line_2;
typedef typename K::Segment_2 Segment_2; typedef typename K::Segment_2 Segment_2;
public: public:
@ -842,48 +843,55 @@ namespace HomogeneousKernelFunctors {
result_type result_type
operator()(const Segment_2& s1, const Segment_2& s2) const operator()(const Segment_2& s1, const Segment_2& s2) const
{
return (*this)(s1.source(), s1.target(),
s2.source(), s2.target());
}
result_type
operator()(const Point_2& s1s, const Point_2& s1t, const Point_2& s2s, const Point_2& s2t) const
{ {
typedef typename K::FT FT; typedef typename K::FT FT;
typename K::Comparison_result cmp_y1 = compare_y(s1.source(), s1.target()); typename K::Comparison_result cmp_y1 = compare_y(s1s, s1t);
if (cmp_y1 == EQUAL) // horizontal if (cmp_y1 == EQUAL) // horizontal
{ {
typename K::Comparison_result cmp_x2 = compare_x(s2.source(), s2.target()); typename K::Comparison_result cmp_x2 = compare_x(s2s, s2t);
if (cmp_x2 == EQUAL) return SMALLER; if (cmp_x2 == EQUAL) return SMALLER;
FT s_hw = s2.source().hw(); FT s_hw = s2s.hw();
FT t_hw = s2.target().hw(); FT t_hw = s2t.hw();
return - CGAL_NTS sign(s2.source().hy()*t_hw - s2.target().hy()*s_hw) * return - CGAL_NTS sign(s2s.hy()*t_hw - s2t.hy()*s_hw) *
CGAL_NTS sign(s2.source().hx()*t_hw - s2.target().hx()*s_hw); CGAL_NTS sign(s2s.hx()*t_hw - s2t.hx()*s_hw);
} }
typename K::Comparison_result cmp_y2 = compare_y(s2.source(), s2.target()); typename K::Comparison_result cmp_y2 = compare_y(s2s, s2t);
if (cmp_y2 == EQUAL) if (cmp_y2 == EQUAL)
{ {
typename K::Comparison_result cmp_x1 = compare_x(s1.source(), s1.target()); typename K::Comparison_result cmp_x1 = compare_x(s1s, s1t);
if (cmp_x1 == EQUAL) return LARGER; if (cmp_x1 == EQUAL) return LARGER;
FT s_hw = s1.source().hw(); FT s_hw = s1s.hw();
FT t_hw = s1.target().hw(); FT t_hw = s1t.hw();
return CGAL_NTS sign(s1.source().hy()*t_hw - s1.target().hy()*s_hw) * return CGAL_NTS sign(s1s.hy()*t_hw - s1t.hy()*s_hw) *
CGAL_NTS sign(s1.source().hx()*t_hw - s1.target().hx()*s_hw); CGAL_NTS sign(s1s.hx()*t_hw - s1t.hx()*s_hw);
} }
typename K::Comparison_result cmp_x1 = compare_x(s1.source(), s1.target()); typename K::Comparison_result cmp_x1 = compare_x(s1s, s1t);
typename K::Comparison_result cmp_x2 = compare_x(s2.source(), s2.target()); typename K::Comparison_result cmp_x2 = compare_x(s2s, s2t);
if (cmp_x1 == EQUAL) if (cmp_x1 == EQUAL)
return cmp_x2 == EQUAL ? EQUAL : LARGER; return cmp_x2 == EQUAL ? EQUAL : LARGER;
if (cmp_x2 == EQUAL) return SMALLER; if (cmp_x2 == EQUAL) return SMALLER;
FT s1_s_hw = s1.source().hw(); FT s1_s_hw = s1s.hw();
FT s1_t_hw = s1.target().hw(); FT s1_t_hw = s1t.hw();
FT s2_s_hw = s2.source().hw(); FT s2_s_hw = s2s.hw();
FT s2_t_hw = s2.target().hw(); FT s2_t_hw = s2t.hw();
FT s1_xdiff = s1.source().hx()*s1_t_hw - s1.target().hx()*s1_s_hw; FT s1_xdiff = s1s.hx()*s1_t_hw - s1t.hx()*s1_s_hw;
FT s1_ydiff = s1.source().hy()*s1_t_hw - s1.target().hy()*s1_s_hw; FT s1_ydiff = s1s.hy()*s1_t_hw - s1t.hy()*s1_s_hw;
FT s2_xdiff = s2.source().hx()*s2_t_hw - s2.target().hx()*s2_s_hw; FT s2_xdiff = s2s.hx()*s2_t_hw - s2t.hx()*s2_s_hw;
FT s2_ydiff = s2.source().hy()*s2_t_hw - s2.target().hy()*s2_s_hw; FT s2_ydiff = s2s.hy()*s2_t_hw - s2t.hy()*s2_s_hw;
typename K::Sign s1_sign = CGAL_NTS sign(s1_ydiff * s1_xdiff); typename K::Sign s1_sign = CGAL_NTS sign(s1_ydiff * s1_xdiff);
typename K::Sign s2_sign = CGAL_NTS sign(s2_ydiff * s2_xdiff); typename K::Sign s2_sign = CGAL_NTS sign(s2_ydiff * s2_xdiff);

View File

@ -852,7 +852,7 @@ const CGAL::Point_3<Kernel>& t);
/// \defgroup compare_slopes_grp CGAL::compare_slopes() /// \defgroup compare_slopes_grp CGAL::compare_slope()
/// \ingroup kernel_global_function /// \ingroup kernel_global_function
/// @{ /// @{
@ -870,9 +870,20 @@ from the left to the right endpoint of the segments.
*/ */
template <typename Kernel> template <typename Kernel>
Comparison_result compare_slope(const CGAL::Segment_2<Kernel> &s1, Comparison_result compare_slope(const CGAL::Segment_2<Kernel> &s1,
const CGAL::Segment_2<Kernel> &s2); const CGAL::Segment_2<Kernel> &s2);
/*! /*!
compares the slopes of the segments `(s1s,s1t)` and `(s2s,s2t)`,
where the slope is the variation of the `y`-coordinate
from the left to the right endpoint of the segments.
*/
template <typename Kernel>
Comparison_result compare_slope(const CGAL::Point_2<Kernel> &s1s,
const CGAL::Point_2<Kernel> &s1t,
const CGAL::Point_2<Kernel> &s2s,
const CGAL::Point_2<Kernel> &s2t);
/*!
compares the slopes of the segments `(p,q)` and `(r,s)`, compares the slopes of the segments `(p,q)` and `(r,s)`,
where the slope is the variation of the `z`-coordinate from the first where the slope is the variation of the `z`-coordinate from the first
to the second point of the segment divided by the length of the segment. to the second point of the segment divided by the length of the segment.

View File

@ -1041,6 +1041,16 @@ public:
Comparison_result operator()(const Kernel::Segment_2& s1, Comparison_result operator()(const Kernel::Segment_2& s1,
const Kernel::Segment_2& s2); const Kernel::Segment_2& s2);
/*!
compares the slopes of the segments `(s1s,s1t)` and `(s2s,s2t)`,
where the slope is the variation of the `y`-coordinate
from the left to the right endpoint of the segments.
*/
Comparison_result operator()(const Kernel::Point_2& s1s,
const Kernel::Point_2& s1t,
const Kernel::Point_2& s2s,
const Kernel::Point_2& s2t));
/// @} /// @}
}; /* end Kernel::CompareSlope_2 */ }; /* end Kernel::CompareSlope_2 */
@ -9712,7 +9722,6 @@ public:
const Kernel::Point_3&s, const Kernel::Point_3&s,
const Kernel::Point_3&t); const Kernel::Point_3&t);
/// @} /// @}
}; /* end Kernel::SideOfOrientedSphere_3 */ }; /* end Kernel::SideOfOrientedSphere_3 */

View File

@ -347,6 +347,15 @@ compare_slope(const Segment_2<K> &s1, const Segment_2<K> &s2)
return internal::compare_slope(s1, s2, K()); return internal::compare_slope(s1, s2, K());
} }
template < class K >
inline
typename K::Comparison_result
compare_slope(const Point_2<K> &s1s, const Point_2<K> &s1t,
const Point_2<K> &s2s, const Point_2<K> &s2t)
{
return internal::compare_slope(s1s, s1t, s2s, s2t, K());
}
#ifndef CGAL_NO_DEPRECATED_CODE #ifndef CGAL_NO_DEPRECATED_CODE
// kept for backward compatibility // kept for backward compatibility

View File

@ -373,6 +373,17 @@ compare_slope(const typename K::Segment_2 &s1,
return k.compare_slope_2_object()(s1, s2); return k.compare_slope_2_object()(s1, s2);
} }
template < class K >
inline
typename K::Comparison_result
compare_slope(const typename K::Point_2 &s1s,
const typename K::Point_2 &s1t,
const typename K::Point_2 &s2s,
const typename K::Point_2 &s2t,const K& k)
{
return k.compare_slope_2_object()(s1s, s1t, s2s, s2t);
}
template < class K > template < class K >
inline inline
typename K::Comparison_result typename K::Comparison_result

View File

@ -477,6 +477,7 @@ test_new_2(const R& rep)
= rep.compare_slope_2_object(); = rep.compare_slope_2_object();
Comparison_result tmp34ee = compare_slope(l1, l2); Comparison_result tmp34ee = compare_slope(l1, l2);
Comparison_result tmp34ff = compare_slope(s1, s2); Comparison_result tmp34ff = compare_slope(s1, s2);
Comparison_result tmp34gg = compare_slope(p3, p5, p2, p3);
typename R::Less_distance_to_point_2 less_distance_to_point typename R::Less_distance_to_point_2 less_distance_to_point
= rep.less_distance_to_point_2_object(); = rep.less_distance_to_point_2_object();
@ -680,7 +681,7 @@ test_new_2(const R& rep)
use(tmp32a); use(tmp31d); use(tmp31c); use(tmp31b); use(tmp31a); use(tmp30); use(tmp32a); use(tmp31d); use(tmp31c); use(tmp31b); use(tmp31a); use(tmp30);
use(tmp26); use(tmp25); use(tmp24); use(tmp26); use(tmp25); use(tmp24);
use(tmp29); use(tmp28); use(tmp33a); use(tmp33b); use(tmp34ab); use(tmp34ac); use(tmp29); use(tmp28); use(tmp33a); use(tmp33b); use(tmp34ab); use(tmp34ac);
use(tmp34ff); use(tmp34ee); use(tmp34dd); use(tmp34cc); use(tmp34bb); use(tmp34ff); use(tmp34gg); use(tmp34ee); use(tmp34dd); use(tmp34cc); use(tmp34bb);
use(tmp34aa); use(tmp34aa);
use(tmp39a); use(tmp36a); use(tmp48c); use(tmp49c); use(tmp50c); use(tmp39a); use(tmp36a); use(tmp48c); use(tmp49c); use(tmp50c);
use(tmp24a); use(tmp24b); use(tmp24c); use(tmp24d); use(tmp24e); use(tmp24f); use(tmp24a); use(tmp24b); use(tmp24c); use(tmp24d); use(tmp24e); use(tmp24f);