Eigne Matrix Class Library
Dependents: Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more
Eigen Matrix Class Library for mbed.
Finally, you can use Eigen on your mbed!!!
src/Core/EigenBase.h@0:13a5d365ba16, 2016-10-13 (annotated)
- Committer:
- ykuroda
- Date:
- Thu Oct 13 04:07:23 2016 +0000
- Revision:
- 0:13a5d365ba16
First commint, Eigne Matrix Class Library
Who changed what in which revision?
User | Revision | Line number | New 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) 2009 Benoit Jacob <jacob.benoit.1@gmail.com> |
ykuroda | 0:13a5d365ba16 | 5 | // Copyright (C) 2009 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_EIGENBASE_H |
ykuroda | 0:13a5d365ba16 | 12 | #define EIGEN_EIGENBASE_H |
ykuroda | 0:13a5d365ba16 | 13 | |
ykuroda | 0:13a5d365ba16 | 14 | namespace Eigen { |
ykuroda | 0:13a5d365ba16 | 15 | |
ykuroda | 0:13a5d365ba16 | 16 | /** Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor MatrixBase(T). |
ykuroda | 0:13a5d365ba16 | 17 | * |
ykuroda | 0:13a5d365ba16 | 18 | * In other words, an EigenBase object is an object that can be copied into a MatrixBase. |
ykuroda | 0:13a5d365ba16 | 19 | * |
ykuroda | 0:13a5d365ba16 | 20 | * Besides MatrixBase-derived classes, this also includes special matrix classes such as diagonal matrices, etc. |
ykuroda | 0:13a5d365ba16 | 21 | * |
ykuroda | 0:13a5d365ba16 | 22 | * Notice that this class is trivial, it is only used to disambiguate overloaded functions. |
ykuroda | 0:13a5d365ba16 | 23 | * |
ykuroda | 0:13a5d365ba16 | 24 | * \sa \ref TopicClassHierarchy |
ykuroda | 0:13a5d365ba16 | 25 | */ |
ykuroda | 0:13a5d365ba16 | 26 | template<typename Derived> struct EigenBase |
ykuroda | 0:13a5d365ba16 | 27 | { |
ykuroda | 0:13a5d365ba16 | 28 | // typedef typename internal::plain_matrix_type<Derived>::type PlainObject; |
ykuroda | 0:13a5d365ba16 | 29 | |
ykuroda | 0:13a5d365ba16 | 30 | typedef typename internal::traits<Derived>::StorageKind StorageKind; |
ykuroda | 0:13a5d365ba16 | 31 | typedef typename internal::traits<Derived>::Index Index; |
ykuroda | 0:13a5d365ba16 | 32 | |
ykuroda | 0:13a5d365ba16 | 33 | /** \returns a reference to the derived object */ |
ykuroda | 0:13a5d365ba16 | 34 | Derived& derived() { return *static_cast<Derived*>(this); } |
ykuroda | 0:13a5d365ba16 | 35 | /** \returns a const reference to the derived object */ |
ykuroda | 0:13a5d365ba16 | 36 | const Derived& derived() const { return *static_cast<const Derived*>(this); } |
ykuroda | 0:13a5d365ba16 | 37 | |
ykuroda | 0:13a5d365ba16 | 38 | inline Derived& const_cast_derived() const |
ykuroda | 0:13a5d365ba16 | 39 | { return *static_cast<Derived*>(const_cast<EigenBase*>(this)); } |
ykuroda | 0:13a5d365ba16 | 40 | inline const Derived& const_derived() const |
ykuroda | 0:13a5d365ba16 | 41 | { return *static_cast<const Derived*>(this); } |
ykuroda | 0:13a5d365ba16 | 42 | |
ykuroda | 0:13a5d365ba16 | 43 | /** \returns the number of rows. \sa cols(), RowsAtCompileTime */ |
ykuroda | 0:13a5d365ba16 | 44 | inline Index rows() const { return derived().rows(); } |
ykuroda | 0:13a5d365ba16 | 45 | /** \returns the number of columns. \sa rows(), ColsAtCompileTime*/ |
ykuroda | 0:13a5d365ba16 | 46 | inline Index cols() const { return derived().cols(); } |
ykuroda | 0:13a5d365ba16 | 47 | /** \returns the number of coefficients, which is rows()*cols(). |
ykuroda | 0:13a5d365ba16 | 48 | * \sa rows(), cols(), SizeAtCompileTime. */ |
ykuroda | 0:13a5d365ba16 | 49 | inline Index size() const { return rows() * cols(); } |
ykuroda | 0:13a5d365ba16 | 50 | |
ykuroda | 0:13a5d365ba16 | 51 | /** \internal Don't use it, but do the equivalent: \code dst = *this; \endcode */ |
ykuroda | 0:13a5d365ba16 | 52 | template<typename Dest> inline void evalTo(Dest& dst) const |
ykuroda | 0:13a5d365ba16 | 53 | { derived().evalTo(dst); } |
ykuroda | 0:13a5d365ba16 | 54 | |
ykuroda | 0:13a5d365ba16 | 55 | /** \internal Don't use it, but do the equivalent: \code dst += *this; \endcode */ |
ykuroda | 0:13a5d365ba16 | 56 | template<typename Dest> inline void addTo(Dest& dst) const |
ykuroda | 0:13a5d365ba16 | 57 | { |
ykuroda | 0:13a5d365ba16 | 58 | // This is the default implementation, |
ykuroda | 0:13a5d365ba16 | 59 | // derived class can reimplement it in a more optimized way. |
ykuroda | 0:13a5d365ba16 | 60 | typename Dest::PlainObject res(rows(),cols()); |
ykuroda | 0:13a5d365ba16 | 61 | evalTo(res); |
ykuroda | 0:13a5d365ba16 | 62 | dst += res; |
ykuroda | 0:13a5d365ba16 | 63 | } |
ykuroda | 0:13a5d365ba16 | 64 | |
ykuroda | 0:13a5d365ba16 | 65 | /** \internal Don't use it, but do the equivalent: \code dst -= *this; \endcode */ |
ykuroda | 0:13a5d365ba16 | 66 | template<typename Dest> inline void subTo(Dest& dst) const |
ykuroda | 0:13a5d365ba16 | 67 | { |
ykuroda | 0:13a5d365ba16 | 68 | // This is the default implementation, |
ykuroda | 0:13a5d365ba16 | 69 | // derived class can reimplement it in a more optimized way. |
ykuroda | 0:13a5d365ba16 | 70 | typename Dest::PlainObject res(rows(),cols()); |
ykuroda | 0:13a5d365ba16 | 71 | evalTo(res); |
ykuroda | 0:13a5d365ba16 | 72 | dst -= res; |
ykuroda | 0:13a5d365ba16 | 73 | } |
ykuroda | 0:13a5d365ba16 | 74 | |
ykuroda | 0:13a5d365ba16 | 75 | /** \internal Don't use it, but do the equivalent: \code dst.applyOnTheRight(*this); \endcode */ |
ykuroda | 0:13a5d365ba16 | 76 | template<typename Dest> inline void applyThisOnTheRight(Dest& dst) const |
ykuroda | 0:13a5d365ba16 | 77 | { |
ykuroda | 0:13a5d365ba16 | 78 | // This is the default implementation, |
ykuroda | 0:13a5d365ba16 | 79 | // derived class can reimplement it in a more optimized way. |
ykuroda | 0:13a5d365ba16 | 80 | dst = dst * this->derived(); |
ykuroda | 0:13a5d365ba16 | 81 | } |
ykuroda | 0:13a5d365ba16 | 82 | |
ykuroda | 0:13a5d365ba16 | 83 | /** \internal Don't use it, but do the equivalent: \code dst.applyOnTheLeft(*this); \endcode */ |
ykuroda | 0:13a5d365ba16 | 84 | template<typename Dest> inline void applyThisOnTheLeft(Dest& dst) const |
ykuroda | 0:13a5d365ba16 | 85 | { |
ykuroda | 0:13a5d365ba16 | 86 | // This is the default implementation, |
ykuroda | 0:13a5d365ba16 | 87 | // derived class can reimplement it in a more optimized way. |
ykuroda | 0:13a5d365ba16 | 88 | dst = this->derived() * dst; |
ykuroda | 0:13a5d365ba16 | 89 | } |
ykuroda | 0:13a5d365ba16 | 90 | |
ykuroda | 0:13a5d365ba16 | 91 | }; |
ykuroda | 0:13a5d365ba16 | 92 | |
ykuroda | 0:13a5d365ba16 | 93 | /*************************************************************************** |
ykuroda | 0:13a5d365ba16 | 94 | * Implementation of matrix base methods |
ykuroda | 0:13a5d365ba16 | 95 | ***************************************************************************/ |
ykuroda | 0:13a5d365ba16 | 96 | |
ykuroda | 0:13a5d365ba16 | 97 | /** \brief Copies the generic expression \a other into *this. |
ykuroda | 0:13a5d365ba16 | 98 | * |
ykuroda | 0:13a5d365ba16 | 99 | * \details The expression must provide a (templated) evalTo(Derived& dst) const |
ykuroda | 0:13a5d365ba16 | 100 | * function which does the actual job. In practice, this allows any user to write |
ykuroda | 0:13a5d365ba16 | 101 | * its own special matrix without having to modify MatrixBase |
ykuroda | 0:13a5d365ba16 | 102 | * |
ykuroda | 0:13a5d365ba16 | 103 | * \returns a reference to *this. |
ykuroda | 0:13a5d365ba16 | 104 | */ |
ykuroda | 0:13a5d365ba16 | 105 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 106 | template<typename OtherDerived> |
ykuroda | 0:13a5d365ba16 | 107 | Derived& DenseBase<Derived>::operator=(const EigenBase<OtherDerived> &other) |
ykuroda | 0:13a5d365ba16 | 108 | { |
ykuroda | 0:13a5d365ba16 | 109 | other.derived().evalTo(derived()); |
ykuroda | 0:13a5d365ba16 | 110 | return derived(); |
ykuroda | 0:13a5d365ba16 | 111 | } |
ykuroda | 0:13a5d365ba16 | 112 | |
ykuroda | 0:13a5d365ba16 | 113 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 114 | template<typename OtherDerived> |
ykuroda | 0:13a5d365ba16 | 115 | Derived& DenseBase<Derived>::operator+=(const EigenBase<OtherDerived> &other) |
ykuroda | 0:13a5d365ba16 | 116 | { |
ykuroda | 0:13a5d365ba16 | 117 | other.derived().addTo(derived()); |
ykuroda | 0:13a5d365ba16 | 118 | return derived(); |
ykuroda | 0:13a5d365ba16 | 119 | } |
ykuroda | 0:13a5d365ba16 | 120 | |
ykuroda | 0:13a5d365ba16 | 121 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 122 | template<typename OtherDerived> |
ykuroda | 0:13a5d365ba16 | 123 | Derived& DenseBase<Derived>::operator-=(const EigenBase<OtherDerived> &other) |
ykuroda | 0:13a5d365ba16 | 124 | { |
ykuroda | 0:13a5d365ba16 | 125 | other.derived().subTo(derived()); |
ykuroda | 0:13a5d365ba16 | 126 | return derived(); |
ykuroda | 0:13a5d365ba16 | 127 | } |
ykuroda | 0:13a5d365ba16 | 128 | |
ykuroda | 0:13a5d365ba16 | 129 | } // end namespace Eigen |
ykuroda | 0:13a5d365ba16 | 130 | |
ykuroda | 0:13a5d365ba16 | 131 | #endif // EIGEN_EIGENBASE_H |