Quaternion.h 32 KB

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