svd_common.h 18 KB

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