Eigen  3.2.91
Rotation2D.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_ROTATION2D_H
11 #define EIGEN_ROTATION2D_H
12 
13 namespace Eigen {
14 
32 namespace internal {
33 
34 template<typename _Scalar> struct traits<Rotation2D<_Scalar> >
35 {
36  typedef _Scalar Scalar;
37 };
38 } // end namespace internal
39 
40 template<typename _Scalar>
41 class Rotation2D : public RotationBase<Rotation2D<_Scalar>,2>
42 {
43  typedef RotationBase<Rotation2D<_Scalar>,2> Base;
44 
45 public:
46 
47  using Base::operator*;
48 
49  enum { Dim = 2 };
51  typedef _Scalar Scalar;
54 
55 protected:
56 
57  Scalar m_angle;
58 
59 public:
60 
62  explicit inline Rotation2D(const Scalar& a) : m_angle(a) {}
63 
66 
68  inline Scalar angle() const { return m_angle; }
69 
71  inline Scalar& angle() { return m_angle; }
72 
74  inline Scalar smallestPositiveAngle() const {
75  Scalar tmp = fmod(m_angle,Scalar(2)*EIGEN_PI);
76  return tmp<Scalar(0) ? tmp + Scalar(2)*EIGEN_PI : tmp;
77  }
78 
80  inline Scalar smallestAngle() const {
81  Scalar tmp = fmod(m_angle,Scalar(2)*EIGEN_PI);
82  if(tmp>Scalar(EIGEN_PI)) tmp -= Scalar(2)*Scalar(EIGEN_PI);
83  else if(tmp<-Scalar(EIGEN_PI)) tmp += Scalar(2)*Scalar(EIGEN_PI);
84  return tmp;
85  }
86 
88  inline Rotation2D inverse() const { return Rotation2D(-m_angle); }
89 
91  inline Rotation2D operator*(const Rotation2D& other) const
92  { return Rotation2D(m_angle + other.m_angle); }
93 
95  inline Rotation2D& operator*=(const Rotation2D& other)
96  { m_angle += other.m_angle; return *this; }
97 
99  Vector2 operator* (const Vector2& vec) const
100  { return toRotationMatrix() * vec; }
101 
102  template<typename Derived>
103  Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
104  Matrix2 toRotationMatrix() const;
105 
109  inline Rotation2D slerp(const Scalar& t, const Rotation2D& other) const
110  {
111  Scalar dist = Rotation2D(other.m_angle-m_angle).smallestAngle();
112  return Rotation2D(m_angle + dist*t);
113  }
114 
120  template<typename NewScalarType>
121  inline typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type cast() const
122  { return typename internal::cast_return_type<Rotation2D,Rotation2D<NewScalarType> >::type(*this); }
123 
125  template<typename OtherScalarType>
126  inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
127  {
128  m_angle = Scalar(other.angle());
129  }
130 
131  static inline Rotation2D Identity() { return Rotation2D(0); }
132 
137  bool isApprox(const Rotation2D& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
138  { return internal::isApprox(m_angle,other.m_angle, prec); }
139 
140 };
141 
148 
153 template<typename Scalar>
154 template<typename Derived>
156 {
157  using std::atan2;
158  EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime==2 && Derived::ColsAtCompileTime==2,YOU_MADE_A_PROGRAMMING_MISTAKE)
159  m_angle = atan2(mat.coeff(1,0), mat.coeff(0,0));
160  return *this;
161 }
162 
165 template<typename Scalar>
168 {
169  using std::sin;
170  using std::cos;
171  Scalar sinA = sin(m_angle);
172  Scalar cosA = cos(m_angle);
173  return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
174 }
175 
176 } // end namespace Eigen
177 
178 #endif // EIGEN_ROTATION2D_H
internal::cast_return_type< Rotation2D, Rotation2D< NewScalarType > >::type cast() const
Definition: Rotation2D.h:121
Matrix2 toRotationMatrix() const
Definition: Rotation2D.h:167
Rotation2D slerp(const Scalar &t, const Rotation2D &other) const
Definition: Rotation2D.h:109
Scalar smallestPositiveAngle() const
Definition: Rotation2D.h:74
Definition: LDLT.h:16
Rotation2D(const Scalar &a)
Definition: Rotation2D.h:62
Rotation2D(const Rotation2D< OtherScalarType > &other)
Definition: Rotation2D.h:126
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:107
Rotation2D()
Definition: Rotation2D.h:65
bool isApprox(const Rotation2D &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Rotation2D.h:137
Scalar smallestAngle() const
Definition: Rotation2D.h:80
Rotation2D operator*(const Rotation2D &other) const
Definition: Rotation2D.h:91
Rotation2D inverse() const
Definition: Rotation2D.h:88
Definition: Eigen_Colamd.h:54
Rotation2D & operator*=(const Rotation2D &other)
Definition: Rotation2D.h:95
Represents a rotation/orientation in a 2 dimensional space.
Definition: ForwardDeclarations.h:265
Rotation2D< double > Rotation2Dd
Definition: Rotation2D.h:147
_Scalar Scalar
Definition: Rotation2D.h:51
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:178
Scalar angle() const
Definition: Rotation2D.h:68
Rotation2D< float > Rotation2Df
Definition: Rotation2D.h:144
Scalar & angle()
Definition: Rotation2D.h:71
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48