2008-06-03 15:32:12 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2008-06-03 15:32:12 +08:00
|
|
|
//
|
2009-02-17 17:53:05 +08:00
|
|
|
// Copyright (C) 2008-2009 Gael Guennebaud <g.gael@free.fr>
|
2008-06-03 15:32:12 +08:00
|
|
|
//
|
|
|
|
// 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/>.
|
|
|
|
|
|
|
|
#include "main.h"
|
|
|
|
#include <Eigen/Geometry>
|
2008-06-16 01:22:41 +08:00
|
|
|
#include <Eigen/LU>
|
2009-01-23 00:39:08 +08:00
|
|
|
#include <Eigen/SVD>
|
2008-06-03 15:32:12 +08:00
|
|
|
|
2009-03-08 19:35:30 +08:00
|
|
|
template<typename Scalar, int Mode> void transformations(void)
|
2008-06-03 15:32:12 +08:00
|
|
|
{
|
|
|
|
/* this test covers the following files:
|
2008-06-15 16:33:44 +08:00
|
|
|
Cross.h Quaternion.h, Transform.cpp
|
2008-06-03 15:32:12 +08:00
|
|
|
*/
|
2008-06-16 01:22:41 +08:00
|
|
|
typedef Matrix<Scalar,2,2> Matrix2;
|
2008-06-03 15:32:12 +08:00
|
|
|
typedef Matrix<Scalar,3,3> Matrix3;
|
|
|
|
typedef Matrix<Scalar,4,4> Matrix4;
|
2008-06-16 01:22:41 +08:00
|
|
|
typedef Matrix<Scalar,2,1> Vector2;
|
2008-06-03 15:32:12 +08:00
|
|
|
typedef Matrix<Scalar,3,1> Vector3;
|
|
|
|
typedef Matrix<Scalar,4,1> Vector4;
|
2008-10-26 07:10:21 +08:00
|
|
|
typedef Quaternion<Scalar> Quaternionx;
|
|
|
|
typedef AngleAxis<Scalar> AngleAxisx;
|
2009-03-08 19:35:30 +08:00
|
|
|
typedef Transform<Scalar,2,Mode> Transform2;
|
|
|
|
typedef Transform<Scalar,3,Mode> Transform3;
|
|
|
|
typedef typename Transform3::MatrixType MatrixType;
|
2009-01-29 00:26:06 +08:00
|
|
|
typedef DiagonalMatrix<Scalar,2> AlignedScaling2;
|
|
|
|
typedef DiagonalMatrix<Scalar,3> AlignedScaling3;
|
2008-08-30 08:08:23 +08:00
|
|
|
typedef Translation<Scalar,2> Translation2;
|
|
|
|
typedef Translation<Scalar,3> Translation3;
|
2008-06-03 15:32:12 +08:00
|
|
|
|
2008-09-11 19:19:34 +08:00
|
|
|
Scalar largeEps = test_precision<Scalar>();
|
|
|
|
if (ei_is_same_type<Scalar,float>::ret)
|
2009-01-28 04:47:12 +08:00
|
|
|
largeEps = 1e-2f;
|
2008-09-11 19:19:34 +08:00
|
|
|
|
2008-09-02 01:31:21 +08:00
|
|
|
Vector3 v0 = Vector3::Random(),
|
|
|
|
v1 = Vector3::Random(),
|
|
|
|
v2 = Vector3::Random();
|
|
|
|
Vector2 u0 = Vector2::Random();
|
2009-02-17 17:53:05 +08:00
|
|
|
Matrix3 matrot1, m;
|
2008-06-03 15:32:12 +08:00
|
|
|
|
2009-01-05 22:47:38 +08:00
|
|
|
Scalar a = ei_random<Scalar>(-Scalar(M_PI), Scalar(M_PI));
|
2008-06-03 21:43:29 +08:00
|
|
|
|
2008-10-26 07:10:21 +08:00
|
|
|
VERIFY_IS_APPROX(v0, AngleAxisx(a, v0.normalized()) * v0);
|
2009-01-05 22:47:38 +08:00
|
|
|
VERIFY_IS_APPROX(-v0, AngleAxisx(Scalar(M_PI), v0.unitOrthogonal()) * v0);
|
2008-11-04 03:14:17 +08:00
|
|
|
VERIFY_IS_APPROX(ei_cos(a)*v0.squaredNorm(), v0.dot(AngleAxisx(a, v0.unitOrthogonal()) * v0));
|
2008-10-26 07:10:21 +08:00
|
|
|
m = AngleAxisx(a, v0.normalized()).toRotationMatrix().adjoint();
|
|
|
|
VERIFY_IS_APPROX(Matrix3::Identity(), m * AngleAxisx(a, v0.normalized()));
|
|
|
|
VERIFY_IS_APPROX(Matrix3::Identity(), AngleAxisx(a, v0.normalized()) * m);
|
2008-08-25 07:16:51 +08:00
|
|
|
|
2009-02-17 17:53:05 +08:00
|
|
|
Quaternionx q1, q2;
|
2008-10-26 07:10:21 +08:00
|
|
|
q1 = AngleAxisx(a, v0.normalized());
|
|
|
|
q2 = AngleAxisx(a, v1.normalized());
|
2008-06-03 15:32:12 +08:00
|
|
|
|
2008-06-03 21:43:29 +08:00
|
|
|
// rotation matrix conversion
|
2009-01-22 04:06:16 +08:00
|
|
|
matrot1 = AngleAxisx(Scalar(0.1), Vector3::UnitX())
|
|
|
|
* AngleAxisx(Scalar(0.2), Vector3::UnitY())
|
|
|
|
* AngleAxisx(Scalar(0.3), Vector3::UnitZ());
|
2008-07-19 21:03:23 +08:00
|
|
|
VERIFY_IS_APPROX(matrot1 * v1,
|
2009-01-22 04:06:16 +08:00
|
|
|
AngleAxisx(Scalar(0.1), Vector3(1,0,0)).toRotationMatrix()
|
|
|
|
* (AngleAxisx(Scalar(0.2), Vector3(0,1,0)).toRotationMatrix()
|
|
|
|
* (AngleAxisx(Scalar(0.3), Vector3(0,0,1)).toRotationMatrix() * v1)));
|
2008-06-16 01:22:41 +08:00
|
|
|
|
|
|
|
// angle-axis conversion
|
2009-11-13 23:41:51 +08:00
|
|
|
AngleAxisx aa = AngleAxisx(q1);
|
2008-10-26 07:10:21 +08:00
|
|
|
VERIFY_IS_APPROX(q1 * v1, Quaternionx(aa) * v1);
|
|
|
|
VERIFY_IS_NOT_APPROX(q1 * v1, Quaternionx(AngleAxisx(aa.angle()*2,aa.axis())) * v1);
|
2008-06-16 01:22:41 +08:00
|
|
|
|
2010-02-02 16:27:41 +08:00
|
|
|
aa.fromRotationMatrix(aa.toRotationMatrix());
|
|
|
|
VERIFY_IS_APPROX(q1 * v1, Quaternionx(aa) * v1);
|
|
|
|
VERIFY_IS_NOT_APPROX(q1 * v1, Quaternionx(AngleAxisx(aa.angle()*2,aa.axis())) * v1);
|
|
|
|
|
2008-06-16 01:22:41 +08:00
|
|
|
// AngleAxis
|
2008-10-26 07:10:21 +08:00
|
|
|
VERIFY_IS_APPROX(AngleAxisx(a,v1.normalized()).toRotationMatrix(),
|
|
|
|
Quaternionx(AngleAxisx(a,v1.normalized())).toRotationMatrix());
|
2008-06-16 01:22:41 +08:00
|
|
|
|
2008-10-26 07:10:21 +08:00
|
|
|
AngleAxisx aa1;
|
2008-06-16 01:22:41 +08:00
|
|
|
m = q1.toRotationMatrix();
|
2008-06-21 23:01:49 +08:00
|
|
|
aa1 = m;
|
2008-10-26 07:10:21 +08:00
|
|
|
VERIFY_IS_APPROX(AngleAxisx(m).toRotationMatrix(),
|
|
|
|
Quaternionx(m).toRotationMatrix());
|
2008-06-16 01:22:41 +08:00
|
|
|
|
2008-06-15 16:33:44 +08:00
|
|
|
// Transform
|
|
|
|
// TODO complete the tests !
|
|
|
|
a = 0;
|
2009-01-22 04:06:16 +08:00
|
|
|
while (ei_abs(a)<Scalar(0.1))
|
|
|
|
a = ei_random<Scalar>(-Scalar(0.4)*Scalar(M_PI), Scalar(0.4)*Scalar(M_PI));
|
2008-10-26 07:10:21 +08:00
|
|
|
q1 = AngleAxisx(a, v0.normalized());
|
2008-06-15 16:33:44 +08:00
|
|
|
Transform3 t0, t1, t2;
|
2009-09-20 07:59:49 +08:00
|
|
|
|
|
|
|
// first test setIdentity() and Identity()
|
2008-06-15 16:33:44 +08:00
|
|
|
t0.setIdentity();
|
2009-09-20 07:59:49 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), Transform3::MatrixType::Identity());
|
|
|
|
t0.matrix().setZero();
|
|
|
|
t0 = Transform3::Identity();
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), Transform3::MatrixType::Identity());
|
2009-11-19 01:15:19 +08:00
|
|
|
|
2008-08-22 01:02:47 +08:00
|
|
|
t0.linear() = q1.toRotationMatrix();
|
2008-06-15 16:33:44 +08:00
|
|
|
t1.setIdentity();
|
2008-08-22 01:02:47 +08:00
|
|
|
t1.linear() = q1.toRotationMatrix();
|
2008-06-15 16:33:44 +08:00
|
|
|
|
2008-09-02 01:31:21 +08:00
|
|
|
v0 << 50, 2, 1;//= ei_random_matrix<Vector3>().cwiseProduct(Vector3(10,2,0.5));
|
2008-06-15 16:33:44 +08:00
|
|
|
t0.scale(v0);
|
|
|
|
t1.prescale(v0);
|
|
|
|
|
2010-01-05 10:24:43 +08:00
|
|
|
VERIFY_IS_APPROX( (t0 * Vector3(1,0,0)).template head<3>().norm(), v0.x());
|
2009-01-28 04:47:12 +08:00
|
|
|
//VERIFY(!ei_isApprox((t1 * Vector3(1,0,0)).norm(), v0.x()));
|
2008-06-15 16:33:44 +08:00
|
|
|
|
|
|
|
t0.setIdentity();
|
|
|
|
t1.setIdentity();
|
|
|
|
v1 << 1, 2, 3;
|
2008-08-22 01:02:47 +08:00
|
|
|
t0.linear() = q1.toRotationMatrix();
|
2008-06-15 16:33:44 +08:00
|
|
|
t0.pretranslate(v0);
|
|
|
|
t0.scale(v1);
|
2008-08-22 01:02:47 +08:00
|
|
|
t1.linear() = q1.conjugate().toRotationMatrix();
|
2009-12-11 05:00:35 +08:00
|
|
|
t1.prescale(v1.cwiseInverse());
|
2008-06-15 16:33:44 +08:00
|
|
|
t1.translate(-v0);
|
|
|
|
|
2009-03-08 19:35:30 +08:00
|
|
|
VERIFY((t0 * t1).matrix().isIdentity(test_precision<Scalar>()));
|
2008-06-15 19:54:18 +08:00
|
|
|
|
|
|
|
t1.fromPositionOrientationScale(v0, q1, v1);
|
|
|
|
VERIFY_IS_APPROX(t1.matrix(), t0.matrix());
|
2008-08-29 05:44:56 +08:00
|
|
|
VERIFY_IS_APPROX(t1*v1, t0*v1);
|
2008-06-16 01:22:41 +08:00
|
|
|
|
2008-08-31 07:10:46 +08:00
|
|
|
t0.setIdentity(); t0.scale(v0).rotate(q1.toRotationMatrix());
|
|
|
|
t1.setIdentity(); t1.scale(v0).rotate(q1);
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
2008-10-26 07:10:21 +08:00
|
|
|
t0.setIdentity(); t0.scale(v0).rotate(AngleAxisx(q1));
|
2008-08-31 07:10:46 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(t0.scale(a).matrix(), t1.scale(Vector3::Constant(a)).matrix());
|
|
|
|
VERIFY_IS_APPROX(t0.prescale(a).matrix(), t1.prescale(Vector3::Constant(a)).matrix());
|
|
|
|
|
2008-12-20 00:16:39 +08:00
|
|
|
// More transform constructors, operator=, operator*=
|
2008-12-20 00:31:22 +08:00
|
|
|
|
2009-01-05 22:26:34 +08:00
|
|
|
Matrix3 mat3 = Matrix3::Random();
|
|
|
|
Matrix4 mat4;
|
|
|
|
mat4 << mat3 , Vector3::Zero() , Vector4::Zero().transpose();
|
|
|
|
Transform3 tmat3(mat3), tmat4(mat4);
|
2009-03-08 19:35:30 +08:00
|
|
|
if(Mode!=int(AffineCompact))
|
|
|
|
tmat4.matrix()(3,3) = Scalar(1);
|
2009-01-05 22:26:34 +08:00
|
|
|
VERIFY_IS_APPROX(tmat3.matrix(), tmat4.matrix());
|
|
|
|
|
2009-01-05 22:47:38 +08:00
|
|
|
Scalar a3 = ei_random<Scalar>(-Scalar(M_PI), Scalar(M_PI));
|
2008-12-19 22:52:03 +08:00
|
|
|
Vector3 v3 = Vector3::Random().normalized();
|
|
|
|
AngleAxisx aa3(a3, v3);
|
|
|
|
Transform3 t3(aa3);
|
|
|
|
Transform3 t4;
|
|
|
|
t4 = aa3;
|
|
|
|
VERIFY_IS_APPROX(t3.matrix(), t4.matrix());
|
2008-12-19 23:55:35 +08:00
|
|
|
t4.rotate(AngleAxisx(-a3,v3));
|
2009-03-08 19:35:30 +08:00
|
|
|
VERIFY_IS_APPROX(t4.matrix(), MatrixType::Identity());
|
2008-12-20 00:16:39 +08:00
|
|
|
t4 *= aa3;
|
|
|
|
VERIFY_IS_APPROX(t3.matrix(), t4.matrix());
|
2008-12-19 22:52:03 +08:00
|
|
|
|
|
|
|
v3 = Vector3::Random();
|
|
|
|
Translation3 tv3(v3);
|
|
|
|
Transform3 t5(tv3);
|
|
|
|
t4 = tv3;
|
|
|
|
VERIFY_IS_APPROX(t5.matrix(), t4.matrix());
|
|
|
|
t4.translate(-v3);
|
2009-03-08 19:35:30 +08:00
|
|
|
VERIFY_IS_APPROX(t4.matrix(), MatrixType::Identity());
|
2008-12-20 00:16:39 +08:00
|
|
|
t4 *= tv3;
|
|
|
|
VERIFY_IS_APPROX(t5.matrix(), t4.matrix());
|
2008-12-19 22:52:03 +08:00
|
|
|
|
2009-01-29 00:26:06 +08:00
|
|
|
AlignedScaling3 sv3(v3);
|
2008-12-19 22:52:03 +08:00
|
|
|
Transform3 t6(sv3);
|
|
|
|
t4 = sv3;
|
|
|
|
VERIFY_IS_APPROX(t6.matrix(), t4.matrix());
|
2009-12-11 05:00:35 +08:00
|
|
|
t4.scale(v3.cwiseInverse());
|
2009-03-08 19:35:30 +08:00
|
|
|
VERIFY_IS_APPROX(t4.matrix(), MatrixType::Identity());
|
2008-12-20 00:16:39 +08:00
|
|
|
t4 *= sv3;
|
|
|
|
VERIFY_IS_APPROX(t6.matrix(), t4.matrix());
|
2008-12-19 22:52:03 +08:00
|
|
|
|
2008-12-20 00:31:22 +08:00
|
|
|
// matrix * transform
|
2009-03-08 19:35:30 +08:00
|
|
|
VERIFY_IS_APPROX((t3.matrix()*t4).matrix(), (t3*t4).matrix());
|
2008-12-20 00:31:22 +08:00
|
|
|
|
2008-12-19 22:52:03 +08:00
|
|
|
// chained Transform product
|
2008-12-20 00:16:39 +08:00
|
|
|
VERIFY_IS_APPROX(((t3*t4)*t5).matrix(), (t3*(t4*t5)).matrix());
|
|
|
|
|
|
|
|
// check that Transform product doesn't have aliasing problems
|
|
|
|
t5 = t4;
|
|
|
|
t5 = t5*t5;
|
|
|
|
VERIFY_IS_APPROX(t5, t4*t4);
|
|
|
|
|
2008-06-16 01:22:41 +08:00
|
|
|
// 2D transformation
|
|
|
|
Transform2 t20, t21;
|
2008-09-02 01:31:21 +08:00
|
|
|
Vector2 v20 = Vector2::Random();
|
|
|
|
Vector2 v21 = Vector2::Random();
|
2008-08-23 23:14:20 +08:00
|
|
|
for (int k=0; k<2; ++k)
|
2009-01-22 04:06:16 +08:00
|
|
|
if (ei_abs(v21[k])<Scalar(1e-3)) v21[k] = Scalar(1e-3);
|
2008-06-16 01:22:41 +08:00
|
|
|
t21.setIdentity();
|
2008-08-22 01:02:47 +08:00
|
|
|
t21.linear() = Rotation2D<Scalar>(a).toRotationMatrix();
|
2008-06-16 01:22:41 +08:00
|
|
|
VERIFY_IS_APPROX(t20.fromPositionOrientationScale(v20,a,v21).matrix(),
|
|
|
|
t21.pretranslate(v20).scale(v21).matrix());
|
|
|
|
|
|
|
|
t21.setIdentity();
|
2008-08-22 01:02:47 +08:00
|
|
|
t21.linear() = Rotation2D<Scalar>(-a).toRotationMatrix();
|
2008-08-23 23:14:20 +08:00
|
|
|
VERIFY( (t20.fromPositionOrientationScale(v20,a,v21)
|
2009-12-11 05:00:35 +08:00
|
|
|
* (t21.prescale(v21.cwiseInverse()).translate(-v20))).matrix().isIdentity(test_precision<Scalar>()) );
|
2008-08-30 08:08:23 +08:00
|
|
|
|
|
|
|
// Transform - new API
|
|
|
|
// 3D
|
|
|
|
t0.setIdentity();
|
|
|
|
t0.rotate(q1).scale(v0).translate(v0);
|
2009-01-29 00:26:06 +08:00
|
|
|
// mat * aligned scaling and mat * translation
|
|
|
|
t1 = (Matrix3(q1) * AlignedScaling3(v0)) * Translation3(v0);
|
2008-08-30 08:08:23 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
2009-01-29 00:26:06 +08:00
|
|
|
t1 = (Matrix3(q1) * Scaling(v0)) * Translation3(v0);
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
t1 = (q1 * Scaling(v0)) * Translation3(v0);
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
// mat * transformation and aligned scaling * translation
|
|
|
|
t1 = Matrix3(q1) * (AlignedScaling3(v0) * Translation3(v0));
|
2008-08-30 08:08:23 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
|
|
|
t0.setIdentity();
|
|
|
|
t0.prerotate(q1).prescale(v0).pretranslate(v0);
|
2009-01-29 00:26:06 +08:00
|
|
|
// translation * aligned scaling and transformation * mat
|
|
|
|
t1 = (Translation3(v0) * AlignedScaling3(v0)) * Matrix3(q1);
|
2008-08-30 08:08:23 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
// scaling * mat and translation * mat
|
2009-01-29 00:26:06 +08:00
|
|
|
t1 = Translation3(v0) * (AlignedScaling3(v0) * Matrix3(q1));
|
2008-08-30 08:08:23 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
2008-10-26 07:10:21 +08:00
|
|
|
|
2008-08-30 08:08:23 +08:00
|
|
|
t0.setIdentity();
|
|
|
|
t0.scale(v0).translate(v0).rotate(q1);
|
2009-01-29 00:26:06 +08:00
|
|
|
// translation * mat and aligned scaling * transformation
|
|
|
|
t1 = AlignedScaling3(v0) * (Translation3(v0) * Matrix3(q1));
|
2008-08-30 08:08:23 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
2009-01-29 00:26:06 +08:00
|
|
|
// transformation * aligned scaling
|
2008-08-30 08:08:23 +08:00
|
|
|
t0.scale(v0);
|
2009-01-29 00:26:06 +08:00
|
|
|
t1 = t1 * AlignedScaling3(v0);
|
2008-08-30 08:08:23 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
// transformation * translation
|
|
|
|
t0.translate(v0);
|
|
|
|
t1 = t1 * Translation3(v0);
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
2008-10-26 07:10:21 +08:00
|
|
|
// translation * transformation
|
2008-08-30 08:08:23 +08:00
|
|
|
t0.pretranslate(v0);
|
|
|
|
t1 = Translation3(v0) * t1;
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
2008-08-31 04:11:04 +08:00
|
|
|
// transform * quaternion
|
|
|
|
t0.rotate(q1);
|
|
|
|
t1 = t1 * q1;
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
|
|
|
// translation * quaternion
|
|
|
|
t0.translate(v1).rotate(q1);
|
|
|
|
t1 = t1 * (Translation3(v1) * q1);
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
2009-01-29 00:26:06 +08:00
|
|
|
// aligned scaling * quaternion
|
2008-08-31 04:11:04 +08:00
|
|
|
t0.scale(v1).rotate(q1);
|
2009-01-29 00:26:06 +08:00
|
|
|
t1 = t1 * (AlignedScaling3(v1) * q1);
|
2008-08-31 04:11:04 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
|
|
|
// quaternion * transform
|
|
|
|
t0.prerotate(q1);
|
|
|
|
t1 = q1 * t1;
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
|
|
|
// quaternion * translation
|
|
|
|
t0.rotate(q1).translate(v1);
|
|
|
|
t1 = t1 * (q1 * Translation3(v1));
|
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
2009-01-29 00:26:06 +08:00
|
|
|
// quaternion * aligned scaling
|
2008-08-31 04:11:04 +08:00
|
|
|
t0.rotate(q1).scale(v1);
|
2009-01-29 00:26:06 +08:00
|
|
|
t1 = t1 * (q1 * AlignedScaling3(v1));
|
2008-08-31 04:11:04 +08:00
|
|
|
VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
|
|
|
|
|
2008-08-30 08:08:23 +08:00
|
|
|
// translation * vector
|
|
|
|
t0.setIdentity();
|
|
|
|
t0.translate(v0);
|
2010-01-05 10:24:43 +08:00
|
|
|
VERIFY_IS_APPROX((t0 * v1).template head<3>(), Translation3(v0) * v1);
|
2008-08-30 08:08:23 +08:00
|
|
|
|
2009-01-29 00:26:06 +08:00
|
|
|
// AlignedScaling * vector
|
2008-08-30 08:08:23 +08:00
|
|
|
t0.setIdentity();
|
|
|
|
t0.scale(v0);
|
2010-01-05 10:24:43 +08:00
|
|
|
VERIFY_IS_APPROX((t0 * v1).template head<3>(), AlignedScaling3(v0) * v1);
|
2008-08-30 20:42:06 +08:00
|
|
|
|
|
|
|
// test transform inversion
|
2009-10-29 06:19:29 +08:00
|
|
|
t0.setIdentity();
|
|
|
|
t0.translate(v0);
|
|
|
|
t0.linear().setRandom();
|
|
|
|
Matrix4 t044 = Matrix4::Zero();
|
|
|
|
t044(3,3) = 1;
|
|
|
|
t044.block(0,0,t0.matrix().rows(),4) = t0.matrix();
|
|
|
|
VERIFY_IS_APPROX(t0.inverse(Affine).matrix(), t044.inverse().block(0,0,t0.matrix().rows(),4));
|
|
|
|
t0.setIdentity();
|
|
|
|
t0.translate(v0).rotate(q1);
|
|
|
|
t044 = Matrix4::Zero();
|
|
|
|
t044(3,3) = 1;
|
|
|
|
t044.block(0,0,t0.matrix().rows(),4) = t0.matrix();
|
|
|
|
VERIFY_IS_APPROX(t0.inverse(Isometry).matrix(), t044.inverse().block(0,0,t0.matrix().rows(),4));
|
2008-08-30 23:22:45 +08:00
|
|
|
|
2009-01-29 00:26:06 +08:00
|
|
|
// test extract rotation and aligned scaling
|
|
|
|
// t0.setIdentity();
|
|
|
|
// t0.translate(v0).rotate(q1).scale(v1);
|
|
|
|
// VERIFY_IS_APPROX(t0.rotation() * v1, Matrix3(q1) * v1);
|
2009-01-23 00:39:08 +08:00
|
|
|
|
|
|
|
Matrix3 mat_rotation, mat_scaling;
|
|
|
|
t0.setIdentity();
|
|
|
|
t0.translate(v0).rotate(q1).scale(v1);
|
|
|
|
t0.computeRotationScaling(&mat_rotation, &mat_scaling);
|
|
|
|
VERIFY_IS_APPROX(t0.linear(), mat_rotation * mat_scaling);
|
|
|
|
VERIFY_IS_APPROX(mat_rotation*mat_rotation.adjoint(), Matrix3::Identity());
|
|
|
|
VERIFY_IS_APPROX(mat_rotation.determinant(), Scalar(1));
|
|
|
|
t0.computeScalingRotation(&mat_scaling, &mat_rotation);
|
|
|
|
VERIFY_IS_APPROX(t0.linear(), mat_scaling * mat_rotation);
|
|
|
|
VERIFY_IS_APPROX(mat_rotation*mat_rotation.adjoint(), Matrix3::Identity());
|
|
|
|
VERIFY_IS_APPROX(mat_rotation.determinant(), Scalar(1));
|
2008-10-26 07:10:21 +08:00
|
|
|
|
|
|
|
// test casting
|
2009-03-08 19:35:30 +08:00
|
|
|
Transform<float,3,Mode> t1f = t1.template cast<float>();
|
2008-10-26 07:10:21 +08:00
|
|
|
VERIFY_IS_APPROX(t1f.template cast<Scalar>(),t1);
|
2009-03-08 19:35:30 +08:00
|
|
|
Transform<double,3,Mode> t1d = t1.template cast<double>();
|
2008-10-26 07:10:21 +08:00
|
|
|
VERIFY_IS_APPROX(t1d.template cast<Scalar>(),t1);
|
|
|
|
|
|
|
|
Translation3 tr1(v0);
|
|
|
|
Translation<float,3> tr1f = tr1.template cast<float>();
|
|
|
|
VERIFY_IS_APPROX(tr1f.template cast<Scalar>(),tr1);
|
|
|
|
Translation<double,3> tr1d = tr1.template cast<double>();
|
|
|
|
VERIFY_IS_APPROX(tr1d.template cast<Scalar>(),tr1);
|
|
|
|
|
|
|
|
AngleAxis<float> aa1f = aa1.template cast<float>();
|
|
|
|
VERIFY_IS_APPROX(aa1f.template cast<Scalar>(),aa1);
|
|
|
|
AngleAxis<double> aa1d = aa1.template cast<double>();
|
|
|
|
VERIFY_IS_APPROX(aa1d.template cast<Scalar>(),aa1);
|
|
|
|
|
|
|
|
Rotation2D<Scalar> r2d1(ei_random<Scalar>());
|
|
|
|
Rotation2D<float> r2d1f = r2d1.template cast<float>();
|
|
|
|
VERIFY_IS_APPROX(r2d1f.template cast<Scalar>(),r2d1);
|
|
|
|
Rotation2D<double> r2d1d = r2d1.template cast<double>();
|
|
|
|
VERIFY_IS_APPROX(r2d1d.template cast<Scalar>(),r2d1);
|
2009-11-19 01:15:19 +08:00
|
|
|
|
2008-06-03 15:32:12 +08:00
|
|
|
}
|
|
|
|
|
2009-02-17 17:53:05 +08:00
|
|
|
void test_geo_transformations()
|
2008-06-03 15:32:12 +08:00
|
|
|
{
|
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_1(( transformations<double,Affine>() ));
|
|
|
|
CALL_SUBTEST_2(( transformations<float,AffineCompact>() ));
|
|
|
|
CALL_SUBTEST_3(( transformations<double,Projective>() ));
|
2008-06-03 15:32:12 +08:00
|
|
|
}
|
|
|
|
}
|