Eigne Matrix Class Library

Dependents:   MPC_current_control HydraulicControlBoard_SW AHRS Test_ekf ... more

Committer:
jsoh91
Date:
Tue Sep 24 00:18:23 2019 +0000
Revision:
1:3b8049da21b8
Parent:
0:13a5d365ba16
ignore and revise some of error parts

Who changed what in which revision?

UserRevisionLine numberNew contents of line
ykuroda 0:13a5d365ba16 1 // This file is part of Eigen, a lightweight C++ template library
ykuroda 0:13a5d365ba16 2 // for linear algebra.
ykuroda 0:13a5d365ba16 3 //
ykuroda 0:13a5d365ba16 4 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
ykuroda 0:13a5d365ba16 5 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
ykuroda 0:13a5d365ba16 6 //
ykuroda 0:13a5d365ba16 7 // This Source Code Form is subject to the terms of the Mozilla
ykuroda 0:13a5d365ba16 8 // Public License v. 2.0. If a copy of the MPL was not distributed
ykuroda 0:13a5d365ba16 9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
ykuroda 0:13a5d365ba16 10
ykuroda 0:13a5d365ba16 11 #ifndef EIGEN_FUZZY_H
ykuroda 0:13a5d365ba16 12 #define EIGEN_FUZZY_H
ykuroda 0:13a5d365ba16 13
ykuroda 0:13a5d365ba16 14 namespace Eigen {
ykuroda 0:13a5d365ba16 15
ykuroda 0:13a5d365ba16 16 namespace internal
ykuroda 0:13a5d365ba16 17 {
ykuroda 0:13a5d365ba16 18
ykuroda 0:13a5d365ba16 19 template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
ykuroda 0:13a5d365ba16 20 struct isApprox_selector
ykuroda 0:13a5d365ba16 21 {
ykuroda 0:13a5d365ba16 22 static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec)
ykuroda 0:13a5d365ba16 23 {
ykuroda 0:13a5d365ba16 24 using std::min;
ykuroda 0:13a5d365ba16 25 typename internal::nested<Derived,2>::type nested(x);
ykuroda 0:13a5d365ba16 26 typename internal::nested<OtherDerived,2>::type otherNested(y);
ykuroda 0:13a5d365ba16 27 return (nested - otherNested).cwiseAbs2().sum() <= prec * prec * (min)(nested.cwiseAbs2().sum(), otherNested.cwiseAbs2().sum());
ykuroda 0:13a5d365ba16 28 }
ykuroda 0:13a5d365ba16 29 };
ykuroda 0:13a5d365ba16 30
ykuroda 0:13a5d365ba16 31 template<typename Derived, typename OtherDerived>
ykuroda 0:13a5d365ba16 32 struct isApprox_selector<Derived, OtherDerived, true>
ykuroda 0:13a5d365ba16 33 {
ykuroda 0:13a5d365ba16 34 static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar&)
ykuroda 0:13a5d365ba16 35 {
ykuroda 0:13a5d365ba16 36 return x.matrix() == y.matrix();
ykuroda 0:13a5d365ba16 37 }
ykuroda 0:13a5d365ba16 38 };
ykuroda 0:13a5d365ba16 39
ykuroda 0:13a5d365ba16 40 template<typename Derived, typename OtherDerived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
ykuroda 0:13a5d365ba16 41 struct isMuchSmallerThan_object_selector
ykuroda 0:13a5d365ba16 42 {
ykuroda 0:13a5d365ba16 43 static bool run(const Derived& x, const OtherDerived& y, const typename Derived::RealScalar& prec)
ykuroda 0:13a5d365ba16 44 {
ykuroda 0:13a5d365ba16 45 return x.cwiseAbs2().sum() <= numext::abs2(prec) * y.cwiseAbs2().sum();
ykuroda 0:13a5d365ba16 46 }
ykuroda 0:13a5d365ba16 47 };
ykuroda 0:13a5d365ba16 48
ykuroda 0:13a5d365ba16 49 template<typename Derived, typename OtherDerived>
ykuroda 0:13a5d365ba16 50 struct isMuchSmallerThan_object_selector<Derived, OtherDerived, true>
ykuroda 0:13a5d365ba16 51 {
ykuroda 0:13a5d365ba16 52 static bool run(const Derived& x, const OtherDerived&, const typename Derived::RealScalar&)
ykuroda 0:13a5d365ba16 53 {
ykuroda 0:13a5d365ba16 54 return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
ykuroda 0:13a5d365ba16 55 }
ykuroda 0:13a5d365ba16 56 };
ykuroda 0:13a5d365ba16 57
ykuroda 0:13a5d365ba16 58 template<typename Derived, bool is_integer = NumTraits<typename Derived::Scalar>::IsInteger>
ykuroda 0:13a5d365ba16 59 struct isMuchSmallerThan_scalar_selector
ykuroda 0:13a5d365ba16 60 {
ykuroda 0:13a5d365ba16 61 static bool run(const Derived& x, const typename Derived::RealScalar& y, const typename Derived::RealScalar& prec)
ykuroda 0:13a5d365ba16 62 {
ykuroda 0:13a5d365ba16 63 return x.cwiseAbs2().sum() <= numext::abs2(prec * y);
ykuroda 0:13a5d365ba16 64 }
ykuroda 0:13a5d365ba16 65 };
ykuroda 0:13a5d365ba16 66
ykuroda 0:13a5d365ba16 67 template<typename Derived>
ykuroda 0:13a5d365ba16 68 struct isMuchSmallerThan_scalar_selector<Derived, true>
ykuroda 0:13a5d365ba16 69 {
ykuroda 0:13a5d365ba16 70 static bool run(const Derived& x, const typename Derived::RealScalar&, const typename Derived::RealScalar&)
ykuroda 0:13a5d365ba16 71 {
ykuroda 0:13a5d365ba16 72 return x.matrix() == Derived::Zero(x.rows(), x.cols()).matrix();
ykuroda 0:13a5d365ba16 73 }
ykuroda 0:13a5d365ba16 74 };
ykuroda 0:13a5d365ba16 75
ykuroda 0:13a5d365ba16 76 } // end namespace internal
ykuroda 0:13a5d365ba16 77
ykuroda 0:13a5d365ba16 78
ykuroda 0:13a5d365ba16 79 /** \returns \c true if \c *this is approximately equal to \a other, within the precision
ykuroda 0:13a5d365ba16 80 * determined by \a prec.
ykuroda 0:13a5d365ba16 81 *
ykuroda 0:13a5d365ba16 82 * \note The fuzzy compares are done multiplicatively. Two vectors \f$ v \f$ and \f$ w \f$
ykuroda 0:13a5d365ba16 83 * are considered to be approximately equal within precision \f$ p \f$ if
ykuroda 0:13a5d365ba16 84 * \f[ \Vert v - w \Vert \leqslant p\,\min(\Vert v\Vert, \Vert w\Vert). \f]
ykuroda 0:13a5d365ba16 85 * For matrices, the comparison is done using the Hilbert-Schmidt norm (aka Frobenius norm
ykuroda 0:13a5d365ba16 86 * L2 norm).
ykuroda 0:13a5d365ba16 87 *
ykuroda 0:13a5d365ba16 88 * \note Because of the multiplicativeness of this comparison, one can't use this function
ykuroda 0:13a5d365ba16 89 * to check whether \c *this is approximately equal to the zero matrix or vector.
ykuroda 0:13a5d365ba16 90 * Indeed, \c isApprox(zero) returns false unless \c *this itself is exactly the zero matrix
ykuroda 0:13a5d365ba16 91 * or vector. If you want to test whether \c *this is zero, use internal::isMuchSmallerThan(const
ykuroda 0:13a5d365ba16 92 * RealScalar&, RealScalar) instead.
ykuroda 0:13a5d365ba16 93 *
ykuroda 0:13a5d365ba16 94 * \sa internal::isMuchSmallerThan(const RealScalar&, RealScalar) const
ykuroda 0:13a5d365ba16 95 */
ykuroda 0:13a5d365ba16 96 template<typename Derived>
ykuroda 0:13a5d365ba16 97 template<typename OtherDerived>
ykuroda 0:13a5d365ba16 98 bool DenseBase<Derived>::isApprox(
ykuroda 0:13a5d365ba16 99 const DenseBase<OtherDerived>& other,
ykuroda 0:13a5d365ba16 100 const RealScalar& prec
ykuroda 0:13a5d365ba16 101 ) const
ykuroda 0:13a5d365ba16 102 {
ykuroda 0:13a5d365ba16 103 return internal::isApprox_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
ykuroda 0:13a5d365ba16 104 }
ykuroda 0:13a5d365ba16 105
ykuroda 0:13a5d365ba16 106 /** \returns \c true if the norm of \c *this is much smaller than \a other,
ykuroda 0:13a5d365ba16 107 * within the precision determined by \a prec.
ykuroda 0:13a5d365ba16 108 *
ykuroda 0:13a5d365ba16 109 * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is
ykuroda 0:13a5d365ba16 110 * considered to be much smaller than \f$ x \f$ within precision \f$ p \f$ if
ykuroda 0:13a5d365ba16 111 * \f[ \Vert v \Vert \leqslant p\,\vert x\vert. \f]
ykuroda 0:13a5d365ba16 112 *
ykuroda 0:13a5d365ba16 113 * For matrices, the comparison is done using the Hilbert-Schmidt norm. For this reason,
ykuroda 0:13a5d365ba16 114 * the value of the reference scalar \a other should come from the Hilbert-Schmidt norm
ykuroda 0:13a5d365ba16 115 * of a reference matrix of same dimensions.
ykuroda 0:13a5d365ba16 116 *
ykuroda 0:13a5d365ba16 117 * \sa isApprox(), isMuchSmallerThan(const DenseBase<OtherDerived>&, RealScalar) const
ykuroda 0:13a5d365ba16 118 */
ykuroda 0:13a5d365ba16 119 template<typename Derived>
ykuroda 0:13a5d365ba16 120 bool DenseBase<Derived>::isMuchSmallerThan(
ykuroda 0:13a5d365ba16 121 const typename NumTraits<Scalar>::Real& other,
ykuroda 0:13a5d365ba16 122 const RealScalar& prec
ykuroda 0:13a5d365ba16 123 ) const
ykuroda 0:13a5d365ba16 124 {
ykuroda 0:13a5d365ba16 125 return internal::isMuchSmallerThan_scalar_selector<Derived>::run(derived(), other, prec);
ykuroda 0:13a5d365ba16 126 }
ykuroda 0:13a5d365ba16 127
ykuroda 0:13a5d365ba16 128 /** \returns \c true if the norm of \c *this is much smaller than the norm of \a other,
ykuroda 0:13a5d365ba16 129 * within the precision determined by \a prec.
ykuroda 0:13a5d365ba16 130 *
ykuroda 0:13a5d365ba16 131 * \note The fuzzy compares are done multiplicatively. A vector \f$ v \f$ is
ykuroda 0:13a5d365ba16 132 * considered to be much smaller than a vector \f$ w \f$ within precision \f$ p \f$ if
ykuroda 0:13a5d365ba16 133 * \f[ \Vert v \Vert \leqslant p\,\Vert w\Vert. \f]
ykuroda 0:13a5d365ba16 134 * For matrices, the comparison is done using the Hilbert-Schmidt norm.
ykuroda 0:13a5d365ba16 135 *
ykuroda 0:13a5d365ba16 136 * \sa isApprox(), isMuchSmallerThan(const RealScalar&, RealScalar) const
ykuroda 0:13a5d365ba16 137 */
ykuroda 0:13a5d365ba16 138 template<typename Derived>
ykuroda 0:13a5d365ba16 139 template<typename OtherDerived>
ykuroda 0:13a5d365ba16 140 bool DenseBase<Derived>::isMuchSmallerThan(
ykuroda 0:13a5d365ba16 141 const DenseBase<OtherDerived>& other,
ykuroda 0:13a5d365ba16 142 const RealScalar& prec
ykuroda 0:13a5d365ba16 143 ) const
ykuroda 0:13a5d365ba16 144 {
ykuroda 0:13a5d365ba16 145 return internal::isMuchSmallerThan_object_selector<Derived, OtherDerived>::run(derived(), other.derived(), prec);
ykuroda 0:13a5d365ba16 146 }
ykuroda 0:13a5d365ba16 147
ykuroda 0:13a5d365ba16 148 } // end namespace Eigen
ykuroda 0:13a5d365ba16 149
ykuroda 0:13a5d365ba16 150 #endif // EIGEN_FUZZY_H