Public Types | Public Member Functions | Static Public Member Functions | Private Attributes | Friends | List of all members
Eigen::EulerAngles< _Scalar, _System > Class Template Reference

Represents a rotation in a 3 dimensional space as three Euler angles. More...

#include <EulerAngles.h>

Inheritance diagram for Eigen::EulerAngles< _Scalar, _System >:
Inheritance graph
[legend]

Public Types

typedef AngleAxis< ScalarAngleAxisType
 
typedef RotationBase< EulerAngles< _Scalar, _System >, 3 > Base
 
typedef Matrix< Scalar, 3, 3 > Matrix3
 
typedef Quaternion< ScalarQuaternionType
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef _Scalar Scalar
 
typedef _System System
 
typedef Matrix< Scalar, 3, 1 > Vector3
 
- Public Types inherited from Eigen::RotationBase< EulerAngles< _Scalar, _System >, 3 >
enum  
 
typedef Matrix< Scalar, Dim, DimRotationMatrixType
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef Matrix< Scalar, Dim, 1 > VectorType
 

Public Member Functions

Scalaralpha ()
 
Scalar alpha () const
 
Vector3angles ()
 
const Vector3angles () const
 
Scalarbeta ()
 
Scalar beta () const
 
template<typename NewScalarType >
EulerAngles< NewScalarType, Systemcast () const
 
 EulerAngles ()
 
template<typename Derived >
 EulerAngles (const MatrixBase< Derived > &other)
 
template<typename Derived >
 EulerAngles (const RotationBase< Derived, 3 > &rot)
 
 EulerAngles (const Scalar &alpha, const Scalar &beta, const Scalar &gamma)
 
 EulerAngles (const Scalar *data)
 
Scalargamma ()
 
Scalar gamma () const
 
EulerAngles inverse () const
 
bool isApprox (const EulerAngles &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
 
 operator QuaternionType () const
 
EulerAngles operator- () const
 
template<class Derived >
EulerAnglesoperator= (const MatrixBase< Derived > &other)
 
template<typename Derived >
EulerAnglesoperator= (const RotationBase< Derived, 3 > &rot)
 
Matrix3 toRotationMatrix () const
 
- Public Member Functions inherited from Eigen::RotationBase< EulerAngles< _Scalar, _System >, 3 >
EIGEN_DEVICE_FUNC VectorType _transformVector (const OtherVectorType &v) const
 
EIGEN_DEVICE_FUNC Derived & derived ()
 
const EIGEN_DEVICE_FUNC Derived & derived () const
 
EIGEN_DEVICE_FUNC Derived inverse () const
 
EIGEN_DEVICE_FUNC RotationMatrixType matrix () 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 Transform< Scalar, Dim, Isometry > operator* (const Translation< Scalar, Dim > &t) const
 
EIGEN_DEVICE_FUNC RotationMatrixType operator* (const UniformScaling< Scalar > &s) const
 
EIGEN_DEVICE_FUNC RotationMatrixType toRotationMatrix () const
 

Static Public Member Functions

static Vector3 AlphaAxisVector ()
 
static Vector3 BetaAxisVector ()
 
static Vector3 GammaAxisVector ()
 

Private Attributes

Vector3 m_angles
 

Friends

std::ostream & operator<< (std::ostream &s, const EulerAngles< Scalar, System > &eulerAngles)
 

Detailed Description

template<typename _Scalar, class _System>
class Eigen::EulerAngles< _Scalar, _System >

Represents a rotation in a 3 dimensional space as three Euler angles.

Euler rotation is a set of three rotation of three angles over three fixed axes, defined by the EulerSystem given as a template parameter.

Here is how intrinsic Euler angles works:

Note
This class support only intrinsic Euler angles for simplicity, see EulerSystem how to easily overcome this for extrinsic systems.

Rotation representation and conversions

It has been proved(see Wikipedia link below) that every rotation can be represented by Euler angles, but there is no single representation (e.g. unlike rotation matrices). Therefore, you can convert from Eigen rotation and to them (including rotation matrices, which is not called "rotations" by Eigen design).

Euler angles usually used for:

However, Euler angles are slow comparing to quaternion or matrices, because their unnatural math definition, although it's simple for human. To overcome this, this class provide easy movement from the math friendly representation to the human friendly representation, and vise-versa.

All the user need to do is a safe simple C++ type conversion, and this class take care for the math. Additionally, some axes related computation is done in compile time.

Euler angles ranges in conversions

Rotations representation as EulerAngles are not single (unlike matrices), and even have infinite EulerAngles representations.
For example, add or subtract 2*PI from either angle of EulerAngles and you'll get the same rotation. This is the general reason for infinite representation, but it's not the only general reason for not having a single representation.

When converting rotation to EulerAngles, this class convert it to specific ranges When converting some rotation to EulerAngles, the rules for ranges are as follow:

See also
EulerAngles(const MatrixBase<Derived>&)
EulerAngles(const RotationBase<Derived, 3>&)

Convenient user typedefs

Convenient typedefs for EulerAngles exist for float and double scalar, in a form of EulerAngles{A}{B}{C}{scalar}, e.g. EulerAnglesXYZd, EulerAnglesZYZf.

Only for positive axes{+x,+y,+z} Euler systems are have convenient typedef. If you need negative axes{-x,-y,-z}, it is recommended to create you own typedef with a word that represent what you need.

Example

// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2015 Tal Hadad <tal_hd@hotmail.com>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
#include "main.h"
#include <unsupported/Eigen/EulerAngles>
using namespace Eigen;
// Unfortunately, we need to specialize it in order to work. (We could add it in main.h test framework)
template <typename Scalar, class System>
{
return verifyIsApprox(a.angles(), b.angles());
}
// Verify that x is in the approxed range [a, b]
#define VERIFY_APPROXED_RANGE(a, x, b) \
do { \
VERIFY_IS_APPROX_OR_LESS_THAN(a, x); \
VERIFY_IS_APPROX_OR_LESS_THAN(x, b); \
} while(0)
const char X = EULER_X;
const char Y = EULER_Y;
const char Z = EULER_Z;
template<typename Scalar, class EulerSystem>
{
typedef EulerAngles<Scalar, EulerSystem> EulerAnglesType;
const Scalar ONE = Scalar(1);
const Scalar HALF_PI = Scalar(EIGEN_PI / 2);
const Scalar PI = Scalar(EIGEN_PI);
// It's very important calc the acceptable precision depending on the distance from the pole.
const Scalar longitudeRadius = std::abs(
std::cos(e.beta()) :
std::sin(e.beta())
);
Scalar precision = test_precision<Scalar>() / longitudeRadius;
Scalar betaRangeStart, betaRangeEnd;
{
betaRangeStart = -HALF_PI;
betaRangeEnd = HALF_PI;
}
else
{
{
betaRangeStart = 0;
betaRangeEnd = PI;
}
else
{
betaRangeStart = -PI;
betaRangeEnd = 0;
}
}
const Vector3 I_ = EulerAnglesType::AlphaAxisVector();
const Vector3 J_ = EulerAnglesType::BetaAxisVector();
const Vector3 K_ = EulerAnglesType::GammaAxisVector();
// Is approx checks
VERIFY(e.isApprox(e));
VERIFY_IS_NOT_APPROX(e, EulerAnglesType(e.alpha() + ONE, e.beta() + ONE, e.gamma() + ONE));
const Matrix3 m(e);
VERIFY_IS_APPROX(Scalar(m.determinant()), ONE);
EulerAnglesType ebis(m);
// When no roll(acting like polar representation), we have the best precision.
// One of those cases is when the Euler angles are on the pole, and because it's singular case,
// the computation returns no roll.
if (ebis.beta() == 0)
precision = test_precision<Scalar>();
// Check that eabis in range
VERIFY_APPROXED_RANGE(-PI, ebis.alpha(), PI);
VERIFY_APPROXED_RANGE(betaRangeStart, ebis.beta(), betaRangeEnd);
VERIFY_APPROXED_RANGE(-PI, ebis.gamma(), PI);
const Matrix3 mbis(AngleAxisType(ebis.alpha(), I_) * AngleAxisType(ebis.beta(), J_) * AngleAxisType(ebis.gamma(), K_));
VERIFY_IS_APPROX(Scalar(mbis.determinant()), ONE);
VERIFY_IS_APPROX(mbis, ebis.toRotationMatrix());
/*std::cout << "===================\n" <<
"e: " << e << std::endl <<
"eabis: " << eabis.transpose() << std::endl <<
"m: " << m << std::endl <<
"mbis: " << mbis << std::endl <<
"X: " << (m * Vector3::UnitX()).transpose() << std::endl <<
"X: " << (mbis * Vector3::UnitX()).transpose() << std::endl;*/
VERIFY(m.isApprox(mbis, precision));
// Test if ea and eabis are the same
// Need to check both singular and non-singular cases
// There are two singular cases.
// 1. When I==K and sin(ea(1)) == 0
// 2. When I!=K and cos(ea(1)) == 0
// TODO: Make this test work well, and use range saturation function.
/*// If I==K, and ea[1]==0, then there no unique solution.
// The remark apply in the case where I!=K, and |ea[1]| is close to +-pi/2.
if( (i!=k || ea[1]!=0) && (i==k || !internal::isApprox(abs(ea[1]),Scalar(EIGEN_PI/2),test_precision<Scalar>())) )
VERIFY_IS_APPROX(ea, eabis);*/
// Quaternions
const QuaternionType q(e);
ebis = q;
const QuaternionType qbis(ebis);
VERIFY(internal::isApprox<Scalar>(std::abs(q.dot(qbis)), ONE, precision));
//VERIFY_IS_APPROX(eabis, eabis2);// Verify that the euler angles are still the same
// A suggestion for simple product test when will be supported.
/*EulerAnglesType e2(PI/2, PI/2, PI/2);
Matrix3 m2(e2);
VERIFY_IS_APPROX(e*e2, m*m2);*/
}
template<signed char A, signed char B, signed char C, typename Scalar>
{
}
template<signed char A, signed char B, signed char C, typename Scalar>
{
verify_euler_vec<+A,+B,+C>(ea);
verify_euler_vec<+A,+B,-C>(ea);
verify_euler_vec<+A,-B,+C>(ea);
verify_euler_vec<+A,-B,-C>(ea);
verify_euler_vec<-A,+B,+C>(ea);
verify_euler_vec<-A,+B,-C>(ea);
verify_euler_vec<-A,-B,+C>(ea);
verify_euler_vec<-A,-B,-C>(ea);
}
template<typename Scalar> void check_all_var(const Matrix<Scalar,3,1>& ea)
{
verify_euler_all_neg<X,Y,Z>(ea);
verify_euler_all_neg<X,Y,X>(ea);
verify_euler_all_neg<X,Z,Y>(ea);
verify_euler_all_neg<X,Z,X>(ea);
verify_euler_all_neg<Y,Z,X>(ea);
verify_euler_all_neg<Y,Z,Y>(ea);
verify_euler_all_neg<Y,X,Z>(ea);
verify_euler_all_neg<Y,X,Y>(ea);
verify_euler_all_neg<Z,X,Y>(ea);
verify_euler_all_neg<Z,X,Z>(ea);
verify_euler_all_neg<Z,Y,X>(ea);
verify_euler_all_neg<Z,Y,Z>(ea);
}
template<typename Scalar> void check_singular_cases(const Scalar& singularBeta)
{
const Scalar PI = Scalar(EIGEN_PI);
{
check_all_var(Vector3(PI/4, singularBeta, PI/3));
check_all_var(Vector3(PI/4, singularBeta - epsilon, PI/3));
check_all_var(Vector3(PI/4, singularBeta - Scalar(1.5)*epsilon, PI/3));
check_all_var(Vector3(PI/4, singularBeta - 2*epsilon, PI/3));
check_all_var(Vector3(PI*Scalar(0.8), singularBeta - epsilon, Scalar(0.9)*PI));
check_all_var(Vector3(PI*Scalar(-0.9), singularBeta + epsilon, PI*Scalar(0.3)));
check_all_var(Vector3(PI*Scalar(-0.6), singularBeta + Scalar(1.5)*epsilon, PI*Scalar(0.3)));
check_all_var(Vector3(PI*Scalar(-0.5), singularBeta + 2*epsilon, PI*Scalar(0.4)));
check_all_var(Vector3(PI*Scalar(0.9), singularBeta + epsilon, Scalar(0.8)*PI));
}
// This one for sanity, it had a problem with near pole cases in float scalar.
check_all_var(Vector3(PI*Scalar(0.8), singularBeta - Scalar(1E-6), Scalar(0.9)*PI));
}
template<typename Scalar> void eulerangles_manual()
{
const Vector3 Zero = Vector3::Zero();
const Scalar PI = Scalar(EIGEN_PI);
// singular cases
// non-singular cases
VectorX alpha = VectorX::LinSpaced(20, Scalar(-0.99) * PI, PI);
VectorX beta = VectorX::LinSpaced(20, Scalar(-0.49) * PI, Scalar(0.49) * PI);
VectorX gamma = VectorX::LinSpaced(20, Scalar(-0.99) * PI, PI);
for (int i = 0; i < alpha.size(); ++i) {
for (int j = 0; j < beta.size(); ++j) {
for (int k = 0; k < gamma.size(); ++k) {
}
}
}
}
template<typename Scalar> void eulerangles_rand()
{
typedef Array<Scalar,3,1> Array3;
typedef Quaternion<Scalar> Quaternionx;
Scalar a = internal::random<Scalar>(-Scalar(EIGEN_PI), Scalar(EIGEN_PI));
Quaternionx q1;
q1 = AngleAxisType(a, Vector3::Random().normalized());
m = q1;
Vector3 ea = m.eulerAngles(0,1,2);
ea = m.eulerAngles(0,1,0);
// Check with purely random Quaternion:
q1.coeffs() = Quaternionx::Coefficients::Random().normalized();
m = q1;
ea = m.eulerAngles(0,1,2);
ea = m.eulerAngles(0,1,0);
// Check with random angles in range [0:pi]x[-pi:pi]x[-pi:pi].
ea = (Array3::Random() + Array3(1,0,0))*Scalar(EIGEN_PI)*Array3(0.5,1,1);
ea[2] = ea[0] = internal::random<Scalar>(0,Scalar(EIGEN_PI));
ea[0] = ea[1] = internal::random<Scalar>(0,Scalar(EIGEN_PI));
ea[1] = 0;
ea.head(2).setZero();
ea.setZero();
}
{
// Simple cast test
EulerAnglesXYZd onesEd(1, 1, 1);
EulerAnglesXYZf onesEf = onesEd.cast<float>();
VERIFY_IS_APPROX(onesEd, onesEf.cast<double>());
// Simple Construction from Vector3 test
VERIFY_IS_APPROX(onesEd, EulerAnglesXYZd(Vector3d::Ones()));
CALL_SUBTEST_1( eulerangles_manual<float>() );
CALL_SUBTEST_2( eulerangles_manual<double>() );
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST_3( eulerangles_rand<float>() );
CALL_SUBTEST_4( eulerangles_rand<double>() );
}
// TODO: Add tests for auto diff
// TODO: Add tests for complex numbers
}

Output:

Additional reading

If you're want to get more idea about how Euler system work in Eigen see EulerSystem.

More information about Euler angles: https://en.wikipedia.org/wiki/Euler_angles

Template Parameters
_Scalarthe scalar type, i.e. the type of the angles.
_Systemthe EulerSystem to use, which represents the axes of rotation.

Definition at line 100 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

Member Typedef Documentation

◆ AngleAxisType

template<typename _Scalar , class _System >
typedef AngleAxis<Scalar> Eigen::EulerAngles< _Scalar, _System >::AngleAxisType

the equivalent angle-axis type

Definition at line 115 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ Base

template<typename _Scalar , class _System >
typedef RotationBase<EulerAngles<_Scalar, _System>, 3> Eigen::EulerAngles< _Scalar, _System >::Base

◆ Matrix3

template<typename _Scalar , class _System >
typedef Matrix<Scalar,3,3> Eigen::EulerAngles< _Scalar, _System >::Matrix3

the equivalent rotation matrix type

Definition at line 112 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ QuaternionType

template<typename _Scalar , class _System >
typedef Quaternion<Scalar> Eigen::EulerAngles< _Scalar, _System >::QuaternionType

the equivalent quaternion type

Definition at line 114 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ RealScalar

template<typename _Scalar , class _System >
typedef NumTraits<Scalar>::Real Eigen::EulerAngles< _Scalar, _System >::RealScalar

◆ Scalar

template<typename _Scalar , class _System >
typedef _Scalar Eigen::EulerAngles< _Scalar, _System >::Scalar

the scalar type of the angles

Definition at line 106 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ System

template<typename _Scalar , class _System >
typedef _System Eigen::EulerAngles< _Scalar, _System >::System

the EulerSystem to use, which represents the axes of rotation.

Definition at line 110 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ Vector3

template<typename _Scalar , class _System >
typedef Matrix<Scalar,3,1> Eigen::EulerAngles< _Scalar, _System >::Vector3

the equivalent 3 dimension vector type

Definition at line 113 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

Constructor & Destructor Documentation

◆ EulerAngles() [1/5]

template<typename _Scalar , class _System >
Eigen::EulerAngles< _Scalar, _System >::EulerAngles ( )
inline

Default constructor without initialization.

Definition at line 140 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ EulerAngles() [2/5]

template<typename _Scalar , class _System >
Eigen::EulerAngles< _Scalar, _System >::EulerAngles ( const Scalar alpha,
const Scalar beta,
const Scalar gamma 
)
inline

Constructs and initialize an EulerAngles (alpha, beta, gamma).

Definition at line 142 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ EulerAngles() [3/5]

template<typename _Scalar , class _System >
Eigen::EulerAngles< _Scalar, _System >::EulerAngles ( const Scalar data)
inlineexplicit

Constructs and initialize an EulerAngles from the array data {alpha, beta, gamma}

Definition at line 147 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ EulerAngles() [4/5]

template<typename _Scalar , class _System >
template<typename Derived >
Eigen::EulerAngles< _Scalar, _System >::EulerAngles ( const MatrixBase< Derived > &  other)
inlineexplicit

Constructs and initializes an EulerAngles from either:

  • a 3x3 rotation matrix expression(i.e. pure orthogonal matrix with determinant of +1),
  • a 3D vector expression representing Euler angles.
Note
If other is a 3x3 rotation matrix, the angles range rules will be as follow:
Alpha and gamma angles will be in the range [-PI, PI].
As for Beta angle:
  • If the system is Tait-Bryan, the beta angle will be in the range [-PI/2, PI/2].
  • otherwise:
    • If the beta axis is positive, the beta angle will be in the range [0, PI]
    • If the beta axis is negative, the beta angle will be in the range [-PI, 0]

Definition at line 162 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ EulerAngles() [5/5]

template<typename _Scalar , class _System >
template<typename Derived >
Eigen::EulerAngles< _Scalar, _System >::EulerAngles ( const RotationBase< Derived, 3 > &  rot)
inline

Constructs and initialize Euler angles from a rotation rot.

Note
If rot is an EulerAngles (even when it represented as RotationBase explicitly), angles ranges are undefined. Otherwise, alpha and gamma angles will be in the range [-PI, PI].
As for Beta angle:
  • If the system is Tait-Bryan, the beta angle will be in the range [-PI/2, PI/2].
  • otherwise:
    • If the beta axis is positive, the beta angle will be in the range [0, PI]
    • If the beta axis is negative, the beta angle will be in the range [-PI, 0]

Definition at line 176 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

Member Function Documentation

◆ alpha() [1/2]

template<typename _Scalar , class _System >
Scalar& Eigen::EulerAngles< _Scalar, _System >::alpha ( )
inline
Returns
A read-write reference to the angle of the first angle.

Definition at line 200 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ alpha() [2/2]

template<typename _Scalar , class _System >
Scalar Eigen::EulerAngles< _Scalar, _System >::alpha ( ) const
inline
Returns
The value of the first angle.

Definition at line 198 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ AlphaAxisVector()

template<typename _Scalar , class _System >
static Vector3 Eigen::EulerAngles< _Scalar, _System >::AlphaAxisVector ( )
inlinestatic
Returns
the axis vector of the first (alpha) rotation

Definition at line 118 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ angles() [1/2]

template<typename _Scalar , class _System >
Vector3& Eigen::EulerAngles< _Scalar, _System >::angles ( )
inline
Returns
A read-write reference to the angle values stored in a vector (alpha, beta, gamma).

Definition at line 195 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ angles() [2/2]

template<typename _Scalar , class _System >
const Vector3& Eigen::EulerAngles< _Scalar, _System >::angles ( ) const
inline
Returns
The angle values stored in a vector (alpha, beta, gamma).

Definition at line 193 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ beta() [1/2]

template<typename _Scalar , class _System >
Scalar& Eigen::EulerAngles< _Scalar, _System >::beta ( )
inline
Returns
A read-write reference to the angle of the second angle.

Definition at line 205 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ beta() [2/2]

template<typename _Scalar , class _System >
Scalar Eigen::EulerAngles< _Scalar, _System >::beta ( ) const
inline
Returns
The value of the second angle.

Definition at line 203 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ BetaAxisVector()

template<typename _Scalar , class _System >
static Vector3 Eigen::EulerAngles< _Scalar, _System >::BetaAxisVector ( )
inlinestatic
Returns
the axis vector of the second (beta) rotation

Definition at line 124 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ cast()

template<typename _Scalar , class _System >
template<typename NewScalarType >
EulerAngles<NewScalarType, System> Eigen::EulerAngles< _Scalar, _System >::cast ( ) const
inline
Returns
*this with scalar type casted to NewScalarType

Definition at line 292 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ gamma() [1/2]

template<typename _Scalar , class _System >
Scalar& Eigen::EulerAngles< _Scalar, _System >::gamma ( )
inline
Returns
A read-write reference to the angle of the third angle.

Definition at line 210 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ gamma() [2/2]

template<typename _Scalar , class _System >
Scalar Eigen::EulerAngles< _Scalar, _System >::gamma ( ) const
inline
Returns
The value of the third angle.

Definition at line 208 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ GammaAxisVector()

template<typename _Scalar , class _System >
static Vector3 Eigen::EulerAngles< _Scalar, _System >::GammaAxisVector ( )
inlinestatic
Returns
the axis vector of the third (gamma) rotation

Definition at line 130 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ inverse()

template<typename _Scalar , class _System >
EulerAngles Eigen::EulerAngles< _Scalar, _System >::inverse ( ) const
inline
Returns
The Euler angles rotation inverse (which is as same as the negative), (-alpha, -beta, -gamma).

Definition at line 215 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ isApprox()

template<typename _Scalar , class _System >
bool Eigen::EulerAngles< _Scalar, _System >::isApprox ( const EulerAngles< _Scalar, _System > &  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 264 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ operator QuaternionType()

template<typename _Scalar , class _System >
Eigen::EulerAngles< _Scalar, _System >::operator QuaternionType ( ) const
inline

Convert the Euler angles to quaternion.

Definition at line 276 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ operator-()

template<typename _Scalar , class _System >
EulerAngles Eigen::EulerAngles< _Scalar, _System >::operator- ( ) const
inline
Returns
The Euler angles rotation negative (which is as same as the inverse), (-alpha, -beta, -gamma).

Definition at line 225 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ operator=() [1/2]

template<typename _Scalar , class _System >
template<class Derived >
EulerAngles& Eigen::EulerAngles< _Scalar, _System >::operator= ( const MatrixBase< Derived > &  other)
inline

Set *this from either:

  • a 3x3 rotation matrix expression(i.e. pure orthogonal matrix with determinant of +1),
  • a 3D vector expression representing Euler angles.

See EulerAngles(const MatrixBase<Derived, 3>&) for more information about angles ranges output.

Definition at line 238 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ operator=() [2/2]

template<typename _Scalar , class _System >
template<typename Derived >
EulerAngles& Eigen::EulerAngles< _Scalar, _System >::operator= ( const RotationBase< Derived, 3 > &  rot)
inline

Set *this from a rotation.

See EulerAngles(const RotationBase<Derived, 3>&) for more information about angles ranges output.

Definition at line 255 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

◆ toRotationMatrix()

template<typename _Scalar , class _System >
Matrix3 Eigen::EulerAngles< _Scalar, _System >::toRotationMatrix ( ) const
inline
Returns
an equivalent 3x3 rotation matrix.

Definition at line 269 of file unsupported/Eigen/src/EulerAngles/EulerAngles.h.

Friends And Related Function Documentation

◆ operator<<

template<typename _Scalar , class _System >
std::ostream& operator<< ( std::ostream &  s,
const EulerAngles< Scalar, System > &  eulerAngles 
)
friend

Member Data Documentation

◆ m_angles

template<typename _Scalar , class _System >
Vector3 Eigen::EulerAngles< _Scalar, _System >::m_angles
private

The documentation for this class was generated from the following file:
Eigen::EULER_X
@ EULER_X
Definition: EulerSystem.h:63
Eigen::EulerAngles
Represents a rotation in a 3 dimensional space as three Euler angles.
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:100
Eigen::EulerAngles::Matrix3
Matrix< Scalar, 3, 3 > Matrix3
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:112
Eigen
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
B
Matrix< SCALARB, Dynamic, Dynamic, opt_B > B
Definition: bench_gemm.cpp:49
Eigen::EULER_Z
@ EULER_Z
Definition: EulerSystem.h:65
EIGEN_PI
#define EIGEN_PI
Definition: Eigen/src/Core/MathFunctions.h:16
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
ceres::sin
Jet< T, N > sin(const Jet< T, N > &f)
Definition: jet.h:439
b
Scalar * b
Definition: benchVecAdd.cpp:17
gtsam::Y
GaussianFactorGraphValuePair Y
Definition: HybridGaussianProductFactor.cpp:29
Eigen::AngleAxis
Represents a 3D rotation as a rotation angle around an arbitrary 3D axis.
Definition: ForwardDeclarations.h:290
Eigen::EulerAngles::QuaternionType
Quaternion< Scalar > QuaternionType
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:114
check_singular_cases
void check_singular_cases(const Scalar &singularBeta)
Definition: test/EulerAngles.cpp:175
VERIFY_APPROXED_RANGE
#define VERIFY_APPROXED_RANGE(a, x, b)
Definition: test/EulerAngles.cpp:24
Eigen::Array
General-purpose arrays with easy API for coefficient-wise operations.
Definition: Array.h:45
X
#define X
Definition: icosphere.cpp:20
ceres::cos
Jet< T, N > cos(const Jet< T, N > &f)
Definition: jet.h:426
Eigen::verifyIsApprox
bool verifyIsApprox(const Type1 &a, const Type2 &b)
Definition: main.h:585
A
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:48
CALL_SUBTEST_4
#define CALL_SUBTEST_4(FUNC)
Definition: split_test_helper.h:22
epsilon
static double epsilon
Definition: testRot3.cpp:37
VERIFY_IS_NOT_APPROX
#define VERIFY_IS_NOT_APPROX(a, b)
Definition: integer_types.cpp:17
CALL_SUBTEST_3
#define CALL_SUBTEST_3(FUNC)
Definition: split_test_helper.h:16
verify_euler_all_neg
void verify_euler_all_neg(const Matrix< Scalar, 3, 1 > &ea)
Definition: test/EulerAngles.cpp:144
CALL_SUBTEST_1
#define CALL_SUBTEST_1(FUNC)
Definition: split_test_helper.h:4
Eigen::EULER_Y
@ EULER_Y
Definition: EulerSystem.h:64
j
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2
Eigen::numext::q
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:1984
check_all_var
void check_all_var(const Matrix< Scalar, 3, 1 > &ea)
Definition: geo_eulerangles.cpp:41
Eigen::EulerAngles::Vector3
Matrix< Scalar, 3, 1 > Vector3
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:113
Eigen::EulerAngles::alpha
Scalar alpha() const
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:198
eulerangles_manual
void eulerangles_manual()
Definition: test/EulerAngles.cpp:197
Eigen::g_repeat
static int g_repeat
Definition: main.h:169
m
Matrix3f m
Definition: AngleAxis_mimic_euler.cpp:1
eulerangles_rand
void eulerangles_rand()
Definition: test/EulerAngles.cpp:229
Eigen::EulerAngles::gamma
Scalar gamma() const
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:208
CALL_SUBTEST_2
#define CALL_SUBTEST_2(FUNC)
Definition: split_test_helper.h:10
Eigen::EulerSystem
Represents a fixed Euler rotation system.
Definition: EulerSystem.h:126
E
DiscreteKey E(5, 2)
VERIFY_IS_APPROX
#define VERIFY_IS_APPROX(a, b)
Definition: integer_types.cpp:15
Eigen::EulerSystem::IsTaitBryan
@ IsTaitBryan
Definition: EulerSystem.h:156
verify_euler
void verify_euler(const Matrix< Scalar, 3, 1 > &ea, int i, int j, int k)
Definition: geo_eulerangles.cpp:17
Eigen::Quaternion
The quaternion class used to represent 3D orientations and rotations.
Definition: ForwardDeclarations.h:293
a
ArrayXXi a
Definition: Array_initializer_list_23_cxx11.cpp:1
C
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:50
verify_euler_vec
void verify_euler_vec(const Matrix< Scalar, 3, 1 > &ea)
Definition: test/EulerAngles.cpp:138
precision
cout precision(2)
main.h
Eigen::EulerAngles::beta
Scalar beta() const
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:203
EIGEN_DECLARE_TEST
#define EIGEN_DECLARE_TEST(X)
Definition: main.h:201
Eigen::EulerAngles::AngleAxisType
AngleAxis< Scalar > AngleAxisType
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:115
Eigen::EulerSystem::IsBetaOpposite
@ IsBetaOpposite
Definition: EulerSystem.h:148
VectorX
Matrix< Scalar, Dynamic, 1 > VectorX
Definition: sparse_lu.cpp:41
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: 3rdparty/Eigen/Eigen/src/Core/Matrix.h:178
Eigen::EulerAngles::Scalar
_Scalar Scalar
Definition: unsupported/Eigen/src/EulerAngles/EulerAngles.h:106
abs
#define abs(x)
Definition: datatypes.h:17
Z
#define Z
Definition: icosphere.cpp:21
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:232
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
VERIFY
#define VERIFY(a)
Definition: main.h:380


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:10:49