Rotation2D.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_ROTATION2D_H
00011 #define EIGEN_ROTATION2D_H
00012 
00013 namespace Eigen { 
00014 
00032 namespace internal {
00033 
00034 template<typename _Scalar> struct traits<Rotation2D<_Scalar> >
00035 {
00036   typedef _Scalar Scalar;
00037 };
00038 } // end namespace internal
00039 
00040 template<typename _Scalar>
00041 class Rotation2D : public RotationBase<Rotation2D<_Scalar>,2>
00042 {
00043   typedef RotationBase<Rotation2D<_Scalar>,2> Base;
00044 
00045 public:
00046 
00047   using Base::operator*;
00048 
00049   enum { Dim = 2 };
00051   typedef _Scalar Scalar;
00052   typedef Matrix<Scalar,2,1> Vector2;
00053   typedef Matrix<Scalar,2,2> Matrix2;
00054 
00055 protected:
00056 
00057   Scalar m_angle;
00058 
00059 public:
00060 
00062   inline Rotation2D(const Scalar& a) : m_angle(a) {}
00063 
00065   inline Scalar angle() const { return m_angle; }
00066 
00068   inline Scalar& angle() { return m_angle; }
00069 
00071   inline Rotation2D inverse() const { return -m_angle; }
00072 
00074   inline Rotation2D operator*(const Rotation2D& other) const
00075   { return m_angle + other.m_angle; }
00076 
00078   inline Rotation2D& operator*=(const Rotation2D& other)
00079   { m_angle += other.m_angle; return *this; }
00080 
00082   Vector2 operator* (const Vector2& vec) const
00083   { return toRotationMatrix() * vec; }
00084 
00085   template<typename Derived>
00086   Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
00087   Matrix2 toRotationMatrix(void) const;
00088 
00092   inline Rotation2D slerp(const Scalar& t, const Rotation2D& other) const
00093   { return m_angle * (1-t) + other.angle() * t; }
00094 
00100   template<typename NewScalarType>
00101   inline typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type cast() const
00102   { return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*this); }
00103 
00105   template<typename OtherScalarType>
00106   inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
00107   {
00108     m_angle = Scalar(other.angle());
00109   }
00110 
00111   static inline Rotation2D Identity() { return Rotation2D(0); }
00112 
00117   bool isApprox(const Rotation2D& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
00118   { return internal::isApprox(m_angle,other.m_angle, prec); }
00119 };
00120 
00123 typedef Rotation2D<float> Rotation2Df;
00126 typedef Rotation2D<double> Rotation2Dd;
00127 
00132 template<typename Scalar>
00133 template<typename Derived>
00134 Rotation2D<Scalar>& Rotation2D<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
00135 {
00136   using std::atan2;
00137   EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
00138   m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
00139   return *this;
00140 }
00141 
00144 template<typename Scalar>
00145 typename Rotation2D<Scalar>::Matrix2
00146 Rotation2D<Scalar>::toRotationMatrix(void) const
00147 {
00148   using std::sin;
00149   using std::cos;
00150   Scalar sinA = sin(m_angle);
00151   Scalar cosA = cos(m_angle);
00152   return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
00153 }
00154 
00155 } // end namespace Eigen
00156 
00157 #endif // EIGEN_ROTATION2D_H


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Sat Jun 8 2019 19:38:49