Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:57:02

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2015 Tal Hadad <tal_hd@hotmail.com>
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_EULERANGLESCLASS_H// TODO: Fix previous "EIGEN_EULERANGLES_H" definition?
0011 #define EIGEN_EULERANGLESCLASS_H
0012 
0013 namespace Eigen
0014 {
0015   /** \class EulerAngles
0016     *
0017     * \ingroup EulerAngles_Module
0018     *
0019     * \brief Represents a rotation in a 3 dimensional space as three Euler angles.
0020     *
0021     * Euler rotation is a set of three rotation of three angles over three fixed axes, defined by the EulerSystem given as a template parameter.
0022     * 
0023     * Here is how intrinsic Euler angles works:
0024     *  - first, rotate the axes system over the alpha axis in angle alpha
0025     *  - then, rotate the axes system over the beta axis(which was rotated in the first stage) in angle beta
0026     *  - then, rotate the axes system over the gamma axis(which was rotated in the two stages above) in angle gamma
0027     *
0028     * \note This class support only intrinsic Euler angles for simplicity,
0029     *  see EulerSystem how to easily overcome this for extrinsic systems.
0030     *
0031     * ### Rotation representation and conversions ###
0032     *
0033     * It has been proved(see Wikipedia link below) that every rotation can be represented
0034     *  by Euler angles, but there is no single representation (e.g. unlike rotation matrices).
0035     * Therefore, you can convert from Eigen rotation and to them
0036     *  (including rotation matrices, which is not called "rotations" by Eigen design).
0037     *
0038     * Euler angles usually used for:
0039     *  - convenient human representation of rotation, especially in interactive GUI.
0040     *  - gimbal systems and robotics
0041     *  - efficient encoding(i.e. 3 floats only) of rotation for network protocols.
0042     *
0043     * However, Euler angles are slow comparing to quaternion or matrices,
0044     *  because their unnatural math definition, although it's simple for human.
0045     * To overcome this, this class provide easy movement from the math friendly representation
0046     *  to the human friendly representation, and vise-versa.
0047     *
0048     * All the user need to do is a safe simple C++ type conversion,
0049     *  and this class take care for the math.
0050     * Additionally, some axes related computation is done in compile time.
0051     *
0052     * #### Euler angles ranges in conversions ####
0053     * Rotations representation as EulerAngles are not single (unlike matrices),
0054     *  and even have infinite EulerAngles representations.<BR>
0055     * For example, add or subtract 2*PI from either angle of EulerAngles
0056     *  and you'll get the same rotation.
0057     * This is the general reason for infinite representation,
0058     *  but it's not the only general reason for not having a single representation.
0059     *
0060     * When converting rotation to EulerAngles, this class convert it to specific ranges
0061     * When converting some rotation to EulerAngles, the rules for ranges are as follow:
0062     * - If the rotation we converting from is an EulerAngles
0063     *  (even when it represented as RotationBase explicitly), angles ranges are __undefined__.
0064     * - otherwise, alpha and gamma angles will be in the range [-PI, PI].<BR>
0065     *   As for Beta angle:
0066     *    - If the system is Tait-Bryan, the beta angle will be in the range [-PI/2, PI/2].
0067     *    - otherwise:
0068     *      - If the beta axis is positive, the beta angle will be in the range [0, PI]
0069     *      - If the beta axis is negative, the beta angle will be in the range [-PI, 0]
0070     *
0071     * \sa EulerAngles(const MatrixBase<Derived>&)
0072     * \sa EulerAngles(const RotationBase<Derived, 3>&)
0073     *
0074     * ### Convenient user typedefs ###
0075     *
0076     * Convenient typedefs for EulerAngles exist for float and double scalar,
0077     *  in a form of EulerAngles{A}{B}{C}{scalar},
0078     *  e.g. \ref EulerAnglesXYZd, \ref EulerAnglesZYZf.
0079     *
0080     * Only for positive axes{+x,+y,+z} Euler systems are have convenient typedef.
0081     * If you need negative axes{-x,-y,-z}, it is recommended to create you own typedef with
0082     *  a word that represent what you need.
0083     *
0084     * ### Example ###
0085     *
0086     * \include EulerAngles.cpp
0087     * Output: \verbinclude EulerAngles.out
0088     *
0089     * ### Additional reading ###
0090     *
0091     * If you're want to get more idea about how Euler system work in Eigen see EulerSystem.
0092     *
0093     * More information about Euler angles: https://en.wikipedia.org/wiki/Euler_angles
0094     *
0095     * \tparam _Scalar the scalar type, i.e. the type of the angles.
0096     *
0097     * \tparam _System the EulerSystem to use, which represents the axes of rotation.
0098     */
0099   template <typename _Scalar, class _System>
0100   class EulerAngles : public RotationBase<EulerAngles<_Scalar, _System>, 3>
0101   {
0102     public:
0103       typedef RotationBase<EulerAngles<_Scalar, _System>, 3> Base;
0104       
0105       /** the scalar type of the angles */
0106       typedef _Scalar Scalar;
0107       typedef typename NumTraits<Scalar>::Real RealScalar;
0108       
0109       /** the EulerSystem to use, which represents the axes of rotation. */
0110       typedef _System System;
0111     
0112       typedef Matrix<Scalar,3,3> Matrix3; /*!< the equivalent rotation matrix type */
0113       typedef Matrix<Scalar,3,1> Vector3; /*!< the equivalent 3 dimension vector type */
0114       typedef Quaternion<Scalar> QuaternionType; /*!< the equivalent quaternion type */
0115       typedef AngleAxis<Scalar> AngleAxisType; /*!< the equivalent angle-axis type */
0116       
0117       /** \returns the axis vector of the first (alpha) rotation */
0118       static Vector3 AlphaAxisVector() {
0119         const Vector3& u = Vector3::Unit(System::AlphaAxisAbs - 1);
0120         return System::IsAlphaOpposite ? -u : u;
0121       }
0122       
0123       /** \returns the axis vector of the second (beta) rotation */
0124       static Vector3 BetaAxisVector() {
0125         const Vector3& u = Vector3::Unit(System::BetaAxisAbs - 1);
0126         return System::IsBetaOpposite ? -u : u;
0127       }
0128       
0129       /** \returns the axis vector of the third (gamma) rotation */
0130       static Vector3 GammaAxisVector() {
0131         const Vector3& u = Vector3::Unit(System::GammaAxisAbs - 1);
0132         return System::IsGammaOpposite ? -u : u;
0133       }
0134 
0135     private:
0136       Vector3 m_angles;
0137 
0138     public:
0139       /** Default constructor without initialization. */
0140       EulerAngles() {}
0141       /** Constructs and initialize an EulerAngles (\p alpha, \p beta, \p gamma). */
0142       EulerAngles(const Scalar& alpha, const Scalar& beta, const Scalar& gamma) :
0143         m_angles(alpha, beta, gamma) {}
0144       
0145       // TODO: Test this constructor
0146       /** Constructs and initialize an EulerAngles from the array data {alpha, beta, gamma} */
0147       explicit EulerAngles(const Scalar* data) : m_angles(data) {}
0148       
0149       /** Constructs and initializes an EulerAngles from either:
0150         *  - a 3x3 rotation matrix expression(i.e. pure orthogonal matrix with determinant of +1),
0151         *  - a 3D vector expression representing Euler angles.
0152         *
0153         * \note If \p other is a 3x3 rotation matrix, the angles range rules will be as follow:<BR>
0154         *  Alpha and gamma angles will be in the range [-PI, PI].<BR>
0155         *  As for Beta angle:
0156         *   - If the system is Tait-Bryan, the beta angle will be in the range [-PI/2, PI/2].
0157         *   - otherwise:
0158         *     - If the beta axis is positive, the beta angle will be in the range [0, PI]
0159         *     - If the beta axis is negative, the beta angle will be in the range [-PI, 0]
0160        */
0161       template<typename Derived>
0162       explicit EulerAngles(const MatrixBase<Derived>& other) { *this = other; }
0163       
0164       /** Constructs and initialize Euler angles from a rotation \p rot.
0165         *
0166         * \note If \p rot is an EulerAngles (even when it represented as RotationBase explicitly),
0167         *  angles ranges are __undefined__.
0168         *  Otherwise, alpha and gamma angles will be in the range [-PI, PI].<BR>
0169         *  As for Beta angle:
0170         *   - If the system is Tait-Bryan, the beta angle will be in the range [-PI/2, PI/2].
0171         *   - otherwise:
0172         *     - If the beta axis is positive, the beta angle will be in the range [0, PI]
0173         *     - If the beta axis is negative, the beta angle will be in the range [-PI, 0]
0174       */
0175       template<typename Derived>
0176       EulerAngles(const RotationBase<Derived, 3>& rot) { System::CalcEulerAngles(*this, rot.toRotationMatrix()); }
0177       
0178       /*EulerAngles(const QuaternionType& q)
0179       {
0180         // TODO: Implement it in a faster way for quaternions
0181         // According to http://www.euclideanspace.com/maths/geometry/rotations/conversions/quaternionToEuler/
0182         //  we can compute only the needed matrix cells and then convert to euler angles. (see ZYX example below)
0183         // Currently we compute all matrix cells from quaternion.
0184 
0185         // Special case only for ZYX
0186         //Scalar y2 = q.y() * q.y();
0187         //m_angles[0] = std::atan2(2*(q.w()*q.z() + q.x()*q.y()), (1 - 2*(y2 + q.z()*q.z())));
0188         //m_angles[1] = std::asin( 2*(q.w()*q.y() - q.z()*q.x()));
0189         //m_angles[2] = std::atan2(2*(q.w()*q.x() + q.y()*q.z()), (1 - 2*(q.x()*q.x() + y2)));
0190       }*/
0191 
0192       /** \returns The angle values stored in a vector (alpha, beta, gamma). */
0193       const Vector3& angles() const { return m_angles; }
0194       /** \returns A read-write reference to the angle values stored in a vector (alpha, beta, gamma). */
0195       Vector3& angles() { return m_angles; }
0196 
0197       /** \returns The value of the first angle. */
0198       Scalar alpha() const { return m_angles[0]; }
0199       /** \returns A read-write reference to the angle of the first angle. */
0200       Scalar& alpha() { return m_angles[0]; }
0201 
0202       /** \returns The value of the second angle. */
0203       Scalar beta() const { return m_angles[1]; }
0204       /** \returns A read-write reference to the angle of the second angle. */
0205       Scalar& beta() { return m_angles[1]; }
0206 
0207       /** \returns The value of the third angle. */
0208       Scalar gamma() const { return m_angles[2]; }
0209       /** \returns A read-write reference to the angle of the third angle. */
0210       Scalar& gamma() { return m_angles[2]; }
0211 
0212       /** \returns The Euler angles rotation inverse (which is as same as the negative),
0213         *  (-alpha, -beta, -gamma).
0214       */
0215       EulerAngles inverse() const
0216       {
0217         EulerAngles res;
0218         res.m_angles = -m_angles;
0219         return res;
0220       }
0221 
0222       /** \returns The Euler angles rotation negative (which is as same as the inverse),
0223         *  (-alpha, -beta, -gamma).
0224       */
0225       EulerAngles operator -() const
0226       {
0227         return inverse();
0228       }
0229       
0230       /** Set \c *this from either:
0231         *  - a 3x3 rotation matrix expression(i.e. pure orthogonal matrix with determinant of +1),
0232         *  - a 3D vector expression representing Euler angles.
0233         *
0234         * See EulerAngles(const MatrixBase<Derived, 3>&) for more information about
0235         *  angles ranges output.
0236       */
0237       template<class Derived>
0238       EulerAngles& operator=(const MatrixBase<Derived>& other)
0239       {
0240         EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename Derived::Scalar>::value),
0241          YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
0242         
0243         internal::eulerangles_assign_impl<System, Derived>::run(*this, other.derived());
0244         return *this;
0245       }
0246 
0247       // TODO: Assign and construct from another EulerAngles (with different system)
0248       
0249       /** Set \c *this from a rotation.
0250         *
0251         * See EulerAngles(const RotationBase<Derived, 3>&) for more information about
0252         *  angles ranges output.
0253       */
0254       template<typename Derived>
0255       EulerAngles& operator=(const RotationBase<Derived, 3>& rot) {
0256         System::CalcEulerAngles(*this, rot.toRotationMatrix());
0257         return *this;
0258       }
0259       
0260       /** \returns \c true if \c *this is approximately equal to \a other, within the precision
0261         * determined by \a prec.
0262         *
0263         * \sa MatrixBase::isApprox() */
0264       bool isApprox(const EulerAngles& other,
0265         const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
0266       { return angles().isApprox(other.angles(), prec); }
0267 
0268       /** \returns an equivalent 3x3 rotation matrix. */
0269       Matrix3 toRotationMatrix() const
0270       {
0271         // TODO: Calc it faster
0272         return static_cast<QuaternionType>(*this).toRotationMatrix();
0273       }
0274 
0275       /** Convert the Euler angles to quaternion. */
0276       operator QuaternionType() const
0277       {
0278         return
0279           AngleAxisType(alpha(), AlphaAxisVector()) *
0280           AngleAxisType(beta(), BetaAxisVector())   *
0281           AngleAxisType(gamma(), GammaAxisVector());
0282       }
0283       
0284       friend std::ostream& operator<<(std::ostream& s, const EulerAngles<Scalar, System>& eulerAngles)
0285       {
0286         s << eulerAngles.angles().transpose();
0287         return s;
0288       }
0289       
0290       /** \returns \c *this with scalar type casted to \a NewScalarType */
0291       template <typename NewScalarType>
0292       EulerAngles<NewScalarType, System> cast() const
0293       {
0294         EulerAngles<NewScalarType, System> e;
0295         e.angles() = angles().template cast<NewScalarType>();
0296         return e;
0297       }
0298   };
0299 
0300 #define EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(AXES, SCALAR_TYPE, SCALAR_POSTFIX) \
0301   /** \ingroup EulerAngles_Module */ \
0302   typedef EulerAngles<SCALAR_TYPE, EulerSystem##AXES> EulerAngles##AXES##SCALAR_POSTFIX;
0303 
0304 #define EIGEN_EULER_ANGLES_TYPEDEFS(SCALAR_TYPE, SCALAR_POSTFIX) \
0305   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XYZ, SCALAR_TYPE, SCALAR_POSTFIX) \
0306   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XYX, SCALAR_TYPE, SCALAR_POSTFIX) \
0307   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XZY, SCALAR_TYPE, SCALAR_POSTFIX) \
0308   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XZX, SCALAR_TYPE, SCALAR_POSTFIX) \
0309  \
0310   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YZX, SCALAR_TYPE, SCALAR_POSTFIX) \
0311   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YZY, SCALAR_TYPE, SCALAR_POSTFIX) \
0312   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YXZ, SCALAR_TYPE, SCALAR_POSTFIX) \
0313   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YXY, SCALAR_TYPE, SCALAR_POSTFIX) \
0314  \
0315   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZXY, SCALAR_TYPE, SCALAR_POSTFIX) \
0316   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZXZ, SCALAR_TYPE, SCALAR_POSTFIX) \
0317   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZYX, SCALAR_TYPE, SCALAR_POSTFIX) \
0318   EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZYZ, SCALAR_TYPE, SCALAR_POSTFIX)
0319 
0320 EIGEN_EULER_ANGLES_TYPEDEFS(float, f)
0321 EIGEN_EULER_ANGLES_TYPEDEFS(double, d)
0322 
0323   namespace internal
0324   {
0325     template<typename _Scalar, class _System>
0326     struct traits<EulerAngles<_Scalar, _System> >
0327     {
0328       typedef _Scalar Scalar;
0329     };
0330     
0331     // set from a rotation matrix
0332     template<class System, class Other>
0333     struct eulerangles_assign_impl<System,Other,3,3>
0334     {
0335       typedef typename Other::Scalar Scalar;
0336       static void run(EulerAngles<Scalar, System>& e, const Other& m)
0337       {
0338         System::CalcEulerAngles(e, m);
0339       }
0340     };
0341     
0342     // set from a vector of Euler angles
0343     template<class System, class Other>
0344     struct eulerangles_assign_impl<System,Other,3,1>
0345     {
0346       typedef typename Other::Scalar Scalar;
0347       static void run(EulerAngles<Scalar, System>& e, const Other& vec)
0348       {
0349         e.angles() = vec;
0350       }
0351     };
0352   }
0353 }
0354 
0355 #endif // EIGEN_EULERANGLESCLASS_H