From a0f9754bbe0abd24d01e9a5cbe8bf3f20e7d121e Mon Sep 17 00:00:00 2001 From: Sylvain Pion Date: Wed, 22 Aug 2007 15:41:11 +0000 Subject: [PATCH] Add radius(Interval_nt). Replace compare_relative_precision() by has_smaller_relative_precision() fixing the case of [-Inf,+Inf]. --- Number_types/include/CGAL/Interval_nt.h | 19 +++++++++++++++---- Number_types/include/CGAL/Lazy_exact_nt.h | 4 ++-- 2 files changed, 17 insertions(+), 6 deletions(-) diff --git a/Number_types/include/CGAL/Interval_nt.h b/Number_types/include/CGAL/Interval_nt.h index fc60cb4ed20..5a92a49cf13 100644 --- a/Number_types/include/CGAL/Interval_nt.h +++ b/Number_types/include/CGAL/Interval_nt.h @@ -433,10 +433,21 @@ width (const Interval_nt & d) // Non-documented template inline -Comparison_result -compare_relative_precision(const Interval_nt & d, double prec) +double +radius (const Interval_nt & d) { - return CGAL::compare(width(d), prec * magnitude(d)); + return width(d)/2; // This could be improved to avoid overflow. +} + +// Non-documented +// This is the relative precision of to_double() (the center of the interval), +// hence we use radius() instead of width(). +template +inline +bool +has_smaller_relative_precision(const Interval_nt & d, double prec) +{ + return magnitude(d) == 0 || radius(d) < prec * magnitude(d); } // Non-documented @@ -446,7 +457,7 @@ relative_precision(const Interval_nt & d) { if (magnitude(d) == 0.0) return 0.0; - return width(d) / magnitude(d); + return radius(d) / magnitude(d); } diff --git a/Number_types/include/CGAL/Lazy_exact_nt.h b/Number_types/include/CGAL/Lazy_exact_nt.h index d59eae0264e..bd244322e45 100644 --- a/Number_types/include/CGAL/Lazy_exact_nt.h +++ b/Number_types/include/CGAL/Lazy_exact_nt.h @@ -1078,8 +1078,8 @@ template < typename ET > class Real_embeddable_traits< Lazy_exact_nt > return r; // If it's precise enough, then OK. - if (compare_relative_precision(app, - Lazy_exact_nt::get_relative_precision_of_to_double()) <= 0) + if (has_smaller_relative_precision(app, + Lazy_exact_nt::get_relative_precision_of_to_double())) return CGAL_NTS to_double(app); CGAL_PROFILER(std::string("failures of : ") + std::string(CGAL_PRETTY_FUNCTION));