AngleAxis.h 8.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247
  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_ANGLEAXIS_H
  10. #define EIGEN_ANGLEAXIS_H
  11. namespace Eigen {
  12. /** \geometry_module \ingroup Geometry_Module
  13. *
  14. * \class AngleAxis
  15. *
  16. * \brief Represents a 3D rotation as a rotation angle around an arbitrary 3D axis
  17. *
  18. * \param _Scalar the scalar type, i.e., the type of the coefficients.
  19. *
  20. * \warning When setting up an AngleAxis object, the axis vector \b must \b be \b normalized.
  21. *
  22. * The following two typedefs are provided for convenience:
  23. * \li \c AngleAxisf for \c float
  24. * \li \c AngleAxisd for \c double
  25. *
  26. * Combined with MatrixBase::Unit{X,Y,Z}, AngleAxis can be used to easily
  27. * mimic Euler-angles. Here is an example:
  28. * \include AngleAxis_mimic_euler.cpp
  29. * Output: \verbinclude AngleAxis_mimic_euler.out
  30. *
  31. * \note This class is not aimed to be used to store a rotation transformation,
  32. * but rather to make easier the creation of other rotation (Quaternion, rotation Matrix)
  33. * and transformation objects.
  34. *
  35. * \sa class Quaternion, class Transform, MatrixBase::UnitX()
  36. */
  37. namespace internal {
  38. template<typename _Scalar> struct traits<AngleAxis<_Scalar> >
  39. {
  40. typedef _Scalar Scalar;
  41. };
  42. }
  43. template<typename _Scalar>
  44. class AngleAxis : public RotationBase<AngleAxis<_Scalar>,3>
  45. {
  46. typedef RotationBase<AngleAxis<_Scalar>,3> Base;
  47. public:
  48. using Base::operator*;
  49. enum { Dim = 3 };
  50. /** the scalar type of the coefficients */
  51. typedef _Scalar Scalar;
  52. typedef Matrix<Scalar,3,3> Matrix3;
  53. typedef Matrix<Scalar,3,1> Vector3;
  54. typedef Quaternion<Scalar> QuaternionType;
  55. protected:
  56. Vector3 m_axis;
  57. Scalar m_angle;
  58. public:
  59. /** Default constructor without initialization. */
  60. EIGEN_DEVICE_FUNC AngleAxis() {}
  61. /** Constructs and initialize the angle-axis rotation from an \a angle in radian
  62. * and an \a axis which \b must \b be \b normalized.
  63. *
  64. * \warning If the \a axis vector is not normalized, then the angle-axis object
  65. * represents an invalid rotation. */
  66. template<typename Derived>
  67. EIGEN_DEVICE_FUNC
  68. inline AngleAxis(const Scalar& angle, const MatrixBase<Derived>& axis) : m_axis(axis), m_angle(angle) {}
  69. /** Constructs and initialize the angle-axis rotation from a quaternion \a q.
  70. * This function implicitly normalizes the quaternion \a q.
  71. */
  72. template<typename QuatDerived>
  73. EIGEN_DEVICE_FUNC inline explicit AngleAxis(const QuaternionBase<QuatDerived>& q) { *this = q; }
  74. /** Constructs and initialize the angle-axis rotation from a 3x3 rotation matrix. */
  75. template<typename Derived>
  76. EIGEN_DEVICE_FUNC inline explicit AngleAxis(const MatrixBase<Derived>& m) { *this = m; }
  77. /** \returns the value of the rotation angle in radian */
  78. EIGEN_DEVICE_FUNC Scalar angle() const { return m_angle; }
  79. /** \returns a read-write reference to the stored angle in radian */
  80. EIGEN_DEVICE_FUNC Scalar& angle() { return m_angle; }
  81. /** \returns the rotation axis */
  82. EIGEN_DEVICE_FUNC const Vector3& axis() const { return m_axis; }
  83. /** \returns a read-write reference to the stored rotation axis.
  84. *
  85. * \warning The rotation axis must remain a \b unit vector.
  86. */
  87. EIGEN_DEVICE_FUNC Vector3& axis() { return m_axis; }
  88. /** Concatenates two rotations */
  89. EIGEN_DEVICE_FUNC inline QuaternionType operator* (const AngleAxis& other) const
  90. { return QuaternionType(*this) * QuaternionType(other); }
  91. /** Concatenates two rotations */
  92. EIGEN_DEVICE_FUNC inline QuaternionType operator* (const QuaternionType& other) const
  93. { return QuaternionType(*this) * other; }
  94. /** Concatenates two rotations */
  95. friend EIGEN_DEVICE_FUNC inline QuaternionType operator* (const QuaternionType& a, const AngleAxis& b)
  96. { return a * QuaternionType(b); }
  97. /** \returns the inverse rotation, i.e., an angle-axis with opposite rotation angle */
  98. EIGEN_DEVICE_FUNC AngleAxis inverse() const
  99. { return AngleAxis(-m_angle, m_axis); }
  100. template<class QuatDerived>
  101. EIGEN_DEVICE_FUNC AngleAxis& operator=(const QuaternionBase<QuatDerived>& q);
  102. template<typename Derived>
  103. EIGEN_DEVICE_FUNC AngleAxis& operator=(const MatrixBase<Derived>& m);
  104. template<typename Derived>
  105. EIGEN_DEVICE_FUNC AngleAxis& fromRotationMatrix(const MatrixBase<Derived>& m);
  106. EIGEN_DEVICE_FUNC Matrix3 toRotationMatrix(void) const;
  107. /** \returns \c *this with scalar type casted to \a NewScalarType
  108. *
  109. * Note that if \a NewScalarType is equal to the current scalar type of \c *this
  110. * then this function smartly returns a const reference to \c *this.
  111. */
  112. template<typename NewScalarType>
  113. EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<AngleAxis,AngleAxis<NewScalarType> >::type cast() const
  114. { return typename internal::cast_return_type<AngleAxis,AngleAxis<NewScalarType> >::type(*this); }
  115. /** Copy constructor with scalar type conversion */
  116. template<typename OtherScalarType>
  117. EIGEN_DEVICE_FUNC inline explicit AngleAxis(const AngleAxis<OtherScalarType>& other)
  118. {
  119. m_axis = other.axis().template cast<Scalar>();
  120. m_angle = Scalar(other.angle());
  121. }
  122. EIGEN_DEVICE_FUNC static inline const AngleAxis Identity() { return AngleAxis(Scalar(0), Vector3::UnitX()); }
  123. /** \returns \c true if \c *this is approximately equal to \a other, within the precision
  124. * determined by \a prec.
  125. *
  126. * \sa MatrixBase::isApprox() */
  127. EIGEN_DEVICE_FUNC bool isApprox(const AngleAxis& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
  128. { return m_axis.isApprox(other.m_axis, prec) && internal::isApprox(m_angle,other.m_angle, prec); }
  129. };
  130. /** \ingroup Geometry_Module
  131. * single precision angle-axis type */
  132. typedef AngleAxis<float> AngleAxisf;
  133. /** \ingroup Geometry_Module
  134. * double precision angle-axis type */
  135. typedef AngleAxis<double> AngleAxisd;
  136. /** Set \c *this from a \b unit quaternion.
  137. *
  138. * The resulting axis is normalized, and the computed angle is in the [0,pi] range.
  139. *
  140. * This function implicitly normalizes the quaternion \a q.
  141. */
  142. template<typename Scalar>
  143. template<typename QuatDerived>
  144. EIGEN_DEVICE_FUNC AngleAxis<Scalar>& AngleAxis<Scalar>::operator=(const QuaternionBase<QuatDerived>& q)
  145. {
  146. EIGEN_USING_STD_MATH(atan2)
  147. EIGEN_USING_STD_MATH(abs)
  148. Scalar n = q.vec().norm();
  149. if(n<NumTraits<Scalar>::epsilon())
  150. n = q.vec().stableNorm();
  151. if (n != Scalar(0))
  152. {
  153. m_angle = Scalar(2)*atan2(n, abs(q.w()));
  154. if(q.w() < Scalar(0))
  155. n = -n;
  156. m_axis = q.vec() / n;
  157. }
  158. else
  159. {
  160. m_angle = Scalar(0);
  161. m_axis << Scalar(1), Scalar(0), Scalar(0);
  162. }
  163. return *this;
  164. }
  165. /** Set \c *this from a 3x3 rotation matrix \a mat.
  166. */
  167. template<typename Scalar>
  168. template<typename Derived>
  169. EIGEN_DEVICE_FUNC AngleAxis<Scalar>& AngleAxis<Scalar>::operator=(const MatrixBase<Derived>& mat)
  170. {
  171. // Since a direct conversion would not be really faster,
  172. // let's use the robust Quaternion implementation:
  173. return *this = QuaternionType(mat);
  174. }
  175. /**
  176. * \brief Sets \c *this from a 3x3 rotation matrix.
  177. **/
  178. template<typename Scalar>
  179. template<typename Derived>
  180. EIGEN_DEVICE_FUNC AngleAxis<Scalar>& AngleAxis<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
  181. {
  182. return *this = QuaternionType(mat);
  183. }
  184. /** Constructs and \returns an equivalent 3x3 rotation matrix.
  185. */
  186. template<typename Scalar>
  187. typename AngleAxis<Scalar>::Matrix3
  188. EIGEN_DEVICE_FUNC AngleAxis<Scalar>::toRotationMatrix(void) const
  189. {
  190. EIGEN_USING_STD_MATH(sin)
  191. EIGEN_USING_STD_MATH(cos)
  192. Matrix3 res;
  193. Vector3 sin_axis = sin(m_angle) * m_axis;
  194. Scalar c = cos(m_angle);
  195. Vector3 cos1_axis = (Scalar(1)-c) * m_axis;
  196. Scalar tmp;
  197. tmp = cos1_axis.x() * m_axis.y();
  198. res.coeffRef(0,1) = tmp - sin_axis.z();
  199. res.coeffRef(1,0) = tmp + sin_axis.z();
  200. tmp = cos1_axis.x() * m_axis.z();
  201. res.coeffRef(0,2) = tmp + sin_axis.y();
  202. res.coeffRef(2,0) = tmp - sin_axis.y();
  203. tmp = cos1_axis.y() * m_axis.z();
  204. res.coeffRef(1,2) = tmp - sin_axis.x();
  205. res.coeffRef(2,1) = tmp + sin_axis.x();
  206. res.diagonal() = (cos1_axis.cwiseProduct(m_axis)).array() + c;
  207. return res;
  208. }
  209. } // end namespace Eigen
  210. #endif // EIGEN_ANGLEAXIS_H