Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-19 09:06:18

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_ROTATION2D_H
0011 #define EIGEN_ROTATION2D_H
0012 
0013 namespace RivetEigen { 
0014 
0015 /** \geometry_module \ingroup Geometry_Module
0016   *
0017   * \class Rotation2D
0018   *
0019   * \brief Represents a rotation/orientation in a 2 dimensional space.
0020   *
0021   * \tparam _Scalar the scalar type, i.e., the type of the coefficients
0022   *
0023   * This class is equivalent to a single scalar representing a counter clock wise rotation
0024   * as a single angle in radian. It provides some additional features such as the automatic
0025   * conversion from/to a 2x2 rotation matrix. Moreover this class aims to provide a similar
0026   * interface to Quaternion in order to facilitate the writing of generic algorithms
0027   * dealing with rotations.
0028   *
0029   * \sa class Quaternion, class Transform
0030   */
0031 
0032 namespace internal {
0033 
0034 template<typename _Scalar> struct traits<Rotation2D<_Scalar> >
0035 {
0036   typedef _Scalar Scalar;
0037 };
0038 } // end namespace internal
0039 
0040 template<typename _Scalar>
0041 class Rotation2D : public RotationBase<Rotation2D<_Scalar>,2>
0042 {
0043   typedef RotationBase<Rotation2D<_Scalar>,2> Base;
0044 
0045 public:
0046 
0047   using Base::operator*;
0048 
0049   enum { Dim = 2 };
0050   /** the scalar type of the coefficients */
0051   typedef _Scalar Scalar;
0052   typedef Matrix<Scalar,2,1> Vector2;
0053   typedef Matrix<Scalar,2,2> Matrix2;
0054 
0055 protected:
0056 
0057   Scalar m_angle;
0058 
0059 public:
0060 
0061   /** Construct a 2D counter clock wise rotation from the angle \a a in radian. */
0062   EIGEN_DEVICE_FUNC explicit inline Rotation2D(const Scalar& a) : m_angle(a) {}
0063   
0064   /** Default constructor wihtout initialization. The represented rotation is undefined. */
0065   EIGEN_DEVICE_FUNC Rotation2D() {}
0066 
0067   /** Construct a 2D rotation from a 2x2 rotation matrix \a mat.
0068     *
0069     * \sa fromRotationMatrix()
0070     */
0071   template<typename Derived>
0072   EIGEN_DEVICE_FUNC explicit Rotation2D(const MatrixBase<Derived>& m)
0073   {
0074     fromRotationMatrix(m.derived());
0075   }
0076 
0077   /** \returns the rotation angle */
0078   EIGEN_DEVICE_FUNC inline Scalar angle() const { return m_angle; }
0079 
0080   /** \returns a read-write reference to the rotation angle */
0081   EIGEN_DEVICE_FUNC inline Scalar& angle() { return m_angle; }
0082   
0083   /** \returns the rotation angle in [0,2pi] */
0084   EIGEN_DEVICE_FUNC inline Scalar smallestPositiveAngle() const {
0085     Scalar tmp = numext::fmod(m_angle,Scalar(2*EIGEN_PI));
0086     return tmp<Scalar(0) ? tmp + Scalar(2*EIGEN_PI) : tmp;
0087   }
0088   
0089   /** \returns the rotation angle in [-pi,pi] */
0090   EIGEN_DEVICE_FUNC inline Scalar smallestAngle() const {
0091     Scalar tmp = numext::fmod(m_angle,Scalar(2*EIGEN_PI));
0092     if(tmp>Scalar(EIGEN_PI))       tmp -= Scalar(2*EIGEN_PI);
0093     else if(tmp<-Scalar(EIGEN_PI)) tmp += Scalar(2*EIGEN_PI);
0094     return tmp;
0095   }
0096 
0097   /** \returns the inverse rotation */
0098   EIGEN_DEVICE_FUNC inline Rotation2D inverse() const { return Rotation2D(-m_angle); }
0099 
0100   /** Concatenates two rotations */
0101   EIGEN_DEVICE_FUNC inline Rotation2D operator*(const Rotation2D& other) const
0102   { return Rotation2D(m_angle + other.m_angle); }
0103 
0104   /** Concatenates two rotations */
0105   EIGEN_DEVICE_FUNC inline Rotation2D& operator*=(const Rotation2D& other)
0106   { m_angle += other.m_angle; return *this; }
0107 
0108   /** Applies the rotation to a 2D vector */
0109   EIGEN_DEVICE_FUNC Vector2 operator* (const Vector2& vec) const
0110   { return toRotationMatrix() * vec; }
0111   
0112   template<typename Derived>
0113   EIGEN_DEVICE_FUNC Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
0114   EIGEN_DEVICE_FUNC Matrix2 toRotationMatrix() const;
0115 
0116   /** Set \c *this from a 2x2 rotation matrix \a mat.
0117     * In other words, this function extract the rotation angle from the rotation matrix.
0118     *
0119     * This method is an alias for fromRotationMatrix()
0120     *
0121     * \sa fromRotationMatrix()
0122     */
0123   template<typename Derived>
0124   EIGEN_DEVICE_FUNC Rotation2D& operator=(const  MatrixBase<Derived>& m)
0125   { return fromRotationMatrix(m.derived()); }
0126 
0127   /** \returns the spherical interpolation between \c *this and \a other using
0128     * parameter \a t. It is in fact equivalent to a linear interpolation.
0129     */
0130   EIGEN_DEVICE_FUNC inline Rotation2D slerp(const Scalar& t, const Rotation2D& other) const
0131   {
0132     Scalar dist = Rotation2D(other.m_angle-m_angle).smallestAngle();
0133     return Rotation2D(m_angle + dist*t);
0134   }
0135 
0136   /** \returns \c *this with scalar type casted to \a NewScalarType
0137     *
0138     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
0139     * then this function smartly returns a const reference to \c *this.
0140     */
0141   template<typename NewScalarType>
0142   EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type cast() const
0143   { return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*this); }
0144 
0145   /** Copy constructor with scalar type conversion */
0146   template<typename OtherScalarType>
0147   EIGEN_DEVICE_FUNC inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
0148   {
0149     m_angle = Scalar(other.angle());
0150   }
0151 
0152   EIGEN_DEVICE_FUNC static inline Rotation2D Identity() { return Rotation2D(0); }
0153 
0154   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
0155     * determined by \a prec.
0156     *
0157     * \sa MatrixBase::isApprox() */
0158   EIGEN_DEVICE_FUNC bool isApprox(const Rotation2D& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
0159   { return internal::isApprox(m_angle,other.m_angle, prec); }
0160   
0161 };
0162 
0163 /** \ingroup Geometry_Module
0164   * single precision 2D rotation type */
0165 typedef Rotation2D<float> Rotation2Df;
0166 /** \ingroup Geometry_Module
0167   * double precision 2D rotation type */
0168 typedef Rotation2D<double> Rotation2Dd;
0169 
0170 /** Set \c *this from a 2x2 rotation matrix \a mat.
0171   * In other words, this function extract the rotation angle
0172   * from the rotation matrix.
0173   */
0174 template<typename Scalar>
0175 template<typename Derived>
0176 EIGEN_DEVICE_FUNC Rotation2D<Scalar>& Rotation2D<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
0177 {
0178   EIGEN_USING_STD(atan2)
0179   EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
0180   m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
0181   return *this;
0182 }
0183 
0184 /** Constructs and \returns an equivalent 2x2 rotation matrix.
0185   */
0186 template<typename Scalar>
0187 typename Rotation2D<Scalar>::Matrix2
0188 EIGEN_DEVICE_FUNC Rotation2D<Scalar>::toRotationMatrix(void) const
0189 {
0190   EIGEN_USING_STD(sin)
0191   EIGEN_USING_STD(cos)
0192   Scalar sinA = sin(m_angle);
0193   Scalar cosA = cos(m_angle);
0194   return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
0195 }
0196 
0197 } // end namespace RivetEigen
0198 
0199 #endif // EIGEN_ROTATION2D_H