Fuzzy.h 5.6 KB

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