Updated fuzzy comparisons to use L2 norm as all my experiments

tends to show L2 norm works very well here.
(the legacy implementation is still available via a preprocessor token
 to allow further experiments if needed...)
This commit is contained in:
Gael Guennebaud 2008-06-06 18:37:53 +00:00
parent 8769bfd9aa
commit a172385720

View File

@ -26,6 +26,78 @@
#ifndef EIGEN_FUZZY_H
#define EIGEN_FUZZY_H
#ifndef EIGEN_LEGACY_COMPARES
/** \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 ei_isMuchSmallerThan(const
* RealScalar&, RealScalar) instead.
*
* \sa ei_isMuchSmallerThan(const RealScalar&, RealScalar) const
*/
template<typename Derived>
template<typename OtherDerived>
bool MatrixBase<Derived>::isApprox(
const MatrixBase<OtherDerived>& other,
typename NumTraits<Scalar>::Real prec
) const
{
const typename ei_nested<Derived,2>::type nested(derived());
const typename ei_nested<OtherDerived,2>::type otherNested(other.derived());
return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * std::min(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum());
}
/** \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.
*
* \sa isApprox(), isMuchSmallerThan(const MatrixBase<OtherDerived>&, RealScalar) const
*/
template<typename Derived>
bool MatrixBase<Derived>::isMuchSmallerThan(
const typename NumTraits<Scalar>::Real& other,
typename NumTraits<Scalar>::Real prec
) const
{
return cwiseAbs2().sum() <= prec * prec * other * other * cols() * rows();
}
/** \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 MatrixBase<Derived>::isMuchSmallerThan(
const MatrixBase<OtherDerived>& other,
typename NumTraits<Scalar>::Real prec
) const
{
return this->cwiseAbs2().sum() <= prec * prec * other.cwiseAbs2().sum();
}
#else
template<typename Derived, typename OtherDerived=Derived, bool IsVector=Derived::IsVectorAtCompileTime>
struct ei_fuzzy_selector;
@ -154,4 +226,6 @@ struct ei_fuzzy_selector<Derived,OtherDerived,false>
}
};
#endif
#endif // EIGEN_FUZZY_H