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.

355 lines
12 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. // discard stack allocation as that too bypasses malloc
  11. #define EIGEN_STACK_ALLOCATION_LIMIT 0
  12. #define EIGEN_RUNTIME_NO_MALLOC
  13. #include "main.h"
  14. #include <Eigen/SVD>
  15. template<typename MatrixType, int QRPreconditioner>
  16. void jacobisvd_check_full(const MatrixType& m, const JacobiSVD<MatrixType, QRPreconditioner>& svd)
  17. {
  18. typedef typename MatrixType::Index Index;
  19. Index rows = m.rows();
  20. Index cols = m.cols();
  21. enum {
  22. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  23. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  24. };
  25. typedef typename MatrixType::Scalar Scalar;
  26. typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixUType;
  27. typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime> MatrixVType;
  28. MatrixType sigma = MatrixType::Zero(rows,cols);
  29. sigma.diagonal() = svd.singularValues().template cast<Scalar>();
  30. MatrixUType u = svd.matrixU();
  31. MatrixVType v = svd.matrixV();
  32. VERIFY_IS_APPROX(m, u * sigma * v.adjoint());
  33. VERIFY_IS_UNITARY(u);
  34. VERIFY_IS_UNITARY(v);
  35. }
  36. template<typename MatrixType, int QRPreconditioner>
  37. void jacobisvd_compare_to_full(const MatrixType& m,
  38. unsigned int computationOptions,
  39. const JacobiSVD<MatrixType, QRPreconditioner>& referenceSvd)
  40. {
  41. typedef typename MatrixType::Index Index;
  42. Index rows = m.rows();
  43. Index cols = m.cols();
  44. Index diagSize = (std::min)(rows, cols);
  45. JacobiSVD<MatrixType, QRPreconditioner> svd(m, computationOptions);
  46. VERIFY_IS_APPROX(svd.singularValues(), referenceSvd.singularValues());
  47. if(computationOptions & ComputeFullU)
  48. VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU());
  49. if(computationOptions & ComputeThinU)
  50. VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU().leftCols(diagSize));
  51. if(computationOptions & ComputeFullV)
  52. VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV());
  53. if(computationOptions & ComputeThinV)
  54. VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV().leftCols(diagSize));
  55. }
  56. template<typename MatrixType, int QRPreconditioner>
  57. void jacobisvd_solve(const MatrixType& m, unsigned int computationOptions)
  58. {
  59. typedef typename MatrixType::Scalar Scalar;
  60. typedef typename MatrixType::Index Index;
  61. Index rows = m.rows();
  62. Index cols = m.cols();
  63. enum {
  64. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  65. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  66. };
  67. typedef Matrix<Scalar, RowsAtCompileTime, Dynamic> RhsType;
  68. typedef Matrix<Scalar, ColsAtCompileTime, Dynamic> SolutionType;
  69. RhsType rhs = RhsType::Random(rows, internal::random<Index>(1, cols));
  70. JacobiSVD<MatrixType, QRPreconditioner> svd(m, computationOptions);
  71. SolutionType x = svd.solve(rhs);
  72. // evaluate normal equation which works also for least-squares solutions
  73. VERIFY_IS_APPROX(m.adjoint()*m*x,m.adjoint()*rhs);
  74. }
  75. template<typename MatrixType, int QRPreconditioner>
  76. void jacobisvd_test_all_computation_options(const MatrixType& m)
  77. {
  78. if (QRPreconditioner == NoQRPreconditioner && m.rows() != m.cols())
  79. return;
  80. JacobiSVD<MatrixType, QRPreconditioner> fullSvd(m, ComputeFullU|ComputeFullV);
  81. jacobisvd_check_full(m, fullSvd);
  82. jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeFullU | ComputeFullV);
  83. #if defined __INTEL_COMPILER
  84. // remark #111: statement is unreachable
  85. #pragma warning disable 111
  86. #endif
  87. if(QRPreconditioner == FullPivHouseholderQRPreconditioner)
  88. return;
  89. jacobisvd_compare_to_full(m, ComputeFullU, fullSvd);
  90. jacobisvd_compare_to_full(m, ComputeFullV, fullSvd);
  91. jacobisvd_compare_to_full(m, 0, fullSvd);
  92. if (MatrixType::ColsAtCompileTime == Dynamic) {
  93. // thin U/V are only available with dynamic number of columns
  94. jacobisvd_compare_to_full(m, ComputeFullU|ComputeThinV, fullSvd);
  95. jacobisvd_compare_to_full(m, ComputeThinV, fullSvd);
  96. jacobisvd_compare_to_full(m, ComputeThinU|ComputeFullV, fullSvd);
  97. jacobisvd_compare_to_full(m, ComputeThinU , fullSvd);
  98. jacobisvd_compare_to_full(m, ComputeThinU|ComputeThinV, fullSvd);
  99. jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeFullU | ComputeThinV);
  100. jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeThinU | ComputeFullV);
  101. jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeThinU | ComputeThinV);
  102. // test reconstruction
  103. typedef typename MatrixType::Index Index;
  104. Index diagSize = (std::min)(m.rows(), m.cols());
  105. JacobiSVD<MatrixType, QRPreconditioner> svd(m, ComputeThinU | ComputeThinV);
  106. VERIFY_IS_APPROX(m, svd.matrixU().leftCols(diagSize) * svd.singularValues().asDiagonal() * svd.matrixV().leftCols(diagSize).adjoint());
  107. }
  108. }
  109. template<typename MatrixType>
  110. void jacobisvd(const MatrixType& a = MatrixType(), bool pickrandom = true)
  111. {
  112. MatrixType m = pickrandom ? MatrixType::Random(a.rows(), a.cols()) : a;
  113. jacobisvd_test_all_computation_options<MatrixType, FullPivHouseholderQRPreconditioner>(m);
  114. jacobisvd_test_all_computation_options<MatrixType, ColPivHouseholderQRPreconditioner>(m);
  115. jacobisvd_test_all_computation_options<MatrixType, HouseholderQRPreconditioner>(m);
  116. jacobisvd_test_all_computation_options<MatrixType, NoQRPreconditioner>(m);
  117. }
  118. template<typename MatrixType> void jacobisvd_verify_assert(const MatrixType& m)
  119. {
  120. typedef typename MatrixType::Scalar Scalar;
  121. typedef typename MatrixType::Index Index;
  122. Index rows = m.rows();
  123. Index cols = m.cols();
  124. enum {
  125. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  126. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  127. };
  128. typedef Matrix<Scalar, RowsAtCompileTime, 1> RhsType;
  129. RhsType rhs(rows);
  130. JacobiSVD<MatrixType> svd;
  131. VERIFY_RAISES_ASSERT(svd.matrixU())
  132. VERIFY_RAISES_ASSERT(svd.singularValues())
  133. VERIFY_RAISES_ASSERT(svd.matrixV())
  134. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  135. MatrixType a = MatrixType::Zero(rows, cols);
  136. a.setZero();
  137. svd.compute(a, 0);
  138. VERIFY_RAISES_ASSERT(svd.matrixU())
  139. VERIFY_RAISES_ASSERT(svd.matrixV())
  140. svd.singularValues();
  141. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  142. if (ColsAtCompileTime == Dynamic)
  143. {
  144. svd.compute(a, ComputeThinU);
  145. svd.matrixU();
  146. VERIFY_RAISES_ASSERT(svd.matrixV())
  147. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  148. svd.compute(a, ComputeThinV);
  149. svd.matrixV();
  150. VERIFY_RAISES_ASSERT(svd.matrixU())
  151. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  152. JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner> svd_fullqr;
  153. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeFullU|ComputeThinV))
  154. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeThinV))
  155. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeFullV))
  156. }
  157. else
  158. {
  159. VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinU))
  160. VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinV))
  161. }
  162. }
  163. template<typename MatrixType>
  164. void jacobisvd_method()
  165. {
  166. enum { Size = MatrixType::RowsAtCompileTime };
  167. typedef typename MatrixType::RealScalar RealScalar;
  168. typedef Matrix<RealScalar, Size, 1> RealVecType;
  169. MatrixType m = MatrixType::Identity();
  170. VERIFY_IS_APPROX(m.jacobiSvd().singularValues(), RealVecType::Ones());
  171. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixU());
  172. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixV());
  173. VERIFY_IS_APPROX(m.jacobiSvd(ComputeFullU|ComputeFullV).solve(m), m);
  174. }
  175. // work around stupid msvc error when constructing at compile time an expression that involves
  176. // a division by zero, even if the numeric type has floating point
  177. template<typename Scalar>
  178. EIGEN_DONT_INLINE Scalar zero() { return Scalar(0); }
  179. // workaround aggressive optimization in ICC
  180. template<typename T> EIGEN_DONT_INLINE T sub(T a, T b) { return a - b; }
  181. template<typename MatrixType>
  182. void jacobisvd_inf_nan()
  183. {
  184. // all this function does is verify we don't iterate infinitely on nan/inf values
  185. JacobiSVD<MatrixType> svd;
  186. typedef typename MatrixType::Scalar Scalar;
  187. Scalar some_inf = Scalar(1) / zero<Scalar>();
  188. VERIFY(sub(some_inf, some_inf) != sub(some_inf, some_inf));
  189. svd.compute(MatrixType::Constant(10,10,some_inf), ComputeFullU | ComputeFullV);
  190. Scalar some_nan = zero<Scalar>() / zero<Scalar>();
  191. VERIFY(some_nan != some_nan);
  192. svd.compute(MatrixType::Constant(10,10,some_nan), ComputeFullU | ComputeFullV);
  193. MatrixType m = MatrixType::Zero(10,10);
  194. m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_inf;
  195. svd.compute(m, ComputeFullU | ComputeFullV);
  196. m = MatrixType::Zero(10,10);
  197. m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_nan;
  198. svd.compute(m, ComputeFullU | ComputeFullV);
  199. }
  200. // Regression test for bug 286: JacobiSVD loops indefinitely with some
  201. // matrices containing denormal numbers.
  202. void jacobisvd_bug286()
  203. {
  204. #if defined __INTEL_COMPILER
  205. // shut up warning #239: floating point underflow
  206. #pragma warning push
  207. #pragma warning disable 239
  208. #endif
  209. Matrix2d M;
  210. M << -7.90884e-313, -4.94e-324,
  211. 0, 5.60844e-313;
  212. #if defined __INTEL_COMPILER
  213. #pragma warning pop
  214. #endif
  215. JacobiSVD<Matrix2d> svd;
  216. svd.compute(M); // just check we don't loop indefinitely
  217. }
  218. void jacobisvd_preallocate()
  219. {
  220. Vector3f v(3.f, 2.f, 1.f);
  221. MatrixXf m = v.asDiagonal();
  222. internal::set_is_malloc_allowed(false);
  223. VERIFY_RAISES_ASSERT(VectorXf tmp(10);)
  224. JacobiSVD<MatrixXf> svd;
  225. internal::set_is_malloc_allowed(true);
  226. svd.compute(m);
  227. VERIFY_IS_APPROX(svd.singularValues(), v);
  228. JacobiSVD<MatrixXf> svd2(3,3);
  229. internal::set_is_malloc_allowed(false);
  230. svd2.compute(m);
  231. internal::set_is_malloc_allowed(true);
  232. VERIFY_IS_APPROX(svd2.singularValues(), v);
  233. VERIFY_RAISES_ASSERT(svd2.matrixU());
  234. VERIFY_RAISES_ASSERT(svd2.matrixV());
  235. svd2.compute(m, ComputeFullU | ComputeFullV);
  236. VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
  237. VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
  238. internal::set_is_malloc_allowed(false);
  239. svd2.compute(m);
  240. internal::set_is_malloc_allowed(true);
  241. JacobiSVD<MatrixXf> svd3(3,3,ComputeFullU|ComputeFullV);
  242. internal::set_is_malloc_allowed(false);
  243. svd2.compute(m);
  244. internal::set_is_malloc_allowed(true);
  245. VERIFY_IS_APPROX(svd2.singularValues(), v);
  246. VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
  247. VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
  248. internal::set_is_malloc_allowed(false);
  249. svd2.compute(m, ComputeFullU|ComputeFullV);
  250. internal::set_is_malloc_allowed(true);
  251. }
  252. void test_jacobisvd()
  253. {
  254. CALL_SUBTEST_3(( jacobisvd_verify_assert(Matrix3f()) ));
  255. CALL_SUBTEST_4(( jacobisvd_verify_assert(Matrix4d()) ));
  256. CALL_SUBTEST_7(( jacobisvd_verify_assert(MatrixXf(10,12)) ));
  257. CALL_SUBTEST_8(( jacobisvd_verify_assert(MatrixXcd(7,5)) ));
  258. for(int i = 0; i < g_repeat; i++) {
  259. Matrix2cd m;
  260. m << 0, 1,
  261. 0, 1;
  262. CALL_SUBTEST_1(( jacobisvd(m, false) ));
  263. m << 1, 0,
  264. 1, 0;
  265. CALL_SUBTEST_1(( jacobisvd(m, false) ));
  266. Matrix2d n;
  267. n << 0, 0,
  268. 0, 0;
  269. CALL_SUBTEST_2(( jacobisvd(n, false) ));
  270. n << 0, 0,
  271. 0, 1;
  272. CALL_SUBTEST_2(( jacobisvd(n, false) ));
  273. CALL_SUBTEST_3(( jacobisvd<Matrix3f>() ));
  274. CALL_SUBTEST_4(( jacobisvd<Matrix4d>() ));
  275. CALL_SUBTEST_5(( jacobisvd<Matrix<float,3,5> >() ));
  276. CALL_SUBTEST_6(( jacobisvd<Matrix<double,Dynamic,2> >(Matrix<double,Dynamic,2>(10,2)) ));
  277. int r = internal::random<int>(1, 30),
  278. c = internal::random<int>(1, 30);
  279. TEST_SET_BUT_UNUSED_VARIABLE(r)
  280. TEST_SET_BUT_UNUSED_VARIABLE(c)
  281. CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(r,c)) ));
  282. CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(r,c)) ));
  283. (void) r;
  284. (void) c;
  285. // Test on inf/nan matrix
  286. CALL_SUBTEST_7( jacobisvd_inf_nan<MatrixXf>() );
  287. }
  288. CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/2), internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/2))) ));
  289. CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/3), internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/3))) ));
  290. // test matrixbase method
  291. CALL_SUBTEST_1(( jacobisvd_method<Matrix2cd>() ));
  292. CALL_SUBTEST_3(( jacobisvd_method<Matrix3f>() ));
  293. // Test problem size constructors
  294. CALL_SUBTEST_7( JacobiSVD<MatrixXf>(10,10) );
  295. // Check that preallocation avoids subsequent mallocs
  296. CALL_SUBTEST_9( jacobisvd_preallocate() );
  297. // Regression check for bug 286
  298. CALL_SUBTEST_2( jacobisvd_bug286() );
  299. }