Rotation2D.h 6.7 KB

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