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@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
13namespace 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  * \tparam _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
32namespace internal {
33
34template<typename _Scalar> struct traits<Rotation2D<_Scalar> >
35{
36  typedef _Scalar Scalar;
37};
38} // end namespace internal
39
40template<typename _Scalar>
41class Rotation2D : public RotationBase<Rotation2D<_Scalar>,2>
42{
43  typedef RotationBase<Rotation2D<_Scalar>,2> Base;
44
45public:
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
55protected:
56
57  Scalar m_angle;
58
59public:
60
61  /** 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  EIGEN_DEVICE_FUNC Rotation2D() {}
66
67  /** Construct a 2D rotation from a 2x2 rotation matrix \a mat.
68    *
69    * \sa fromRotationMatrix()
70    */
71  template<typename Derived>
72  EIGEN_DEVICE_FUNC explicit Rotation2D(const MatrixBase<Derived>& m)
73  {
74    fromRotationMatrix(m.derived());
75  }
76
77  /** \returns the rotation angle */
78  EIGEN_DEVICE_FUNC inline Scalar angle() const { return m_angle; }
79
80  /** \returns a read-write reference to the rotation angle */
81  EIGEN_DEVICE_FUNC inline Scalar& angle() { return m_angle; }
82
83  /** \returns the rotation angle in [0,2pi] */
84  EIGEN_DEVICE_FUNC inline Scalar smallestPositiveAngle() const {
85    Scalar tmp = numext::fmod(m_angle,Scalar(2*EIGEN_PI));
86    return tmp<Scalar(0) ? tmp + Scalar(2*EIGEN_PI) : tmp;
87  }
88
89  /** \returns the rotation angle in [-pi,pi] */
90  EIGEN_DEVICE_FUNC inline Scalar smallestAngle() const {
91    Scalar tmp = numext::fmod(m_angle,Scalar(2*EIGEN_PI));
92    if(tmp>Scalar(EIGEN_PI))       tmp -= Scalar(2*EIGEN_PI);
93    else if(tmp<-Scalar(EIGEN_PI)) tmp += Scalar(2*EIGEN_PI);
94    return tmp;
95  }
96
97  /** \returns the inverse rotation */
98  EIGEN_DEVICE_FUNC inline Rotation2D inverse() const { return Rotation2D(-m_angle); }
99
100  /** Concatenates two rotations */
101  EIGEN_DEVICE_FUNC inline Rotation2D operator*(const Rotation2D& other) const
102  { return Rotation2D(m_angle + other.m_angle); }
103
104  /** Concatenates two rotations */
105  EIGEN_DEVICE_FUNC inline Rotation2D& operator*=(const Rotation2D& other)
106  { m_angle += other.m_angle; return *this; }
107
108  /** Applies the rotation to a 2D vector */
109  EIGEN_DEVICE_FUNC Vector2 operator* (const Vector2& vec) const
110  { return toRotationMatrix() * vec; }
111
112  template<typename Derived>
113  EIGEN_DEVICE_FUNC Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
114  EIGEN_DEVICE_FUNC Matrix2 toRotationMatrix() const;
115
116  /** Set \c *this from a 2x2 rotation matrix \a mat.
117    * In other words, this function extract the rotation angle from the rotation matrix.
118    *
119    * This method is an alias for fromRotationMatrix()
120    *
121    * \sa fromRotationMatrix()
122    */
123  template<typename Derived>
124  EIGEN_DEVICE_FUNC Rotation2D& operator=(const  MatrixBase<Derived>& m)
125  { return fromRotationMatrix(m.derived()); }
126
127  /** \returns the spherical interpolation between \c *this and \a other using
128    * parameter \a t. It is in fact equivalent to a linear interpolation.
129    */
130  EIGEN_DEVICE_FUNC inline Rotation2D slerp(const Scalar& t, const Rotation2D& other) const
131  {
132    Scalar dist = Rotation2D(other.m_angle-m_angle).smallestAngle();
133    return Rotation2D(m_angle + dist*t);
134  }
135
136  /** \returns \c *this with scalar type casted to \a NewScalarType
137    *
138    * Note that if \a NewScalarType is equal to the current scalar type of \c *this
139    * then this function smartly returns a const reference to \c *this.
140    */
141  template<typename NewScalarType>
142  EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type cast() const
143  { return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*this); }
144
145  /** Copy constructor with scalar type conversion */
146  template<typename OtherScalarType>
147  EIGEN_DEVICE_FUNC inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
148  {
149    m_angle = Scalar(other.angle());
150  }
151
152  EIGEN_DEVICE_FUNC static inline Rotation2D Identity() { return Rotation2D(0); }
153
154  /** \returns \c true if \c *this is approximately equal to \a other, within the precision
155    * determined by \a prec.
156    *
157    * \sa MatrixBase::isApprox() */
158  EIGEN_DEVICE_FUNC bool isApprox(const Rotation2D& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
159  { return internal::isApprox(m_angle,other.m_angle, prec); }
160
161};
162
163/** \ingroup Geometry_Module
164  * single precision 2D rotation type */
165typedef Rotation2D<float> Rotation2Df;
166/** \ingroup Geometry_Module
167  * double precision 2D rotation type */
168typedef Rotation2D<double> Rotation2Dd;
169
170/** Set \c *this from a 2x2 rotation matrix \a mat.
171  * In other words, this function extract the rotation angle
172  * from the rotation matrix.
173  */
174template<typename Scalar>
175template<typename Derived>
176EIGEN_DEVICE_FUNC Rotation2D<Scalar>& Rotation2D<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
177{
178  EIGEN_USING_STD_MATH(atan2)
179  EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
180  m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
181  return *this;
182}
183
184/** Constructs and \returns an equivalent 2x2 rotation matrix.
185  */
186template<typename Scalar>
187typename Rotation2D<Scalar>::Matrix2
188EIGEN_DEVICE_FUNC Rotation2D<Scalar>::toRotationMatrix(void) const
189{
190  EIGEN_USING_STD_MATH(sin)
191  EIGEN_USING_STD_MATH(cos)
192  Scalar sinA = sin(m_angle);
193  Scalar cosA = cos(m_angle);
194  return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
195}
196
197} // end namespace Eigen
198
199#endif // EIGEN_ROTATION2D_H
200