mirror of https://github.com/CGAL/cgal
Add some useful functions for Interval_nt:
- width() - magnitude() - compare_relative_precision(Interval, double) And use the last one in to_double(Lazy_exact_nt).
This commit is contained in:
parent
0389b256f5
commit
20a5fb5c95
|
|
@ -280,6 +280,7 @@ Uncertain<bool>
|
|||
operator!=(const Interval_nt<Protected> &a, int b)
|
||||
{ return ! (a == b); }
|
||||
|
||||
// Non-documented
|
||||
// Returns true if the interval is a unique representable double.
|
||||
template <bool Protected>
|
||||
inline
|
||||
|
|
@ -292,6 +293,7 @@ fit_in_double (const Interval_nt<Protected> & d, double &r)
|
|||
return b;
|
||||
}
|
||||
|
||||
// Non-documented
|
||||
template <bool Protected>
|
||||
inline
|
||||
bool
|
||||
|
|
@ -300,6 +302,34 @@ is_singleton (const Interval_nt<Protected> & d)
|
|||
return d.is_point();
|
||||
}
|
||||
|
||||
// Non-documented
|
||||
template <bool Protected>
|
||||
inline
|
||||
double
|
||||
magnitude (const Interval_nt<Protected> & d)
|
||||
{
|
||||
return (std::max)(std::fabs(d.inf()), std::fabs(d.sup()));
|
||||
}
|
||||
|
||||
// Non-documented
|
||||
template <bool Protected>
|
||||
inline
|
||||
double
|
||||
width (const Interval_nt<Protected> & d)
|
||||
{
|
||||
return d.sup() - d.inf();
|
||||
}
|
||||
|
||||
// Non-documented
|
||||
template <bool Protected>
|
||||
inline
|
||||
Comparison_result
|
||||
compare_relative_precision(const Interval_nt<Protected> & d, double prec)
|
||||
{
|
||||
return CGAL::compare(width(d), prec * magnitude(d));
|
||||
}
|
||||
|
||||
|
||||
template< bool Protected >
|
||||
class Is_valid< Interval_nt<Protected> >
|
||||
: public Unary_function< Interval_nt<Protected>, bool > {
|
||||
|
|
|
|||
|
|
@ -1039,9 +1039,8 @@ template < typename ET > class Real_embeddable_traits< Lazy_exact_nt<ET> >
|
|||
return r;
|
||||
|
||||
// If it's precise enough, then OK.
|
||||
if ((app.sup() - app.inf())
|
||||
< Lazy_exact_nt<ET>::get_relative_precision_of_to_double()
|
||||
* (std::max)(std::fabs(app.inf()), std::fabs(app.sup())) )
|
||||
if (compare_relative_precision(app,
|
||||
Lazy_exact_nt<ET>::get_relative_precision_of_to_double()) <= 0)
|
||||
return CGAL_NTS to_double(app);
|
||||
|
||||
CGAL_PROFILER(std::string("failures of : ") + std::string(CGAL_PRETTY_FUNCTION));
|
||||
|
|
|
|||
Loading…
Reference in New Issue