svd_common.h 19 KB

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