You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

778 lines
28 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2009 Mathieu Gautier <mathieu.gautier@cea.fr>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #ifndef EIGEN_QUATERNION_H
  11. #define EIGEN_QUATERNION_H
  12. namespace Eigen {
  13. /***************************************************************************
  14. * Definition of QuaternionBase<Derived>
  15. * The implementation is at the end of the file
  16. ***************************************************************************/
  17. namespace internal {
  18. template<typename Other,
  19. int OtherRows=Other::RowsAtCompileTime,
  20. int OtherCols=Other::ColsAtCompileTime>
  21. struct quaternionbase_assign_impl;
  22. }
  23. /** \geometry_module \ingroup Geometry_Module
  24. * \class QuaternionBase
  25. * \brief Base class for quaternion expressions
  26. * \tparam Derived derived type (CRTP)
  27. * \sa class Quaternion
  28. */
  29. template<class Derived>
  30. class QuaternionBase : public RotationBase<Derived, 3>
  31. {
  32. typedef RotationBase<Derived, 3> Base;
  33. public:
  34. using Base::operator*;
  35. using Base::derived;
  36. typedef typename internal::traits<Derived>::Scalar Scalar;
  37. typedef typename NumTraits<Scalar>::Real RealScalar;
  38. typedef typename internal::traits<Derived>::Coefficients Coefficients;
  39. enum {
  40. Flags = Eigen::internal::traits<Derived>::Flags
  41. };
  42. // typedef typename Matrix<Scalar,4,1> Coefficients;
  43. /** the type of a 3D vector */
  44. typedef Matrix<Scalar,3,1> Vector3;
  45. /** the equivalent rotation matrix type */
  46. typedef Matrix<Scalar,3,3> Matrix3;
  47. /** the equivalent angle-axis type */
  48. typedef AngleAxis<Scalar> AngleAxisType;
  49. /** \returns the \c x coefficient */
  50. inline Scalar x() const { return this->derived().coeffs().coeff(0); }
  51. /** \returns the \c y coefficient */
  52. inline Scalar y() const { return this->derived().coeffs().coeff(1); }
  53. /** \returns the \c z coefficient */
  54. inline Scalar z() const { return this->derived().coeffs().coeff(2); }
  55. /** \returns the \c w coefficient */
  56. inline Scalar w() const { return this->derived().coeffs().coeff(3); }
  57. /** \returns a reference to the \c x coefficient */
  58. inline Scalar& x() { return this->derived().coeffs().coeffRef(0); }
  59. /** \returns a reference to the \c y coefficient */
  60. inline Scalar& y() { return this->derived().coeffs().coeffRef(1); }
  61. /** \returns a reference to the \c z coefficient */
  62. inline Scalar& z() { return this->derived().coeffs().coeffRef(2); }
  63. /** \returns a reference to the \c w coefficient */
  64. inline Scalar& w() { return this->derived().coeffs().coeffRef(3); }
  65. /** \returns a read-only vector expression of the imaginary part (x,y,z) */
  66. inline const VectorBlock<const Coefficients,3> vec() const { return coeffs().template head<3>(); }
  67. /** \returns a vector expression of the imaginary part (x,y,z) */
  68. inline VectorBlock<Coefficients,3> vec() { return coeffs().template head<3>(); }
  69. /** \returns a read-only vector expression of the coefficients (x,y,z,w) */
  70. inline const typename internal::traits<Derived>::Coefficients& coeffs() const { return derived().coeffs(); }
  71. /** \returns a vector expression of the coefficients (x,y,z,w) */
  72. inline typename internal::traits<Derived>::Coefficients& coeffs() { return derived().coeffs(); }
  73. EIGEN_STRONG_INLINE QuaternionBase<Derived>& operator=(const QuaternionBase<Derived>& other);
  74. template<class OtherDerived> EIGEN_STRONG_INLINE Derived& operator=(const QuaternionBase<OtherDerived>& other);
  75. // disabled this copy operator as it is giving very strange compilation errors when compiling
  76. // test_stdvector with GCC 4.4.2. This looks like a GCC bug though, so feel free to re-enable it if it's
  77. // useful; however notice that we already have the templated operator= above and e.g. in MatrixBase
  78. // we didn't have to add, in addition to templated operator=, such a non-templated copy operator.
  79. // Derived& operator=(const QuaternionBase& other)
  80. // { return operator=<Derived>(other); }
  81. Derived& operator=(const AngleAxisType& aa);
  82. template<class OtherDerived> Derived& operator=(const MatrixBase<OtherDerived>& m);
  83. /** \returns a quaternion representing an identity rotation
  84. * \sa MatrixBase::Identity()
  85. */
  86. static inline Quaternion<Scalar> Identity() { return Quaternion<Scalar>(1, 0, 0, 0); }
  87. /** \sa QuaternionBase::Identity(), MatrixBase::setIdentity()
  88. */
  89. inline QuaternionBase& setIdentity() { coeffs() << 0, 0, 0, 1; return *this; }
  90. /** \returns the squared norm of the quaternion's coefficients
  91. * \sa QuaternionBase::norm(), MatrixBase::squaredNorm()
  92. */
  93. inline Scalar squaredNorm() const { return coeffs().squaredNorm(); }
  94. /** \returns the norm of the quaternion's coefficients
  95. * \sa QuaternionBase::squaredNorm(), MatrixBase::norm()
  96. */
  97. inline Scalar norm() const { return coeffs().norm(); }
  98. /** Normalizes the quaternion \c *this
  99. * \sa normalized(), MatrixBase::normalize() */
  100. inline void normalize() { coeffs().normalize(); }
  101. /** \returns a normalized copy of \c *this
  102. * \sa normalize(), MatrixBase::normalized() */
  103. inline Quaternion<Scalar> normalized() const { return Quaternion<Scalar>(coeffs().normalized()); }
  104. /** \returns the dot product of \c *this and \a other
  105. * Geometrically speaking, the dot product of two unit quaternions
  106. * corresponds to the cosine of half the angle between the two rotations.
  107. * \sa angularDistance()
  108. */
  109. template<class OtherDerived> inline Scalar dot(const QuaternionBase<OtherDerived>& other) const { return coeffs().dot(other.coeffs()); }
  110. template<class OtherDerived> Scalar angularDistance(const QuaternionBase<OtherDerived>& other) const;
  111. /** \returns an equivalent 3x3 rotation matrix */
  112. Matrix3 toRotationMatrix() const;
  113. /** \returns the quaternion which transform \a a into \a b through a rotation */
  114. template<typename Derived1, typename Derived2>
  115. Derived& setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
  116. template<class OtherDerived> EIGEN_STRONG_INLINE Quaternion<Scalar> operator* (const QuaternionBase<OtherDerived>& q) const;
  117. template<class OtherDerived> EIGEN_STRONG_INLINE Derived& operator*= (const QuaternionBase<OtherDerived>& q);
  118. /** \returns the quaternion describing the inverse rotation */
  119. Quaternion<Scalar> inverse() const;
  120. /** \returns the conjugated quaternion */
  121. Quaternion<Scalar> conjugate() const;
  122. /** \returns an interpolation for a constant motion between \a other and \c *this
  123. * \a t in [0;1]
  124. * see http://en.wikipedia.org/wiki/Slerp
  125. */
  126. template<class OtherDerived> Quaternion<Scalar> slerp(Scalar t, const QuaternionBase<OtherDerived>& other) const;
  127. /** \returns \c true if \c *this is approximately equal to \a other, within the precision
  128. * determined by \a prec.
  129. *
  130. * \sa MatrixBase::isApprox() */
  131. template<class OtherDerived>
  132. bool isApprox(const QuaternionBase<OtherDerived>& other, RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
  133. { return coeffs().isApprox(other.coeffs(), prec); }
  134. /** return the result vector of \a v through the rotation*/
  135. EIGEN_STRONG_INLINE Vector3 _transformVector(Vector3 v) const;
  136. /** \returns \c *this with scalar type casted to \a NewScalarType
  137. *
  138. * Note that if \a NewScalarType is equal to the current scalar type of \c *this
  139. * then this function smartly returns a const reference to \c *this.
  140. */
  141. template<typename NewScalarType>
  142. inline typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type cast() const
  143. {
  144. return typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type(derived());
  145. }
  146. #ifdef EIGEN_QUATERNIONBASE_PLUGIN
  147. # include EIGEN_QUATERNIONBASE_PLUGIN
  148. #endif
  149. };
  150. /***************************************************************************
  151. * Definition/implementation of Quaternion<Scalar>
  152. ***************************************************************************/
  153. /** \geometry_module \ingroup Geometry_Module
  154. *
  155. * \class Quaternion
  156. *
  157. * \brief The quaternion class used to represent 3D orientations and rotations
  158. *
  159. * \param _Scalar the scalar type, i.e., the type of the coefficients
  160. *
  161. * This class represents a quaternion \f$ w+xi+yj+zk \f$ that is a convenient representation of
  162. * orientations and rotations of objects in three dimensions. Compared to other representations
  163. * like Euler angles or 3x3 matrices, quatertions offer the following advantages:
  164. * \li \b compact storage (4 scalars)
  165. * \li \b efficient to compose (28 flops),
  166. * \li \b stable spherical interpolation
  167. *
  168. * The following two typedefs are provided for convenience:
  169. * \li \c Quaternionf for \c float
  170. * \li \c Quaterniond for \c double
  171. *
  172. * \sa class AngleAxis, class Transform
  173. */
  174. namespace internal {
  175. template<typename _Scalar,int _Options>
  176. struct traits<Quaternion<_Scalar,_Options> >
  177. {
  178. typedef Quaternion<_Scalar,_Options> PlainObject;
  179. typedef _Scalar Scalar;
  180. typedef Matrix<_Scalar,4,1,_Options> Coefficients;
  181. enum{
  182. IsAligned = internal::traits<Coefficients>::Flags & AlignedBit,
  183. Flags = IsAligned ? (AlignedBit | LvalueBit) : LvalueBit
  184. };
  185. };
  186. }
  187. template<typename _Scalar, int _Options>
  188. class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
  189. {
  190. typedef QuaternionBase<Quaternion<_Scalar,_Options> > Base;
  191. enum { IsAligned = internal::traits<Quaternion>::IsAligned };
  192. public:
  193. typedef _Scalar Scalar;
  194. EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Quaternion)
  195. using Base::operator*=;
  196. typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
  197. typedef typename Base::AngleAxisType AngleAxisType;
  198. /** Default constructor leaving the quaternion uninitialized. */
  199. inline Quaternion() {}
  200. /** Constructs and initializes the quaternion \f$ w+xi+yj+zk \f$ from
  201. * its four coefficients \a w, \a x, \a y and \a z.
  202. *
  203. * \warning Note the order of the arguments: the real \a w coefficient first,
  204. * while internally the coefficients are stored in the following order:
  205. * [\c x, \c y, \c z, \c w]
  206. */
  207. inline Quaternion(Scalar w, Scalar x, Scalar y, Scalar z) : m_coeffs(x, y, z, w){}
  208. /** Constructs and initialize a quaternion from the array data */
  209. inline Quaternion(const Scalar* data) : m_coeffs(data) {}
  210. /** Copy constructor */
  211. template<class Derived> EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) { this->Base::operator=(other); }
  212. /** Constructs and initializes a quaternion from the angle-axis \a aa */
  213. explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
  214. /** Constructs and initializes a quaternion from either:
  215. * - a rotation matrix expression,
  216. * - a 4D vector expression representing quaternion coefficients.
  217. */
  218. template<typename Derived>
  219. explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; }
  220. /** Explicit copy constructor with scalar conversion */
  221. template<typename OtherScalar, int OtherOptions>
  222. explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other)
  223. { m_coeffs = other.coeffs().template cast<Scalar>(); }
  224. template<typename Derived1, typename Derived2>
  225. static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
  226. inline Coefficients& coeffs() { return m_coeffs;}
  227. inline const Coefficients& coeffs() const { return m_coeffs;}
  228. EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(IsAligned)
  229. protected:
  230. Coefficients m_coeffs;
  231. #ifndef EIGEN_PARSED_BY_DOXYGEN
  232. static EIGEN_STRONG_INLINE void _check_template_params()
  233. {
  234. EIGEN_STATIC_ASSERT( (_Options & DontAlign) == _Options,
  235. INVALID_MATRIX_TEMPLATE_PARAMETERS)
  236. }
  237. #endif
  238. };
  239. /** \ingroup Geometry_Module
  240. * single precision quaternion type */
  241. typedef Quaternion<float> Quaternionf;
  242. /** \ingroup Geometry_Module
  243. * double precision quaternion type */
  244. typedef Quaternion<double> Quaterniond;
  245. /***************************************************************************
  246. * Specialization of Map<Quaternion<Scalar>>
  247. ***************************************************************************/
  248. namespace internal {
  249. template<typename _Scalar, int _Options>
  250. struct traits<Map<Quaternion<_Scalar>, _Options> >:
  251. traits<Quaternion<_Scalar, _Options> >
  252. {
  253. typedef _Scalar Scalar;
  254. typedef Map<Matrix<_Scalar,4,1>, _Options> Coefficients;
  255. typedef traits<Quaternion<_Scalar, _Options> > TraitsBase;
  256. enum {
  257. IsAligned = TraitsBase::IsAligned,
  258. Flags = TraitsBase::Flags
  259. };
  260. };
  261. }
  262. namespace internal {
  263. template<typename _Scalar, int _Options>
  264. struct traits<Map<const Quaternion<_Scalar>, _Options> >:
  265. traits<Quaternion<_Scalar> >
  266. {
  267. typedef _Scalar Scalar;
  268. typedef Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;
  269. typedef traits<Quaternion<_Scalar, _Options> > TraitsBase;
  270. enum {
  271. IsAligned = TraitsBase::IsAligned,
  272. Flags = TraitsBase::Flags & ~LvalueBit
  273. };
  274. };
  275. }
  276. /** \brief Quaternion expression mapping a constant memory buffer
  277. *
  278. * \param _Scalar the type of the Quaternion coefficients
  279. * \param _Options see class Map
  280. *
  281. * This is a specialization of class Map for Quaternion. This class allows to view
  282. * a 4 scalar memory buffer as an Eigen's Quaternion object.
  283. *
  284. * \sa class Map, class Quaternion, class QuaternionBase
  285. */
  286. template<typename _Scalar, int _Options>
  287. class Map<const Quaternion<_Scalar>, _Options >
  288. : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
  289. {
  290. typedef QuaternionBase<Map<const Quaternion<_Scalar>, _Options> > Base;
  291. public:
  292. typedef _Scalar Scalar;
  293. typedef typename internal::traits<Map>::Coefficients Coefficients;
  294. EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Map)
  295. using Base::operator*=;
  296. /** Constructs a Mapped Quaternion object from the pointer \a coeffs
  297. *
  298. * The pointer \a coeffs must reference the four coeffecients of Quaternion in the following order:
  299. * \code *coeffs == {x, y, z, w} \endcode
  300. *
  301. * If the template parameter _Options is set to #Aligned, then the pointer coeffs must be aligned. */
  302. EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
  303. inline const Coefficients& coeffs() const { return m_coeffs;}
  304. protected:
  305. const Coefficients m_coeffs;
  306. };
  307. /** \brief Expression of a quaternion from a memory buffer
  308. *
  309. * \param _Scalar the type of the Quaternion coefficients
  310. * \param _Options see class Map
  311. *
  312. * This is a specialization of class Map for Quaternion. This class allows to view
  313. * a 4 scalar memory buffer as an Eigen's Quaternion object.
  314. *
  315. * \sa class Map, class Quaternion, class QuaternionBase
  316. */
  317. template<typename _Scalar, int _Options>
  318. class Map<Quaternion<_Scalar>, _Options >
  319. : public QuaternionBase<Map<Quaternion<_Scalar>, _Options> >
  320. {
  321. typedef QuaternionBase<Map<Quaternion<_Scalar>, _Options> > Base;
  322. public:
  323. typedef _Scalar Scalar;
  324. typedef typename internal::traits<Map>::Coefficients Coefficients;
  325. EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Map)
  326. using Base::operator*=;
  327. /** Constructs a Mapped Quaternion object from the pointer \a coeffs
  328. *
  329. * The pointer \a coeffs must reference the four coeffecients of Quaternion in the following order:
  330. * \code *coeffs == {x, y, z, w} \endcode
  331. *
  332. * If the template parameter _Options is set to #Aligned, then the pointer coeffs must be aligned. */
  333. EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
  334. inline Coefficients& coeffs() { return m_coeffs; }
  335. inline const Coefficients& coeffs() const { return m_coeffs; }
  336. protected:
  337. Coefficients m_coeffs;
  338. };
  339. /** \ingroup Geometry_Module
  340. * Map an unaligned array of single precision scalar as a quaternion */
  341. typedef Map<Quaternion<float>, 0> QuaternionMapf;
  342. /** \ingroup Geometry_Module
  343. * Map an unaligned array of double precision scalar as a quaternion */
  344. typedef Map<Quaternion<double>, 0> QuaternionMapd;
  345. /** \ingroup Geometry_Module
  346. * Map a 16-bits aligned array of double precision scalars as a quaternion */
  347. typedef Map<Quaternion<float>, Aligned> QuaternionMapAlignedf;
  348. /** \ingroup Geometry_Module
  349. * Map a 16-bits aligned array of double precision scalars as a quaternion */
  350. typedef Map<Quaternion<double>, Aligned> QuaternionMapAlignedd;
  351. /***************************************************************************
  352. * Implementation of QuaternionBase methods
  353. ***************************************************************************/
  354. // Generic Quaternion * Quaternion product
  355. // This product can be specialized for a given architecture via the Arch template argument.
  356. namespace internal {
  357. template<int Arch, class Derived1, class Derived2, typename Scalar, int _Options> struct quat_product
  358. {
  359. static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a, const QuaternionBase<Derived2>& b){
  360. return Quaternion<Scalar>
  361. (
  362. a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
  363. a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
  364. a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
  365. a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()
  366. );
  367. }
  368. };
  369. }
  370. /** \returns the concatenation of two rotations as a quaternion-quaternion product */
  371. template <class Derived>
  372. template <class OtherDerived>
  373. EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar>
  374. QuaternionBase<Derived>::operator* (const QuaternionBase<OtherDerived>& other) const
  375. {
  376. EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
  377. YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
  378. return internal::quat_product<Architecture::Target, Derived, OtherDerived,
  379. typename internal::traits<Derived>::Scalar,
  380. internal::traits<Derived>::IsAligned && internal::traits<OtherDerived>::IsAligned>::run(*this, other);
  381. }
  382. /** \sa operator*(Quaternion) */
  383. template <class Derived>
  384. template <class OtherDerived>
  385. EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*= (const QuaternionBase<OtherDerived>& other)
  386. {
  387. derived() = derived() * other.derived();
  388. return derived();
  389. }
  390. /** Rotation of a vector by a quaternion.
  391. * \remarks If the quaternion is used to rotate several points (>1)
  392. * then it is much more efficient to first convert it to a 3x3 Matrix.
  393. * Comparison of the operation cost for n transformations:
  394. * - Quaternion2: 30n
  395. * - Via a Matrix3: 24 + 15n
  396. */
  397. template <class Derived>
  398. EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
  399. QuaternionBase<Derived>::_transformVector(Vector3 v) const
  400. {
  401. // Note that this algorithm comes from the optimization by hand
  402. // of the conversion to a Matrix followed by a Matrix/Vector product.
  403. // It appears to be much faster than the common algorithm found
  404. // in the litterature (30 versus 39 flops). It also requires two
  405. // Vector3 as temporaries.
  406. Vector3 uv = this->vec().cross(v);
  407. uv += uv;
  408. return v + this->w() * uv + this->vec().cross(uv);
  409. }
  410. template<class Derived>
  411. EIGEN_STRONG_INLINE QuaternionBase<Derived>& QuaternionBase<Derived>::operator=(const QuaternionBase<Derived>& other)
  412. {
  413. coeffs() = other.coeffs();
  414. return derived();
  415. }
  416. template<class Derived>
  417. template<class OtherDerived>
  418. EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
  419. {
  420. coeffs() = other.coeffs();
  421. return derived();
  422. }
  423. /** Set \c *this from an angle-axis \a aa and returns a reference to \c *this
  424. */
  425. template<class Derived>
  426. EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa)
  427. {
  428. Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
  429. this->w() = internal::cos(ha);
  430. this->vec() = internal::sin(ha) * aa.axis();
  431. return derived();
  432. }
  433. /** Set \c *this from the expression \a xpr:
  434. * - if \a xpr is a 4x1 vector, then \a xpr is assumed to be a quaternion
  435. * - if \a xpr is a 3x3 matrix, then \a xpr is assumed to be rotation matrix
  436. * and \a xpr is converted to a quaternion
  437. */
  438. template<class Derived>
  439. template<class MatrixDerived>
  440. inline Derived& QuaternionBase<Derived>::operator=(const MatrixBase<MatrixDerived>& xpr)
  441. {
  442. EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
  443. YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
  444. internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
  445. return derived();
  446. }
  447. /** Convert the quaternion to a 3x3 rotation matrix. The quaternion is required to
  448. * be normalized, otherwise the result is undefined.
  449. */
  450. template<class Derived>
  451. inline typename QuaternionBase<Derived>::Matrix3
  452. QuaternionBase<Derived>::toRotationMatrix(void) const
  453. {
  454. // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
  455. // if not inlined then the cost of the return by value is huge ~ +35%,
  456. // however, not inlining this function is an order of magnitude slower, so
  457. // it has to be inlined, and so the return by value is not an issue
  458. Matrix3 res;
  459. const Scalar tx = Scalar(2)*this->x();
  460. const Scalar ty = Scalar(2)*this->y();
  461. const Scalar tz = Scalar(2)*this->z();
  462. const Scalar twx = tx*this->w();
  463. const Scalar twy = ty*this->w();
  464. const Scalar twz = tz*this->w();
  465. const Scalar txx = tx*this->x();
  466. const Scalar txy = ty*this->x();
  467. const Scalar txz = tz*this->x();
  468. const Scalar tyy = ty*this->y();
  469. const Scalar tyz = tz*this->y();
  470. const Scalar tzz = tz*this->z();
  471. res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
  472. res.coeffRef(0,1) = txy-twz;
  473. res.coeffRef(0,2) = txz+twy;
  474. res.coeffRef(1,0) = txy+twz;
  475. res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
  476. res.coeffRef(1,2) = tyz-twx;
  477. res.coeffRef(2,0) = txz-twy;
  478. res.coeffRef(2,1) = tyz+twx;
  479. res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
  480. return res;
  481. }
  482. /** Sets \c *this to be a quaternion representing a rotation between
  483. * the two arbitrary vectors \a a and \a b. In other words, the built
  484. * rotation represent a rotation sending the line of direction \a a
  485. * to the line of direction \a b, both lines passing through the origin.
  486. *
  487. * \returns a reference to \c *this.
  488. *
  489. * Note that the two input vectors do \b not have to be normalized, and
  490. * do not need to have the same norm.
  491. */
  492. template<class Derived>
  493. template<typename Derived1, typename Derived2>
  494. inline Derived& QuaternionBase<Derived>::setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b)
  495. {
  496. using std::max;
  497. Vector3 v0 = a.normalized();
  498. Vector3 v1 = b.normalized();
  499. Scalar c = v1.dot(v0);
  500. // if dot == -1, vectors are nearly opposites
  501. // => accuraletly compute the rotation axis by computing the
  502. // intersection of the two planes. This is done by solving:
  503. // x^T v0 = 0
  504. // x^T v1 = 0
  505. // under the constraint:
  506. // ||x|| = 1
  507. // which yields a singular value problem
  508. if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
  509. {
  510. c = max<Scalar>(c,-1);
  511. Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
  512. JacobiSVD<Matrix<Scalar,2,3> > svd(m, ComputeFullV);
  513. Vector3 axis = svd.matrixV().col(2);
  514. Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
  515. this->w() = internal::sqrt(w2);
  516. this->vec() = axis * internal::sqrt(Scalar(1) - w2);
  517. return derived();
  518. }
  519. Vector3 axis = v0.cross(v1);
  520. Scalar s = internal::sqrt((Scalar(1)+c)*Scalar(2));
  521. Scalar invs = Scalar(1)/s;
  522. this->vec() = axis * invs;
  523. this->w() = s * Scalar(0.5);
  524. return derived();
  525. }
  526. /** Returns a quaternion representing a rotation between
  527. * the two arbitrary vectors \a a and \a b. In other words, the built
  528. * rotation represent a rotation sending the line of direction \a a
  529. * to the line of direction \a b, both lines passing through the origin.
  530. *
  531. * \returns resulting quaternion
  532. *
  533. * Note that the two input vectors do \b not have to be normalized, and
  534. * do not need to have the same norm.
  535. */
  536. template<typename Scalar, int Options>
  537. template<typename Derived1, typename Derived2>
  538. Quaternion<Scalar,Options> Quaternion<Scalar,Options>::FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b)
  539. {
  540. Quaternion quat;
  541. quat.setFromTwoVectors(a, b);
  542. return quat;
  543. }
  544. /** \returns the multiplicative inverse of \c *this
  545. * Note that in most cases, i.e., if you simply want the opposite rotation,
  546. * and/or the quaternion is normalized, then it is enough to use the conjugate.
  547. *
  548. * \sa QuaternionBase::conjugate()
  549. */
  550. template <class Derived>
  551. inline Quaternion<typename internal::traits<Derived>::Scalar> QuaternionBase<Derived>::inverse() const
  552. {
  553. // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
  554. Scalar n2 = this->squaredNorm();
  555. if (n2 > 0)
  556. return Quaternion<Scalar>(conjugate().coeffs() / n2);
  557. else
  558. {
  559. // return an invalid result to flag the error
  560. return Quaternion<Scalar>(Coefficients::Zero());
  561. }
  562. }
  563. /** \returns the conjugate of the \c *this which is equal to the multiplicative inverse
  564. * if the quaternion is normalized.
  565. * The conjugate of a quaternion represents the opposite rotation.
  566. *
  567. * \sa Quaternion2::inverse()
  568. */
  569. template <class Derived>
  570. inline Quaternion<typename internal::traits<Derived>::Scalar>
  571. QuaternionBase<Derived>::conjugate() const
  572. {
  573. return Quaternion<Scalar>(this->w(),-this->x(),-this->y(),-this->z());
  574. }
  575. /** \returns the angle (in radian) between two rotations
  576. * \sa dot()
  577. */
  578. template <class Derived>
  579. template <class OtherDerived>
  580. inline typename internal::traits<Derived>::Scalar
  581. QuaternionBase<Derived>::angularDistance(const QuaternionBase<OtherDerived>& other) const
  582. {
  583. using std::acos;
  584. double d = internal::abs(this->dot(other));
  585. if (d>=1.0)
  586. return Scalar(0);
  587. return static_cast<Scalar>(2 * acos(d));
  588. }
  589. /** \returns the spherical linear interpolation between the two quaternions
  590. * \c *this and \a other at the parameter \a t
  591. */
  592. template <class Derived>
  593. template <class OtherDerived>
  594. Quaternion<typename internal::traits<Derived>::Scalar>
  595. QuaternionBase<Derived>::slerp(Scalar t, const QuaternionBase<OtherDerived>& other) const
  596. {
  597. using std::acos;
  598. static const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
  599. Scalar d = this->dot(other);
  600. Scalar absD = internal::abs(d);
  601. Scalar scale0;
  602. Scalar scale1;
  603. if(absD>=one)
  604. {
  605. scale0 = Scalar(1) - t;
  606. scale1 = t;
  607. }
  608. else
  609. {
  610. // theta is the angle between the 2 quaternions
  611. Scalar theta = acos(absD);
  612. Scalar sinTheta = internal::sin(theta);
  613. scale0 = internal::sin( ( Scalar(1) - t ) * theta) / sinTheta;
  614. scale1 = internal::sin( ( t * theta) ) / sinTheta;
  615. }
  616. if(d<0) scale1 = -scale1;
  617. return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
  618. }
  619. namespace internal {
  620. // set from a rotation matrix
  621. template<typename Other>
  622. struct quaternionbase_assign_impl<Other,3,3>
  623. {
  624. typedef typename Other::Scalar Scalar;
  625. typedef DenseIndex Index;
  626. template<class Derived> static inline void run(QuaternionBase<Derived>& q, const Other& mat)
  627. {
  628. // This algorithm comes from "Quaternion Calculus and Fast Animation",
  629. // Ken Shoemake, 1987 SIGGRAPH course notes
  630. Scalar t = mat.trace();
  631. if (t > Scalar(0))
  632. {
  633. t = sqrt(t + Scalar(1.0));
  634. q.w() = Scalar(0.5)*t;
  635. t = Scalar(0.5)/t;
  636. q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t;
  637. q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t;
  638. q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t;
  639. }
  640. else
  641. {
  642. DenseIndex i = 0;
  643. if (mat.coeff(1,1) > mat.coeff(0,0))
  644. i = 1;
  645. if (mat.coeff(2,2) > mat.coeff(i,i))
  646. i = 2;
  647. DenseIndex j = (i+1)%3;
  648. DenseIndex k = (j+1)%3;
  649. t = sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0));
  650. q.coeffs().coeffRef(i) = Scalar(0.5) * t;
  651. t = Scalar(0.5)/t;
  652. q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t;
  653. q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t;
  654. q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t;
  655. }
  656. }
  657. };
  658. // set from a vector of coefficients assumed to be a quaternion
  659. template<typename Other>
  660. struct quaternionbase_assign_impl<Other,4,1>
  661. {
  662. typedef typename Other::Scalar Scalar;
  663. template<class Derived> static inline void run(QuaternionBase<Derived>& q, const Other& vec)
  664. {
  665. q.coeffs() = vec;
  666. }
  667. };
  668. } // end namespace internal
  669. } // end namespace Eigen
  670. #endif // EIGEN_QUATERNION_H