Rotation2D.h
Go to the documentation of this file.
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 {
44 
45 public:
46 
47  using Base::operator*;
48 
49  enum { Dim = 2 };
51  typedef _Scalar Scalar;
54 
55 protected:
56 
58 
59 public:
60 
62  EIGEN_DEVICE_FUNC explicit inline Rotation2D(const Scalar& a) : m_angle(a) {}
63 
65  EIGEN_DEVICE_FUNC Rotation2D() {}
66 
71  template<typename Derived>
72  EIGEN_DEVICE_FUNC explicit Rotation2D(const MatrixBase<Derived>& m)
73  {
74  fromRotationMatrix(m.derived());
75  }
76 
78  EIGEN_DEVICE_FUNC inline Scalar angle() const { return m_angle; }
79 
81  EIGEN_DEVICE_FUNC inline Scalar& angle() { return m_angle; }
82 
84  EIGEN_DEVICE_FUNC inline Scalar smallestPositiveAngle() const {
86  return tmp<Scalar(0) ? tmp + Scalar(2*EIGEN_PI) : tmp;
87  }
88 
90  EIGEN_DEVICE_FUNC inline Scalar smallestAngle() const {
92  if(tmp>Scalar(EIGEN_PI)) tmp -= Scalar(2*EIGEN_PI);
93  else if(tmp<-Scalar(EIGEN_PI)) tmp += Scalar(2*EIGEN_PI);
94  return tmp;
95  }
96 
98  EIGEN_DEVICE_FUNC inline Rotation2D inverse() const { return Rotation2D(-m_angle); }
99 
101  EIGEN_DEVICE_FUNC inline Rotation2D operator*(const Rotation2D& other) const
102  { return Rotation2D(m_angle + other.m_angle); }
103 
105  EIGEN_DEVICE_FUNC inline Rotation2D& operator*=(const Rotation2D& other)
106  { m_angle += other.m_angle; return *this; }
107 
109  EIGEN_DEVICE_FUNC Vector2 operator* (const Vector2& vec) const
110  { return toRotationMatrix() * vec; }
111 
112  template<typename Derived>
113  EIGEN_DEVICE_FUNC Rotation2D& fromRotationMatrix(const MatrixBase<Derived>& m);
114  EIGEN_DEVICE_FUNC Matrix2 toRotationMatrix() const;
115 
123  template<typename Derived>
124  EIGEN_DEVICE_FUNC Rotation2D& operator=(const MatrixBase<Derived>& m)
125  { return fromRotationMatrix(m.derived()); }
126 
130  EIGEN_DEVICE_FUNC 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>
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); }
144 
146  template<typename OtherScalarType>
147  EIGEN_DEVICE_FUNC inline explicit Rotation2D(const Rotation2D<OtherScalarType>& other)
148  {
149  m_angle = Scalar(other.angle());
150  }
151 
152  EIGEN_DEVICE_FUNC static inline Rotation2D Identity() { return Rotation2D(0); }
153 
158  EIGEN_DEVICE_FUNC 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 
174 template<typename Scalar>
175 template<typename Derived>
177 {
178  EIGEN_USING_STD_MATH(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 
186 template<typename Scalar>
188 EIGEN_DEVICE_FUNC Rotation2D<Scalar>::toRotationMatrix(void) const
189 {
190  EIGEN_USING_STD_MATH(sin)
191  EIGEN_USING_STD_MATH(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
sin
const EIGEN_DEVICE_FUNC SinReturnType sin() const
Definition: ArrayCwiseUnaryOps.h:220
Eigen
Definition: common.h:73
Eigen::internal::cast_return_type
Definition: XprHelper.h:489
Eigen::atan2
const AutoDiffScalar< Matrix< typename internal::traits< typename internal::remove_all< DerTypeA >::type >::Scalar, Dynamic, 1 > > atan2(const AutoDiffScalar< DerTypeA > &a, const AutoDiffScalar< DerTypeB > &b)
Definition: AutoDiffScalar.h:623
EIGEN_PI
#define EIGEN_PI
Definition: Eigen/src/Core/MathFunctions.h:15
Eigen::Rotation2D::Rotation2D
EIGEN_DEVICE_FUNC Rotation2D()
Definition: Rotation2D.h:65
Eigen::Rotation2D::Rotation2D
EIGEN_DEVICE_FUNC Rotation2D(const Rotation2D< OtherScalarType > &other)
Definition: Rotation2D.h:147
Eigen::Rotation2D::operator=
EIGEN_DEVICE_FUNC Rotation2D & operator=(const MatrixBase< Derived > &m)
Definition: Rotation2D.h:124
Eigen::internal::isApprox
EIGEN_DEVICE_FUNC bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Definition: Eigen/src/Core/MathFunctions.h:1361
Eigen::internal::traits< Rotation2D< _Scalar > >::Scalar
_Scalar Scalar
Definition: Rotation2D.h:36
Eigen::Rotation2Df
Rotation2D< float > Rotation2Df
Definition: Rotation2D.h:165
Eigen::Rotation2D::fromRotationMatrix
EIGEN_DEVICE_FUNC Rotation2D & fromRotationMatrix(const MatrixBase< Derived > &m)
Eigen::Rotation2D::Vector2
Matrix< Scalar, 2, 1 > Vector2
Definition: Rotation2D.h:52
Eigen::Rotation2D::inverse
EIGEN_DEVICE_FUNC Rotation2D inverse() const
Definition: Rotation2D.h:98
Eigen::Rotation2D::Identity
static EIGEN_DEVICE_FUNC Rotation2D Identity()
Definition: Rotation2D.h:152
Eigen::Rotation2D::operator*
EIGEN_DEVICE_FUNC Rotation2D operator*(const Rotation2D &other) const
Definition: Rotation2D.h:101
Eigen::Rotation2D::smallestAngle
EIGEN_DEVICE_FUNC Scalar smallestAngle() const
Definition: Rotation2D.h:90
Eigen::Rotation2D::angle
EIGEN_DEVICE_FUNC Scalar angle() const
Definition: Rotation2D.h:78
Eigen::RotationBase
Common base class for compact rotation representations.
Definition: ForwardDeclarations.h:266
Eigen::Rotation2Dd
Rotation2D< double > Rotation2Dd
Definition: Rotation2D.h:168
Eigen::Rotation2D::Rotation2D
EIGEN_DEVICE_FUNC Rotation2D(const Scalar &a)
Definition: Rotation2D.h:62
Eigen::numext::fmod
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T fmod(const T &a, const T &b)
Definition: Eigen/src/Core/MathFunctions.h:1244
Eigen::Rotation2D::Rotation2D
EIGEN_DEVICE_FUNC Rotation2D(const MatrixBase< Derived > &m)
Definition: Rotation2D.h:72
Eigen::Rotation2D::isApprox
EIGEN_DEVICE_FUNC bool isApprox(const Rotation2D &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Rotation2D.h:158
Eigen::Rotation2D::Scalar
_Scalar Scalar
Definition: Rotation2D.h:51
Eigen::Rotation2D
Represents a rotation/orientation in a 2 dimensional space.
Definition: ForwardDeclarations.h:269
Eigen::Rotation2D::cast
EIGEN_DEVICE_FUNC internal::cast_return_type< Rotation2D, Rotation2D< NewScalarType > >::type cast() const
Definition: Rotation2D.h:142
Eigen::internal::traits
Definition: ForwardDeclarations.h:17
EIGEN_STATIC_ASSERT
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:124
Eigen::Rotation2D::smallestPositiveAngle
EIGEN_DEVICE_FUNC Scalar smallestPositiveAngle() const
Definition: Rotation2D.h:84
Eigen::Rotation2D::Matrix2
Matrix< Scalar, 2, 2 > Matrix2
Definition: Rotation2D.h:53
a
Scalar * a
Definition: cholesky.cpp:26
Eigen::Rotation2D::Dim
@ Dim
Definition: Rotation2D.h:49
Eigen::Rotation2D::operator*=
EIGEN_DEVICE_FUNC Rotation2D & operator*=(const Rotation2D &other)
Definition: Rotation2D.h:105
Eigen::Rotation2D::Base
RotationBase< Rotation2D< _Scalar >, 2 > Base
Definition: Rotation2D.h:43
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:178
Eigen::Rotation2D::toRotationMatrix
EIGEN_DEVICE_FUNC Matrix2 toRotationMatrix() const
Definition: Rotation2D.h:188
internal
Definition: BandTriangularSolver.h:13
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Eigen::Rotation2D::slerp
EIGEN_DEVICE_FUNC Rotation2D slerp(const Scalar &t, const Rotation2D &other) const
Definition: Rotation2D.h:130
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:150
cos
const EIGEN_DEVICE_FUNC CosReturnType cos() const
Definition: ArrayCwiseUnaryOps.h:202
Eigen::Rotation2D::m_angle
Scalar m_angle
Definition: Rotation2D.h:57
Eigen::Rotation2D::angle
EIGEN_DEVICE_FUNC Scalar & angle()
Definition: Rotation2D.h:81
mat
else mat
Definition: eigenvalues.cpp:43


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Wed Mar 2 2022 00:06:10