Eigen  3.4.90 (git rev a4098ac676528a83cfb73d4d26ce1b42ec05f47c)
Eigen::Quaternion< Scalar_, Options_ > Class Template Reference

Detailed Description

template<typename Scalar_, int Options_>
class Eigen::Quaternion< Scalar_, Options_ >

The quaternion class used to represent 3D orientations and rotations.

This is defined in the Geometry module.

#include <Eigen/Geometry>
Template Parameters
Scalar_the scalar type, i.e., the type of the coefficients
Options_controls the memory alignment of the coefficients. Can be # AutoAlign or # DontAlign. Default is AutoAlign.

This class represents a quaternion \( w+xi+yj+zk \) that is a convenient representation of orientations and rotations of objects in three dimensions. Compared to other representations like Euler angles or 3x3 matrices, quaternions offer the following advantages:

  • compact storage (4 scalars)
  • efficient to compose (28 flops),
  • stable spherical interpolation

The following two typedefs are provided for convenience:

  • Quaternionf for float
  • Quaterniond for double
Warning
Operations interpreting the quaternion as rotation have undefined behavior if the quaternion is not normalized.
See also
class AngleAxis, class Transform
+ Inheritance diagram for Eigen::Quaternion< Scalar_, Options_ >:

Public Member Functions

template<typename Derived1 , typename Derived2 >
Quaternion< Scalar, Options > FromTwoVectors (const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
 
Quaternionoperator= (Quaternion &&other) EIGEN_NOEXCEPT_IF(std
 
 Quaternion ()
 
 Quaternion (const AngleAxisType &aa)
 
template<typename Derived >
 Quaternion (const MatrixBase< Derived > &other)
 
template<typename OtherScalar , int OtherOptions>
 Quaternion (const Quaternion< OtherScalar, OtherOptions > &other)
 
template<class Derived >
 Quaternion (const QuaternionBase< Derived > &other)
 
 Quaternion (const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)
 
 Quaternion (const Scalar *data)
 
 Quaternion (Quaternion &&other) EIGEN_NOEXCEPT_IF(std
 
- Public Member Functions inherited from Eigen::QuaternionBase< Quaternion< Scalar_, Options_ > >
Vector3 _transformVector (const Vector3 &v) const
 
internal::traits< Quaternion< Scalar_, Options_ > >::Scalar angularDistance (const QuaternionBase< OtherDerived > &other) const
 
internal::cast_return_type< Quaternion< Scalar_, Options_ >, Quaternion< NewScalarType > >::type cast () const
 
internal::traits< Quaternion< Scalar_, Options_ > >::Coefficients & coeffs ()
 
const internal::traits< Quaternion< Scalar_, Options_ > >::Coefficients & coeffs () const
 
Quaternion< Scalarconjugate () const
 
Scalar dot (const QuaternionBase< OtherDerived > &other) const
 
Quaternion< Scalarinverse () const
 
bool isApprox (const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
 
Scalar norm () const
 
void normalize ()
 
Quaternion< Scalarnormalized () const
 
bool operator!= (const QuaternionBase< OtherDerived > &other) const
 
Quaternion< typename internal::traits< Quaternion< Scalar_, Options_ > >::Scalaroperator* (const QuaternionBase< OtherDerived > &other) const
 
Quaternion< Scalar_, Options_ > & operator*= (const QuaternionBase< OtherDerived > &q)
 
Quaternion< Scalar_, Options_ > & operator= (const AngleAxisType &aa)
 
Quaternion< Scalar_, Options_ > & operator= (const MatrixBase< MatrixDerived > &xpr)
 
bool operator== (const QuaternionBase< OtherDerived > &other) const
 
Quaternion< Scalar_, Options_ > & setFromTwoVectors (const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
 
QuaternionBasesetIdentity ()
 
Quaternion< typename internal::traits< Quaternion< Scalar_, Options_ > >::Scalarslerp (const Scalar &t, const QuaternionBase< OtherDerived > &other) const
 
Scalar squaredNorm () const
 
Matrix3 toRotationMatrix () const
 
VectorBlock< Coefficients, 3 > vec ()
 
const VectorBlock< const Coefficients, 3 > vec () const
 
NonConstCoeffReturnType w ()
 
CoeffReturnType w () const
 
NonConstCoeffReturnType x ()
 
CoeffReturnType x () const
 
NonConstCoeffReturnType y ()
 
CoeffReturnType y () const
 
NonConstCoeffReturnType z ()
 
CoeffReturnType z () const
 
- Public Member Functions inherited from Eigen::RotationBase< Derived, Dim_ >
Derived inverse () const
 
RotationMatrixType matrix () const
 
template<typename OtherDerived >
internal::rotation_base_generic_product_selector< Derived, OtherDerived, OtherDerived::IsVectorAtCompileTime >::ReturnType operator* (const EigenBase< OtherDerived > &e) const
 
template<int Mode, int Options>
Transform< Scalar, Dim, Mode > operator* (const Transform< Scalar, Dim, Mode, Options > &t) const
 
Transform< Scalar, Dim, Isometryoperator* (const Translation< Scalar, Dim > &t) const
 
RotationMatrixType operator* (const UniformScaling< Scalar > &s) const
 
RotationMatrixType toRotationMatrix () const
 

Static Public Member Functions

static Quaternion UnitRandom ()
 
- Static Public Member Functions inherited from Eigen::QuaternionBase< Quaternion< Scalar_, Options_ > >
static Quaternion< ScalarIdentity ()
 

Additional Inherited Members

- Public Types inherited from Eigen::QuaternionBase< Quaternion< Scalar_, Options_ > >
typedef AngleAxis< ScalarAngleAxisType
 
typedef Matrix< Scalar, 3, 3 > Matrix3
 
typedef Matrix< Scalar, 3, 1 > Vector3
 
- Public Types inherited from Eigen::RotationBase< Derived, Dim_ >
typedef Matrix< Scalar, Dim, Dim > RotationMatrixType
 
typedef internal::traits< Derived >::Scalar Scalar
 

Constructor & Destructor Documentation

◆ Quaternion() [1/8]

template<typename Scalar_ , int Options_>
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( )
inline

Default constructor leaving the quaternion uninitialized.

◆ Quaternion() [2/8]

template<typename Scalar_ , int Options_>
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( const Scalar &  w,
const Scalar &  x,
const Scalar &  y,
const Scalar &  z 
)
inline

Constructs and initializes the quaternion \( w+xi+yj+zk \) from its four coefficients w, x, y and z.

Warning
Note the order of the arguments: the real w coefficient first, while internally the coefficients are stored in the following order: [x, y, z, w]

◆ Quaternion() [3/8]

template<typename Scalar_ , int Options_>
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( const Scalar *  data)
inlineexplicit

Constructs and initialize a quaternion from the array data

◆ Quaternion() [4/8]

template<typename Scalar_ , int Options_>
template<class Derived >
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( const QuaternionBase< Derived > &  other)
inline

Copy constructor

◆ Quaternion() [5/8]

template<typename Scalar_ , int Options_>
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( const AngleAxisType aa)
inlineexplicit

Constructs and initializes a quaternion from the angle-axis aa

◆ Quaternion() [6/8]

template<typename Scalar_ , int Options_>
template<typename Derived >
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( const MatrixBase< Derived > &  other)
inlineexplicit

Constructs and initializes a quaternion from either:

  • a rotation matrix expression,
  • a 4D vector expression representing quaternion coefficients in the order [x, y, z, w].

◆ Quaternion() [7/8]

template<typename Scalar_ , int Options_>
template<typename OtherScalar , int OtherOptions>
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( const Quaternion< OtherScalar, OtherOptions > &  other)
inlineexplicit

Explicit copy constructor with scalar conversion

◆ Quaternion() [8/8]

template<typename Scalar_ , int Options_>
Eigen::Quaternion< Scalar_, Options_ >::Quaternion ( Quaternion< Scalar_, Options_ > &&  other)
inline

Default move constructor

Member Function Documentation

◆ FromTwoVectors()

template<typename Scalar_ , int Options_>
template<typename Derived1 , typename Derived2 >
Quaternion< Scalar, Options > Eigen::Quaternion< Scalar_, Options_ >::FromTwoVectors ( const MatrixBase< Derived1 > &  a,
const MatrixBase< Derived2 > &  b 
)

Returns a quaternion representing a rotation between the two arbitrary vectors a and b. In other words, the built rotation represent a rotation sending the line of direction a to the line of direction b, both lines passing through the origin.

Returns
resulting quaternion

Note that the two input vectors do not have to be normalized, and do not need to have the same norm.

◆ operator=()

template<typename Scalar_ , int Options_>
Quaternion & Eigen::Quaternion< Scalar_, Options_ >::operator= ( Quaternion< Scalar_, Options_ > &&  other)
inline

Default move assignment operator

◆ UnitRandom()

template<typename Scalar , int Options>
Quaternion< Scalar, Options > Eigen::Quaternion< Scalar, Options >::UnitRandom
static
Returns
a random unit quaternion following a uniform distribution law on SO(3)
Note
The implementation is based on http://planning.cs.uiuc.edu/node198.html

The documentation for this class was generated from the following file: