mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-06-04 18:54:00 +08:00
Geometry/EulerAngles: introduce canonicalEulerAngles
This commit is contained in:
parent
7d9bb90f15
commit
c18f94e3b0
@ -399,9 +399,12 @@ template<typename Derived> class MatrixBase
|
||||
EIGEN_DEVICE_FUNC
|
||||
inline PlainObject unitOrthogonal(void) const;
|
||||
|
||||
EIGEN_DEVICE_FUNC
|
||||
EIGEN_DEPRECATED EIGEN_DEVICE_FUNC
|
||||
inline Matrix<Scalar,3,1> eulerAngles(Index a0, Index a1, Index a2) const;
|
||||
|
||||
EIGEN_DEVICE_FUNC
|
||||
inline Matrix<Scalar,3,1> canonicalEulerAngles(Index a0, Index a1, Index a2) const;
|
||||
|
||||
// put this as separate enum value to work around possible GCC 4.3 bug (?)
|
||||
enum { HomogeneousReturnTypeDirection = ColsAtCompileTime==1&&RowsAtCompileTime==1 ? ((internal::traits<Derived>::Flags&RowMajorBit)==RowMajorBit ? Horizontal : Vertical)
|
||||
: ColsAtCompileTime==1 ? Vertical : Horizontal };
|
||||
|
@ -2,6 +2,7 @@
|
||||
// for linear algebra.
|
||||
//
|
||||
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
|
||||
// Copyright (C) 2023 Juraj Oršulić, University of Zagreb <juraj.orsulic@fer.hr>
|
||||
//
|
||||
// This Source Code Form is subject to the terms of the Mozilla
|
||||
// Public License v. 2.0. If a copy of the MPL was not distributed
|
||||
@ -12,12 +13,12 @@
|
||||
|
||||
#include "./InternalHeaderCheck.h"
|
||||
|
||||
namespace Eigen {
|
||||
namespace Eigen {
|
||||
|
||||
/** \geometry_module \ingroup Geometry_Module
|
||||
*
|
||||
*
|
||||
* \returns the Euler-angles of the rotation matrix \c *this using the convention defined by the triplet (\a a0,\a a1,\a a2)
|
||||
* \returns the canonical 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:
|
||||
@ -29,85 +30,188 @@ namespace Eigen {
|
||||
* * AngleAxisf(ea[1], Vector3f::UnitX())
|
||||
* * AngleAxisf(ea[2], Vector3f::UnitZ()); \endcode
|
||||
* This corresponds to the right-multiply conventions (with right hand side frames).
|
||||
*
|
||||
* The returned angles are in the ranges [0:pi]x[-pi:pi]x[-pi:pi].
|
||||
*
|
||||
*
|
||||
* For Tait-Bryan angle configurations (a0 != a2), the returned angles are in the ranges [-pi:pi]x[-pi/2:pi/2]x[-pi:pi].
|
||||
* For proper Euler angle configurations (a0 == a2), the returned angles are in the ranges [-pi:pi]x[0:pi]x[-pi:pi].
|
||||
*
|
||||
* \sa class AngleAxis
|
||||
*/
|
||||
template<typename Derived>
|
||||
EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
|
||||
MatrixBase<Derived>::eulerAngles(Index a0, Index a1, Index a2) const
|
||||
MatrixBase<Derived>::canonicalEulerAngles(Index a0, Index a1, Index a2) const
|
||||
{
|
||||
EIGEN_USING_STD(atan2)
|
||||
EIGEN_USING_STD(sin)
|
||||
EIGEN_USING_STD(cos)
|
||||
/* Implemented from Graphics Gems IV */
|
||||
EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived,3,3)
|
||||
EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived, 3, 3)
|
||||
|
||||
Matrix<Scalar,3,1> res;
|
||||
typedef Matrix<typename Derived::Scalar,2,1> Vector2;
|
||||
Matrix<Scalar, 3, 1> res;
|
||||
|
||||
const Index odd = ((a0+1)%3 == a1) ? 0 : 1;
|
||||
const Index odd = ((a0 + 1) % 3 == a1) ? 0 : 1;
|
||||
const Index i = a0;
|
||||
const Index j = (a0 + 1 + odd)%3;
|
||||
const Index k = (a0 + 2 - odd)%3;
|
||||
|
||||
if (a0==a2)
|
||||
const Index j = (a0 + 1 + odd) % 3;
|
||||
const Index k = (a0 + 2 - odd) % 3;
|
||||
|
||||
if (a0 == a2)
|
||||
{
|
||||
res[0] = atan2(coeff(j,i), coeff(k,i));
|
||||
if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0)))
|
||||
// Proper Euler angles (same first and last axis).
|
||||
// The i, j, k indices enable addressing the input matrix as the XYX archetype matrix (see Graphics Gems IV),
|
||||
// where e.g. coeff(k, i) means third column, first row in the XYX archetype matrix:
|
||||
// c2 s2s1 s2c1
|
||||
// s2s3 -c2s1s3 + c1c3 -c2c1s3 - s1c3
|
||||
// -s2c3 c2s1c3 + c1s3 c2c1c3 - s1s3
|
||||
|
||||
// Note: s2 is always positive.
|
||||
Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
|
||||
if (odd)
|
||||
{
|
||||
if(res[0] > Scalar(0)) {
|
||||
res[0] -= Scalar(EIGEN_PI);
|
||||
}
|
||||
else {
|
||||
res[0] += Scalar(EIGEN_PI);
|
||||
}
|
||||
Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
|
||||
res[1] = -atan2(s2, coeff(i,i));
|
||||
res[0] = numext::atan2(coeff(j, i), coeff(k, i));
|
||||
// s2 is always positive, so res[1] will be within the canonical [0, pi] range
|
||||
res[1] = numext::atan2(s2, coeff(i, i));
|
||||
}
|
||||
else
|
||||
{
|
||||
Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
|
||||
res[1] = atan2(s2, coeff(i,i));
|
||||
// In the !odd case, signs of all three angles are flipped at the very end. To keep the solution within the canonical range,
|
||||
// we flip the solution and make res[1] always negative here (since s2 is always positive, -atan2(s2, c2) will always be negative).
|
||||
// The final flip at the end due to !odd will thus make res[1] positive and canonical.
|
||||
// NB: in the general case, there are two correct solutions, but only one is canonical. For proper Euler angles,
|
||||
// flipping from one solution to the other involves flipping the sign of the second angle res[1] and adding/subtracting pi
|
||||
// to the first and third angles. The addition/subtraction of pi to the first angle res[0] is handled here by flipping
|
||||
// the signs of arguments to atan2, while the calculation of the third angle does not need special adjustment since
|
||||
// it uses the adjusted res[0] as the input and produces a correct result.
|
||||
res[0] = numext::atan2(-coeff(j, i), -coeff(k, i));
|
||||
res[1] = -numext::atan2(s2, coeff(i, i));
|
||||
}
|
||||
|
||||
|
||||
// With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
|
||||
// we can compute their respective rotation, and apply its inverse to M. Since the result must
|
||||
// be a rotation around x, we have:
|
||||
//
|
||||
// c2 s1.s2 c1.s2 1 0 0
|
||||
// c2 s1.s2 c1.s2 1 0 0
|
||||
// 0 c1 -s1 * M = 0 c3 s3
|
||||
// -s2 s1.c2 c1.c2 0 -s3 c3
|
||||
//
|
||||
// Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
|
||||
|
||||
Scalar s1 = sin(res[0]);
|
||||
Scalar c1 = cos(res[0]);
|
||||
res[2] = atan2(c1*coeff(j,k)-s1*coeff(k,k), c1*coeff(j,j) - s1 * coeff(k,j));
|
||||
}
|
||||
|
||||
Scalar s1 = numext::sin(res[0]);
|
||||
Scalar c1 = numext::cos(res[0]);
|
||||
res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
|
||||
}
|
||||
else
|
||||
{
|
||||
res[0] = atan2(coeff(j,k), coeff(k,k));
|
||||
Scalar c2 = Vector2(coeff(i,i), coeff(i,j)).norm();
|
||||
if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0))) {
|
||||
if(res[0] > Scalar(0)) {
|
||||
res[0] -= Scalar(EIGEN_PI);
|
||||
}
|
||||
else {
|
||||
res[0] += Scalar(EIGEN_PI);
|
||||
}
|
||||
res[1] = atan2(-coeff(i,k), -c2);
|
||||
}
|
||||
else
|
||||
res[1] = atan2(-coeff(i,k), c2);
|
||||
Scalar s1 = sin(res[0]);
|
||||
Scalar c1 = cos(res[0]);
|
||||
res[2] = atan2(s1*coeff(k,i)-c1*coeff(j,i), c1*coeff(j,j) - s1 * coeff(k,j));
|
||||
// Tait-Bryan angles (all three axes are different; typically used for yaw-pitch-roll calculations).
|
||||
// The i, j, k indices enable addressing the input matrix as the XYZ archetype matrix (see Graphics Gems IV),
|
||||
// where e.g. coeff(k, i) means third column, first row in the XYZ archetype matrix:
|
||||
// c2c3 s2s1c3 - c1s3 s2c1c3 + s1s3
|
||||
// c2s3 s2s1s3 + c1c3 s2c1s3 - s1c3
|
||||
// -s2 c2s1 c2c1
|
||||
|
||||
res[0] = numext::atan2(coeff(j, k), coeff(k, k));
|
||||
|
||||
Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
|
||||
// c2 is always positive, so the following atan2 will always return a result in the correct canonical middle angle range [-pi/2, pi/2]
|
||||
res[1] = numext::atan2(-coeff(i, k), c2);
|
||||
|
||||
Scalar s1 = numext::sin(res[0]);
|
||||
Scalar c1 = numext::cos(res[0]);
|
||||
res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
|
||||
}
|
||||
if (!odd)
|
||||
{
|
||||
res = -res;
|
||||
|
||||
}
|
||||
|
||||
return res;
|
||||
}
|
||||
|
||||
/** \geometry_module \ingroup Geometry_Module
|
||||
*
|
||||
*
|
||||
* \returns the Euler-angles of the rotation matrix \c *this using the convention defined by the triplet (\a a0,\a a1,\a a2)
|
||||
*
|
||||
* NB: The returned angles are in non-canonical ranges [0:pi]x[-pi:pi]x[-pi:pi]. For canonical Tait-Bryan/proper Euler ranges, use canonicalEulerAngles.
|
||||
*
|
||||
* \sa MatrixBase::canonicalEulerAngles
|
||||
* \sa class AngleAxis
|
||||
*/
|
||||
template<typename Derived>
|
||||
EIGEN_DEPRECATED EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
|
||||
MatrixBase<Derived>::eulerAngles(Index a0, Index a1, Index a2) const
|
||||
{
|
||||
/* Implemented from Graphics Gems IV */
|
||||
EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived, 3, 3)
|
||||
|
||||
Matrix<Scalar, 3, 1> res;
|
||||
|
||||
const Index odd = ((a0 + 1) % 3 == a1) ? 0 : 1;
|
||||
const Index i = a0;
|
||||
const Index j = (a0 + 1 + odd) % 3;
|
||||
const Index k = (a0 + 2 - odd) % 3;
|
||||
|
||||
if (a0 == a2)
|
||||
{
|
||||
res[0] = numext::atan2(coeff(j, i), coeff(k, i));
|
||||
if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0)))
|
||||
{
|
||||
if (res[0] > Scalar(0))
|
||||
{
|
||||
res[0] -= Scalar(EIGEN_PI);
|
||||
}
|
||||
else
|
||||
{
|
||||
res[0] += Scalar(EIGEN_PI);
|
||||
}
|
||||
|
||||
Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
|
||||
res[1] = -numext::atan2(s2, coeff(i, i));
|
||||
}
|
||||
else
|
||||
{
|
||||
Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
|
||||
res[1] = numext::atan2(s2, coeff(i, i));
|
||||
}
|
||||
|
||||
// With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
|
||||
// we can compute their respective rotation, and apply its inverse to M. Since the result must
|
||||
// be a rotation around x, we have:
|
||||
//
|
||||
// c2 s1.s2 c1.s2 1 0 0
|
||||
// 0 c1 -s1 * M = 0 c3 s3
|
||||
// -s2 s1.c2 c1.c2 0 -s3 c3
|
||||
//
|
||||
// Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
|
||||
|
||||
Scalar s1 = numext::sin(res[0]);
|
||||
Scalar c1 = numext::cos(res[0]);
|
||||
res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
|
||||
}
|
||||
else
|
||||
{
|
||||
res[0] = numext::atan2(coeff(j, k), coeff(k, k));
|
||||
Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
|
||||
if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0)))
|
||||
{
|
||||
if (res[0] > Scalar(0))
|
||||
{
|
||||
res[0] -= Scalar(EIGEN_PI);
|
||||
}
|
||||
else
|
||||
{
|
||||
res[0] += Scalar(EIGEN_PI);
|
||||
}
|
||||
res[1] = numext::atan2(-coeff(i, k), -c2);
|
||||
}
|
||||
else
|
||||
{
|
||||
res[1] = numext::atan2(-coeff(i, k), c2);
|
||||
}
|
||||
Scalar s1 = numext::sin(res[0]);
|
||||
Scalar c1 = numext::cos(res[0]);
|
||||
res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
|
||||
}
|
||||
if (!odd)
|
||||
{
|
||||
res = -res;
|
||||
}
|
||||
|
||||
return res;
|
||||
}
|
||||
|
||||
|
@ -2,11 +2,15 @@
|
||||
// for linear algebra.
|
||||
//
|
||||
// Copyright (C) 2008-2012 Gael Guennebaud <gael.guennebaud@inria.fr>
|
||||
// Copyright (C) 2023 Juraj Oršulić, University of Zagreb <juraj.orsulic@fer.hr>
|
||||
//
|
||||
// This Source Code Form is subject to the terms of the Mozilla
|
||||
// Public License v. 2.0. If a copy of the MPL was not distributed
|
||||
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
||||
|
||||
// Silence warnings about using the deprecated non-canonical .eulerAngles(), which are still being tested.
|
||||
#define EIGEN_NO_DEPRECATED_WARNING
|
||||
|
||||
#include "main.h"
|
||||
#include <Eigen/Geometry>
|
||||
#include <Eigen/LU>
|
||||
@ -14,53 +18,89 @@
|
||||
|
||||
|
||||
template<typename Scalar>
|
||||
void verify_euler(const Matrix<Scalar,3,1>& ea, int i, int j, int k)
|
||||
void verify_euler(const Matrix<Scalar, 3, 1>& ea, int i, int j, int k)
|
||||
{
|
||||
typedef Matrix<Scalar,3,3> Matrix3;
|
||||
typedef Matrix<Scalar,3,1> Vector3;
|
||||
typedef Matrix<Scalar, 3, 3> Matrix3;
|
||||
typedef Matrix<Scalar, 3, 1> Vector3;
|
||||
typedef AngleAxis<Scalar> AngleAxisx;
|
||||
using std::abs;
|
||||
Matrix3 m(AngleAxisx(ea[0], Vector3::Unit(i)) * AngleAxisx(ea[1], Vector3::Unit(j)) * AngleAxisx(ea[2], Vector3::Unit(k)));
|
||||
Vector3 eabis = m.eulerAngles(i, j, k);
|
||||
Matrix3 mbis(AngleAxisx(eabis[0], Vector3::Unit(i)) * AngleAxisx(eabis[1], Vector3::Unit(j)) * AngleAxisx(eabis[2], Vector3::Unit(k)));
|
||||
VERIFY_IS_APPROX(m, mbis);
|
||||
/* If I==K, and ea[1]==0, then there no unique solution. */
|
||||
/* The remark apply in the case where I!=K, and |ea[1]| is close to pi/2. */
|
||||
if((i!=k || !numext::is_exactly_zero(ea[1])) && (i == k || !internal::isApprox(abs(ea[1]), Scalar(EIGEN_PI / 2), test_precision<Scalar>())) )
|
||||
VERIFY((ea-eabis).norm() <= test_precision<Scalar>());
|
||||
|
||||
// approx_or_less_than does not work for 0
|
||||
VERIFY(0 < eabis[0] || test_isMuchSmallerThan(eabis[0], Scalar(1)));
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[0], Scalar(EIGEN_PI));
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(-Scalar(EIGEN_PI), eabis[1]);
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[1], Scalar(EIGEN_PI));
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(-Scalar(EIGEN_PI), eabis[2]);
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[2], Scalar(EIGEN_PI));
|
||||
const Matrix3 m(AngleAxisx(ea[0], Vector3::Unit(i)) * AngleAxisx(ea[1], Vector3::Unit(j)) * AngleAxisx(ea[2], Vector3::Unit(k)));
|
||||
|
||||
// Test non-canonical eulerAngles
|
||||
{
|
||||
Vector3 eabis = m.eulerAngles(i, j, k);
|
||||
Matrix3 mbis(AngleAxisx(eabis[0], Vector3::Unit(i)) * AngleAxisx(eabis[1], Vector3::Unit(j)) * AngleAxisx(eabis[2], Vector3::Unit(k)));
|
||||
VERIFY_IS_APPROX(m, mbis);
|
||||
|
||||
// approx_or_less_than does not work for 0
|
||||
VERIFY(0 < eabis[0] || test_isMuchSmallerThan(eabis[0], Scalar(1)));
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[0], Scalar(EIGEN_PI));
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(-Scalar(EIGEN_PI), eabis[1]);
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[1], Scalar(EIGEN_PI));
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(-Scalar(EIGEN_PI), eabis[2]);
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[2], Scalar(EIGEN_PI));
|
||||
}
|
||||
|
||||
// Test canonicalEulerAngles
|
||||
{
|
||||
Vector3 eabis = m.canonicalEulerAngles(i, j, k);
|
||||
Matrix3 mbis(AngleAxisx(eabis[0], Vector3::Unit(i)) * AngleAxisx(eabis[1], Vector3::Unit(j)) * AngleAxisx(eabis[2], Vector3::Unit(k)));
|
||||
VERIFY_IS_APPROX(m, mbis);
|
||||
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(-Scalar(EIGEN_PI), eabis[0]);
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[0], Scalar(EIGEN_PI));
|
||||
if (i != k)
|
||||
{
|
||||
// Tait-Bryan sequence
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(-Scalar(EIGEN_PI / 2), eabis[1]);
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[1], Scalar(EIGEN_PI / 2));
|
||||
}
|
||||
else
|
||||
{
|
||||
// Proper Euler sequence
|
||||
// approx_or_less_than does not work for 0
|
||||
VERIFY(0 < eabis[1] || test_isMuchSmallerThan(eabis[1], Scalar(1)));
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[1], Scalar(EIGEN_PI));
|
||||
}
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(-Scalar(EIGEN_PI), eabis[2]);
|
||||
VERIFY_IS_APPROX_OR_LESS_THAN(eabis[2], Scalar(EIGEN_PI));
|
||||
}
|
||||
}
|
||||
|
||||
template<typename Scalar> void check_all_var(const Matrix<Scalar,3,1>& ea)
|
||||
template<typename Scalar> void check_all_var(const Matrix<Scalar, 3, 1>& ea)
|
||||
{
|
||||
verify_euler(ea, 0,1,2);
|
||||
verify_euler(ea, 0,1,0);
|
||||
verify_euler(ea, 0,2,1);
|
||||
verify_euler(ea, 0,2,0);
|
||||
auto verify_permutation = [](const Matrix<Scalar, 3, 1>& eap)
|
||||
{
|
||||
verify_euler(eap, 0, 1, 2);
|
||||
verify_euler(eap, 0, 1, 0);
|
||||
verify_euler(eap, 0, 2, 1);
|
||||
verify_euler(eap, 0, 2, 0);
|
||||
|
||||
verify_euler(ea, 1,2,0);
|
||||
verify_euler(ea, 1,2,1);
|
||||
verify_euler(ea, 1,0,2);
|
||||
verify_euler(ea, 1,0,1);
|
||||
verify_euler(eap, 1, 2, 0);
|
||||
verify_euler(eap, 1, 2, 1);
|
||||
verify_euler(eap, 1, 0, 2);
|
||||
verify_euler(eap, 1, 0, 1);
|
||||
|
||||
verify_euler(ea, 2,0,1);
|
||||
verify_euler(ea, 2,0,2);
|
||||
verify_euler(ea, 2,1,0);
|
||||
verify_euler(ea, 2,1,2);
|
||||
verify_euler(eap, 2, 0, 1);
|
||||
verify_euler(eap, 2, 0, 2);
|
||||
verify_euler(eap, 2, 1, 0);
|
||||
verify_euler(eap, 2, 1, 2);
|
||||
};
|
||||
|
||||
int i, j, k;
|
||||
for (i = 0; i < 3; i++)
|
||||
for (j = 0; j < 3; j++)
|
||||
for (k = 0; k < 3; k++)
|
||||
{
|
||||
Matrix<Scalar,3,1> eap(ea(i), ea(j), ea(k));
|
||||
verify_permutation(eap);
|
||||
}
|
||||
}
|
||||
|
||||
template<typename Scalar> void eulerangles()
|
||||
{
|
||||
typedef Matrix<Scalar,3,3> Matrix3;
|
||||
typedef Matrix<Scalar,3,1> Vector3;
|
||||
typedef Array<Scalar,3,1> Array3;
|
||||
typedef Matrix<Scalar, 3, 3> Matrix3;
|
||||
typedef Matrix<Scalar, 3, 1> Vector3;
|
||||
typedef Array<Scalar, 3, 1> Array3;
|
||||
typedef Quaternion<Scalar> Quaternionx;
|
||||
typedef AngleAxis<Scalar> AngleAxisx;
|
||||
|
||||
@ -69,43 +109,97 @@ template<typename Scalar> void eulerangles()
|
||||
q1 = AngleAxisx(a, Vector3::Random().normalized());
|
||||
Matrix3 m;
|
||||
m = q1;
|
||||
|
||||
Vector3 ea = m.eulerAngles(0,1,2);
|
||||
|
||||
Vector3 ea = m.eulerAngles(0, 1, 2);
|
||||
check_all_var(ea);
|
||||
ea = m.eulerAngles(0,1,0);
|
||||
ea = m.eulerAngles(0, 1, 0);
|
||||
check_all_var(ea);
|
||||
|
||||
|
||||
// Check with purely random Quaternion:
|
||||
q1.coeffs() = Quaternionx::Coefficients::Random().normalized();
|
||||
m = q1;
|
||||
ea = m.eulerAngles(0,1,2);
|
||||
ea = m.eulerAngles(0, 1, 2);
|
||||
check_all_var(ea);
|
||||
ea = m.eulerAngles(0,1,0);
|
||||
ea = m.eulerAngles(0, 1, 0);
|
||||
check_all_var(ea);
|
||||
|
||||
// Check with random angles in range [0:pi]x[-pi:pi]x[-pi:pi].
|
||||
ea = (Array3::Random() + Array3(1,0,0))*Scalar(EIGEN_PI)*Array3(0.5,1,1);
|
||||
|
||||
// Check with random angles in range [-pi:pi]x[-pi:pi]x[-pi:pi].
|
||||
ea = Array3::Random() * Scalar(EIGEN_PI);
|
||||
check_all_var(ea);
|
||||
|
||||
ea[2] = ea[0] = internal::random<Scalar>(0,Scalar(EIGEN_PI));
|
||||
|
||||
auto test_with_some_zeros = [](const Vector3& eaz)
|
||||
{
|
||||
check_all_var(eaz);
|
||||
Vector3 ea_glz = eaz;
|
||||
ea_glz[0] = Scalar(0);
|
||||
check_all_var(ea_glz);
|
||||
ea_glz[0] = internal::random<Scalar>(-0.001, 0.001);
|
||||
check_all_var(ea_glz);
|
||||
ea_glz[2] = Scalar(0);
|
||||
check_all_var(ea_glz);
|
||||
ea_glz[2] = internal::random<Scalar>(-0.001, 0.001);
|
||||
check_all_var(ea_glz);
|
||||
};
|
||||
// Check gimbal lock configurations and a bit noisy gimbal locks
|
||||
Vector3 ea_gl = ea;
|
||||
ea_gl[1] = EIGEN_PI/2;
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_gl[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_gl[1] = -EIGEN_PI/2;
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_gl[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_gl[1] = EIGEN_PI/2;
|
||||
ea_gl[2] = ea_gl[0];
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_gl[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_gl[1] = -EIGEN_PI/2;
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_gl[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
|
||||
// Similar to above, but with pi instead of pi/2
|
||||
Vector3 ea_pi = ea;
|
||||
ea_pi[1] = EIGEN_PI;
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_pi[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_pi[1] = -EIGEN_PI;
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_pi[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_pi[1] = EIGEN_PI;
|
||||
ea_pi[2] = ea_pi[0];
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_pi[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_pi[1] = -EIGEN_PI;
|
||||
test_with_some_zeros(ea_gl);
|
||||
ea_pi[1] += internal::random<Scalar>(-0.001, 0.001);
|
||||
test_with_some_zeros(ea_gl);
|
||||
|
||||
ea[2] = ea[0] = internal::random<Scalar>(0, Scalar(EIGEN_PI));
|
||||
check_all_var(ea);
|
||||
|
||||
ea[0] = ea[1] = internal::random<Scalar>(0,Scalar(EIGEN_PI));
|
||||
|
||||
ea[0] = ea[1] = internal::random<Scalar>(0, Scalar(EIGEN_PI));
|
||||
check_all_var(ea);
|
||||
|
||||
|
||||
ea[1] = 0;
|
||||
check_all_var(ea);
|
||||
|
||||
|
||||
ea.head(2).setZero();
|
||||
check_all_var(ea);
|
||||
|
||||
|
||||
ea.setZero();
|
||||
check_all_var(ea);
|
||||
}
|
||||
|
||||
EIGEN_DECLARE_TEST(geo_eulerangles)
|
||||
{
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
for(int i = 0; i < g_repeat; i++)
|
||||
{
|
||||
CALL_SUBTEST_1( eulerangles<float>() );
|
||||
CALL_SUBTEST_2( eulerangles<double>() );
|
||||
}
|
||||
|
Loading…
x
Reference in New Issue
Block a user