Fuzzy.h 5.46 KB
Newer Older
LM's avatar
LM committed
1 2 3 4 5 6
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
//
Don Gagne's avatar
Don Gagne committed
7 8 9
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
LM's avatar
LM committed
10 11 12 13

#ifndef EIGEN_FUZZY_H
#define EIGEN_FUZZY_H

Don Gagne's avatar
Don Gagne committed
14 15
namespace Eigen { 

LM's avatar
LM committed
16 17 18 19 20 21
namespace internal
{

template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
struct isApprox_selector
{
Don Gagne's avatar
Don Gagne committed
22
  static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec)
LM's avatar
LM committed
23 24
  {
    using std::min;
Don Gagne's avatar
Don Gagne committed
25 26 27
    typename internal::nested<Derived,2>::type nested(x);
    typename internal::nested<OtherDerived,2>::type otherNested(y);
    return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * (min)(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum());
LM's avatar
LM committed
28 29 30 31 32 33
  }
};

template<typename Derived, typename OtherDerived>
struct isApprox_selector<Derived, OtherDerived, true>
{
Don Gagne's avatar
Don Gagne committed
34
  static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar&)
LM's avatar
LM committed
35 36 37 38 39 40 41 42
  {
    return x.matrix() == y.matrix();
  }
};

template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
struct isMuchSmallerThan_object_selector
{
Don Gagne's avatar
Don Gagne committed
43
  static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec)
LM's avatar
LM committed
44
  {
Don Gagne's avatar
Don Gagne committed
45
    return x.cwiseAbs2().sum() <= numext::abs2(prec) * y.cwiseAbs2().sum();
LM's avatar
LM committed
46 47 48 49 50 51
  }
};

template<typename Derived, typename OtherDerived>
struct isMuchSmallerThan_object_selector<Derived, OtherDerived, true>
{
Don Gagne's avatar
Don Gagne committed
52
  static bool run(const Derived& x, const OtherDerived&, const typename Derived::RealScalar&)
LM's avatar
LM committed
53 54 55 56 57 58 59 60
  {
    return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
  }
};

template<typename Derived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
struct isMuchSmallerThan_scalar_selector
{
Don Gagne's avatar
Don Gagne committed
61
  static bool run(const Derived& x, const typename Derived::RealScalar& y, const typename Derived::RealScalar& prec)
LM's avatar
LM committed
62
  {
Don Gagne's avatar
Don Gagne committed
63
    return x.cwiseAbs2().sum() <= numext::abs2(prec * y);
LM's avatar
LM committed
64 65 66 67 68 69
  }
};

template<typename Derived>
struct isMuchSmallerThan_scalar_selector<Derived, true>
{
Don Gagne's avatar
Don Gagne committed
70
  static bool run(const Derived& x, const typename Derived::RealScalar&, const typename Derived::RealScalar&)
LM's avatar
LM committed
71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99
  {
    return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
  }
};

} // end namespace internal


/** \returns \c true if \c *this is approximately equal to \a other, within the precision
  * determined by \a prec.
  *
  * \note The fuzzy compares are done multiplicatively. Two vectors \f$ v \f$ and \f$ w \f$
  * are considered to be approximately equal within precision \f$ p \f$ if
  * \f[ \Vert v - w \Vert \leqslant p\,\min(\Vert v\Vert, \Vert w\Vert). \f]
  * For matrices, the comparison is done using the Hilbert-Schmidt norm (aka Frobenius norm
  * L2 norm).
  *
  * \note Because of the multiplicativeness of this comparison, one can't use this function
  * to check whether \c *this is approximately equal to the zero matrix or vector.
  * Indeed, \c isApprox(zero) returns false unless \c *this itself is exactly the zero matrix
  * or vector. If you want to test whether \c *this is zero, use internal::isMuchSmallerThan(const
  * RealScalar&, RealScalar) instead.
  *
  * \sa internal::isMuchSmallerThan(const RealScalar&, RealScalar) const
  */
template<typename Derived>
template<typename OtherDerived>
bool DenseBase<Derived>::isApprox(
  const DenseBase<OtherDerived>& other,
Don Gagne's avatar
Don Gagne committed
100
  const RealScalar& prec
LM's avatar
LM committed
101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121
) const
{
  return internal::isApprox_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
}

/** \returns \c true if the norm of \c *this is much smaller than \a other,
  * within the precision determined by \a prec.
  *
  * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is
  * considered to be much smaller than \f$ x \f$ within precision \f$ p \f$ if
  * \f[ \Vert v \Vert \leqslant p\,\vert x\vert. \f]
  *
  * For matrices, the comparison is done using the Hilbert-Schmidt norm. For this reason,
  * the value of the reference scalar \a other should come from the Hilbert-Schmidt norm
  * of a reference matrix of same dimensions.
  *
  * \sa isApprox(), isMuchSmallerThan(const DenseBase<OtherDerived>&, RealScalar) const
  */
template<typename Derived>
bool DenseBase<Derived>::isMuchSmallerThan(
  const typename NumTraits<Scalar>::Real& other,
Don Gagne's avatar
Don Gagne committed
122
  const RealScalar& prec
LM's avatar
LM committed
123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141
) const
{
  return internal::isMuchSmallerThan_scalar_selector<Derived>::run(derived(), other, prec);
}

/** \returns \c true if the norm of \c *this is much smaller than the norm of \a other,
  * within the precision determined by \a prec.
  *
  * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is
  * considered to be much smaller than a vector \f$ w \f$ within precision \f$ p \f$ if
  * \f[ \Vert v \Vert \leqslant p\,\Vert w\Vert. \f]
  * For matrices, the comparison is done using the Hilbert-Schmidt norm.
  *
  * \sa isApprox(), isMuchSmallerThan(const RealScalar&, RealScalar) const
  */
template<typename Derived>
template<typename OtherDerived>
bool DenseBase<Derived>::isMuchSmallerThan(
  const DenseBase<OtherDerived>& other,
Don Gagne's avatar
Don Gagne committed
142
  const RealScalar& prec
LM's avatar
LM committed
143 144 145 146 147
) const
{
  return internal::isMuchSmallerThan_object_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
}

Don Gagne's avatar
Don Gagne committed
148 149
} // end namespace Eigen

LM's avatar
LM committed
150
#endif // EIGEN_FUZZY_H