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.

478 lines
18 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2014 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. #ifndef SVD_DEFAULT
  11. #error a macro SVD_DEFAULT(MatrixType) must be defined prior to including svd_common.h
  12. #endif
  13. #ifndef SVD_FOR_MIN_NORM
  14. #error a macro SVD_FOR_MIN_NORM(MatrixType) must be defined prior to including svd_common.h
  15. #endif
  16. #include "svd_fill.h"
  17. // Check that the matrix m is properly reconstructed and that the U and V factors are unitary
  18. // The SVD must have already been computed.
  19. template<typename SvdType, typename MatrixType>
  20. void svd_check_full(const MatrixType& m, const SvdType& svd)
  21. {
  22. typedef typename MatrixType::Index Index;
  23. Index rows = m.rows();
  24. Index cols = m.cols();
  25. enum {
  26. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  27. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  28. };
  29. typedef typename MatrixType::Scalar Scalar;
  30. typedef typename MatrixType::RealScalar RealScalar;
  31. typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixUType;
  32. typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime> MatrixVType;
  33. MatrixType sigma = MatrixType::Zero(rows,cols);
  34. sigma.diagonal() = svd.singularValues().template cast<Scalar>();
  35. MatrixUType u = svd.matrixU();
  36. MatrixVType v = svd.matrixV();
  37. RealScalar scaling = m.cwiseAbs().maxCoeff();
  38. if(scaling<=(std::numeric_limits<RealScalar>::min)())
  39. scaling = RealScalar(1);
  40. VERIFY_IS_APPROX(m/scaling, u * (sigma/scaling) * v.adjoint());
  41. VERIFY_IS_UNITARY(u);
  42. VERIFY_IS_UNITARY(v);
  43. }
  44. // Compare partial SVD defined by computationOptions to a full SVD referenceSvd
  45. template<typename SvdType, typename MatrixType>
  46. void svd_compare_to_full(const MatrixType& m,
  47. unsigned int computationOptions,
  48. const SvdType& referenceSvd)
  49. {
  50. typedef typename MatrixType::RealScalar RealScalar;
  51. Index rows = m.rows();
  52. Index cols = m.cols();
  53. Index diagSize = (std::min)(rows, cols);
  54. RealScalar prec = test_precision<RealScalar>();
  55. SvdType svd(m, computationOptions);
  56. VERIFY_IS_APPROX(svd.singularValues(), referenceSvd.singularValues());
  57. if(computationOptions & (ComputeFullV|ComputeThinV))
  58. {
  59. VERIFY( (svd.matrixV().adjoint()*svd.matrixV()).isIdentity(prec) );
  60. VERIFY_IS_APPROX( svd.matrixV().leftCols(diagSize) * svd.singularValues().asDiagonal() * svd.matrixV().leftCols(diagSize).adjoint(),
  61. referenceSvd.matrixV().leftCols(diagSize) * referenceSvd.singularValues().asDiagonal() * referenceSvd.matrixV().leftCols(diagSize).adjoint());
  62. }
  63. if(computationOptions & (ComputeFullU|ComputeThinU))
  64. {
  65. VERIFY( (svd.matrixU().adjoint()*svd.matrixU()).isIdentity(prec) );
  66. VERIFY_IS_APPROX( svd.matrixU().leftCols(diagSize) * svd.singularValues().cwiseAbs2().asDiagonal() * svd.matrixU().leftCols(diagSize).adjoint(),
  67. referenceSvd.matrixU().leftCols(diagSize) * referenceSvd.singularValues().cwiseAbs2().asDiagonal() * referenceSvd.matrixU().leftCols(diagSize).adjoint());
  68. }
  69. // The following checks are not critical.
  70. // For instance, with Dived&Conquer SVD, if only the factor 'V' is computedt then different matrix-matrix product implementation will be used
  71. // and the resulting 'V' factor might be significantly different when the SVD decomposition is not unique, especially with single precision float.
  72. ++g_test_level;
  73. if(computationOptions & ComputeFullU) VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU());
  74. if(computationOptions & ComputeThinU) VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU().leftCols(diagSize));
  75. if(computationOptions & ComputeFullV) VERIFY_IS_APPROX(svd.matrixV().cwiseAbs(), referenceSvd.matrixV().cwiseAbs());
  76. if(computationOptions & ComputeThinV) VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV().leftCols(diagSize));
  77. --g_test_level;
  78. }
  79. //
  80. template<typename SvdType, typename MatrixType>
  81. void svd_least_square(const MatrixType& m, unsigned int computationOptions)
  82. {
  83. typedef typename MatrixType::Scalar Scalar;
  84. typedef typename MatrixType::RealScalar RealScalar;
  85. typedef typename MatrixType::Index Index;
  86. Index rows = m.rows();
  87. Index cols = m.cols();
  88. enum {
  89. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  90. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  91. };
  92. typedef Matrix<Scalar, RowsAtCompileTime, Dynamic> RhsType;
  93. typedef Matrix<Scalar, ColsAtCompileTime, Dynamic> SolutionType;
  94. RhsType rhs = RhsType::Random(rows, internal::random<Index>(1, cols));
  95. SvdType svd(m, computationOptions);
  96. if(internal::is_same<RealScalar,double>::value) svd.setThreshold(1e-8);
  97. else if(internal::is_same<RealScalar,float>::value) svd.setThreshold(2e-4);
  98. SolutionType x = svd.solve(rhs);
  99. RealScalar residual = (m*x-rhs).norm();
  100. RealScalar rhs_norm = rhs.norm();
  101. if(!test_isMuchSmallerThan(residual,rhs.norm()))
  102. {
  103. // ^^^ If the residual is very small, then we have an exact solution, so we are already good.
  104. // evaluate normal equation which works also for least-squares solutions
  105. if(internal::is_same<RealScalar,double>::value || svd.rank()==m.diagonal().size())
  106. {
  107. using std::sqrt;
  108. // This test is not stable with single precision.
  109. // This is probably because squaring m signicantly affects the precision.
  110. if(internal::is_same<RealScalar,float>::value) ++g_test_level;
  111. VERIFY_IS_APPROX(m.adjoint()*(m*x),m.adjoint()*rhs);
  112. if(internal::is_same<RealScalar,float>::value) --g_test_level;
  113. }
  114. // Check that there is no significantly better solution in the neighborhood of x
  115. for(Index k=0;k<x.rows();++k)
  116. {
  117. using std::abs;
  118. SolutionType y(x);
  119. y.row(k) = (1.+2*NumTraits<RealScalar>::epsilon())*x.row(k);
  120. RealScalar residual_y = (m*y-rhs).norm();
  121. VERIFY( test_isMuchSmallerThan(abs(residual_y-residual), rhs_norm) || residual < residual_y );
  122. if(internal::is_same<RealScalar,float>::value) ++g_test_level;
  123. VERIFY( test_isApprox(residual_y,residual) || residual < residual_y );
  124. if(internal::is_same<RealScalar,float>::value) --g_test_level;
  125. y.row(k) = (1.-2*NumTraits<RealScalar>::epsilon())*x.row(k);
  126. residual_y = (m*y-rhs).norm();
  127. VERIFY( test_isMuchSmallerThan(abs(residual_y-residual), rhs_norm) || residual < residual_y );
  128. if(internal::is_same<RealScalar,float>::value) ++g_test_level;
  129. VERIFY( test_isApprox(residual_y,residual) || residual < residual_y );
  130. if(internal::is_same<RealScalar,float>::value) --g_test_level;
  131. }
  132. }
  133. }
  134. // check minimal norm solutions, the inoput matrix m is only used to recover problem size
  135. template<typename MatrixType>
  136. void svd_min_norm(const MatrixType& m, unsigned int computationOptions)
  137. {
  138. typedef typename MatrixType::Scalar Scalar;
  139. typedef typename MatrixType::Index Index;
  140. Index cols = m.cols();
  141. enum {
  142. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  143. };
  144. typedef Matrix<Scalar, ColsAtCompileTime, Dynamic> SolutionType;
  145. // generate a full-rank m x n problem with m<n
  146. enum {
  147. RankAtCompileTime2 = ColsAtCompileTime==Dynamic ? Dynamic : (ColsAtCompileTime)/2+1,
  148. RowsAtCompileTime3 = ColsAtCompileTime==Dynamic ? Dynamic : ColsAtCompileTime+1
  149. };
  150. typedef Matrix<Scalar, RankAtCompileTime2, ColsAtCompileTime> MatrixType2;
  151. typedef Matrix<Scalar, RankAtCompileTime2, 1> RhsType2;
  152. typedef Matrix<Scalar, ColsAtCompileTime, RankAtCompileTime2> MatrixType2T;
  153. Index rank = RankAtCompileTime2==Dynamic ? internal::random<Index>(1,cols) : Index(RankAtCompileTime2);
  154. MatrixType2 m2(rank,cols);
  155. int guard = 0;
  156. do {
  157. m2.setRandom();
  158. } while(SVD_FOR_MIN_NORM(MatrixType2)(m2).setThreshold(test_precision<Scalar>()).rank()!=rank && (++guard)<10);
  159. VERIFY(guard<10);
  160. RhsType2 rhs2 = RhsType2::Random(rank);
  161. // use QR to find a reference minimal norm solution
  162. HouseholderQR<MatrixType2T> qr(m2.adjoint());
  163. Matrix<Scalar,Dynamic,1> tmp = qr.matrixQR().topLeftCorner(rank,rank).template triangularView<Upper>().adjoint().solve(rhs2);
  164. tmp.conservativeResize(cols);
  165. tmp.tail(cols-rank).setZero();
  166. SolutionType x21 = qr.householderQ() * tmp;
  167. // now check with SVD
  168. SVD_FOR_MIN_NORM(MatrixType2) svd2(m2, computationOptions);
  169. SolutionType x22 = svd2.solve(rhs2);
  170. VERIFY_IS_APPROX(m2*x21, rhs2);
  171. VERIFY_IS_APPROX(m2*x22, rhs2);
  172. VERIFY_IS_APPROX(x21, x22);
  173. // Now check with a rank deficient matrix
  174. typedef Matrix<Scalar, RowsAtCompileTime3, ColsAtCompileTime> MatrixType3;
  175. typedef Matrix<Scalar, RowsAtCompileTime3, 1> RhsType3;
  176. Index rows3 = RowsAtCompileTime3==Dynamic ? internal::random<Index>(rank+1,2*cols) : Index(RowsAtCompileTime3);
  177. Matrix<Scalar,RowsAtCompileTime3,Dynamic> C = Matrix<Scalar,RowsAtCompileTime3,Dynamic>::Random(rows3,rank);
  178. MatrixType3 m3 = C * m2;
  179. RhsType3 rhs3 = C * rhs2;
  180. SVD_FOR_MIN_NORM(MatrixType3) svd3(m3, computationOptions);
  181. SolutionType x3 = svd3.solve(rhs3);
  182. VERIFY_IS_APPROX(m3*x3, rhs3);
  183. VERIFY_IS_APPROX(m3*x21, rhs3);
  184. VERIFY_IS_APPROX(m2*x3, rhs2);
  185. VERIFY_IS_APPROX(x21, x3);
  186. }
  187. // Check full, compare_to_full, least_square, and min_norm for all possible compute-options
  188. template<typename SvdType, typename MatrixType>
  189. void svd_test_all_computation_options(const MatrixType& m, bool full_only)
  190. {
  191. // if (QRPreconditioner == NoQRPreconditioner && m.rows() != m.cols())
  192. // return;
  193. SvdType fullSvd(m, ComputeFullU|ComputeFullV);
  194. CALL_SUBTEST(( svd_check_full(m, fullSvd) ));
  195. CALL_SUBTEST(( svd_least_square<SvdType>(m, ComputeFullU | ComputeFullV) ));
  196. CALL_SUBTEST(( svd_min_norm(m, ComputeFullU | ComputeFullV) ));
  197. #if defined __INTEL_COMPILER
  198. // remark #111: statement is unreachable
  199. #pragma warning disable 111
  200. #endif
  201. if(full_only)
  202. return;
  203. CALL_SUBTEST(( svd_compare_to_full(m, ComputeFullU, fullSvd) ));
  204. CALL_SUBTEST(( svd_compare_to_full(m, ComputeFullV, fullSvd) ));
  205. CALL_SUBTEST(( svd_compare_to_full(m, 0, fullSvd) ));
  206. if (MatrixType::ColsAtCompileTime == Dynamic) {
  207. // thin U/V are only available with dynamic number of columns
  208. CALL_SUBTEST(( svd_compare_to_full(m, ComputeFullU|ComputeThinV, fullSvd) ));
  209. CALL_SUBTEST(( svd_compare_to_full(m, ComputeThinV, fullSvd) ));
  210. CALL_SUBTEST(( svd_compare_to_full(m, ComputeThinU|ComputeFullV, fullSvd) ));
  211. CALL_SUBTEST(( svd_compare_to_full(m, ComputeThinU , fullSvd) ));
  212. CALL_SUBTEST(( svd_compare_to_full(m, ComputeThinU|ComputeThinV, fullSvd) ));
  213. CALL_SUBTEST(( svd_least_square<SvdType>(m, ComputeFullU | ComputeThinV) ));
  214. CALL_SUBTEST(( svd_least_square<SvdType>(m, ComputeThinU | ComputeFullV) ));
  215. CALL_SUBTEST(( svd_least_square<SvdType>(m, ComputeThinU | ComputeThinV) ));
  216. CALL_SUBTEST(( svd_min_norm(m, ComputeFullU | ComputeThinV) ));
  217. CALL_SUBTEST(( svd_min_norm(m, ComputeThinU | ComputeFullV) ));
  218. CALL_SUBTEST(( svd_min_norm(m, ComputeThinU | ComputeThinV) ));
  219. // test reconstruction
  220. typedef typename MatrixType::Index Index;
  221. Index diagSize = (std::min)(m.rows(), m.cols());
  222. SvdType svd(m, ComputeThinU | ComputeThinV);
  223. VERIFY_IS_APPROX(m, svd.matrixU().leftCols(diagSize) * svd.singularValues().asDiagonal() * svd.matrixV().leftCols(diagSize).adjoint());
  224. }
  225. }
  226. // work around stupid msvc error when constructing at compile time an expression that involves
  227. // a division by zero, even if the numeric type has floating point
  228. template<typename Scalar>
  229. EIGEN_DONT_INLINE Scalar zero() { return Scalar(0); }
  230. // workaround aggressive optimization in ICC
  231. template<typename T> EIGEN_DONT_INLINE T sub(T a, T b) { return a - b; }
  232. // all this function does is verify we don't iterate infinitely on nan/inf values
  233. template<typename SvdType, typename MatrixType>
  234. void svd_inf_nan()
  235. {
  236. SvdType svd;
  237. typedef typename MatrixType::Scalar Scalar;
  238. Scalar some_inf = Scalar(1) / zero<Scalar>();
  239. VERIFY(sub(some_inf, some_inf) != sub(some_inf, some_inf));
  240. svd.compute(MatrixType::Constant(10,10,some_inf), ComputeFullU | ComputeFullV);
  241. Scalar nan = std::numeric_limits<Scalar>::quiet_NaN();
  242. VERIFY(nan != nan);
  243. svd.compute(MatrixType::Constant(10,10,nan), ComputeFullU | ComputeFullV);
  244. MatrixType m = MatrixType::Zero(10,10);
  245. m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_inf;
  246. svd.compute(m, ComputeFullU | ComputeFullV);
  247. m = MatrixType::Zero(10,10);
  248. m(internal::random<int>(0,9), internal::random<int>(0,9)) = nan;
  249. svd.compute(m, ComputeFullU | ComputeFullV);
  250. // regression test for bug 791
  251. m.resize(3,3);
  252. m << 0, 2*NumTraits<Scalar>::epsilon(), 0.5,
  253. 0, -0.5, 0,
  254. nan, 0, 0;
  255. svd.compute(m, ComputeFullU | ComputeFullV);
  256. m.resize(4,4);
  257. m << 1, 0, 0, 0,
  258. 0, 3, 1, 2e-308,
  259. 1, 0, 1, nan,
  260. 0, nan, nan, 0;
  261. svd.compute(m, ComputeFullU | ComputeFullV);
  262. }
  263. // Regression test for bug 286: JacobiSVD loops indefinitely with some
  264. // matrices containing denormal numbers.
  265. template<typename>
  266. void svd_underoverflow()
  267. {
  268. #if defined __INTEL_COMPILER
  269. // shut up warning #239: floating point underflow
  270. #pragma warning push
  271. #pragma warning disable 239
  272. #endif
  273. Matrix2d M;
  274. M << -7.90884e-313, -4.94e-324,
  275. 0, 5.60844e-313;
  276. SVD_DEFAULT(Matrix2d) svd;
  277. svd.compute(M,ComputeFullU|ComputeFullV);
  278. CALL_SUBTEST( svd_check_full(M,svd) );
  279. // Check all 2x2 matrices made with the following coefficients:
  280. VectorXd value_set(9);
  281. value_set << 0, 1, -1, 5.60844e-313, -5.60844e-313, 4.94e-324, -4.94e-324, -4.94e-223, 4.94e-223;
  282. Array4i id(0,0,0,0);
  283. int k = 0;
  284. do
  285. {
  286. M << value_set(id(0)), value_set(id(1)), value_set(id(2)), value_set(id(3));
  287. svd.compute(M,ComputeFullU|ComputeFullV);
  288. CALL_SUBTEST( svd_check_full(M,svd) );
  289. id(k)++;
  290. if(id(k)>=value_set.size())
  291. {
  292. while(k<3 && id(k)>=value_set.size()) id(++k)++;
  293. id.head(k).setZero();
  294. k=0;
  295. }
  296. } while((id<int(value_set.size())).all());
  297. #if defined __INTEL_COMPILER
  298. #pragma warning pop
  299. #endif
  300. // Check for overflow:
  301. Matrix3d M3;
  302. M3 << 4.4331978442502944e+307, -5.8585363752028680e+307, 6.4527017443412964e+307,
  303. 3.7841695601406358e+307, 2.4331702789740617e+306, -3.5235707140272905e+307,
  304. -8.7190887618028355e+307, -7.3453213709232193e+307, -2.4367363684472105e+307;
  305. SVD_DEFAULT(Matrix3d) svd3;
  306. svd3.compute(M3,ComputeFullU|ComputeFullV); // just check we don't loop indefinitely
  307. CALL_SUBTEST( svd_check_full(M3,svd3) );
  308. }
  309. // void jacobisvd(const MatrixType& a = MatrixType(), bool pickrandom = true)
  310. template<typename MatrixType>
  311. void svd_all_trivial_2x2( void (*cb)(const MatrixType&,bool) )
  312. {
  313. MatrixType M;
  314. VectorXd value_set(3);
  315. value_set << 0, 1, -1;
  316. Array4i id(0,0,0,0);
  317. int k = 0;
  318. do
  319. {
  320. M << value_set(id(0)), value_set(id(1)), value_set(id(2)), value_set(id(3));
  321. cb(M,false);
  322. id(k)++;
  323. if(id(k)>=value_set.size())
  324. {
  325. while(k<3 && id(k)>=value_set.size()) id(++k)++;
  326. id.head(k).setZero();
  327. k=0;
  328. }
  329. } while((id<int(value_set.size())).all());
  330. }
  331. template<typename>
  332. void svd_preallocate()
  333. {
  334. Vector3f v(3.f, 2.f, 1.f);
  335. MatrixXf m = v.asDiagonal();
  336. internal::set_is_malloc_allowed(false);
  337. VERIFY_RAISES_ASSERT(VectorXf tmp(10);)
  338. SVD_DEFAULT(MatrixXf) svd;
  339. internal::set_is_malloc_allowed(true);
  340. svd.compute(m);
  341. VERIFY_IS_APPROX(svd.singularValues(), v);
  342. SVD_DEFAULT(MatrixXf) svd2(3,3);
  343. internal::set_is_malloc_allowed(false);
  344. svd2.compute(m);
  345. internal::set_is_malloc_allowed(true);
  346. VERIFY_IS_APPROX(svd2.singularValues(), v);
  347. VERIFY_RAISES_ASSERT(svd2.matrixU());
  348. VERIFY_RAISES_ASSERT(svd2.matrixV());
  349. svd2.compute(m, ComputeFullU | ComputeFullV);
  350. VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
  351. VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
  352. internal::set_is_malloc_allowed(false);
  353. svd2.compute(m);
  354. internal::set_is_malloc_allowed(true);
  355. SVD_DEFAULT(MatrixXf) svd3(3,3,ComputeFullU|ComputeFullV);
  356. internal::set_is_malloc_allowed(false);
  357. svd2.compute(m);
  358. internal::set_is_malloc_allowed(true);
  359. VERIFY_IS_APPROX(svd2.singularValues(), v);
  360. VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
  361. VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
  362. internal::set_is_malloc_allowed(false);
  363. svd2.compute(m, ComputeFullU|ComputeFullV);
  364. internal::set_is_malloc_allowed(true);
  365. }
  366. template<typename SvdType,typename MatrixType>
  367. void svd_verify_assert(const MatrixType& m)
  368. {
  369. typedef typename MatrixType::Scalar Scalar;
  370. typedef typename MatrixType::Index Index;
  371. Index rows = m.rows();
  372. Index cols = m.cols();
  373. enum {
  374. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  375. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  376. };
  377. typedef Matrix<Scalar, RowsAtCompileTime, 1> RhsType;
  378. RhsType rhs(rows);
  379. SvdType svd;
  380. VERIFY_RAISES_ASSERT(svd.matrixU())
  381. VERIFY_RAISES_ASSERT(svd.singularValues())
  382. VERIFY_RAISES_ASSERT(svd.matrixV())
  383. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  384. MatrixType a = MatrixType::Zero(rows, cols);
  385. a.setZero();
  386. svd.compute(a, 0);
  387. VERIFY_RAISES_ASSERT(svd.matrixU())
  388. VERIFY_RAISES_ASSERT(svd.matrixV())
  389. svd.singularValues();
  390. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  391. if (ColsAtCompileTime == Dynamic)
  392. {
  393. svd.compute(a, ComputeThinU);
  394. svd.matrixU();
  395. VERIFY_RAISES_ASSERT(svd.matrixV())
  396. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  397. svd.compute(a, ComputeThinV);
  398. svd.matrixV();
  399. VERIFY_RAISES_ASSERT(svd.matrixU())
  400. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  401. }
  402. else
  403. {
  404. VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinU))
  405. VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinV))
  406. }
  407. }
  408. #undef SVD_DEFAULT
  409. #undef SVD_FOR_MIN_NORM