Eigne Matrix Class Library
Dependents: Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more
Fuzzy.h
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com> 00005 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr> 00006 // 00007 // This Source Code Form is subject to the terms of the Mozilla 00008 // Public License v. 2.0. If a copy of the MPL was not distributed 00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 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 } // end namespace internal 00077 00078 00079 /** \returns \c true if \c *this is approximately equal to \a other, within the precision 00080 * determined by \a prec. 00081 * 00082 * \note The fuzzy compares are done multiplicatively. Two vectors \f$ v \f$ and \f$ w \f$ 00083 * are considered to be approximately equal within precision \f$ p \f$ if 00084 * \f[ \Vert v - w \Vert \leqslant p\,\min(\Vert v\Vert, \Vert w\Vert). \f] 00085 * For matrices, the comparison is done using the Hilbert-Schmidt norm (aka Frobenius norm 00086 * L2 norm). 00087 * 00088 * \note Because of the multiplicativeness of this comparison, one can't use this function 00089 * to check whether \c *this is approximately equal to the zero matrix or vector. 00090 * Indeed, \c isApprox(zero) returns false unless \c *this itself is exactly the zero matrix 00091 * or vector. If you want to test whether \c *this is zero, use internal::isMuchSmallerThan(const 00092 * RealScalar&, RealScalar) instead. 00093 * 00094 * \sa internal::isMuchSmallerThan(const RealScalar&, RealScalar) const 00095 */ 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 00106 /** \returns \c true if the norm of \c *this is much smaller than \a other, 00107 * within the precision determined by \a prec. 00108 * 00109 * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is 00110 * considered to be much smaller than \f$ x \f$ within precision \f$ p \f$ if 00111 * \f[ \Vert v \Vert \leqslant p\,\vert x\vert. \f] 00112 * 00113 * For matrices, the comparison is done using the Hilbert-Schmidt norm. For this reason, 00114 * the value of the reference scalar \a other should come from the Hilbert-Schmidt norm 00115 * of a reference matrix of same dimensions. 00116 * 00117 * \sa isApprox(), isMuchSmallerThan(const DenseBase<OtherDerived>&, RealScalar) const 00118 */ 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 00128 /** \returns \c true if the norm of \c *this is much smaller than the norm of \a other, 00129 * within the precision determined by \a prec. 00130 * 00131 * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is 00132 * considered to be much smaller than a vector \f$ w \f$ within precision \f$ p \f$ if 00133 * \f[ \Vert v \Vert \leqslant p\,\Vert w\Vert. \f] 00134 * For matrices, the comparison is done using the Hilbert-Schmidt norm. 00135 * 00136 * \sa isApprox(), isMuchSmallerThan(const RealScalar&, RealScalar) const 00137 */ 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 } // end namespace Eigen 00149 00150 #endif // EIGEN_FUZZY_H
Generated on Tue Jul 12 2022 17:46:53 by 1.7.2