Home | History | Annotate | Download | only in Geometry
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 //
      6 // This Source Code Form is subject to the terms of the Mozilla
      7 // Public License v. 2.0. If a copy of the MPL was not distributed
      8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
      9 
     10 #ifndef EIGEN_ROTATION2D_H
     11 #define EIGEN_ROTATION2D_H
     12 
     13 namespace Eigen {
     14 
     15 /** \geometry_module \ingroup Geometry_Module
     16   *
     17   * \class Rotation2D
     18   *
     19   * \brief Represents a rotation/orientation in a 2 dimensional space.
     20   *
     21   * \param _Scalar the scalar type, i.e., the type of the coefficients
     22   *
     23   * This class is equivalent to a single scalar representing a counter clock wise rotation
     24   * as a single angle in radian. It provides some additional features such as the automatic
     25   * conversion from/to a 2x2 rotation matrix. Moreover this class aims to provide a similar
     26   * interface to Quaternion in order to facilitate the writing of generic algorithms
     27   * dealing with rotations.
     28   *
     29   * \sa class Quaternion, class Transform
     30   */
     31 
     32 namespace internal {
     33 
     34 template<typename _Scalar> struct traits<Rotation2D<_Scalar> >
     35 {
     36   typedef _Scalar Scalar;
     37 };
     38 } // end namespace internal
     39 
     40 template<typename _Scalar>
     41 class Rotation2D : public RotationBase<Rotation2D<_Scalar>,2>
     42 {
     43   typedef RotationBase<Rotation2D<_Scalar>,2> Base;
     44 
     45 public:
     46 
     47   using Base::operator*;
     48 
     49   enum { Dim = 2 };
     50   /** the scalar type of the coefficients */
     51   typedef _Scalar Scalar;
     52   typedef Matrix<Scalar,2,1> Vector2;
     53   typedef Matrix<Scalar,2,2> Matrix2;
     54 
     55 protected:
     56 
     57   Scalar m_angle;
     58 
     59 public:
     60 
     61   /** Construct a 2D counter clock wise rotation from the angle \a a in radian. */
     62   inline Rotation2D(const Scalar& a) : m_angle(a) {}
     63 
     64   /** Default constructor wihtout initialization. The represented rotation is undefined. */
     65   Rotation2D() {}
     66 
     67   /** \returns the rotation angle */
     68   inline Scalar angle() const { return m_angle; }
     69 
     70   /** \returns a read-write reference to the rotation angle */
     71   inline Scalar& angle() { return m_angle; }
     72 
     73   /** \returns the inverse rotation */
     74   inline Rotation2D inverse() const { return -m_angle; }
     75 
     76   /** Concatenates two rotations */
     77   inline Rotation2D operator*(const Rotation2D& other) const
     78   { return m_angle + other.m_angle; }
     79 
     80   /** Concatenates two rotations */
     81   inline Rotation2D& operator*=(const Rotation2D& other)
     82   { m_angle += other.m_angle; return *this; }
     83 
     84   /** Applies the rotation to a 2D vector */
     85   Vector2 operator* (const Vector2& vec) const
     86   { return toRotationMatrix() * vec; }
     87 
     88   template<typename Derived>
     89   Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
     90   Matrix2 toRotationMatrix() const;
     91 
     92   /** \returns the spherical interpolation between \c *this and \a other using
     93     * parameter \a t. It is in fact equivalent to a linear interpolation.
     94     */
     95   inline Rotation2D slerp(const Scalar& t, const Rotation2D& other) const
     96   { return m_angle * (1-t) + other.angle() * t; }
     97 
     98   /** \returns \c *this with scalar type casted to \a NewScalarType
     99     *
    100     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
    101     * then this function smartly returns a const reference to \c *this.
    102     */
    103   template<typename NewScalarType>
    104   inline typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type cast() const
    105   { return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*this); }
    106 
    107   /** Copy constructor with scalar type conversion */
    108   template<typename OtherScalarType>
    109   inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
    110   {
    111     m_angle = Scalar(other.angle());
    112   }
    113 
    114   static inline Rotation2D Identity() { return Rotation2D(0); }
    115 
    116   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
    117     * determined by \a prec.
    118     *
    119     * \sa MatrixBase::isApprox() */
    120   bool isApprox(const Rotation2D& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
    121   { return internal::isApprox(m_angle,other.m_angle, prec); }
    122 };
    123 
    124 /** \ingroup Geometry_Module
    125   * single precision 2D rotation type */
    126 typedef Rotation2D<float> Rotation2Df;
    127 /** \ingroup Geometry_Module
    128   * double precision 2D rotation type */
    129 typedef Rotation2D<double> Rotation2Dd;
    130 
    131 /** Set \c *this from a 2x2 rotation matrix \a mat.
    132   * In other words, this function extract the rotation angle
    133   * from the rotation matrix.
    134   */
    135 template<typename Scalar>
    136 template<typename Derived>
    137 Rotation2D<Scalar>& Rotation2D<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
    138 {
    139   using std::atan2;
    140   EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
    141   m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
    142   return *this;
    143 }
    144 
    145 /** Constructs and \returns an equivalent 2x2 rotation matrix.
    146   */
    147 template<typename Scalar>
    148 typename Rotation2D<Scalar>::Matrix2
    149 Rotation2D<Scalar>::toRotationMatrix(void) const
    150 {
    151   using std::sin;
    152   using std::cos;
    153   Scalar sinA = sin(m_angle);
    154   Scalar cosA = cos(m_angle);
    155   return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
    156 }
    157 
    158 } // end namespace Eigen
    159 
    160 #endif // EIGEN_ROTATION2D_H
    161