Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #ifndef EIGEN_FUZZY_H
00012 #define EIGEN_FUZZY_H
00013
00014 namespace Eigen {
00015
00016 namespace internal
00017 {
00018
00019 template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
00020 struct isApprox_selector
00021 {
00022 static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec)
00023 {
00024 using std::min;
00025 typename internal::nested<Derived,2>::type nested(x);
00026 typename internal::nested<OtherDerived,2>::type otherNested(y);
00027 return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * (min)(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum());
00028 }
00029 };
00030
00031 template<typename Derived, typename OtherDerived>
00032 struct isApprox_selector<Derived, OtherDerived, true>
00033 {
00034 static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar&)
00035 {
00036 return x.matrix() == y.matrix();
00037 }
00038 };
00039
00040 template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
00041 struct isMuchSmallerThan_object_selector
00042 {
00043 static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec)
00044 {
00045 return x.cwiseAbs2().sum() <= numext::abs2(prec) * y.cwiseAbs2().sum();
00046 }
00047 };
00048
00049 template<typename Derived, typename OtherDerived>
00050 struct isMuchSmallerThan_object_selector<Derived, OtherDerived, true>
00051 {
00052 static bool run(const Derived& x, const OtherDerived&, const typename Derived::RealScalar&)
00053 {
00054 return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
00055 }
00056 };
00057
00058 template<typename Derived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
00059 struct isMuchSmallerThan_scalar_selector
00060 {
00061 static bool run(const Derived& x, const typename Derived::RealScalar& y, const typename Derived::RealScalar& prec)
00062 {
00063 return x.cwiseAbs2().sum() <= numext::abs2(prec * y);
00064 }
00065 };
00066
00067 template<typename Derived>
00068 struct isMuchSmallerThan_scalar_selector<Derived, true>
00069 {
00070 static bool run(const Derived& x, const typename Derived::RealScalar&, const typename Derived::RealScalar&)
00071 {
00072 return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
00073 }
00074 };
00075
00076 }
00077
00078
00096 template<typename Derived>
00097 template<typename OtherDerived>
00098 bool DenseBase<Derived>::isApprox(
00099 const DenseBase<OtherDerived>& other,
00100 const RealScalar& prec
00101 ) const
00102 {
00103 return internal::isApprox_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
00104 }
00105
00119 template<typename Derived>
00120 bool DenseBase<Derived>::isMuchSmallerThan(
00121 const typename NumTraits<Scalar>::Real& other,
00122 const RealScalar& prec
00123 ) const
00124 {
00125 return internal::isMuchSmallerThan_scalar_selector<Derived>::run(derived(), other, prec);
00126 }
00127
00138 template<typename Derived>
00139 template<typename OtherDerived>
00140 bool DenseBase<Derived>::isMuchSmallerThan(
00141 const DenseBase<OtherDerived>& other,
00142 const RealScalar& prec
00143 ) const
00144 {
00145 return internal::isMuchSmallerThan_object_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
00146 }
00147
00148 }
00149
00150 #endif // EIGEN_FUZZY_H