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