Eigen::Rotation2D< Scalar_ > Class Template Reference

Represents a rotation/orientation in a 2 dimensional space. More...

+ Inheritance diagram for Eigen::Rotation2D< Scalar_ >:

Public Types

enum  { Dim }
 
typedef Matrix< Scalar, 2, 2 > Matrix2
 
typedef Scalar_ Scalar
 
typedef Matrix< Scalar, 2, 1 > Vector2
 
- Public Types inherited from Eigen::RotationBase< Rotation2D< Scalar_ >, 2 >
enum  
 
typedef Matrix< Scalar, Dim, DimRotationMatrixType
 
typedef internal::traits< Rotation2D< Scalar_ > >::Scalar Scalar
 
typedef Matrix< Scalar, Dim, 1 > VectorType
 

Public Member Functions

Scalarangle ()
 
Scalar angle () const
 
template<typename NewScalarType >
internal::cast_return_type< Rotation2D, Rotation2D< NewScalarType > >::type cast () const
 
template<typename Derived >
Rotation2DfromRotationMatrix (const MatrixBase< Derived > &m)
 
template<typename Derived >
Rotation2D< Scalar > & fromRotationMatrix (const MatrixBase< Derived > &mat)
 
Rotation2D inverse () const
 
bool isApprox (const Rotation2D &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
 
Rotation2D operator* (const Rotation2D &other) const
 
Vector2 operator* (const Vector2 &vec) const
 
Rotation2Doperator*= (const Rotation2D &other)
 
template<typename Derived >
Rotation2Doperator= (const MatrixBase< Derived > &m)
 
 Rotation2D ()
 
template<typename Derived >
 Rotation2D (const MatrixBase< Derived > &m)
 
template<typename OtherScalarType >
 Rotation2D (const Rotation2D< OtherScalarType > &other)
 
 Rotation2D (const Scalar &a)
 
Rotation2D slerp (const Scalar &t, const Rotation2D &other) const
 
Scalar smallestAngle () const
 
Scalar smallestPositiveAngle () const
 
Matrix2 toRotationMatrix () const
 
- Public Member Functions inherited from Eigen::RotationBase< Rotation2D< Scalar_ >, 2 >
VectorType _transformVector (const OtherVectorType &v) const
 
Rotation2D< Scalar_ > & derived ()
 
const Rotation2D< Scalar_ > & derived () const
 
Rotation2D< Scalar_ > inverse () const
 
RotationMatrixType matrix () const
 
internal::rotation_base_generic_product_selector< Rotation2D< Scalar_ >, 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 Rotation2D Identity ()
 

Protected Attributes

Scalar m_angle
 

Private Types

typedef RotationBase< Rotation2D< Scalar_ >, 2 > Base
 

Detailed Description

template<typename Scalar_>
class Eigen::Rotation2D< Scalar_ >

Represents a rotation/orientation in a 2 dimensional space.

This is defined in the Geometry module.

Template Parameters
Scalar_the scalar type, i.e., the type of the coefficients

This class is equivalent to a single scalar representing a counter clock wise rotation as a single angle in radian. It provides some additional features such as the automatic conversion from/to a 2x2 rotation matrix. Moreover this class aims to provide a similar interface to Quaternion in order to facilitate the writing of generic algorithms dealing with rotations.

See also
class Quaternion, class Transform

Definition at line 43 of file Rotation2D.h.

Member Typedef Documentation

◆ Base

template<typename Scalar_ >
typedef RotationBase<Rotation2D<Scalar_>,2> Eigen::Rotation2D< Scalar_ >::Base
private

Definition at line 45 of file Rotation2D.h.

◆ Matrix2

template<typename Scalar_ >
typedef Matrix<Scalar,2,2> Eigen::Rotation2D< Scalar_ >::Matrix2

Definition at line 55 of file Rotation2D.h.

◆ Scalar

template<typename Scalar_ >
typedef Scalar_ Eigen::Rotation2D< Scalar_ >::Scalar

the scalar type of the coefficients

Definition at line 53 of file Rotation2D.h.

◆ Vector2

template<typename Scalar_ >
typedef Matrix<Scalar,2,1> Eigen::Rotation2D< Scalar_ >::Vector2

Definition at line 54 of file Rotation2D.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar_ >
anonymous enum
Enumerator
Dim 

Definition at line 51 of file Rotation2D.h.

51 { Dim = 2 };

Constructor & Destructor Documentation

◆ Rotation2D() [1/4]

template<typename Scalar_ >
Eigen::Rotation2D< Scalar_ >::Rotation2D ( const Scalar a)
inlineexplicit

Construct a 2D counter clock wise rotation from the angle a in radian.

Definition at line 64 of file Rotation2D.h.

64 : m_angle(a) {}

◆ Rotation2D() [2/4]

template<typename Scalar_ >
Eigen::Rotation2D< Scalar_ >::Rotation2D ( )
inline

Default constructor wihtout initialization. The represented rotation is undefined.

Definition at line 67 of file Rotation2D.h.

67 {}

◆ Rotation2D() [3/4]

template<typename Scalar_ >
template<typename Derived >
Eigen::Rotation2D< Scalar_ >::Rotation2D ( const MatrixBase< Derived > &  m)
inlineexplicit

Construct a 2D rotation from a 2x2 rotation matrix mat.

See also
fromRotationMatrix()

Definition at line 74 of file Rotation2D.h.

75  {
76  fromRotationMatrix(m.derived());
77  }
Matrix3f m
Rotation2D & fromRotationMatrix(const MatrixBase< Derived > &m)

◆ Rotation2D() [4/4]

template<typename Scalar_ >
template<typename OtherScalarType >
Eigen::Rotation2D< Scalar_ >::Rotation2D ( const Rotation2D< OtherScalarType > &  other)
inlineexplicit

Copy constructor with scalar type conversion

Definition at line 149 of file Rotation2D.h.

150  {
151  m_angle = Scalar(other.angle());
152  }

Member Function Documentation

◆ angle() [1/2]

template<typename Scalar_ >
Scalar& Eigen::Rotation2D< Scalar_ >::angle ( )
inline
Returns
a read-write reference to the rotation angle

Definition at line 83 of file Rotation2D.h.

83 { return m_angle; }

◆ angle() [2/2]

template<typename Scalar_ >
Scalar Eigen::Rotation2D< Scalar_ >::angle ( ) const
inline
Returns
the rotation angle

Definition at line 80 of file Rotation2D.h.

80 { return m_angle; }

◆ cast()

template<typename Scalar_ >
template<typename NewScalarType >
internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type Eigen::Rotation2D< Scalar_ >::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.

Definition at line 144 of file Rotation2D.h.

145  { return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*this); }

◆ fromRotationMatrix() [1/2]

template<typename Scalar_ >
template<typename Derived >
Rotation2D& Eigen::Rotation2D< Scalar_ >::fromRotationMatrix ( const MatrixBase< Derived > &  m)

◆ fromRotationMatrix() [2/2]

template<typename Scalar_ >
template<typename Derived >
Rotation2D<Scalar>& Eigen::Rotation2D< Scalar_ >::fromRotationMatrix ( const MatrixBase< Derived > &  mat)

Set *this from a 2x2 rotation matrix mat. In other words, this function extract the rotation angle from the rotation matrix.

Definition at line 178 of file Rotation2D.h.

179 {
181  EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
182  m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
183  return *this;
184 }
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
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26

◆ Identity()

template<typename Scalar_ >
static Rotation2D Eigen::Rotation2D< Scalar_ >::Identity ( )
inlinestatic

Definition at line 154 of file Rotation2D.h.

154 { return Rotation2D(0); }

◆ inverse()

template<typename Scalar_ >
Rotation2D Eigen::Rotation2D< Scalar_ >::inverse ( ) const
inline
Returns
the inverse rotation

Definition at line 100 of file Rotation2D.h.

100 { return Rotation2D(-m_angle); }

◆ isApprox()

template<typename Scalar_ >
bool Eigen::Rotation2D< Scalar_ >::isApprox ( const Rotation2D< Scalar_ > &  other,
const typename NumTraits< Scalar >::Real &  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 160 of file Rotation2D.h.

161  { return internal::isApprox(m_angle,other.m_angle, prec); }
bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())

◆ operator*() [1/2]

template<typename Scalar_ >
Rotation2D Eigen::Rotation2D< Scalar_ >::operator* ( const Rotation2D< Scalar_ > &  other) const
inline

Concatenates two rotations

Definition at line 103 of file Rotation2D.h.

104  { return Rotation2D(m_angle + other.m_angle); }

◆ operator*() [2/2]

template<typename Scalar_ >
Vector2 Eigen::Rotation2D< Scalar_ >::operator* ( const Vector2 vec) const
inline

Applies the rotation to a 2D vector

Definition at line 111 of file Rotation2D.h.

112  { return toRotationMatrix() * vec; }
Matrix2 toRotationMatrix() const
Definition: Rotation2D.h:190

◆ operator*=()

template<typename Scalar_ >
Rotation2D& Eigen::Rotation2D< Scalar_ >::operator*= ( const Rotation2D< Scalar_ > &  other)
inline

Concatenates two rotations

Definition at line 107 of file Rotation2D.h.

108  { m_angle += other.m_angle; return *this; }

◆ operator=()

template<typename Scalar_ >
template<typename Derived >
Rotation2D& Eigen::Rotation2D< Scalar_ >::operator= ( const MatrixBase< Derived > &  m)
inline

Set *this from a 2x2 rotation matrix mat. In other words, this function extract the rotation angle from the rotation matrix.

This method is an alias for fromRotationMatrix()

See also
fromRotationMatrix()

Definition at line 126 of file Rotation2D.h.

127  { return fromRotationMatrix(m.derived()); }

◆ slerp()

template<typename Scalar_ >
Rotation2D Eigen::Rotation2D< Scalar_ >::slerp ( const Scalar t,
const Rotation2D< Scalar_ > &  other 
) const
inline
Returns
the spherical interpolation between *this and other using parameter t. It is in fact equivalent to a linear interpolation.

Definition at line 132 of file Rotation2D.h.

133  {
134  Scalar dist = Rotation2D(other.m_angle-m_angle).smallestAngle();
135  return Rotation2D(m_angle + dist*t);
136  }

◆ smallestAngle()

template<typename Scalar_ >
Scalar Eigen::Rotation2D< Scalar_ >::smallestAngle ( ) const
inline
Returns
the rotation angle in [-pi,pi]

Definition at line 92 of file Rotation2D.h.

92  {
94  if(tmp>Scalar(EIGEN_PI)) tmp -= Scalar(2*EIGEN_PI);
95  else if(tmp<-Scalar(EIGEN_PI)) tmp += Scalar(2*EIGEN_PI);
96  return tmp;
97  }
#define EIGEN_PI
Definition: MathFunctions.h:16
EIGEN_ALWAYS_INLINE T fmod(const T &a, const T &b)

◆ smallestPositiveAngle()

template<typename Scalar_ >
Scalar Eigen::Rotation2D< Scalar_ >::smallestPositiveAngle ( ) const
inline
Returns
the rotation angle in [0,2pi]

Definition at line 86 of file Rotation2D.h.

86  {
88  return tmp<Scalar(0) ? tmp + Scalar(2*EIGEN_PI) : tmp;
89  }

◆ toRotationMatrix()

template<typename Scalar >
Rotation2D< Scalar >::Matrix2 Eigen::Rotation2D< Scalar >::toRotationMatrix ( void  ) const

Constructs and

Returns
an equivalent 2x2 rotation matrix.

Definition at line 190 of file Rotation2D.h.

191 {
194  Scalar sinA = sin(m_angle);
195  Scalar cosA = cos(m_angle);
196  return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
197 }
Matrix< Scalar, 2, 2 > Matrix2
Definition: Rotation2D.h:55
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)

Member Data Documentation

◆ m_angle

template<typename Scalar_ >
Scalar Eigen::Rotation2D< Scalar_ >::m_angle
protected

Definition at line 59 of file Rotation2D.h.


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