Eigen::QuaternionBase< Derived > Class Template Reference

Base class for quaternion expressions. More...

#include <Quaternion.h>

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

Public Types

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

Public Member Functions

EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType x () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType y () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType z () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType w () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType x ()
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType y ()
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType z ()
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType w ()
 
EIGEN_DEVICE_FUNC const VectorBlock< const Coefficients, 3 > vec () const
 
EIGEN_DEVICE_FUNC VectorBlock< Coefficients, 3 > vec ()
 
EIGEN_DEVICE_FUNC const internal::traits< Derived >::Coefficientscoeffs () const
 
EIGEN_DEVICE_FUNC internal::traits< Derived >::Coefficientscoeffs ()
 
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase< Derived > & operator= (const QuaternionBase< Derived > &other)
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived & operator= (const QuaternionBase< OtherDerived > &other)
 
EIGEN_DEVICE_FUNC Derived & operator= (const AngleAxisType &aa)
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC Derived & operator= (const MatrixBase< OtherDerived > &m)
 
EIGEN_DEVICE_FUNC QuaternionBasesetIdentity ()
 
EIGEN_DEVICE_FUNC Scalar squaredNorm () const
 
EIGEN_DEVICE_FUNC Scalar norm () const
 
EIGEN_DEVICE_FUNC void normalize ()
 
EIGEN_DEVICE_FUNC Quaternion< Scalarnormalized () const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC Scalar dot (const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC Scalar angularDistance (const QuaternionBase< OtherDerived > &other) const
 
EIGEN_DEVICE_FUNC Matrix3 toRotationMatrix () const
 
template<typename Derived1 , typename Derived2 >
EIGEN_DEVICE_FUNC Derived & setFromTwoVectors (const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion< Scalaroperator* (const QuaternionBase< OtherDerived > &q) const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived & operator*= (const QuaternionBase< OtherDerived > &q)
 
EIGEN_DEVICE_FUNC Quaternion< Scalarinverse () const
 
EIGEN_DEVICE_FUNC Quaternion< Scalarconjugate () const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC Quaternion< Scalarslerp (const Scalar &t, const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC bool operator== (const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC bool operator!= (const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC bool isApprox (const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
 
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Vector3 _transformVector (const Vector3 &v) const
 
template<typename NewScalarType >
EIGEN_DEVICE_FUNC std::enable_if_t< internal::is_same< Scalar, NewScalarType >::value, const Derived & > cast () const
 
template<typename NewScalarType >
EIGEN_DEVICE_FUNC std::enable_if_t<!internal::is_same< Scalar, NewScalarType >::value, Quaternion< NewScalarType > > cast () const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion< typename internal::traits< Derived >::Scalaroperator* (const QuaternionBase< OtherDerived > &other) const
 
template<class MatrixDerived >
EIGEN_DEVICE_FUNC Derived & operator= (const MatrixBase< MatrixDerived > &xpr)
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC internal::traits< Derived >::Scalar angularDistance (const QuaternionBase< OtherDerived > &other) const
 
template<class OtherDerived >
EIGEN_DEVICE_FUNC Quaternion< typename internal::traits< Derived >::Scalarslerp (const Scalar &t, const QuaternionBase< OtherDerived > &other) const
 
- Public Member Functions inherited from Eigen::RotationBase< Derived, 3 >
EIGEN_DEVICE_FUNC const Derived & derived () const
 
EIGEN_DEVICE_FUNC Derived & derived ()
 
EIGEN_DEVICE_FUNC RotationMatrixType toRotationMatrix () const
 
EIGEN_DEVICE_FUNC RotationMatrixType matrix () const
 
EIGEN_DEVICE_FUNC Derived inverse () const
 
EIGEN_DEVICE_FUNC Transform< Scalar, Dim, Isometry > operator* (const Translation< Scalar, Dim > &t) const
 
EIGEN_DEVICE_FUNC RotationMatrixType operator* (const UniformScaling< Scalar > &s) const
 
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE internal::rotation_base_generic_product_selector< Derived, OtherDerived, OtherDerived::IsVectorAtCompileTime >::ReturnType operator* (const EigenBase< OtherDerived > &e) const
 
EIGEN_DEVICE_FUNC Transform< Scalar, Dim, Mode > operator* (const Transform< Scalar, Dim, Mode, Options > &t) const
 
EIGEN_DEVICE_FUNC VectorType _transformVector (const OtherVectorType &v) const
 

Static Public Member Functions

static EIGEN_DEVICE_FUNC Quaternion< ScalarIdentity ()
 

Friends

std::ostream & operator<< (std::ostream &s, const QuaternionBase< Derived > &q)
 

Detailed Description

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

Base class for quaternion expressions.

\geometry_module

Template Parameters
Derivedderived type (CRTP)
See also
class Quaternion

Member Typedef Documentation

◆ AngleAxisType

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

the equivalent angle-axis type

◆ Base

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

◆ Coefficients

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

◆ CoeffReturnType

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

◆ Matrix3

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

the equivalent rotation matrix type

◆ NonConstCoeffReturnType

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

◆ RealScalar

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

◆ Scalar

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

◆ Vector3

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

the type of a 3D vector

Member Enumeration Documentation

◆ anonymous enum

template<class Derived >
anonymous enum
Enumerator
Flags 
Extend namespace for flags.
Definition: fsi_chan_precond_driver.cc:56
Definition: ForwardDeclarations.h:21

Member Function Documentation

◆ _transformVector()

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

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
544  {
545  // Note that this algorithm comes from the optimization by hand
546  // of the conversion to a Matrix followed by a Matrix/Vector product.
547  // It appears to be much faster than the common algorithm found
548  // in the literature (30 versus 39 flops). It also requires two
549  // Vector3 as temporaries.
550  Vector3 uv = this->vec().cross(v);
551  uv += uv;
552  return v + this->w() * uv + this->vec().cross(uv);
553 }
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
Matrix< Scalar, 3, 1 > Vector3
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:53
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType w() const
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:66
EIGEN_DEVICE_FUNC const VectorBlock< const Coefficients, 3 > vec() const
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:78

References v, and w.

◆ angularDistance() [1/2]

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

◆ angularDistance() [2/2]

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC 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()
765  {
767  Quaternion<Scalar> d = (*this) * other.conjugate();
768  return Scalar(2) * atan2(d.vec().norm(), numext::abs(d.w()));
769 }
#define EIGEN_USING_STD(FUNC)
Definition: Macros.h:1090
internal::traits< Derived >::Scalar Scalar
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:42
This class contains the 4 components of a quaternion and the standard operators and functions needed ...
Definition: Kernel/Math/Quaternion.h:42
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
Definition: MathFunctions.h:1355
AutoDiffScalar< Matrix< typename internal::traits< internal::remove_all_t< DerTypeA > >::Scalar, Dynamic, 1 > > atan2(const AutoDiffScalar< DerTypeA > &a, const AutoDiffScalar< DerTypeB > &b)
Definition: AutoDiffScalar.h:558

References Eigen::numext::abs(), Eigen::atan2(), Eigen::QuaternionBase< Derived >::conjugate(), EIGEN_USING_STD, Eigen::QuaternionBase< Derived >::vec(), and Eigen::QuaternionBase< Derived >::w().

◆ cast() [1/2]

template<class Derived >
template<typename NewScalarType >
EIGEN_DEVICE_FUNC std::enable_if_t<internal::is_same<Scalar, NewScalarType>::value, const Derived&> Eigen::QuaternionBase< Derived >::cast ( ) const
inline
215  {
216  return derived();
217  }
EIGEN_DEVICE_FUNC const Derived & derived() const
Definition: RotationBase.h:43

References Eigen::RotationBase< Derived, 3 >::derived().

◆ cast() [2/2]

template<class Derived >
template<typename NewScalarType >
EIGEN_DEVICE_FUNC std::enable_if_t<!internal::is_same<Scalar, NewScalarType>::value, Quaternion<NewScalarType> > Eigen::QuaternionBase< Derived >::cast ( ) const
inline
222  {
223  return Quaternion<NewScalarType>(coeffs().template cast<NewScalarType>());
224  }
EIGEN_DEVICE_FUNC const internal::traits< Derived >::Coefficients & coeffs() const
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:84

References Eigen::QuaternionBase< Derived >::coeffs().

◆ coeffs() [1/2]

template<class Derived >
EIGEN_DEVICE_FUNC internal::traits<Derived>::Coefficients& Eigen::QuaternionBase< Derived >::coeffs ( )
inline
Returns
a vector expression of the coefficients (x,y,z,w)
89 { return derived().coeffs(); }

References Eigen::RotationBase< Derived, 3 >::derived().

◆ coeffs() [2/2]

◆ conjugate()

template<class Derived >
EIGEN_DEVICE_FUNC 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()
755  {
757 }
SCALAR Scalar
Definition: bench_gemm.cpp:45
auto run(Kernel kernel, Args &&... args) -> decltype(kernel(args...))
Definition: gpu_test_helper.h:414

References Eigen::run().

Referenced by Eigen::QuaternionBase< Derived >::angularDistance().

◆ dot()

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC 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()
143  {
144  return coeffs().dot(other.coeffs());
145  }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ Identity()

template<class Derived >
static EIGEN_DEVICE_FUNC Quaternion<Scalar> Eigen::QuaternionBase< Derived >::Identity ( )
inlinestatic
Returns
a quaternion representing an identity rotation
See also
MatrixBase::Identity()
109  {
110  return Quaternion<Scalar>(Scalar(1), Scalar(0), Scalar(0), Scalar(0));
111  }

◆ inverse()

template<class Derived >
EIGEN_DEVICE_FUNC 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()
726  {
727  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
728  Scalar n2 = this->squaredNorm();
729  if (n2 > Scalar(0))
730  return Quaternion<Scalar>(conjugate().coeffs() / n2);
731  else {
732  // return an invalid result to flag the error
734  }
735 }
EIGEN_DEVICE_FUNC Scalar squaredNorm() const
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:123
EIGEN_DEVICE_FUNC Quaternion< Scalar > conjugate() const
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:754
double Zero
Definition: pseudosolid_node_update_elements.cc:35

References oomph::PseudoSolidHelper::Zero.

◆ isApprox()

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC 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()
195  {
196  return coeffs().isApprox(other.coeffs(), prec);
197  }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ norm()

template<class Derived >
EIGEN_DEVICE_FUNC Scalar Eigen::QuaternionBase< Derived >::norm ( ) const
inline
Returns
the norm of the quaternion's coefficients
See also
QuaternionBase::squaredNorm(), MatrixBase::norm()
128 { return coeffs().norm(); }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ normalize()

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

Normalizes the quaternion *this

See also
normalized(), MatrixBase::normalize()
132 { coeffs().normalize(); }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ normalized()

template<class Derived >
EIGEN_DEVICE_FUNC Quaternion<Scalar> Eigen::QuaternionBase< Derived >::normalized ( ) const
inline
Returns
a normalized copy of *this
See also
normalize(), MatrixBase::normalized()
135 { return Quaternion<Scalar>(coeffs().normalized()); }
EIGEN_DEVICE_FUNC Quaternion< Scalar > normalized() const
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:135

References Eigen::QuaternionBase< Derived >::coeffs(), and Eigen::QuaternionBase< Derived >::normalized().

Referenced by Eigen::QuaternionBase< Derived >::normalized().

◆ operator!=()

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC 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==
185  {
186  return coeffs() != other.coeffs();
187  }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ operator*() [1/2]

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar> Eigen::QuaternionBase< Derived >::operator* ( const QuaternionBase< OtherDerived > &  other) const
Returns
the concatenation of two rotations as a quaternion-quaternion product
518  {
521  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
522  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
523  typename internal::traits<Derived>::Scalar>::run(*this, other);
524 }
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
@ Target
Definition: Constants.h:495
@ value
Definition: Meta.h:206

References EIGEN_STATIC_ASSERT, Eigen::run(), and Eigen::Architecture::Target.

◆ operator*() [2/2]

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

◆ operator*=()

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived & Eigen::QuaternionBase< Derived >::operator*= ( const QuaternionBase< OtherDerived > &  other)
See also
operator*(Quaternion)
530  {
531  derived() = derived() * other.derived();
532  return derived();
533 }

References Eigen::RotationBase< Derived, Dim_ >::derived().

◆ operator=() [1/5]

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

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

573  {
576  Scalar ha = Scalar(0.5) * aa.angle(); // Scalar(0.5) to suppress precision loss warnings
577  this->w() = cos(ha);
578  this->vec() = sin(ha) * aa.axis();
579  return derived();
580 }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137

References Eigen::AngleAxis< Scalar_ >::angle(), Eigen::AngleAxis< Scalar_ >::axis(), cos(), EIGEN_USING_STD, sin(), and w.

◆ operator=() [2/5]

template<class Derived >
template<class MatrixDerived >
EIGEN_DEVICE_FUNC 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
590  {
593  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
595  return derived();
596 }
void run(const string &dir_name, LinearSolver *linear_solver_pt, const unsigned nel_1d, bool mess_up_order)
Definition: two_d_poisson_compare_solvers.cc:317

References EIGEN_STATIC_ASSERT, and run().

◆ operator=() [3/5]

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

◆ operator=() [4/5]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase< Derived > & Eigen::QuaternionBase< Derived >::operator= ( const QuaternionBase< Derived > &  other)
557  {
558  coeffs() = other.coeffs();
559  return derived();
560 }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ operator=() [5/5]

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived & Eigen::QuaternionBase< Derived >::operator= ( const QuaternionBase< OtherDerived > &  other)
565  {
566  coeffs() = other.coeffs();
567  return derived();
568 }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ operator==()

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC 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!=
176  {
177  return coeffs() == other.coeffs();
178  }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ setFromTwoVectors()

template<class Derived >
template<typename Derived1 , typename Derived2 >
EIGEN_DEVICE_FUNC 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.

649  {
651  Vector3 v0 = a.normalized();
652  Vector3 v1 = b.normalized();
653  Scalar c = v1.dot(v0);
654 
655  // if dot == -1, vectors are nearly opposites
656  // => accurately compute the rotation axis by computing the
657  // intersection of the two planes. This is done by solving:
658  // x^T v0 = 0
659  // x^T v1 = 0
660  // under the constraint:
661  // ||x|| = 1
662  // which yields a singular value problem
663  if (c < Scalar(-1) + NumTraits<Scalar>::dummy_precision()) {
664  c = numext::maxi(c, Scalar(-1));
665  Matrix<Scalar, 2, 3> m;
666  m << v0.transpose(), v1.transpose();
667  JacobiSVD<Matrix<Scalar, 2, 3>, ComputeFullV> svd(m);
668  Vector3 axis = svd.matrixV().col(2);
669 
670  Scalar w2 = (Scalar(1) + c) * Scalar(0.5);
671  this->w() = sqrt(w2);
672  this->vec() = axis * sqrt(Scalar(1) - w2);
673  return derived();
674  }
675  Vector3 axis = v0.cross(v1);
676  Scalar s = sqrt((Scalar(1) + c) * Scalar(2));
677  Scalar invs = Scalar(1) / s;
678  this->vec() = axis * invs;
679  this->w() = s * Scalar(0.5);
680 
681  return derived();
682 }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
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())
Scalar * b
Definition: benchVecAdd.cpp:17
@ ComputeFullV
Definition: Constants.h:393
RealScalar s
Definition: level1_cplx_impl.h:130
const Scalar * a
Definition: level2_cplx_impl.h:32
int * m
Definition: level2_cplx_impl.h:294
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:926
int c
Definition: calibrate.py:100

References a, b, calibrate::c, Eigen::ComputeFullV, EIGEN_USING_STD, m, Eigen::numext::maxi(), s, sqrt(), svd(), v1(), and w.

Referenced by Eigen::Quaternion< Scalar_, Options_ >::FromTwoVectors().

◆ setIdentity()

template<class Derived >
EIGEN_DEVICE_FUNC QuaternionBase& Eigen::QuaternionBase< Derived >::setIdentity ( )
inline
See also
QuaternionBase::Identity(), MatrixBase::setIdentity()
115  {
116  coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(1);
117  return *this;
118  }

References Eigen::QuaternionBase< Derived >::coeffs().

Referenced by main().

◆ slerp() [1/2]

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

◆ slerp() [2/2]

template<class Derived >
template<class OtherDerived >
EIGEN_DEVICE_FUNC 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.

780  {
783  const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
784  Scalar d = this->dot(other);
785  Scalar absD = numext::abs(d);
786 
787  Scalar scale0;
788  Scalar scale1;
789 
790  if (absD >= one) {
791  scale0 = Scalar(1) - t;
792  scale1 = t;
793  } else {
794  // theta is the angle between the 2 quaternions
795  Scalar theta = acos(absD);
796  Scalar sinTheta = sin(theta);
797 
798  scale0 = sin((Scalar(1) - t) * theta) / sinTheta;
799  scale1 = sin((t * theta)) / sinTheta;
800  }
801  if (d < Scalar(0)) scale1 = -scale1;
802 
803  return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
804 }
AnnoyingScalar acos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:138
EIGEN_DEVICE_FUNC Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:143
double theta
Definition: two_d_biharmonic.cc:236
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43
t
Definition: plotPSD.py:36

References Eigen::numext::abs(), acos(), Eigen::QuaternionBase< Derived >::coeffs(), dot(), EIGEN_USING_STD, oomph::SarahBL::epsilon, sin(), plotPSD::t, and BiharmonicTestFunctions2::theta.

◆ squaredNorm()

template<class Derived >
EIGEN_DEVICE_FUNC Scalar Eigen::QuaternionBase< Derived >::squaredNorm ( ) const
inline
Returns
the squared norm of the quaternion's coefficients
See also
QuaternionBase::norm(), MatrixBase::squaredNorm()
123 { return coeffs().squaredNorm(); }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ toRotationMatrix()

template<class Derived >
EIGEN_DEVICE_FUNC 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.

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

References res, w, plotDoE::x, and y.

Referenced by Eigen::EulerAngles< Scalar_, _System >::toRotationMatrix().

◆ vec() [1/2]

template<class Derived >
EIGEN_DEVICE_FUNC VectorBlock<Coefficients, 3> Eigen::QuaternionBase< Derived >::vec ( )
inline
Returns
a vector expression of the imaginary part (x,y,z)
81 { return coeffs().template head<3>(); }

References Eigen::QuaternionBase< Derived >::coeffs().

◆ vec() [2/2]

template<class Derived >
EIGEN_DEVICE_FUNC const VectorBlock<const Coefficients, 3> Eigen::QuaternionBase< Derived >::vec ( ) const
inline
Returns
a read-only vector expression of the imaginary part (x,y,z)
78 { return coeffs().template head<3>(); }

References Eigen::QuaternionBase< Derived >::coeffs().

Referenced by Eigen::QuaternionBase< Derived >::angularDistance().

◆ w() [1/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::w ( )
inline
Returns
a reference to the w coefficient (if Derived is a non-const lvalue)
75 { return this->derived().coeffs().w(); }

References Eigen::RotationBase< Derived, 3 >::derived().

◆ w() [2/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType Eigen::QuaternionBase< Derived >::w ( ) const
inline
Returns
the w coefficient
66 { return this->derived().coeffs().coeff(3); }

References Eigen::RotationBase< Derived, 3 >::derived().

Referenced by Eigen::QuaternionBase< Derived >::angularDistance().

◆ x() [1/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::x ( )
inline
Returns
a reference to the x coefficient (if Derived is a non-const lvalue)
69 { return this->derived().coeffs().x(); }

References Eigen::RotationBase< Derived, 3 >::derived().

◆ x() [2/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType Eigen::QuaternionBase< Derived >::x ( ) const
inline
Returns
the x coefficient
60 { return this->derived().coeffs().coeff(0); }

References Eigen::RotationBase< Derived, 3 >::derived().

◆ y() [1/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::y ( )
inline
Returns
a reference to the y coefficient (if Derived is a non-const lvalue)
71 { return this->derived().coeffs().y(); }

References Eigen::RotationBase< Derived, 3 >::derived().

◆ y() [2/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType Eigen::QuaternionBase< Derived >::y ( ) const
inline
Returns
the y coefficient
62 { return this->derived().coeffs().coeff(1); }

References Eigen::RotationBase< Derived, 3 >::derived().

◆ z() [1/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR NonConstCoeffReturnType Eigen::QuaternionBase< Derived >::z ( )
inline
Returns
a reference to the z coefficient (if Derived is a non-const lvalue)
73 { return this->derived().coeffs().z(); }

References Eigen::RotationBase< Derived, 3 >::derived().

◆ z() [2/2]

template<class Derived >
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR CoeffReturnType Eigen::QuaternionBase< Derived >::z ( ) const
inline
Returns
the z coefficient
64 { return this->derived().coeffs().coeff(2); }

References Eigen::RotationBase< Derived, 3 >::derived().

Friends And Related Function Documentation

◆ operator<<

template<class Derived >
std::ostream& operator<< ( std::ostream &  s,
const QuaternionBase< Derived > &  q 
)
friend
228  {
229  s << q.x() << "i + " << q.y() << "j + " << q.z() << "k"
230  << " + " << q.w();
231  return s;
232  }
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019

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