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.

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