Fuzzy.h revision c981c48f5bc9aefeffc0bcb0cc3934c2fae179dd
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
5// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11#ifndef EIGEN_FUZZY_H
12#define EIGEN_FUZZY_H
13
14namespace Eigen {
15
16namespace internal
17{
18
19template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
20struct isApprox_selector
21{
22  static bool run(const Derived& x, const OtherDerived& y, typename Derived::RealScalar prec)
23  {
24    using std::min;
25    typename internal::nested<Derived,2>::type nested(x);
26    typename internal::nested<OtherDerived,2>::type otherNested(y);
27    return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * (min)(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum());
28  }
29};
30
31template<typename Derived, typename OtherDerived>
32struct isApprox_selector<Derived, OtherDerived, true>
33{
34  static bool run(const Derived& x, const OtherDerived& y, typename Derived::RealScalar)
35  {
36    return x.matrix() == y.matrix();
37  }
38};
39
40template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
41struct isMuchSmallerThan_object_selector
42{
43  static bool run(const Derived& x, const OtherDerived& y, typename Derived::RealScalar prec)
44  {
45    return x.cwiseAbs2().sum() <= abs2(prec) * y.cwiseAbs2().sum();
46  }
47};
48
49template<typename Derived, typename OtherDerived>
50struct isMuchSmallerThan_object_selector<Derived, OtherDerived, true>
51{
52  static bool run(const Derived& x, const OtherDerived&, typename Derived::RealScalar)
53  {
54    return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
55  }
56};
57
58template<typename Derived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
59struct isMuchSmallerThan_scalar_selector
60{
61  static bool run(const Derived& x, const typename Derived::RealScalar& y, typename Derived::RealScalar prec)
62  {
63    return x.cwiseAbs2().sum() <= abs2(prec * y);
64  }
65};
66
67template<typename Derived>
68struct isMuchSmallerThan_scalar_selector<Derived, true>
69{
70  static bool run(const Derived& x, const typename Derived::RealScalar&, typename Derived::RealScalar)
71  {
72    return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
73  }
74};
75
76} // end namespace internal
77
78
79/** \returns \c true if \c *this is approximately equal to \a other, within the precision
80  * determined by \a prec.
81  *
82  * \note The fuzzy compares are done multiplicatively. Two vectors \f$ v \f$ and \f$ w \f$
83  * are considered to be approximately equal within precision \f$ p \f$ if
84  * \f[ \Vert v - w \Vert \leqslant p\,\min(\Vert v\Vert, \Vert w\Vert). \f]
85  * For matrices, the comparison is done using the Hilbert-Schmidt norm (aka Frobenius norm
86  * L2 norm).
87  *
88  * \note Because of the multiplicativeness of this comparison, one can't use this function
89  * to check whether \c *this is approximately equal to the zero matrix or vector.
90  * Indeed, \c isApprox(zero) returns false unless \c *this itself is exactly the zero matrix
91  * or vector. If you want to test whether \c *this is zero, use internal::isMuchSmallerThan(const
92  * RealScalar&, RealScalar) instead.
93  *
94  * \sa internal::isMuchSmallerThan(const RealScalar&, RealScalar) const
95  */
96template<typename Derived>
97template<typename OtherDerived>
98bool DenseBase<Derived>::isApprox(
99  const DenseBase<OtherDerived>& other,
100  RealScalar prec
101) const
102{
103  return internal::isApprox_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
104}
105
106/** \returns \c true if the norm of \c *this is much smaller than \a other,
107  * within the precision determined by \a prec.
108  *
109  * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is
110  * considered to be much smaller than \f$ x \f$ within precision \f$ p \f$ if
111  * \f[ \Vert v \Vert \leqslant p\,\vert x\vert. \f]
112  *
113  * For matrices, the comparison is done using the Hilbert-Schmidt norm. For this reason,
114  * the value of the reference scalar \a other should come from the Hilbert-Schmidt norm
115  * of a reference matrix of same dimensions.
116  *
117  * \sa isApprox(), isMuchSmallerThan(const DenseBase<OtherDerived>&, RealScalar) const
118  */
119template<typename Derived>
120bool DenseBase<Derived>::isMuchSmallerThan(
121  const typename NumTraits<Scalar>::Real& other,
122  RealScalar prec
123) const
124{
125  return internal::isMuchSmallerThan_scalar_selector<Derived>::run(derived(), other, prec);
126}
127
128/** \returns \c true if the norm of \c *this is much smaller than the norm of \a other,
129  * within the precision determined by \a prec.
130  *
131  * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is
132  * considered to be much smaller than a vector \f$ w \f$ within precision \f$ p \f$ if
133  * \f[ \Vert v \Vert \leqslant p\,\Vert w\Vert. \f]
134  * For matrices, the comparison is done using the Hilbert-Schmidt norm.
135  *
136  * \sa isApprox(), isMuchSmallerThan(const RealScalar&, RealScalar) const
137  */
138template<typename Derived>
139template<typename OtherDerived>
140bool DenseBase<Derived>::isMuchSmallerThan(
141  const DenseBase<OtherDerived>& other,
142  RealScalar prec
143) const
144{
145  return internal::isMuchSmallerThan_object_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
146}
147
148} // end namespace Eigen
149
150#endif // EIGEN_FUZZY_H
151