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 | |
14 | namespace Eigen { |
15 | |
16 | namespace internal |
17 | { |
18 | |
19 | template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger> |
20 | struct isApprox_selector |
21 | { |
22 | EIGEN_DEVICE_FUNC |
23 | static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec) |
24 | { |
25 | typename internal::nested_eval<Derived,2>::type nested(x); |
26 | typename internal::nested_eval<OtherDerived,2>::type otherNested(y); |
27 | return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * numext::mini(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum()); |
28 | } |
29 | }; |
30 | |
31 | template<typename Derived, typename OtherDerived> |
32 | struct isApprox_selector<Derived, OtherDerived, true> |
33 | { |
34 | EIGEN_DEVICE_FUNC |
35 | static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar&) |
36 | { |
37 | return x.matrix() == y.matrix(); |
38 | } |
39 | }; |
40 | |
41 | template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger> |
42 | struct isMuchSmallerThan_object_selector |
43 | { |
44 | EIGEN_DEVICE_FUNC |
45 | static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec) |
46 | { |
47 | return x.cwiseAbs2().sum() <= numext::abs2(prec) * y.cwiseAbs2().sum(); |
48 | } |
49 | }; |
50 | |
51 | template<typename Derived, typename OtherDerived> |
52 | struct isMuchSmallerThan_object_selector<Derived, OtherDerived, true> |
53 | { |
54 | EIGEN_DEVICE_FUNC |
55 | static bool run(const Derived& x, const OtherDerived&, const typename Derived::RealScalar&) |
56 | { |
57 | return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix(); |
58 | } |
59 | }; |
60 | |
61 | template<typename Derived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger> |
62 | struct isMuchSmallerThan_scalar_selector |
63 | { |
64 | EIGEN_DEVICE_FUNC |
65 | static bool run(const Derived& x, const typename Derived::RealScalar& y, const typename Derived::RealScalar& prec) |
66 | { |
67 | return x.cwiseAbs2().sum() <= numext::abs2(prec * y); |
68 | } |
69 | }; |
70 | |
71 | template<typename Derived> |
72 | struct isMuchSmallerThan_scalar_selector<Derived, true> |
73 | { |
74 | EIGEN_DEVICE_FUNC |
75 | static bool run(const Derived& x, const typename Derived::RealScalar&, const typename Derived::RealScalar&) |
76 | { |
77 | return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix(); |
78 | } |
79 | }; |
80 | |
81 | } // end namespace internal |
82 | |
83 | |
84 | /** \returns \c true if \c *this is approximately equal to \a other, within the precision |
85 | * determined by \a prec. |
86 | * |
87 | * \note The fuzzy compares are done multiplicatively. Two vectors \f$ v \f$ and \f$ w \f$ |
88 | * are considered to be approximately equal within precision \f$ p \f$ if |
89 | * \f[ \Vert v - w \Vert \leqslant p\,\min(\Vert v\Vert, \Vert w\Vert). \f] |
90 | * For matrices, the comparison is done using the Hilbert-Schmidt norm (aka Frobenius norm |
91 | * L2 norm). |
92 | * |
93 | * \note Because of the multiplicativeness of this comparison, one can't use this function |
94 | * to check whether \c *this is approximately equal to the zero matrix or vector. |
95 | * Indeed, \c isApprox(zero) returns false unless \c *this itself is exactly the zero matrix |
96 | * or vector. If you want to test whether \c *this is zero, use internal::isMuchSmallerThan(const |
97 | * RealScalar&, RealScalar) instead. |
98 | * |
99 | * \sa internal::isMuchSmallerThan(const RealScalar&, RealScalar) const |
100 | */ |
101 | template<typename Derived> |
102 | template<typename OtherDerived> |
103 | bool DenseBase<Derived>::isApprox( |
104 | const DenseBase<OtherDerived>& other, |
105 | const RealScalar& prec |
106 | ) const |
107 | { |
108 | return internal::isApprox_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec); |
109 | } |
110 | |
111 | /** \returns \c true if the norm of \c *this is much smaller than \a other, |
112 | * within the precision determined by \a prec. |
113 | * |
114 | * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is |
115 | * considered to be much smaller than \f$ x \f$ within precision \f$ p \f$ if |
116 | * \f[ \Vert v \Vert \leqslant p\,\vert x\vert. \f] |
117 | * |
118 | * For matrices, the comparison is done using the Hilbert-Schmidt norm. For this reason, |
119 | * the value of the reference scalar \a other should come from the Hilbert-Schmidt norm |
120 | * of a reference matrix of same dimensions. |
121 | * |
122 | * \sa isApprox(), isMuchSmallerThan(const DenseBase<OtherDerived>&, RealScalar) const |
123 | */ |
124 | template<typename Derived> |
125 | bool DenseBase<Derived>::isMuchSmallerThan( |
126 | const typename NumTraits<Scalar>::Real& other, |
127 | const RealScalar& prec |
128 | ) const |
129 | { |
130 | return internal::isMuchSmallerThan_scalar_selector<Derived>::run(derived(), other, prec); |
131 | } |
132 | |
133 | /** \returns \c true if the norm of \c *this is much smaller than the norm of \a other, |
134 | * within the precision determined by \a prec. |
135 | * |
136 | * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is |
137 | * considered to be much smaller than a vector \f$ w \f$ within precision \f$ p \f$ if |
138 | * \f[ \Vert v \Vert \leqslant p\,\Vert w\Vert. \f] |
139 | * For matrices, the comparison is done using the Hilbert-Schmidt norm. |
140 | * |
141 | * \sa isApprox(), isMuchSmallerThan(const RealScalar&, RealScalar) const |
142 | */ |
143 | template<typename Derived> |
144 | template<typename OtherDerived> |
145 | bool DenseBase<Derived>::isMuchSmallerThan( |
146 | const DenseBase<OtherDerived>& other, |
147 | const RealScalar& prec |
148 | ) const |
149 | { |
150 | return internal::isMuchSmallerThan_object_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec); |
151 | } |
152 | |
153 | } // end namespace Eigen |
154 | |
155 | #endif // EIGEN_FUZZY_H |
156 | |