Rotation2D.h 6.72 KB
Newer Older
LM's avatar
LM committed
1 2 3 4 5
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
//
Don Gagne's avatar
Don Gagne committed
6 7 8
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
LM's avatar
LM committed
9 10 11 12

#ifndef EIGEN_ROTATION2D_H
#define EIGEN_ROTATION2D_H

Don Gagne's avatar
Don Gagne committed
13 14
namespace Eigen { 

LM's avatar
LM committed
15 16 17 18 19 20
/** \geometry_module \ingroup Geometry_Module
  *
  * \class Rotation2D
  *
  * \brief Represents a rotation/orientation in a 2 dimensional space.
  *
21
  * \tparam _Scalar the scalar type, i.e., the type of the coefficients
LM's avatar
LM committed
22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61
  *
  * This class is equivalent to a single scalar representing a counter clock wise rotation
  * as a single angle in radian. It provides some additional features such as the automatic
  * conversion from/to a 2x2 rotation matrix. Moreover this class aims to provide a similar
  * interface to Quaternion in order to facilitate the writing of generic algorithms
  * dealing with rotations.
  *
  * \sa class Quaternion, class Transform
  */

namespace internal {

template<typename _Scalar> struct traits<Rotation2D<_Scalar> >
{
  typedef _Scalar Scalar;
};
} // end namespace internal

template<typename _Scalar>
class Rotation2D : public RotationBase<Rotation2D<_Scalar>,2>
{
  typedef RotationBase<Rotation2D<_Scalar>,2> Base;

public:

  using Base::operator*;

  enum { Dim = 2 };
  /** the scalar type of the coefficients */
  typedef _Scalar Scalar;
  typedef Matrix<Scalar,2,1> Vector2;
  typedef Matrix<Scalar,2,2> Matrix2;

protected:

  Scalar m_angle;

public:

  /** Construct a 2D counter clock wise rotation from the angle \a a in radian. */
62
  EIGEN_DEVICE_FUNC explicit inline Rotation2D(const Scalar& a) : m_angle(a) {}
63 64
  
  /** Default constructor wihtout initialization. The represented rotation is undefined. */
65 66 67 68 69 70 71 72 73 74 75
  EIGEN_DEVICE_FUNC Rotation2D() {}

  /** Construct a 2D rotation from a 2x2 rotation matrix \a mat.
    *
    * \sa fromRotationMatrix()
    */
  template<typename Derived>
  EIGEN_DEVICE_FUNC explicit Rotation2D(const MatrixBase<Derived>& m)
  {
    fromRotationMatrix(m.derived());
  }
LM's avatar
LM committed
76 77

  /** \returns the rotation angle */
78
  EIGEN_DEVICE_FUNC inline Scalar angle() const { return m_angle; }
LM's avatar
LM committed
79 80

  /** \returns a read-write reference to the rotation angle */
81 82 83 84 85 86 87 88 89 90 91 92 93 94 95
  EIGEN_DEVICE_FUNC inline Scalar& angle() { return m_angle; }
  
  /** \returns the rotation angle in [0,2pi] */
  EIGEN_DEVICE_FUNC inline Scalar smallestPositiveAngle() const {
    Scalar tmp = numext::fmod(m_angle,Scalar(2*EIGEN_PI));
    return tmp<Scalar(0) ? tmp + Scalar(2*EIGEN_PI) : tmp;
  }
  
  /** \returns the rotation angle in [-pi,pi] */
  EIGEN_DEVICE_FUNC inline Scalar smallestAngle() const {
    Scalar tmp = numext::fmod(m_angle,Scalar(2*EIGEN_PI));
    if(tmp>Scalar(EIGEN_PI))       tmp -= Scalar(2*EIGEN_PI);
    else if(tmp<-Scalar(EIGEN_PI)) tmp += Scalar(2*EIGEN_PI);
    return tmp;
  }
LM's avatar
LM committed
96 97

  /** \returns the inverse rotation */
98
  EIGEN_DEVICE_FUNC inline Rotation2D inverse() const { return Rotation2D(-m_angle); }
LM's avatar
LM committed
99 100

  /** Concatenates two rotations */
101 102
  EIGEN_DEVICE_FUNC inline Rotation2D operator*(const Rotation2D& other) const
  { return Rotation2D(m_angle + other.m_angle); }
LM's avatar
LM committed
103 104

  /** Concatenates two rotations */
105
  EIGEN_DEVICE_FUNC inline Rotation2D& operator*=(const Rotation2D& other)
Don Gagne's avatar
Don Gagne committed
106
  { m_angle += other.m_angle; return *this; }
LM's avatar
LM committed
107 108

  /** Applies the rotation to a 2D vector */
109
  EIGEN_DEVICE_FUNC Vector2 operator* (const Vector2& vec) const
LM's avatar
LM committed
110
  { return toRotationMatrix() * vec; }
111
  
LM's avatar
LM committed
112
  template<typename Derived>
113 114 115 116 117 118 119 120 121 122 123 124 125
  EIGEN_DEVICE_FUNC Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
  EIGEN_DEVICE_FUNC Matrix2 toRotationMatrix() const;

  /** Set \c *this from a 2x2 rotation matrix \a mat.
    * In other words, this function extract the rotation angle from the rotation matrix.
    *
    * This method is an alias for fromRotationMatrix()
    *
    * \sa fromRotationMatrix()
    */
  template<typename Derived>
  EIGEN_DEVICE_FUNC Rotation2D& operator=(const  MatrixBase<Derived>& m)
  { return fromRotationMatrix(m.derived()); }
LM's avatar
LM committed
126 127 128 129

  /** \returns the spherical interpolation between \c *this and \a other using
    * parameter \a t. It is in fact equivalent to a linear interpolation.
    */
130 131 132 133 134
  EIGEN_DEVICE_FUNC inline Rotation2D slerp(const Scalar& t, const Rotation2D& other) const
  {
    Scalar dist = Rotation2D(other.m_angle-m_angle).smallestAngle();
    return Rotation2D(m_angle + dist*t);
  }
LM's avatar
LM committed
135 136 137 138 139 140 141

  /** \returns \c *this with scalar type casted to \a NewScalarType
    *
    * Note that if \a NewScalarType is equal to the current scalar type of \c *this
    * then this function smartly returns a const reference to \c *this.
    */
  template<typename NewScalarType>
142
  EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type cast() const
LM's avatar
LM committed
143 144 145 146
  { return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*this); }

  /** Copy constructor with scalar type conversion */
  template<typename OtherScalarType>
147
  EIGEN_DEVICE_FUNC inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
LM's avatar
LM committed
148 149 150 151
  {
    m_angle = Scalar(other.angle());
  }

152
  EIGEN_DEVICE_FUNC static inline Rotation2D Identity() { return Rotation2D(0); }
LM's avatar
LM committed
153 154 155 156 157

  /** \returns \c true if \c *this is approximately equal to \a other, within the precision
    * determined by \a prec.
    *
    * \sa MatrixBase::isApprox() */
158
  EIGEN_DEVICE_FUNC bool isApprox(const Rotation2D& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
LM's avatar
LM committed
159
  { return internal::isApprox(m_angle,other.m_angle, prec); }
160
  
LM's avatar
LM committed
161 162 163 164 165 166 167 168 169 170 171 172 173 174 175
};

/** \ingroup Geometry_Module
  * single precision 2D rotation type */
typedef Rotation2D<float> Rotation2Df;
/** \ingroup Geometry_Module
  * double precision 2D rotation type */
typedef Rotation2D<double> Rotation2Dd;

/** Set \c *this from a 2x2 rotation matrix \a mat.
  * In other words, this function extract the rotation angle
  * from the rotation matrix.
  */
template<typename Scalar>
template<typename Derived>
176
EIGEN_DEVICE_FUNC Rotation2D<Scalar>& Rotation2D<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
LM's avatar
LM committed
177
{
178
  EIGEN_USING_STD_MATH(atan2)
LM's avatar
LM committed
179
  EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
Don Gagne's avatar
Don Gagne committed
180
  m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
LM's avatar
LM committed
181 182 183 184 185 186 187
  return *this;
}

/** Constructs and \returns an equivalent 2x2 rotation matrix.
  */
template<typename Scalar>
typename Rotation2D<Scalar>::Matrix2
188
EIGEN_DEVICE_FUNC Rotation2D<Scalar>::toRotationMatrix(void) const
LM's avatar
LM committed
189
{
190 191
  EIGEN_USING_STD_MATH(sin)
  EIGEN_USING_STD_MATH(cos)
Don Gagne's avatar
Don Gagne committed
192 193
  Scalar sinA = sin(m_angle);
  Scalar cosA = cos(m_angle);
LM's avatar
LM committed
194 195 196
  return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
}

Don Gagne's avatar
Don Gagne committed
197 198
} // end namespace Eigen

LM's avatar
LM committed
199
#endif // EIGEN_ROTATION2D_H