10#ifndef EIGEN_ROTATION2D_H
11#define EIGEN_ROTATION2D_H
34template<
typename _Scalar>
struct traits<Rotation2D<_Scalar> >
36 typedef _Scalar Scalar;
40template<
typename _Scalar>
43 typedef RotationBase<Rotation2D<_Scalar>,2>
Base;
47 using Base::operator*;
52 typedef Matrix<Scalar,2,1>
Vector2;
53 typedef Matrix<Scalar,2,2> Matrix2;
71 template<
typename Derived>
72 EIGEN_DEVICE_FUNC
explicit Rotation2D(
const MatrixBase<Derived>& m)
74 fromRotationMatrix(m.derived());
78 EIGEN_DEVICE_FUNC
inline Scalar angle()
const {
return m_angle; }
93 else if(tmp<-
Scalar(EIGEN_PI)) tmp +=
Scalar(2*EIGEN_PI);
102 {
return Rotation2D(m_angle + other.m_angle); }
106 { m_angle += other.m_angle;
return *
this; }
112 template<
typename Derived>
113 EIGEN_DEVICE_FUNC
Rotation2D& fromRotationMatrix(
const MatrixBase<Derived>& m);
123 template<
typename Derived>
125 {
return fromRotationMatrix(m.derived()); }
141 template<
typename NewScalarType>
142 EIGEN_DEVICE_FUNC
inline typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type
cast()
const
143 {
return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*
this); }
146 template<
typename OtherScalarType>
147 EIGEN_DEVICE_FUNC
inline explicit Rotation2D(
const Rotation2D<OtherScalarType>& other)
149 m_angle =
Scalar(other.angle());
159 {
return internal::isApprox(m_angle,other.m_angle, prec); }
174template<
typename Scalar>
175template<
typename Derived>
178 EIGEN_USING_STD(atan2)
179 EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
180 m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
186template<
typename Scalar>
187typename Rotation2D<Scalar>::Matrix2
192 Scalar sinA = sin(m_angle);
193 Scalar cosA = cos(m_angle);
194 return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
Represents a rotation/orientation in a 2 dimensional space.
Definition Rotation2D.h:42
Scalar & angle()
Definition Rotation2D.h:81
Scalar smallestPositiveAngle() const
Definition Rotation2D.h:84
Rotation2D(const Rotation2D< OtherScalarType > &other)
Definition Rotation2D.h:147
Rotation2D()
Definition Rotation2D.h:65
Rotation2D & operator=(const MatrixBase< Derived > &m)
Definition Rotation2D.h:124
Rotation2D(const MatrixBase< Derived > &m)
Definition Rotation2D.h:72
Rotation2D inverse() const
Definition Rotation2D.h:98
Rotation2D & operator*=(const Rotation2D &other)
Definition Rotation2D.h:105
Matrix2 toRotationMatrix() const
Definition Rotation2D.h:188
bool isApprox(const Rotation2D &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition Rotation2D.h:158
Scalar angle() const
Definition Rotation2D.h:78
internal::cast_return_type< Rotation2D, Rotation2D< NewScalarType > >::type cast() const
Definition Rotation2D.h:142
Rotation2D slerp(const Scalar &t, const Rotation2D &other) const
Definition Rotation2D.h:130
_Scalar Scalar
Definition Rotation2D.h:51
Scalar smallestAngle() const
Definition Rotation2D.h:90
Rotation2D operator*(const Rotation2D &other) const
Definition Rotation2D.h:101
Rotation2D(const Scalar &a)
Definition Rotation2D.h:62
Common base class for compact rotation representations.
Definition RotationBase.h:30
friend RotationMatrixType operator*(const EigenBase< OtherDerived > &l, const Rotation2D< _Scalar > &r)
Definition RotationBase.h:76
Rotation2D< float > Rotation2Df
Definition Rotation2D.h:165
Rotation2D< double > Rotation2Dd
Definition Rotation2D.h:168
Matrix< Type, 2, 2 > Matrix2
[c++11]
Definition Matrix.h:540
Namespace containing all symbols from the Eigen library.
Definition Core:141
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:233