aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/Geometry1
-rw-r--r--Eigen/src/Core/MatrixBase.h1
-rw-r--r--Eigen/src/Geometry/EulerAngles.h100
-rw-r--r--test/geometry.cpp22
4 files changed, 124 insertions, 0 deletions
diff --git a/Eigen/Geometry b/Eigen/Geometry
index 19be347fb..95a3e9bed 100644
--- a/Eigen/Geometry
+++ b/Eigen/Geometry
@@ -29,6 +29,7 @@ namespace Eigen {
#include "src/Geometry/Rotation2D.h"
#include "src/Geometry/Quaternion.h"
#include "src/Geometry/AngleAxis.h"
+#include "src/Geometry/EulerAngles.h"
#include "src/Geometry/Transform.h"
#include "src/Geometry/Translation.h"
#include "src/Geometry/Scaling.h"
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index 53090b984..19ab4040c 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -593,6 +593,7 @@ template<typename Derived> class MatrixBase
template<typename OtherDerived>
EvalType cross(const MatrixBase<OtherDerived>& other) const;
EvalType unitOrthogonal(void) const;
+ Matrix<Scalar,3,1> eulerAngles(int a0, int a1, int a2) const;
#ifdef EIGEN_MATRIXBASE_PLUGIN
#include EIGEN_MATRIXBASE_PLUGIN
diff --git a/Eigen/src/Geometry/EulerAngles.h b/Eigen/src/Geometry/EulerAngles.h
new file mode 100644
index 000000000..6087a1433
--- /dev/null
+++ b/Eigen/src/Geometry/EulerAngles.h
@@ -0,0 +1,100 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra. Eigen itself is part of the KDE project.
+//
+// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#ifndef EIGEN_EULERANGLES_H
+#define EIGEN_EULERANGLES_H
+
+/** \geometry_module \ingroup GeometryModule
+ * \nonstableyet
+ *
+ * \returns the Euler-angles of the rotation matrix \c *this using the convention defined by the triplet (\a a0,\a a1,\a a2)
+ *
+ * Each of the three parameters \a a0,\a a1,\a a2 represents the respective rotation axis as an integer in {0,1,2}.
+ * For instance, in:
+ * \code Vector3f ea = mat.eulerAngles(2, 0, 2); \endcode
+ * "2" represents the z axis and "0" the x axis, etc. The returned angles are such that
+ * we have the following equality:
+ * \code
+ * mat == AngleAxisf(ea[0], Vector3f::UnitZ())
+ * * AngleAxisf(ea[1], Vector3f::UnitX())
+ * * AngleAxisf(ea[2], Vector3f::UnitZ()); \endcode
+ * This corresponds to the right-multiply conventions (with right hand side frames).
+ */
+// FIXME perhaps the triplet could be template parameters
+// and/or packed into constants: EulerXYZ, EulerXYX, etc....
+// FIXME should we support the reversed conventions ? (left multiply)
+template<typename Derived>
+inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
+MatrixBase<Derived>::eulerAngles(int a0, int a1, int a2) const
+{
+ EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived,3,3)
+
+ Matrix<Scalar,3,1> res;
+ typedef Matrix<typename Derived::Scalar,2,1> Vector2;
+ const Scalar epsilon = precision<Scalar>();
+
+ const int odd = ((a0+1)%3 == a1) ? 0 : 1;
+ const int i = a0;
+ const int j = (a0 + 1 + odd)%3;
+ const int k = (a0 + 2 - odd)%3;
+
+ if (a0==a2)
+ {
+ Scalar s = Vector2(coeff(j,i) , coeff(k,i)).norm();
+
+ res[1] = std::atan2(s, coeff(i,i));
+ if (s > epsilon)
+ {
+ res[0] = std::atan2(coeff(j,i), coeff(k,i));
+ res[2] = std::atan2(coeff(i,j),-coeff(i,k));
+ }
+ else
+ {
+ res[0] = Scalar(0);
+ res[2] = std::atan2(-coeff(k,j), coeff(j,j));
+ }
+ }
+ else
+ {
+ Scalar c = Vector2(coeff(i,i) , coeff(i,j)).norm();
+
+ res[1] = std::atan2(-coeff(i,k), c);
+ if (c > epsilon)
+ {
+ res[0] = std::atan2(coeff(j,k), coeff(k,k));
+ res[2] = std::atan2(coeff(i,j), coeff(i,i));
+ }
+ else
+ {
+ res[0] = Scalar(0);
+ res[2] = -std::atan2(-coeff(k,j), coeff(j,j));
+ }
+ }
+ if (!odd)
+ res = -res;
+ return res;
+}
+
+
+#endif // EIGEN_EULERANGLES_H
diff --git a/test/geometry.cpp b/test/geometry.cpp
index b6899501c..af23bce82 100644
--- a/test/geometry.cpp
+++ b/test/geometry.cpp
@@ -318,6 +318,28 @@ template<typename Scalar> void geometry(void)
VERIFY_IS_APPROX(r2d1f.template cast<Scalar>(),r2d1);
Rotation2D<double> r2d1d = r2d1.template cast<double>();
VERIFY_IS_APPROX(r2d1d.template cast<Scalar>(),r2d1);
+
+ m = q1;
+
+ #define VERIFY_EULER(I,J,K, X,Y,Z) { \
+ Vector3 ea = m.eulerAngles(I,J,K); \
+ Matrix3 m1 = Matrix3(AngleAxisx(ea[0], Vector3::Unit##X()) * AngleAxisx(ea[1], Vector3::Unit##Y()) * AngleAxisx(ea[2], Vector3::Unit##Z())); \
+ VERIFY_IS_APPROX(m, Matrix3(AngleAxisx(ea[0], Vector3::Unit##X()) * AngleAxisx(ea[1], Vector3::Unit##Y()) * AngleAxisx(ea[2], Vector3::Unit##Z()))); \
+ }
+ VERIFY_EULER(0,1,2, X,Y,Z);
+ VERIFY_EULER(0,1,0, X,Y,X);
+ VERIFY_EULER(0,2,1, X,Z,Y);
+ VERIFY_EULER(0,2,0, X,Z,X);
+
+ VERIFY_EULER(1,2,0, Y,Z,X);
+ VERIFY_EULER(1,2,1, Y,Z,Y);
+ VERIFY_EULER(1,0,2, Y,X,Z);
+ VERIFY_EULER(1,0,1, Y,X,Y);
+
+ VERIFY_EULER(2,0,1, Z,X,Y);
+ VERIFY_EULER(2,0,2, Z,X,Z);
+ VERIFY_EULER(2,1,0, Z,Y,X);
+ VERIFY_EULER(2,1,2, Z,Y,Z);
}
void test_geometry()