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.

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