You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

487 lines
17 KiB

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