Medial Code Documentation
Loading...
Searching...
No Matches
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
13namespace Eigen {
14
32namespace internal {
33
34template<typename _Scalar> struct traits<Rotation2D<_Scalar> >
35{
36 typedef _Scalar Scalar;
37};
38} // end namespace internal
39
40template<typename _Scalar>
41class Rotation2D : public RotationBase<Rotation2D<_Scalar>,2>
42{
44
45public:
46
47 using Base::operator*;
48
49 enum { Dim = 2 };
51 typedef _Scalar Scalar;
54
55protected:
56
57 Scalar m_angle;
58
59public:
60
62 explicit inline Rotation2D(const Scalar& a) : m_angle(a) {}
63
66
71 template<typename Derived>
73 {
74 fromRotationMatrix(m.derived());
75 }
76
78 inline Scalar angle() const { return m_angle; }
79
81 inline Scalar& angle() { return m_angle; }
82
85 Scalar tmp = fmod(m_angle,Scalar(2)*EIGEN_PI);
86 return tmp<Scalar(0) ? tmp + Scalar(2)*EIGEN_PI : tmp;
87 }
88
90 inline Scalar smallestAngle() const {
91 Scalar tmp = fmod(m_angle,Scalar(2)*EIGEN_PI);
92 if(tmp>Scalar(EIGEN_PI)) tmp -= Scalar(2)*Scalar(EIGEN_PI);
93 else if(tmp<-Scalar(EIGEN_PI)) tmp += Scalar(2)*Scalar(EIGEN_PI);
94 return tmp;
95 }
96
98 inline Rotation2D inverse() const { return Rotation2D(-m_angle); }
99
101 inline Rotation2D operator*(const Rotation2D& other) const
102 { return Rotation2D(m_angle + other.m_angle); }
103
105 inline Rotation2D& operator*=(const Rotation2D& other)
106 { m_angle += other.m_angle; return *this; }
107
109 Vector2 operator* (const Vector2& vec) const
110 { return toRotationMatrix() * vec; }
111
112 template<typename Derived>
113 Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
114 Matrix2 toRotationMatrix() const;
115
123 template<typename Derived>
125 { return fromRotationMatrix(m.derived()); }
126
130 inline Rotation2D slerp(const Scalar& t, const Rotation2D& other) const
131 {
132 Scalar dist = Rotation2D(other.m_angle-m_angle).smallestAngle();
133 return Rotation2D(m_angle + dist*t);
134 }
135
141 template<typename NewScalarType>
144
146 template<typename OtherScalarType>
147 inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
148 {
149 m_angle = Scalar(other.angle());
150 }
151
152 static inline Rotation2D Identity() { return Rotation2D(0); }
153
158 bool isApprox(const Rotation2D& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
159 { return internal::isApprox(m_angle,other.m_angle, prec); }
160
161};
162
169
174template<typename Scalar>
175template<typename Derived>
177{
178 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));
181 return *this;
182}
183
186template<typename Scalar>
189{
190 using std::sin;
191 using std::cos;
192 Scalar sinA = sin(m_angle);
193 Scalar cosA = cos(m_angle);
194 return (Matrix2() << cosA, -sinA, sinA, cosA).finished();
195}
196
197} // end namespace Eigen
198
199#endif // EIGEN_ROTATION2D_H
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
\geometry_module
Definition Rotation2D.h:42
Scalar & angle()
Definition Rotation2D.h:81
Scalar smallestPositiveAngle() const
Definition Rotation2D.h:84
Rotation2D(const Rotation2D< OtherScalarType > &other)
Copy constructor with scalar type conversion.
Definition Rotation2D.h:147
Rotation2D()
Default constructor wihtout initialization.
Definition Rotation2D.h:65
Rotation2D & operator=(const MatrixBase< Derived > &m)
Set *this from a 2x2 rotation matrix mat.
Definition Rotation2D.h:124
Rotation2D(const MatrixBase< Derived > &m)
Construct a 2D rotation from a 2x2 rotation matrix mat.
Definition Rotation2D.h:72
Rotation2D inverse() const
Definition Rotation2D.h:98
Rotation2D & operator*=(const Rotation2D &other)
Concatenates two rotations.
Definition Rotation2D.h:105
Matrix2 toRotationMatrix() const
Constructs and.
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
the scalar type of the coefficients
Definition Rotation2D.h:51
Scalar smallestAngle() const
Definition Rotation2D.h:90
Rotation2D operator*(const Rotation2D &other) const
Concatenates two rotations.
Definition Rotation2D.h:101
Rotation2D(const Scalar &a)
Construct a 2D counter clock wise rotation from the angle a in radian.
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
Pseudo expression representing a solving operation.
Definition Solve.h:63
Holds information about the various numeric (i.e.
Definition NumTraits.h:108
Definition XprHelper.h:500
Definition ForwardDeclarations.h:17