10#ifndef EIGEN_ROTATION2D_H
11#define EIGEN_ROTATION2D_H
14#include "./InternalHeaderCheck.h"
37template <
typename Scalar_>
38struct traits<Rotation2D<Scalar_> > {
39 typedef Scalar_ Scalar;
43template <
typename Scalar_>
48 using Base::operator*;
70 template <
typename Derived>
72 fromRotationMatrix(m.
derived());
76 EIGEN_DEVICE_FUNC
inline Scalar angle()
const {
return m_angle; }
83 Scalar tmp = numext::fmod(m_angle,
Scalar(2 * EIGEN_PI));
84 return tmp <
Scalar(0) ? tmp +
Scalar(2 * EIGEN_PI) : tmp;
89 Scalar tmp = numext::fmod(m_angle,
Scalar(2 * EIGEN_PI));
90 if (tmp >
Scalar(EIGEN_PI))
91 tmp -=
Scalar(2 * EIGEN_PI);
92 else if (tmp < -
Scalar(EIGEN_PI))
93 tmp +=
Scalar(2 * EIGEN_PI);
107 m_angle += other.m_angle;
114 template <
typename Derived>
125 template <
typename Derived>
127 return fromRotationMatrix(m.
derived());
143 template <
typename NewScalarType>
144 EIGEN_DEVICE_FUNC
inline typename internal::cast_return_type<Rotation2D, Rotation2D<NewScalarType> >::type
cast()
146 return typename internal::cast_return_type<Rotation2D, Rotation2D<NewScalarType> >::type(*
this);
150 template <
typename OtherScalarType>
161 EIGEN_DEVICE_FUNC
bool isApprox(
const Rotation2D& other,
const typename NumTraits<Scalar>::Real& prec =
163 return internal::isApprox(m_angle, other.m_angle, prec);
178template <
typename Scalar>
179template <
typename Derived>
181 EIGEN_USING_STD(atan2)
182 EIGEN_STATIC_ASSERT(Derived::RowsAtCompileTime == 2 && Derived::ColsAtCompileTime == 2,
183 YOU_MADE_A_PROGRAMMING_MISTAKE)
184 m_angle = atan2(mat.
coeff(1, 0), mat.
coeff(0, 0));
190template <
typename Scalar>
194 Scalar sinA = sin(m_angle);
195 Scalar cosA = cos(m_angle);
196 return (
Matrix2() << cosA, -sinA, sinA, cosA).finished();
Derived & derived()
Definition EigenBase.h:49
EIGEN_CONSTEXPR CoeffReturnType coeff(Index row, Index col) const
Definition DenseCoeffsBase.h:92
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Represents a rotation/orientation in a 2 dimensional space.
Definition Rotation2D.h:44
Rotation2D inverse() const
Definition Rotation2D.h:98
Scalar smallestPositiveAngle() const
Definition Rotation2D.h:82
Scalar smallestAngle() const
Definition Rotation2D.h:88
Rotation2D(const Scalar &a)
Definition Rotation2D.h:61
Rotation2D & operator=(const MatrixBase< Derived > &m)
Definition Rotation2D.h:126
Vector2 operator*(const Vector2 &vec) const
Definition Rotation2D.h:112
Scalar angle() const
Definition Rotation2D.h:76
Matrix2 toRotationMatrix() const
Definition Rotation2D.h:191
Rotation2D operator*(const Rotation2D &other) const
Definition Rotation2D.h:101
internal::cast_return_type< Rotation2D, Rotation2D< NewScalarType > >::type cast() const
Definition Rotation2D.h:144
Rotation2D(const Rotation2D< OtherScalarType > &other)
Definition Rotation2D.h:151
Rotation2D()
Definition Rotation2D.h:64
Rotation2D & operator*=(const Rotation2D &other)
Definition Rotation2D.h:106
Scalar & angle()
Definition Rotation2D.h:79
Rotation2D slerp(const Scalar &t, const Rotation2D &other) const
Definition Rotation2D.h:133
bool isApprox(const Rotation2D &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition Rotation2D.h:161
Scalar_ Scalar
Definition Rotation2D.h:52
Rotation2D(const MatrixBase< Derived > &m)
Definition Rotation2D.h:71
Common base class for compact rotation representations.
Definition RotationBase.h:32
Namespace containing all symbols from the Eigen library.
Definition Core:137
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523