Eigen::QuaternionBase< Derived > Class Template Reference

Base class for quaternion expressions. More...

+ Inheritance diagram for Eigen::QuaternionBase< Derived >:

Public Types

enum  { Flags }
 
typedef AngleAxis< ScalarAngleAxisType
 
typedef RotationBase< Derived, 3 > Base
 
typedef internal::traits< Derived >::Coefficients Coefficients
 
typedef Coefficients::CoeffReturnType CoeffReturnType
 
typedef Matrix< Scalar, 3, 3 > Matrix3
 
typedef std::conditional_t< bool(internal::traits< Derived >::Flags &LvalueBit), Scalar &, CoeffReturnTypeNonConstCoeffReturnType
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef Matrix< Scalar, 3, 1 > Vector3
 
- Public Types inherited from Eigen::RotationBase< Derived, 3 >
enum  
 
typedef Matrix< Scalar, Dim, DimRotationMatrixType
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef Matrix< Scalar, Dim, 1 > VectorType
 

Public Member Functions

Vector3 _transformVector (const Vector3 &v) const
 
template<class OtherDerived >
Scalar angularDistance (const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
internal::traits< Derived >::Scalar angularDistance (const QuaternionBase< OtherDerived > &other) const
 
template<typename NewScalarType >
internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast () const
 
internal::traits< Derived >::Coefficientscoeffs ()
 
const internal::traits< Derived >::Coefficientscoeffs () const
 
Quaternion< Scalarconjugate () const
 
template<class OtherDerived >
Scalar dot (const QuaternionBase< OtherDerived > &other) const
 
Quaternion< Scalarinverse () const
 
template<class OtherDerived >
bool isApprox (const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
 
Scalar norm () const
 
void normalize ()
 
Quaternion< Scalarnormalized () const
 
template<class OtherDerived >
bool operator!= (const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
Quaternion< typename internal::traits< Derived >::Scalaroperator* (const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
Quaternion< Scalaroperator* (const QuaternionBase< OtherDerived > &q) const
 
template<class OtherDerived >
Derived & operator*= (const QuaternionBase< OtherDerived > &q)
 
Derived & operator= (const AngleAxisType &aa)
 
template<class MatrixDerived >
Derived & operator= (const MatrixBase< MatrixDerived > &xpr)
 
template<class OtherDerived >
Derived & operator= (const MatrixBase< OtherDerived > &m)
 
QuaternionBase< Derived > & operator= (const QuaternionBase< Derived > &other)
 
template<class OtherDerived >
Derived & operator= (const QuaternionBase< OtherDerived > &other)
 
template<class OtherDerived >
bool operator== (const QuaternionBase< OtherDerived > &other) const
 
template<typename Derived1 , typename Derived2 >
Derived & setFromTwoVectors (const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
 
QuaternionBasesetIdentity ()
 
template<class OtherDerived >
Quaternion< Scalarslerp (const Scalar &t, const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
Quaternion< typename internal::traits< Derived >::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, 3 >
VectorType _transformVector (const OtherVectorType &v) const
 
Derived & derived ()
 
const Derived & derived () const
 
Derived inverse () const
 
RotationMatrixType matrix () const
 
internal::rotation_base_generic_product_selector< Derived, OtherDerived, OtherDerived::IsVectorAtCompileTime >::ReturnType operator* (const EigenBase< OtherDerived > &e) const
 
Transform< Scalar, Dim, Mode > operator* (const Transform< Scalar, Dim, Mode, Options > &t) const
 
Transform< Scalar, Dim, Isometry > operator* (const Translation< Scalar, Dim > &t) const
 
RotationMatrixType operator* (const UniformScaling< Scalar > &s) const
 
RotationMatrixType toRotationMatrix () const
 

Static Public Member Functions

static Quaternion< ScalarIdentity ()
 

Detailed Description

template<class Derived>
class Eigen::QuaternionBase< Derived >

Base class for quaternion expressions.

This is defined in the Geometry module.

Template Parameters
Derivedderived type (CRTP)
See also
class Quaternion

Definition at line 37 of file Quaternion.h.

Member Typedef Documentation

◆ AngleAxisType

template<class Derived >
typedef AngleAxis<Scalar> Eigen::QuaternionBase< Derived >::AngleAxisType

the equivalent angle-axis type

Definition at line 63 of file Quaternion.h.

◆ Base

template<class Derived >
typedef RotationBase<Derived, 3> Eigen::QuaternionBase< Derived >::Base

Definition at line 40 of file Quaternion.h.

◆ Coefficients

template<class Derived >
typedef internal::traits<Derived>::Coefficients Eigen::QuaternionBase< Derived >::Coefficients

Definition at line 47 of file Quaternion.h.

◆ CoeffReturnType

template<class Derived >
typedef Coefficients::CoeffReturnType Eigen::QuaternionBase< Derived >::CoeffReturnType

Definition at line 48 of file Quaternion.h.

◆ Matrix3

template<class Derived >
typedef Matrix<Scalar,3,3> Eigen::QuaternionBase< Derived >::Matrix3

the equivalent rotation matrix type

Definition at line 61 of file Quaternion.h.

◆ NonConstCoeffReturnType

template<class Derived >
typedef std::conditional_t<bool(internal::traits<Derived>::Flags&LvalueBit), Scalar&, CoeffReturnType> Eigen::QuaternionBase< Derived >::NonConstCoeffReturnType

Definition at line 50 of file Quaternion.h.

◆ RealScalar

template<class Derived >
typedef NumTraits<Scalar>::Real Eigen::QuaternionBase< Derived >::RealScalar

Definition at line 46 of file Quaternion.h.

◆ Scalar

template<class Derived >
typedef internal::traits<Derived>::Scalar Eigen::QuaternionBase< Derived >::Scalar

Definition at line 45 of file Quaternion.h.

◆ Vector3

template<class Derived >
typedef Matrix<Scalar,3,1> Eigen::QuaternionBase< Derived >::Vector3

the type of a 3D vector

Definition at line 59 of file Quaternion.h.

Member Enumeration Documentation

◆ anonymous enum

template<class Derived >
anonymous enum
Enumerator
Flags 

Definition at line 53 of file Quaternion.h.

53  {
54  Flags = Eigen::internal::traits<Derived>::Flags
55  };

Member Function Documentation

◆ _transformVector()

template<class Derived >
QuaternionBase< Derived >::Vector3 Eigen::QuaternionBase< Derived >::_transformVector ( const Vector3 v) const
inline

return the result vector of v through the rotation

Rotation of a vector by a quaternion.

Remarks
If the quaternion is used to rotate several points (>1) then it is much more efficient to first convert it to a 3x3 Matrix. Comparison of the operation cost for n transformations:
  • Quaternion2: 30n
  • Via a Matrix3: 24 + 15n

Definition at line 537 of file Quaternion.h.

538 {
539  // Note that this algorithm comes from the optimization by hand
540  // of the conversion to a Matrix followed by a Matrix/Vector product.
541  // It appears to be much faster than the common algorithm found
542  // in the literature (30 versus 39 flops). It also requires two
543  // Vector3 as temporaries.
544  Vector3 uv = this->vec().cross(v);
545  uv += uv;
546  return v + this->w() * uv + this->vec().cross(uv);
547 }
Array< int, Dynamic, 1 > v
Matrix< Scalar, 3, 1 > Vector3
Definition: Quaternion.h:59
CoeffReturnType w() const
Definition: Quaternion.h:74
const VectorBlock< const Coefficients, 3 > vec() const
Definition: Quaternion.h:86

◆ angularDistance() [1/2]

template<class Derived >
template<class OtherDerived >
Scalar Eigen::QuaternionBase< Derived >::angularDistance ( const QuaternionBase< OtherDerived > &  other) const

◆ angularDistance() [2/2]

template<class Derived >
template<class OtherDerived >
internal::traits<Derived>::Scalar Eigen::QuaternionBase< Derived >::angularDistance ( const QuaternionBase< OtherDerived > &  other) const
inline
Returns
the angle (in radian) between two rotations
See also
dot()

Definition at line 770 of file Quaternion.h.

771 {
773  Quaternion<Scalar> d = (*this) * other.conjugate();
774  return Scalar(2) * atan2( d.vec().norm(), numext::abs(d.w()) );
775 }
const CwiseBinaryOp< atan2< Scalar >, const Derived, const OtherDerived > atan2(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const
#define EIGEN_USING_STD(FUNC)
Definition: Macros.h:1080
internal::traits< Derived >::Scalar Scalar
Definition: Quaternion.h:45
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)

◆ cast()

template<class Derived >
template<typename NewScalarType >
internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type Eigen::QuaternionBase< Derived >::cast ( ) const
inline
Returns
*this with scalar type casted to NewScalarType

Note that if NewScalarType is equal to the current scalar type of *this then this function smartly returns a const reference to *this.

◆ coeffs() [1/2]

template<class Derived >
internal::traits<Derived>::Coefficients& Eigen::QuaternionBase< Derived >::coeffs ( )
inline
Returns
a vector expression of the coefficients (x,y,z,w)

Definition at line 95 of file Quaternion.h.

95 { return derived().coeffs(); }
const Derived & derived() const
Definition: RotationBase.h:43

◆ coeffs() [2/2]

template<class Derived >
const internal::traits<Derived>::Coefficients& Eigen::QuaternionBase< Derived >::coeffs ( ) const
inline
Returns
a read-only vector expression of the coefficients (x,y,z,w)

Definition at line 92 of file Quaternion.h.

92 { return derived().coeffs(); }

◆ conjugate()

template<class Derived >
Quaternion< typename internal::traits< Derived >::Scalar > Eigen::QuaternionBase< Derived >::conjugate ( void  ) const
inline
Returns
the conjugated quaternion
the conjugate of the *this which is equal to the multiplicative inverse if the quaternion is normalized. The conjugate of a quaternion represents the opposite rotation.
See also
Quaternion2::inverse()

Definition at line 757 of file Quaternion.h.

758 {
759  return internal::quat_conj<Architecture::Target, Derived,
760  typename internal::traits<Derived>::Scalar>::run(*this);
761 
762 }

◆ dot()

template<class Derived >
template<class OtherDerived >
Scalar Eigen::QuaternionBase< Derived >::dot ( const QuaternionBase< OtherDerived > &  other) const
inline
Returns
the dot product of *this and other Geometrically speaking, the dot product of two unit quaternions corresponds to the cosine of half the angle between the two rotations.
See also
angularDistance()

Definition at line 141 of file Quaternion.h.

141 { return coeffs().dot(other.coeffs()); }
const internal::traits< Derived >::Coefficients & coeffs() const
Definition: Quaternion.h:92

◆ Identity()

template<class Derived >
static Quaternion<Scalar> Eigen::QuaternionBase< Derived >::Identity ( )
inlinestatic
Returns
a quaternion representing an identity rotation
See also
MatrixBase::Identity()

Definition at line 113 of file Quaternion.h.

113 { return Quaternion<Scalar>(Scalar(1), Scalar(0), Scalar(0), Scalar(0)); }

◆ inverse()

template<class Derived >
Quaternion< typename internal::traits< Derived >::Scalar > Eigen::QuaternionBase< Derived >::inverse
inline
Returns
the quaternion describing the inverse rotation
the multiplicative inverse of *this Note that in most cases, i.e., if you simply want the opposite rotation, and/or the quaternion is normalized, then it is enough to use the conjugate.
See also
QuaternionBase::conjugate()

Definition at line 726 of file Quaternion.h.

727 {
728  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
729  Scalar n2 = this->squaredNorm();
730  if (n2 > Scalar(0))
731  return Quaternion<Scalar>(conjugate().coeffs() / n2);
732  else
733  {
734  // return an invalid result to flag the error
735  return Quaternion<Scalar>(Coefficients::Zero());
736  }
737 }
Scalar squaredNorm() const
Definition: Quaternion.h:122
Quaternion< Scalar > conjugate() const
Definition: Quaternion.h:757

◆ isApprox()

template<class Derived >
template<class OtherDerived >
bool Eigen::QuaternionBase< Derived >::isApprox ( const QuaternionBase< OtherDerived > &  other,
const RealScalar prec = NumTraits<Scalar>::dummy_precision() 
) const
inline
Returns
true if *this is approximately equal to other, within the precision determined by prec.
See also
MatrixBase::isApprox()

Definition at line 184 of file Quaternion.h.

185  { return coeffs().isApprox(other.coeffs(), prec); }

◆ norm()

template<class Derived >
Scalar Eigen::QuaternionBase< Derived >::norm ( ) const
inline
Returns
the norm of the quaternion's coefficients
See also
QuaternionBase::squaredNorm(), MatrixBase::norm()

Definition at line 127 of file Quaternion.h.

127 { return coeffs().norm(); }

◆ normalize()

template<class Derived >
void Eigen::QuaternionBase< Derived >::normalize ( void  )
inline

Normalizes the quaternion *this

See also
normalized(), MatrixBase::normalize()

Definition at line 131 of file Quaternion.h.

131 { coeffs().normalize(); }

◆ normalized()

template<class Derived >
Quaternion<Scalar> Eigen::QuaternionBase< Derived >::normalized ( ) const
inline
Returns
a normalized copy of *this
See also
normalize(), MatrixBase::normalized()

Definition at line 134 of file Quaternion.h.

134 { return Quaternion<Scalar>(coeffs().normalized()); }
Quaternion< Scalar > normalized() const
Definition: Quaternion.h:134

◆ operator!=()

template<class Derived >
template<class OtherDerived >
bool Eigen::QuaternionBase< Derived >::operator!= ( const QuaternionBase< OtherDerived > &  other) const
inline
Returns
true if at least one pair of coefficients of *this and other are not exactly equal to each other.
Warning
When using floating point scalar values you probably should rather use a fuzzy comparison such as isApprox()
See also
isApprox(), operator==

Definition at line 176 of file Quaternion.h.

177  { return coeffs() != other.coeffs(); }

◆ operator*() [1/2]

template<class Derived >
template<class OtherDerived >
Quaternion<typename internal::traits<Derived>::Scalar> Eigen::QuaternionBase< Derived >::operator* ( const QuaternionBase< OtherDerived > &  other) const
inline
Returns
the concatenation of two rotations as a quaternion-quaternion product

Definition at line 511 of file Quaternion.h.

512 {
513  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
514  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
515  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
516  typename internal::traits<Derived>::Scalar>::run(*this, other);
517 }
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26

◆ operator*() [2/2]

template<class Derived >
template<class OtherDerived >
Quaternion<Scalar> Eigen::QuaternionBase< Derived >::operator* ( const QuaternionBase< OtherDerived > &  q) const
inline

◆ operator*=()

template<class Derived >
template<class OtherDerived >
Derived & Eigen::QuaternionBase< Derived >::operator*= ( const QuaternionBase< OtherDerived > &  other)
inline
See also
operator*(Quaternion)

Definition at line 522 of file Quaternion.h.

523 {
524  derived() = derived() * other.derived();
525  return derived();
526 }

◆ operator=() [1/5]

template<class Derived >
Derived & Eigen::QuaternionBase< Derived >::operator= ( const AngleAxisType aa)
inline

Set *this from an angle-axis aa and returns a reference to *this

Definition at line 567 of file Quaternion.h.

568 {
571  Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
572  this->w() = cos(ha);
573  this->vec() = sin(ha) * aa.axis();
574  return derived();
575 }
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_cos_op< typename Derived::Scalar >, const Derived > cos(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)

◆ operator=() [2/5]

template<class Derived >
template<class MatrixDerived >
Derived& Eigen::QuaternionBase< Derived >::operator= ( const MatrixBase< MatrixDerived > &  xpr)
inline

Set *this from the expression xpr:

  • if xpr is a 4x1 vector, then xpr is assumed to be a quaternion
  • if xpr is a 3x3 matrix, then xpr is assumed to be rotation matrix and xpr is converted to a quaternion

Definition at line 585 of file Quaternion.h.

586 {
587  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
588  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
589  internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
590  return derived();
591 }

◆ operator=() [3/5]

template<class Derived >
template<class OtherDerived >
Derived& Eigen::QuaternionBase< Derived >::operator= ( const MatrixBase< OtherDerived > &  m)

◆ operator=() [4/5]

template<class Derived >
QuaternionBase< Derived > & Eigen::QuaternionBase< Derived >::operator= ( const QuaternionBase< Derived > &  other)
inline

Definition at line 550 of file Quaternion.h.

551 {
552  coeffs() = other.coeffs();
553  return derived();
554 }

◆ operator=() [5/5]

template<class Derived >
template<class OtherDerived >
Derived & Eigen::QuaternionBase< Derived >::operator= ( const QuaternionBase< OtherDerived > &  other)
inline

Definition at line 558 of file Quaternion.h.

559 {
560  coeffs() = other.coeffs();
561  return derived();
562 }

◆ operator==()

template<class Derived >
template<class OtherDerived >
bool Eigen::QuaternionBase< Derived >::operator== ( const QuaternionBase< OtherDerived > &  other) const
inline
Returns
true if each coefficients of *this and other are all exactly equal.
Warning
When using floating point scalar values you probably should rather use a fuzzy comparison such as isApprox()
See also
isApprox(), operator!=

Definition at line 168 of file Quaternion.h.

169  { return coeffs() == other.coeffs(); }

◆ setFromTwoVectors()

template<class Derived >
template<typename Derived1 , typename Derived2 >
Derived & Eigen::QuaternionBase< Derived >::setFromTwoVectors ( const MatrixBase< Derived1 > &  a,
const MatrixBase< Derived2 > &  b 
)
inline
Returns
the quaternion which transform a into b through a rotation

Sets *this to be 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
a reference to *this.

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

Definition at line 644 of file Quaternion.h.

645 {
647  Vector3 v0 = a.normalized();
648  Vector3 v1 = b.normalized();
649  Scalar c = v1.dot(v0);
650 
651  // if dot == -1, vectors are nearly opposites
652  // => accurately compute the rotation axis by computing the
653  // intersection of the two planes. This is done by solving:
654  // x^T v0 = 0
655  // x^T v1 = 0
656  // under the constraint:
657  // ||x|| = 1
658  // which yields a singular value problem
659  if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
660  {
661  c = numext::maxi(c,Scalar(-1));
662  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
663  JacobiSVD<Matrix<Scalar,2,3>, ComputeFullV> svd(m);
664  Vector3 axis = svd.matrixV().col(2);
665 
666  Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
667  this->w() = sqrt(w2);
668  this->vec() = axis * sqrt(Scalar(1) - w2);
669  return derived();
670  }
671  Vector3 axis = v0.cross(v1);
672  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
673  Scalar invs = Scalar(1)/s;
674  this->vec() = axis * invs;
675  this->w() = s * Scalar(0.5);
676 
677  return derived();
678 }
Matrix3f m
Array< int, 3, 1 > b
Array33i c
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
M1<< 1, 2, 3, 4, 5, 6, 7, 8, 9;Map< RowVectorXf > v1(M1.data(), M1.size())
@ ComputeFullV
Definition: Constants.h:399
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)

◆ setIdentity()

template<class Derived >
QuaternionBase& Eigen::QuaternionBase< Derived >::setIdentity ( )
inline
See also
QuaternionBase::Identity(), MatrixBase::setIdentity()

Definition at line 117 of file Quaternion.h.

117 { coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(1); return *this; }

◆ slerp() [1/2]

template<class Derived >
template<class OtherDerived >
Quaternion<Scalar> Eigen::QuaternionBase< Derived >::slerp ( const Scalar t,
const QuaternionBase< OtherDerived > &  other 
) const

◆ slerp() [2/2]

template<class Derived >
template<class OtherDerived >
Quaternion<typename internal::traits<Derived>::Scalar> Eigen::QuaternionBase< Derived >::slerp ( const Scalar t,
const QuaternionBase< OtherDerived > &  other 
) const
Returns
the spherical linear interpolation between the two quaternions *this and other at the parameter t in [0;1].

This represents an interpolation for a constant motion between *this and other, see also http://en.wikipedia.org/wiki/Slerp.

Definition at line 788 of file Quaternion.h.

789 {
792  const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
793  Scalar d = this->dot(other);
794  Scalar absD = numext::abs(d);
795 
796  Scalar scale0;
797  Scalar scale1;
798 
799  if(absD>=one)
800  {
801  scale0 = Scalar(1) - t;
802  scale1 = t;
803  }
804  else
805  {
806  // theta is the angle between the 2 quaternions
807  Scalar theta = acos(absD);
808  Scalar sinTheta = sin(theta);
809 
810  scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;
811  scale1 = sin( ( t * theta) ) / sinTheta;
812  }
813  if(d<Scalar(0)) scale1 = -scale1;
814 
815  return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
816 }
Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition: Quaternion.h:141
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_acos_op< typename Derived::Scalar >, const Derived > acos(const Eigen::ArrayBase< Derived > &x)

◆ squaredNorm()

template<class Derived >
Scalar Eigen::QuaternionBase< Derived >::squaredNorm ( ) const
inline
Returns
the squared norm of the quaternion's coefficients
See also
QuaternionBase::norm(), MatrixBase::squaredNorm()

Definition at line 122 of file Quaternion.h.

122 { return coeffs().squaredNorm(); }

◆ toRotationMatrix()

template<class Derived >
QuaternionBase< Derived >::Matrix3 Eigen::QuaternionBase< Derived >::toRotationMatrix ( void  ) const
inline
Returns
an equivalent 3x3 rotation matrix

Convert the quaternion to a 3x3 rotation matrix. The quaternion is required to be normalized, otherwise the result is undefined.

Definition at line 598 of file Quaternion.h.

599 {
600  // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
601  // if not inlined then the cost of the return by value is huge ~ +35%,
602  // however, not inlining this function is an order of magnitude slower, so
603  // it has to be inlined, and so the return by value is not an issue
604  Matrix3 res;
605 
606  const Scalar tx = Scalar(2)*this->x();
607  const Scalar ty = Scalar(2)*this->y();
608  const Scalar tz = Scalar(2)*this->z();
609  const Scalar twx = tx*this->w();
610  const Scalar twy = ty*this->w();
611  const Scalar twz = tz*this->w();
612  const Scalar txx = tx*this->x();
613  const Scalar txy = ty*this->x();
614  const Scalar txz = tz*this->x();
615  const Scalar tyy = ty*this->y();
616  const Scalar tyz = tz*this->y();
617  const Scalar tzz = tz*this->z();
618 
619  res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
620  res.coeffRef(0,1) = txy-twz;
621  res.coeffRef(0,2) = txz+twy;
622  res.coeffRef(1,0) = txy+twz;
623  res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
624  res.coeffRef(1,2) = tyz-twx;
625  res.coeffRef(2,0) = txz-twy;
626  res.coeffRef(2,1) = tyz+twx;
627  res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
628 
629  return res;
630 }
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
CoeffReturnType z() const
Definition: Quaternion.h:72
CoeffReturnType y() const
Definition: Quaternion.h:70
Matrix< Scalar, 3, 3 > Matrix3
Definition: Quaternion.h:61
CoeffReturnType x() const
Definition: Quaternion.h:68

◆ vec() [1/2]

template<class Derived >
VectorBlock<Coefficients,3> Eigen::QuaternionBase< Derived >::vec ( )
inline
Returns
a vector expression of the imaginary part (x,y,z)

Definition at line 89 of file Quaternion.h.

89 { return coeffs().template head<3>(); }

◆ vec() [2/2]

template<class Derived >
const VectorBlock<const Coefficients,3> Eigen::QuaternionBase< Derived >::vec ( ) const
inline
Returns
a read-only vector expression of the imaginary part (x,y,z)

Definition at line 86 of file Quaternion.h.

86 { return coeffs().template head<3>(); }

◆ w() [1/2]

template<class Derived >
NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::w ( )
inline
Returns
a reference to the w coefficient (if Derived is a non-const lvalue)

Definition at line 83 of file Quaternion.h.

83 { return this->derived().coeffs().w(); }

◆ w() [2/2]

template<class Derived >
CoeffReturnType Eigen::QuaternionBase< Derived >::w ( ) const
inline
Returns
the w coefficient

Definition at line 74 of file Quaternion.h.

74 { return this->derived().coeffs().coeff(3); }

◆ x() [1/2]

template<class Derived >
NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::x ( )
inline
Returns
a reference to the x coefficient (if Derived is a non-const lvalue)

Definition at line 77 of file Quaternion.h.

77 { return this->derived().coeffs().x(); }

◆ x() [2/2]

template<class Derived >
CoeffReturnType Eigen::QuaternionBase< Derived >::x ( ) const
inline
Returns
the x coefficient

Definition at line 68 of file Quaternion.h.

68 { return this->derived().coeffs().coeff(0); }

◆ y() [1/2]

template<class Derived >
NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::y ( )
inline
Returns
a reference to the y coefficient (if Derived is a non-const lvalue)

Definition at line 79 of file Quaternion.h.

79 { return this->derived().coeffs().y(); }

◆ y() [2/2]

template<class Derived >
CoeffReturnType Eigen::QuaternionBase< Derived >::y ( ) const
inline
Returns
the y coefficient

Definition at line 70 of file Quaternion.h.

70 { return this->derived().coeffs().coeff(1); }

◆ z() [1/2]

template<class Derived >
NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::z ( )
inline
Returns
a reference to the z coefficient (if Derived is a non-const lvalue)

Definition at line 81 of file Quaternion.h.

81 { return this->derived().coeffs().z(); }

◆ z() [2/2]

template<class Derived >
CoeffReturnType Eigen::QuaternionBase< Derived >::z ( ) const
inline
Returns
the z coefficient

Definition at line 72 of file Quaternion.h.

72 { return this->derived().coeffs().coeff(2); }

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