Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
AngleAxis.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008 Gael Guennebaud <[email protected]>
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_ANGLEAXIS_H
11#define EIGEN_ANGLEAXIS_H
12
13// IWYU pragma: private
14#include "./InternalHeaderCheck.h"
15
16namespace Eigen {
17
44namespace internal {
45template <typename Scalar_>
46struct traits<AngleAxis<Scalar_> > {
47 typedef Scalar_ Scalar;
48};
49} // namespace internal
50
51template <typename Scalar_>
52class AngleAxis : public RotationBase<AngleAxis<Scalar_>, 3> {
54
55 public:
56 using Base::operator*;
57
58 enum { Dim = 3 };
60 typedef Scalar_ Scalar;
64
65 protected:
66 Vector3 m_axis;
67 Scalar m_angle;
68
69 public:
71 EIGEN_DEVICE_FUNC AngleAxis() {}
77 template <typename Derived>
78 EIGEN_DEVICE_FUNC inline AngleAxis(const Scalar& angle, const MatrixBase<Derived>& axis)
79 : m_axis(axis), m_angle(angle) {}
83 template <typename QuatDerived>
84 EIGEN_DEVICE_FUNC inline explicit AngleAxis(const QuaternionBase<QuatDerived>& q) {
85 *this = q;
86 }
88 template <typename Derived>
89 EIGEN_DEVICE_FUNC inline explicit AngleAxis(const MatrixBase<Derived>& m) {
90 *this = m;
91 }
92
94 EIGEN_DEVICE_FUNC Scalar angle() const { return m_angle; }
96 EIGEN_DEVICE_FUNC Scalar& angle() { return m_angle; }
97
99 EIGEN_DEVICE_FUNC const Vector3& axis() const { return m_axis; }
104 EIGEN_DEVICE_FUNC Vector3& axis() { return m_axis; }
105
107 EIGEN_DEVICE_FUNC inline QuaternionType operator*(const AngleAxis& other) const {
108 return QuaternionType(*this) * QuaternionType(other);
109 }
110
112 EIGEN_DEVICE_FUNC inline QuaternionType operator*(const QuaternionType& other) const {
113 return QuaternionType(*this) * other;
114 }
115
117 friend EIGEN_DEVICE_FUNC inline QuaternionType operator*(const QuaternionType& a, const AngleAxis& b) {
118 return a * QuaternionType(b);
119 }
120
122 EIGEN_DEVICE_FUNC AngleAxis inverse() const { return AngleAxis(-m_angle, m_axis); }
123
124 template <class QuatDerived>
125 EIGEN_DEVICE_FUNC AngleAxis& operator=(const QuaternionBase<QuatDerived>& q);
126 template <typename Derived>
127 EIGEN_DEVICE_FUNC AngleAxis& operator=(const MatrixBase<Derived>& m);
128
129 template <typename Derived>
130 EIGEN_DEVICE_FUNC AngleAxis& fromRotationMatrix(const MatrixBase<Derived>& m);
131 EIGEN_DEVICE_FUNC Matrix3 toRotationMatrix(void) const;
132
138 template <typename NewScalarType>
139 EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<AngleAxis, AngleAxis<NewScalarType> >::type cast()
140 const {
141 return typename internal::cast_return_type<AngleAxis, AngleAxis<NewScalarType> >::type(*this);
142 }
143
145 template <typename OtherScalarType>
146 EIGEN_DEVICE_FUNC inline explicit AngleAxis(const AngleAxis<OtherScalarType>& other) {
147 m_axis = other.axis().template cast<Scalar>();
148 m_angle = Scalar(other.angle());
149 }
150
151 EIGEN_DEVICE_FUNC static inline const AngleAxis Identity() { return AngleAxis(Scalar(0), Vector3::UnitX()); }
152
157 EIGEN_DEVICE_FUNC bool isApprox(const AngleAxis& other, const typename NumTraits<Scalar>::Real& prec =
159 return m_axis.isApprox(other.m_axis, prec) && internal::isApprox(m_angle, other.m_angle, prec);
160 }
161};
162
169
176template <typename Scalar>
177template <typename QuatDerived>
179 EIGEN_USING_STD(atan2)
180 EIGEN_USING_STD(abs)
181 Scalar n = q.vec().norm();
182 if (n < NumTraits<Scalar>::epsilon()) n = q.vec().stableNorm();
183
184 if (n != Scalar(0)) {
185 m_angle = Scalar(2) * atan2(n, abs(q.w()));
186 if (q.w() < Scalar(0)) n = -n;
187 m_axis = q.vec() / n;
188 } else {
189 m_angle = Scalar(0);
190 m_axis << Scalar(1), Scalar(0), Scalar(0);
191 }
192 return *this;
193}
194
197template <typename Scalar>
198template <typename Derived>
200 // Since a direct conversion would not be really faster,
201 // let's use the robust Quaternion implementation:
202 return *this = QuaternionType(mat);
203}
204
208template <typename Scalar>
209template <typename Derived>
211 return *this = QuaternionType(mat);
212}
213
216template <typename Scalar>
218 EIGEN_USING_STD(sin)
219 EIGEN_USING_STD(cos)
220 Matrix3 res;
221 Vector3 sin_axis = sin(m_angle) * m_axis;
222 Scalar c = cos(m_angle);
223 Vector3 cos1_axis = (Scalar(1) - c) * m_axis;
224
225 Scalar tmp;
226 tmp = cos1_axis.x() * m_axis.y();
227 res.coeffRef(0, 1) = tmp - sin_axis.z();
228 res.coeffRef(1, 0) = tmp + sin_axis.z();
229
230 tmp = cos1_axis.x() * m_axis.z();
231 res.coeffRef(0, 2) = tmp + sin_axis.y();
232 res.coeffRef(2, 0) = tmp - sin_axis.y();
233
234 tmp = cos1_axis.y() * m_axis.z();
235 res.coeffRef(1, 2) = tmp - sin_axis.x();
236 res.coeffRef(2, 1) = tmp + sin_axis.x();
237
238 res.diagonal() = (cos1_axis.cwiseProduct(m_axis)).array() + c;
239
240 return res;
241}
242
243} // end namespace Eigen
244
245#endif // EIGEN_ANGLEAXIS_H
Represents a 3D rotation as a rotation angle around an arbitrary 3D axis.
Definition AngleAxis.h:52
QuaternionType operator*(const QuaternionType &other) const
Definition AngleAxis.h:112
AngleAxis(const AngleAxis< OtherScalarType > &other)
Definition AngleAxis.h:146
Scalar_ Scalar
Definition AngleAxis.h:60
QuaternionType operator*(const AngleAxis &other) const
Definition AngleAxis.h:107
Scalar & angle()
Definition AngleAxis.h:96
Vector3 & axis()
Definition AngleAxis.h:104
AngleAxis(const QuaternionBase< QuatDerived > &q)
Definition AngleAxis.h:84
internal::cast_return_type< AngleAxis, AngleAxis< NewScalarType > >::type cast() const
Definition AngleAxis.h:139
AngleAxis(const Scalar &angle, const MatrixBase< Derived > &axis)
Definition AngleAxis.h:78
AngleAxis(const MatrixBase< Derived > &m)
Definition AngleAxis.h:89
AngleAxis inverse() const
Definition AngleAxis.h:122
const Vector3 & axis() const
Definition AngleAxis.h:99
AngleAxis()
Definition AngleAxis.h:71
bool isApprox(const AngleAxis &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
Definition AngleAxis.h:157
friend QuaternionType operator*(const QuaternionType &a, const AngleAxis &b)
Definition AngleAxis.h:117
Matrix3 toRotationMatrix(void) const
Definition AngleAxis.h:217
Scalar angle() const
Definition AngleAxis.h:94
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
static const BasisReturnType UnitX()
Definition CwiseNullaryOp.h:895
DiagonalReturnType diagonal()
Definition Diagonal.h:162
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr Scalar & coeffRef(Index rowId, Index colId)
Definition PlainObjectBase.h:217
Base class for quaternion expressions.
Definition Quaternion.h:35
const VectorBlock< const Coefficients, 3 > vec() const
Definition Quaternion.h:78
EIGEN_CONSTEXPR CoeffReturnType w() const
Definition Quaternion.h:66
The quaternion class used to represent 3D orientations and rotations.
Definition Quaternion.h:285
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