Home | History | Annotate | Download | only in test
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      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 #include "main.h"
     11 #include <Eigen/Geometry>
     12 #include <Eigen/LU>
     13 #include <Eigen/SVD>
     14 
     15 template<typename Scalar, int Mode, int Options> void non_projective_only()
     16 {
     17     /* this test covers the following files:
     18      Cross.h Quaternion.h, Transform.cpp
     19   */
     20   typedef Matrix<Scalar,3,1> Vector3;
     21   typedef Quaternion<Scalar> Quaternionx;
     22   typedef AngleAxis<Scalar> AngleAxisx;
     23   typedef Transform<Scalar,3,Mode,Options> Transform3;
     24   typedef DiagonalMatrix<Scalar,3> AlignedScaling3;
     25   typedef Translation<Scalar,3> Translation3;
     26 
     27   Vector3 v0 = Vector3::Random(),
     28           v1 = Vector3::Random();
     29 
     30   Transform3 t0, t1, t2;
     31 
     32   Scalar a = internal::random<Scalar>(-Scalar(M_PI), Scalar(M_PI));
     33 
     34   Quaternionx q1, q2;
     35 
     36   q1 = AngleAxisx(a, v0.normalized());
     37 
     38   t0 = Transform3::Identity();
     39   VERIFY_IS_APPROX(t0.matrix(), Transform3::MatrixType::Identity());
     40 
     41   t0.linear() = q1.toRotationMatrix();
     42 
     43   v0 << 50, 2, 1;
     44   t0.scale(v0);
     45 
     46   VERIFY_IS_APPROX( (t0 * Vector3(1,0,0)).template head<3>().norm(), v0.x());
     47 
     48   t0.setIdentity();
     49   t1.setIdentity();
     50   v1 << 1, 2, 3;
     51   t0.linear() = q1.toRotationMatrix();
     52   t0.pretranslate(v0);
     53   t0.scale(v1);
     54   t1.linear() = q1.conjugate().toRotationMatrix();
     55   t1.prescale(v1.cwiseInverse());
     56   t1.translate(-v0);
     57 
     58   VERIFY((t0 * t1).matrix().isIdentity(test_precision<Scalar>()));
     59 
     60   t1.fromPositionOrientationScale(v0, q1, v1);
     61   VERIFY_IS_APPROX(t1.matrix(), t0.matrix());
     62   VERIFY_IS_APPROX(t1*v1, t0*v1);
     63 
     64   // translation * vector
     65   t0.setIdentity();
     66   t0.translate(v0);
     67   VERIFY_IS_APPROX((t0 * v1).template head<3>(), Translation3(v0) * v1);
     68 
     69   // AlignedScaling * vector
     70   t0.setIdentity();
     71   t0.scale(v0);
     72   VERIFY_IS_APPROX((t0 * v1).template head<3>(), AlignedScaling3(v0) * v1);
     73 }
     74 
     75 template<typename Scalar, int Mode, int Options> void transformations()
     76 {
     77   /* this test covers the following files:
     78      Cross.h Quaternion.h, Transform.cpp
     79   */
     80   using std::cos;
     81   using std::abs;
     82   typedef Matrix<Scalar,3,3> Matrix3;
     83   typedef Matrix<Scalar,4,4> Matrix4;
     84   typedef Matrix<Scalar,2,1> Vector2;
     85   typedef Matrix<Scalar,3,1> Vector3;
     86   typedef Matrix<Scalar,4,1> Vector4;
     87   typedef Quaternion<Scalar> Quaternionx;
     88   typedef AngleAxis<Scalar> AngleAxisx;
     89   typedef Transform<Scalar,2,Mode,Options> Transform2;
     90   typedef Transform<Scalar,3,Mode,Options> Transform3;
     91   typedef typename Transform3::MatrixType MatrixType;
     92   typedef DiagonalMatrix<Scalar,3> AlignedScaling3;
     93   typedef Translation<Scalar,2> Translation2;
     94   typedef Translation<Scalar,3> Translation3;
     95 
     96   Vector3 v0 = Vector3::Random(),
     97           v1 = Vector3::Random();
     98   Matrix3 matrot1, m;
     99 
    100   Scalar a = internal::random<Scalar>(-Scalar(M_PI), Scalar(M_PI));
    101   Scalar s0 = internal::random<Scalar>();
    102 
    103   VERIFY_IS_APPROX(v0, AngleAxisx(a, v0.normalized()) * v0);
    104   VERIFY_IS_APPROX(-v0, AngleAxisx(Scalar(M_PI), v0.unitOrthogonal()) * v0);
    105   VERIFY_IS_APPROX(cos(a)*v0.squaredNorm(), v0.dot(AngleAxisx(a, v0.unitOrthogonal()) * v0));
    106   m = AngleAxisx(a, v0.normalized()).toRotationMatrix().adjoint();
    107   VERIFY_IS_APPROX(Matrix3::Identity(), m * AngleAxisx(a, v0.normalized()));
    108   VERIFY_IS_APPROX(Matrix3::Identity(), AngleAxisx(a, v0.normalized()) * m);
    109 
    110   Quaternionx q1, q2;
    111   q1 = AngleAxisx(a, v0.normalized());
    112   q2 = AngleAxisx(a, v1.normalized());
    113 
    114   // rotation matrix conversion
    115   matrot1 = AngleAxisx(Scalar(0.1), Vector3::UnitX())
    116           * AngleAxisx(Scalar(0.2), Vector3::UnitY())
    117           * AngleAxisx(Scalar(0.3), Vector3::UnitZ());
    118   VERIFY_IS_APPROX(matrot1 * v1,
    119        AngleAxisx(Scalar(0.1), Vector3(1,0,0)).toRotationMatrix()
    120     * (AngleAxisx(Scalar(0.2), Vector3(0,1,0)).toRotationMatrix()
    121     * (AngleAxisx(Scalar(0.3), Vector3(0,0,1)).toRotationMatrix() * v1)));
    122 
    123   // angle-axis conversion
    124   AngleAxisx aa = AngleAxisx(q1);
    125   VERIFY_IS_APPROX(q1 * v1, Quaternionx(aa) * v1);
    126   VERIFY_IS_NOT_APPROX(q1 * v1, Quaternionx(AngleAxisx(aa.angle()*2,aa.axis())) * v1);
    127 
    128   aa.fromRotationMatrix(aa.toRotationMatrix());
    129   VERIFY_IS_APPROX(q1 * v1, Quaternionx(aa) * v1);
    130   VERIFY_IS_NOT_APPROX(q1 * v1, Quaternionx(AngleAxisx(aa.angle()*2,aa.axis())) * v1);
    131 
    132   // AngleAxis
    133   VERIFY_IS_APPROX(AngleAxisx(a,v1.normalized()).toRotationMatrix(),
    134     Quaternionx(AngleAxisx(a,v1.normalized())).toRotationMatrix());
    135 
    136   AngleAxisx aa1;
    137   m = q1.toRotationMatrix();
    138   aa1 = m;
    139   VERIFY_IS_APPROX(AngleAxisx(m).toRotationMatrix(),
    140     Quaternionx(m).toRotationMatrix());
    141 
    142   // Transform
    143   // TODO complete the tests !
    144   a = 0;
    145   while (abs(a)<Scalar(0.1))
    146     a = internal::random<Scalar>(-Scalar(0.4)*Scalar(M_PI), Scalar(0.4)*Scalar(M_PI));
    147   q1 = AngleAxisx(a, v0.normalized());
    148   Transform3 t0, t1, t2;
    149 
    150   // first test setIdentity() and Identity()
    151   t0.setIdentity();
    152   VERIFY_IS_APPROX(t0.matrix(), Transform3::MatrixType::Identity());
    153   t0.matrix().setZero();
    154   t0 = Transform3::Identity();
    155   VERIFY_IS_APPROX(t0.matrix(), Transform3::MatrixType::Identity());
    156 
    157   t0.setIdentity();
    158   t1.setIdentity();
    159   v1 << 1, 2, 3;
    160   t0.linear() = q1.toRotationMatrix();
    161   t0.pretranslate(v0);
    162   t0.scale(v1);
    163   t1.linear() = q1.conjugate().toRotationMatrix();
    164   t1.prescale(v1.cwiseInverse());
    165   t1.translate(-v0);
    166 
    167   VERIFY((t0 * t1).matrix().isIdentity(test_precision<Scalar>()));
    168 
    169   t1.fromPositionOrientationScale(v0, q1, v1);
    170   VERIFY_IS_APPROX(t1.matrix(), t0.matrix());
    171 
    172   t0.setIdentity(); t0.scale(v0).rotate(q1.toRotationMatrix());
    173   t1.setIdentity(); t1.scale(v0).rotate(q1);
    174   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    175 
    176   t0.setIdentity(); t0.scale(v0).rotate(AngleAxisx(q1));
    177   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    178 
    179   VERIFY_IS_APPROX(t0.scale(a).matrix(), t1.scale(Vector3::Constant(a)).matrix());
    180   VERIFY_IS_APPROX(t0.prescale(a).matrix(), t1.prescale(Vector3::Constant(a)).matrix());
    181 
    182   // More transform constructors, operator=, operator*=
    183 
    184   Matrix3 mat3 = Matrix3::Random();
    185   Matrix4 mat4;
    186   mat4 << mat3 , Vector3::Zero() , Vector4::Zero().transpose();
    187   Transform3 tmat3(mat3), tmat4(mat4);
    188   if(Mode!=int(AffineCompact))
    189     tmat4.matrix()(3,3) = Scalar(1);
    190   VERIFY_IS_APPROX(tmat3.matrix(), tmat4.matrix());
    191 
    192   Scalar a3 = internal::random<Scalar>(-Scalar(M_PI), Scalar(M_PI));
    193   Vector3 v3 = Vector3::Random().normalized();
    194   AngleAxisx aa3(a3, v3);
    195   Transform3 t3(aa3);
    196   Transform3 t4;
    197   t4 = aa3;
    198   VERIFY_IS_APPROX(t3.matrix(), t4.matrix());
    199   t4.rotate(AngleAxisx(-a3,v3));
    200   VERIFY_IS_APPROX(t4.matrix(), MatrixType::Identity());
    201   t4 *= aa3;
    202   VERIFY_IS_APPROX(t3.matrix(), t4.matrix());
    203 
    204   v3 = Vector3::Random();
    205   Translation3 tv3(v3);
    206   Transform3 t5(tv3);
    207   t4 = tv3;
    208   VERIFY_IS_APPROX(t5.matrix(), t4.matrix());
    209   t4.translate(-v3);
    210   VERIFY_IS_APPROX(t4.matrix(), MatrixType::Identity());
    211   t4 *= tv3;
    212   VERIFY_IS_APPROX(t5.matrix(), t4.matrix());
    213 
    214   AlignedScaling3 sv3(v3);
    215   Transform3 t6(sv3);
    216   t4 = sv3;
    217   VERIFY_IS_APPROX(t6.matrix(), t4.matrix());
    218   t4.scale(v3.cwiseInverse());
    219   VERIFY_IS_APPROX(t4.matrix(), MatrixType::Identity());
    220   t4 *= sv3;
    221   VERIFY_IS_APPROX(t6.matrix(), t4.matrix());
    222 
    223   // matrix * transform
    224   VERIFY_IS_APPROX((t3.matrix()*t4).matrix(), (t3*t4).matrix());
    225 
    226   // chained Transform product
    227   VERIFY_IS_APPROX(((t3*t4)*t5).matrix(), (t3*(t4*t5)).matrix());
    228 
    229   // check that Transform product doesn't have aliasing problems
    230   t5 = t4;
    231   t5 = t5*t5;
    232   VERIFY_IS_APPROX(t5, t4*t4);
    233 
    234   // 2D transformation
    235   Transform2 t20, t21;
    236   Vector2 v20 = Vector2::Random();
    237   Vector2 v21 = Vector2::Random();
    238   for (int k=0; k<2; ++k)
    239     if (abs(v21[k])<Scalar(1e-3)) v21[k] = Scalar(1e-3);
    240   t21.setIdentity();
    241   t21.linear() = Rotation2D<Scalar>(a).toRotationMatrix();
    242   VERIFY_IS_APPROX(t20.fromPositionOrientationScale(v20,a,v21).matrix(),
    243     t21.pretranslate(v20).scale(v21).matrix());
    244 
    245   t21.setIdentity();
    246   t21.linear() = Rotation2D<Scalar>(-a).toRotationMatrix();
    247   VERIFY( (t20.fromPositionOrientationScale(v20,a,v21)
    248         * (t21.prescale(v21.cwiseInverse()).translate(-v20))).matrix().isIdentity(test_precision<Scalar>()) );
    249 
    250   // Transform - new API
    251   // 3D
    252   t0.setIdentity();
    253   t0.rotate(q1).scale(v0).translate(v0);
    254   // mat * aligned scaling and mat * translation
    255   t1 = (Matrix3(q1) * AlignedScaling3(v0)) * Translation3(v0);
    256   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    257   t1 = (Matrix3(q1) * Eigen::Scaling(v0)) * Translation3(v0);
    258   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    259   t1 = (q1 * Eigen::Scaling(v0)) * Translation3(v0);
    260   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    261   // mat * transformation and aligned scaling * translation
    262   t1 = Matrix3(q1) * (AlignedScaling3(v0) * Translation3(v0));
    263   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    264 
    265 
    266   t0.setIdentity();
    267   t0.scale(s0).translate(v0);
    268   t1 = Eigen::Scaling(s0) * Translation3(v0);
    269   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    270   t0.prescale(s0);
    271   t1 = Eigen::Scaling(s0) * t1;
    272   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    273 
    274   t0 = t3;
    275   t0.scale(s0);
    276   t1 = t3 * Eigen::Scaling(s0,s0,s0);
    277   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    278   t0.prescale(s0);
    279   t1 = Eigen::Scaling(s0,s0,s0) * t1;
    280   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    281 
    282   t0 = t3;
    283   t0.scale(s0);
    284   t1 = t3 * Eigen::Scaling(s0);
    285   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    286   t0.prescale(s0);
    287   t1 = Eigen::Scaling(s0) * t1;
    288   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    289 
    290   t0.setIdentity();
    291   t0.prerotate(q1).prescale(v0).pretranslate(v0);
    292   // translation * aligned scaling and transformation * mat
    293   t1 = (Translation3(v0) * AlignedScaling3(v0)) * Transform3(q1);
    294   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    295   // scaling * mat and translation * mat
    296   t1 = Translation3(v0) * (AlignedScaling3(v0) * Transform3(q1));
    297   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    298 
    299   t0.setIdentity();
    300   t0.scale(v0).translate(v0).rotate(q1);
    301   // translation * mat and aligned scaling * transformation
    302   t1 = AlignedScaling3(v0) * (Translation3(v0) * Transform3(q1));
    303   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    304   // transformation * aligned scaling
    305   t0.scale(v0);
    306   t1 *= AlignedScaling3(v0);
    307   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    308   // transformation * translation
    309   t0.translate(v0);
    310   t1 = t1 * Translation3(v0);
    311   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    312   // translation * transformation
    313   t0.pretranslate(v0);
    314   t1 = Translation3(v0) * t1;
    315   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    316 
    317   // transform * quaternion
    318   t0.rotate(q1);
    319   t1 = t1 * q1;
    320   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    321 
    322   // translation * quaternion
    323   t0.translate(v1).rotate(q1);
    324   t1 = t1 * (Translation3(v1) * q1);
    325   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    326 
    327   // aligned scaling * quaternion
    328   t0.scale(v1).rotate(q1);
    329   t1 = t1 * (AlignedScaling3(v1) * q1);
    330   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    331 
    332   // quaternion * transform
    333   t0.prerotate(q1);
    334   t1 = q1 * t1;
    335   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    336 
    337   // quaternion * translation
    338   t0.rotate(q1).translate(v1);
    339   t1 = t1 * (q1 * Translation3(v1));
    340   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    341 
    342   // quaternion * aligned scaling
    343   t0.rotate(q1).scale(v1);
    344   t1 = t1 * (q1 * AlignedScaling3(v1));
    345   VERIFY_IS_APPROX(t0.matrix(), t1.matrix());
    346 
    347   // test transform inversion
    348   t0.setIdentity();
    349   t0.translate(v0);
    350   t0.linear().setRandom();
    351   Matrix4 t044 = Matrix4::Zero();
    352   t044(3,3) = 1;
    353   t044.block(0,0,t0.matrix().rows(),4) = t0.matrix();
    354   VERIFY_IS_APPROX(t0.inverse(Affine).matrix(), t044.inverse().block(0,0,t0.matrix().rows(),4));
    355   t0.setIdentity();
    356   t0.translate(v0).rotate(q1);
    357   t044 = Matrix4::Zero();
    358   t044(3,3) = 1;
    359   t044.block(0,0,t0.matrix().rows(),4) = t0.matrix();
    360   VERIFY_IS_APPROX(t0.inverse(Isometry).matrix(), t044.inverse().block(0,0,t0.matrix().rows(),4));
    361 
    362   Matrix3 mat_rotation, mat_scaling;
    363   t0.setIdentity();
    364   t0.translate(v0).rotate(q1).scale(v1);
    365   t0.computeRotationScaling(&mat_rotation, &mat_scaling);
    366   VERIFY_IS_APPROX(t0.linear(), mat_rotation * mat_scaling);
    367   VERIFY_IS_APPROX(mat_rotation*mat_rotation.adjoint(), Matrix3::Identity());
    368   VERIFY_IS_APPROX(mat_rotation.determinant(), Scalar(1));
    369   t0.computeScalingRotation(&mat_scaling, &mat_rotation);
    370   VERIFY_IS_APPROX(t0.linear(), mat_scaling * mat_rotation);
    371   VERIFY_IS_APPROX(mat_rotation*mat_rotation.adjoint(), Matrix3::Identity());
    372   VERIFY_IS_APPROX(mat_rotation.determinant(), Scalar(1));
    373 
    374   // test casting
    375   Transform<float,3,Mode> t1f = t1.template cast<float>();
    376   VERIFY_IS_APPROX(t1f.template cast<Scalar>(),t1);
    377   Transform<double,3,Mode> t1d = t1.template cast<double>();
    378   VERIFY_IS_APPROX(t1d.template cast<Scalar>(),t1);
    379 
    380   Translation3 tr1(v0);
    381   Translation<float,3> tr1f = tr1.template cast<float>();
    382   VERIFY_IS_APPROX(tr1f.template cast<Scalar>(),tr1);
    383   Translation<double,3> tr1d = tr1.template cast<double>();
    384   VERIFY_IS_APPROX(tr1d.template cast<Scalar>(),tr1);
    385 
    386   AngleAxis<float> aa1f = aa1.template cast<float>();
    387   VERIFY_IS_APPROX(aa1f.template cast<Scalar>(),aa1);
    388   AngleAxis<double> aa1d = aa1.template cast<double>();
    389   VERIFY_IS_APPROX(aa1d.template cast<Scalar>(),aa1);
    390 
    391   Rotation2D<Scalar> r2d1(internal::random<Scalar>());
    392   Rotation2D<float> r2d1f = r2d1.template cast<float>();
    393   VERIFY_IS_APPROX(r2d1f.template cast<Scalar>(),r2d1);
    394   Rotation2D<double> r2d1d = r2d1.template cast<double>();
    395   VERIFY_IS_APPROX(r2d1d.template cast<Scalar>(),r2d1);
    396 
    397   t20 = Translation2(v20) * (Rotation2D<Scalar>(s0) * Eigen::Scaling(s0));
    398   t21 = Translation2(v20) * Rotation2D<Scalar>(s0) * Eigen::Scaling(s0);
    399   VERIFY_IS_APPROX(t20,t21);
    400 }
    401 
    402 template<typename Scalar> void transform_alignment()
    403 {
    404   typedef Transform<Scalar,3,Projective,AutoAlign> Projective3a;
    405   typedef Transform<Scalar,3,Projective,DontAlign> Projective3u;
    406 
    407   EIGEN_ALIGN16 Scalar array1[16];
    408   EIGEN_ALIGN16 Scalar array2[16];
    409   EIGEN_ALIGN16 Scalar array3[16+1];
    410   Scalar* array3u = array3+1;
    411 
    412   Projective3a *p1 = ::new(reinterpret_cast<void*>(array1)) Projective3a;
    413   Projective3u *p2 = ::new(reinterpret_cast<void*>(array2)) Projective3u;
    414   Projective3u *p3 = ::new(reinterpret_cast<void*>(array3u)) Projective3u;
    415 
    416   p1->matrix().setRandom();
    417   *p2 = *p1;
    418   *p3 = *p1;
    419 
    420   VERIFY_IS_APPROX(p1->matrix(), p2->matrix());
    421   VERIFY_IS_APPROX(p1->matrix(), p3->matrix());
    422 
    423   VERIFY_IS_APPROX( (*p1) * (*p1), (*p2)*(*p3));
    424 
    425   #if defined(EIGEN_VECTORIZE) && EIGEN_ALIGN_STATICALLY
    426   if(internal::packet_traits<Scalar>::Vectorizable)
    427     VERIFY_RAISES_ASSERT((::new(reinterpret_cast<void*>(array3u)) Projective3a));
    428   #endif
    429 }
    430 
    431 template<typename Scalar, int Dim, int Options> void transform_products()
    432 {
    433   typedef Matrix<Scalar,Dim+1,Dim+1> Mat;
    434   typedef Transform<Scalar,Dim,Projective,Options> Proj;
    435   typedef Transform<Scalar,Dim,Affine,Options> Aff;
    436   typedef Transform<Scalar,Dim,AffineCompact,Options> AffC;
    437 
    438   Proj p; p.matrix().setRandom();
    439   Aff a; a.linear().setRandom(); a.translation().setRandom();
    440   AffC ac = a;
    441 
    442   Mat p_m(p.matrix()), a_m(a.matrix());
    443 
    444   VERIFY_IS_APPROX((p*p).matrix(), p_m*p_m);
    445   VERIFY_IS_APPROX((a*a).matrix(), a_m*a_m);
    446   VERIFY_IS_APPROX((p*a).matrix(), p_m*a_m);
    447   VERIFY_IS_APPROX((a*p).matrix(), a_m*p_m);
    448   VERIFY_IS_APPROX((ac*a).matrix(), a_m*a_m);
    449   VERIFY_IS_APPROX((a*ac).matrix(), a_m*a_m);
    450   VERIFY_IS_APPROX((p*ac).matrix(), p_m*a_m);
    451   VERIFY_IS_APPROX((ac*p).matrix(), a_m*p_m);
    452 }
    453 
    454 void test_geo_transformations()
    455 {
    456   for(int i = 0; i < g_repeat; i++) {
    457     CALL_SUBTEST_1(( transformations<double,Affine,AutoAlign>() ));
    458     CALL_SUBTEST_1(( non_projective_only<double,Affine,AutoAlign>() ));
    459 
    460     CALL_SUBTEST_2(( transformations<float,AffineCompact,AutoAlign>() ));
    461     CALL_SUBTEST_2(( non_projective_only<float,AffineCompact,AutoAlign>() ));
    462     CALL_SUBTEST_2(( transform_alignment<float>() ));
    463 
    464     CALL_SUBTEST_3(( transformations<double,Projective,AutoAlign>() ));
    465     CALL_SUBTEST_3(( transformations<double,Projective,DontAlign>() ));
    466     CALL_SUBTEST_3(( transform_alignment<double>() ));
    467 
    468     CALL_SUBTEST_4(( transformations<float,Affine,RowMajor|AutoAlign>() ));
    469     CALL_SUBTEST_4(( non_projective_only<float,Affine,RowMajor>() ));
    470 
    471     CALL_SUBTEST_5(( transformations<double,AffineCompact,RowMajor|AutoAlign>() ));
    472     CALL_SUBTEST_5(( non_projective_only<double,AffineCompact,RowMajor>() ));
    473 
    474     CALL_SUBTEST_6(( transformations<double,Projective,RowMajor|AutoAlign>() ));
    475     CALL_SUBTEST_6(( transformations<double,Projective,RowMajor|DontAlign>() ));
    476 
    477 
    478     CALL_SUBTEST_7(( transform_products<double,3,RowMajor|AutoAlign>() ));
    479     CALL_SUBTEST_7(( transform_products<float,2,AutoAlign>() ));
    480   }
    481 }
    482