cholesky.cpp 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529
  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. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #ifndef EIGEN_NO_ASSERTION_CHECKING
  10. #define EIGEN_NO_ASSERTION_CHECKING
  11. #endif
  12. #define TEST_ENABLE_TEMPORARY_TRACKING
  13. #include "main.h"
  14. #include <Eigen/Cholesky>
  15. #include <Eigen/QR>
  16. template<typename MatrixType, int UpLo>
  17. typename MatrixType::RealScalar matrix_l1_norm(const MatrixType& m) {
  18. if(m.cols()==0) return typename MatrixType::RealScalar(0);
  19. MatrixType symm = m.template selfadjointView<UpLo>();
  20. return symm.cwiseAbs().colwise().sum().maxCoeff();
  21. }
  22. template<typename MatrixType,template <typename,int> class CholType> void test_chol_update(const MatrixType& symm)
  23. {
  24. typedef typename MatrixType::Scalar Scalar;
  25. typedef typename MatrixType::RealScalar RealScalar;
  26. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  27. MatrixType symmLo = symm.template triangularView<Lower>();
  28. MatrixType symmUp = symm.template triangularView<Upper>();
  29. MatrixType symmCpy = symm;
  30. CholType<MatrixType,Lower> chollo(symmLo);
  31. CholType<MatrixType,Upper> cholup(symmUp);
  32. for (int k=0; k<10; ++k)
  33. {
  34. VectorType vec = VectorType::Random(symm.rows());
  35. RealScalar sigma = internal::random<RealScalar>();
  36. symmCpy += sigma * vec * vec.adjoint();
  37. // we are doing some downdates, so it might be the case that the matrix is not SPD anymore
  38. CholType<MatrixType,Lower> chol(symmCpy);
  39. if(chol.info()!=Success)
  40. break;
  41. chollo.rankUpdate(vec, sigma);
  42. VERIFY_IS_APPROX(symmCpy, chollo.reconstructedMatrix());
  43. cholup.rankUpdate(vec, sigma);
  44. VERIFY_IS_APPROX(symmCpy, cholup.reconstructedMatrix());
  45. }
  46. }
  47. template<typename MatrixType> void cholesky(const MatrixType& m)
  48. {
  49. /* this test covers the following files:
  50. LLT.h LDLT.h
  51. */
  52. Index rows = m.rows();
  53. Index cols = m.cols();
  54. typedef typename MatrixType::Scalar Scalar;
  55. typedef typename NumTraits<Scalar>::Real RealScalar;
  56. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
  57. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  58. MatrixType a0 = MatrixType::Random(rows,cols);
  59. VectorType vecB = VectorType::Random(rows), vecX(rows);
  60. MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
  61. SquareMatrixType symm = a0 * a0.adjoint();
  62. // let's make sure the matrix is not singular or near singular
  63. for (int k=0; k<3; ++k)
  64. {
  65. MatrixType a1 = MatrixType::Random(rows,cols);
  66. symm += a1 * a1.adjoint();
  67. }
  68. {
  69. SquareMatrixType symmUp = symm.template triangularView<Upper>();
  70. SquareMatrixType symmLo = symm.template triangularView<Lower>();
  71. LLT<SquareMatrixType,Lower> chollo(symmLo);
  72. VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
  73. vecX = chollo.solve(vecB);
  74. VERIFY_IS_APPROX(symm * vecX, vecB);
  75. matX = chollo.solve(matB);
  76. VERIFY_IS_APPROX(symm * matX, matB);
  77. const MatrixType symmLo_inverse = chollo.solve(MatrixType::Identity(rows,cols));
  78. RealScalar rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Lower>(symmLo)) /
  79. matrix_l1_norm<MatrixType, Lower>(symmLo_inverse);
  80. RealScalar rcond_est = chollo.rcond();
  81. // Verify that the estimated condition number is within a factor of 10 of the
  82. // truth.
  83. VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
  84. // test the upper mode
  85. LLT<SquareMatrixType,Upper> cholup(symmUp);
  86. VERIFY_IS_APPROX(symm, cholup.reconstructedMatrix());
  87. vecX = cholup.solve(vecB);
  88. VERIFY_IS_APPROX(symm * vecX, vecB);
  89. matX = cholup.solve(matB);
  90. VERIFY_IS_APPROX(symm * matX, matB);
  91. // Verify that the estimated condition number is within a factor of 10 of the
  92. // truth.
  93. const MatrixType symmUp_inverse = cholup.solve(MatrixType::Identity(rows,cols));
  94. rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Upper>(symmUp)) /
  95. matrix_l1_norm<MatrixType, Upper>(symmUp_inverse);
  96. rcond_est = cholup.rcond();
  97. VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
  98. MatrixType neg = -symmLo;
  99. chollo.compute(neg);
  100. VERIFY(neg.size()==0 || chollo.info()==NumericalIssue);
  101. VERIFY_IS_APPROX(MatrixType(chollo.matrixL().transpose().conjugate()), MatrixType(chollo.matrixU()));
  102. VERIFY_IS_APPROX(MatrixType(chollo.matrixU().transpose().conjugate()), MatrixType(chollo.matrixL()));
  103. VERIFY_IS_APPROX(MatrixType(cholup.matrixL().transpose().conjugate()), MatrixType(cholup.matrixU()));
  104. VERIFY_IS_APPROX(MatrixType(cholup.matrixU().transpose().conjugate()), MatrixType(cholup.matrixL()));
  105. // test some special use cases of SelfCwiseBinaryOp:
  106. MatrixType m1 = MatrixType::Random(rows,cols), m2(rows,cols);
  107. m2 = m1;
  108. m2 += symmLo.template selfadjointView<Lower>().llt().solve(matB);
  109. VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
  110. m2 = m1;
  111. m2 -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
  112. VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
  113. m2 = m1;
  114. m2.noalias() += symmLo.template selfadjointView<Lower>().llt().solve(matB);
  115. VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
  116. m2 = m1;
  117. m2.noalias() -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
  118. VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
  119. }
  120. // LDLT
  121. {
  122. int sign = internal::random<int>()%2 ? 1 : -1;
  123. if(sign == -1)
  124. {
  125. symm = -symm; // test a negative matrix
  126. }
  127. SquareMatrixType symmUp = symm.template triangularView<Upper>();
  128. SquareMatrixType symmLo = symm.template triangularView<Lower>();
  129. LDLT<SquareMatrixType,Lower> ldltlo(symmLo);
  130. VERIFY(ldltlo.info()==Success);
  131. VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
  132. vecX = ldltlo.solve(vecB);
  133. VERIFY_IS_APPROX(symm * vecX, vecB);
  134. matX = ldltlo.solve(matB);
  135. VERIFY_IS_APPROX(symm * matX, matB);
  136. const MatrixType symmLo_inverse = ldltlo.solve(MatrixType::Identity(rows,cols));
  137. RealScalar rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Lower>(symmLo)) /
  138. matrix_l1_norm<MatrixType, Lower>(symmLo_inverse);
  139. RealScalar rcond_est = ldltlo.rcond();
  140. // Verify that the estimated condition number is within a factor of 10 of the
  141. // truth.
  142. VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
  143. LDLT<SquareMatrixType,Upper> ldltup(symmUp);
  144. VERIFY(ldltup.info()==Success);
  145. VERIFY_IS_APPROX(symm, ldltup.reconstructedMatrix());
  146. vecX = ldltup.solve(vecB);
  147. VERIFY_IS_APPROX(symm * vecX, vecB);
  148. matX = ldltup.solve(matB);
  149. VERIFY_IS_APPROX(symm * matX, matB);
  150. // Verify that the estimated condition number is within a factor of 10 of the
  151. // truth.
  152. const MatrixType symmUp_inverse = ldltup.solve(MatrixType::Identity(rows,cols));
  153. rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Upper>(symmUp)) /
  154. matrix_l1_norm<MatrixType, Upper>(symmUp_inverse);
  155. rcond_est = ldltup.rcond();
  156. VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
  157. VERIFY_IS_APPROX(MatrixType(ldltlo.matrixL().transpose().conjugate()), MatrixType(ldltlo.matrixU()));
  158. VERIFY_IS_APPROX(MatrixType(ldltlo.matrixU().transpose().conjugate()), MatrixType(ldltlo.matrixL()));
  159. VERIFY_IS_APPROX(MatrixType(ldltup.matrixL().transpose().conjugate()), MatrixType(ldltup.matrixU()));
  160. VERIFY_IS_APPROX(MatrixType(ldltup.matrixU().transpose().conjugate()), MatrixType(ldltup.matrixL()));
  161. if(MatrixType::RowsAtCompileTime==Dynamic)
  162. {
  163. // note : each inplace permutation requires a small temporary vector (mask)
  164. // check inplace solve
  165. matX = matB;
  166. VERIFY_EVALUATION_COUNT(matX = ldltlo.solve(matX), 0);
  167. VERIFY_IS_APPROX(matX, ldltlo.solve(matB).eval());
  168. matX = matB;
  169. VERIFY_EVALUATION_COUNT(matX = ldltup.solve(matX), 0);
  170. VERIFY_IS_APPROX(matX, ldltup.solve(matB).eval());
  171. }
  172. // restore
  173. if(sign == -1)
  174. symm = -symm;
  175. // check matrices coming from linear constraints with Lagrange multipliers
  176. if(rows>=3)
  177. {
  178. SquareMatrixType A = symm;
  179. Index c = internal::random<Index>(0,rows-2);
  180. A.bottomRightCorner(c,c).setZero();
  181. // Make sure a solution exists:
  182. vecX.setRandom();
  183. vecB = A * vecX;
  184. vecX.setZero();
  185. ldltlo.compute(A);
  186. VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
  187. vecX = ldltlo.solve(vecB);
  188. VERIFY_IS_APPROX(A * vecX, vecB);
  189. }
  190. // check non-full rank matrices
  191. if(rows>=3)
  192. {
  193. Index r = internal::random<Index>(1,rows-1);
  194. Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,r);
  195. SquareMatrixType A = a * a.adjoint();
  196. // Make sure a solution exists:
  197. vecX.setRandom();
  198. vecB = A * vecX;
  199. vecX.setZero();
  200. ldltlo.compute(A);
  201. VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
  202. vecX = ldltlo.solve(vecB);
  203. VERIFY_IS_APPROX(A * vecX, vecB);
  204. }
  205. // check matrices with a wide spectrum
  206. if(rows>=3)
  207. {
  208. using std::pow;
  209. using std::sqrt;
  210. RealScalar s = (std::min)(16,std::numeric_limits<RealScalar>::max_exponent10/8);
  211. Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,rows);
  212. Matrix<RealScalar,Dynamic,1> d = Matrix<RealScalar,Dynamic,1>::Random(rows);
  213. for(Index k=0; k<rows; ++k)
  214. d(k) = d(k)*pow(RealScalar(10),internal::random<RealScalar>(-s,s));
  215. SquareMatrixType A = a * d.asDiagonal() * a.adjoint();
  216. // Make sure a solution exists:
  217. vecX.setRandom();
  218. vecB = A * vecX;
  219. vecX.setZero();
  220. ldltlo.compute(A);
  221. VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
  222. vecX = ldltlo.solve(vecB);
  223. if(ldltlo.vectorD().real().cwiseAbs().minCoeff()>RealScalar(0))
  224. {
  225. VERIFY_IS_APPROX(A * vecX,vecB);
  226. }
  227. else
  228. {
  229. RealScalar large_tol = sqrt(test_precision<RealScalar>());
  230. VERIFY((A * vecX).isApprox(vecB, large_tol));
  231. ++g_test_level;
  232. VERIFY_IS_APPROX(A * vecX,vecB);
  233. --g_test_level;
  234. }
  235. }
  236. }
  237. // update/downdate
  238. CALL_SUBTEST(( test_chol_update<SquareMatrixType,LLT>(symm) ));
  239. CALL_SUBTEST(( test_chol_update<SquareMatrixType,LDLT>(symm) ));
  240. }
  241. template<typename MatrixType> void cholesky_cplx(const MatrixType& m)
  242. {
  243. // classic test
  244. cholesky(m);
  245. // test mixing real/scalar types
  246. Index rows = m.rows();
  247. Index cols = m.cols();
  248. typedef typename MatrixType::Scalar Scalar;
  249. typedef typename NumTraits<Scalar>::Real RealScalar;
  250. typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RealMatrixType;
  251. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  252. RealMatrixType a0 = RealMatrixType::Random(rows,cols);
  253. VectorType vecB = VectorType::Random(rows), vecX(rows);
  254. MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
  255. RealMatrixType symm = a0 * a0.adjoint();
  256. // let's make sure the matrix is not singular or near singular
  257. for (int k=0; k<3; ++k)
  258. {
  259. RealMatrixType a1 = RealMatrixType::Random(rows,cols);
  260. symm += a1 * a1.adjoint();
  261. }
  262. {
  263. RealMatrixType symmLo = symm.template triangularView<Lower>();
  264. LLT<RealMatrixType,Lower> chollo(symmLo);
  265. VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
  266. vecX = chollo.solve(vecB);
  267. VERIFY_IS_APPROX(symm * vecX, vecB);
  268. // matX = chollo.solve(matB);
  269. // VERIFY_IS_APPROX(symm * matX, matB);
  270. }
  271. // LDLT
  272. {
  273. int sign = internal::random<int>()%2 ? 1 : -1;
  274. if(sign == -1)
  275. {
  276. symm = -symm; // test a negative matrix
  277. }
  278. RealMatrixType symmLo = symm.template triangularView<Lower>();
  279. LDLT<RealMatrixType,Lower> ldltlo(symmLo);
  280. VERIFY(ldltlo.info()==Success);
  281. VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
  282. vecX = ldltlo.solve(vecB);
  283. VERIFY_IS_APPROX(symm * vecX, vecB);
  284. // matX = ldltlo.solve(matB);
  285. // VERIFY_IS_APPROX(symm * matX, matB);
  286. }
  287. }
  288. // regression test for bug 241
  289. template<typename MatrixType> void cholesky_bug241(const MatrixType& m)
  290. {
  291. eigen_assert(m.rows() == 2 && m.cols() == 2);
  292. typedef typename MatrixType::Scalar Scalar;
  293. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  294. MatrixType matA;
  295. matA << 1, 1, 1, 1;
  296. VectorType vecB;
  297. vecB << 1, 1;
  298. VectorType vecX = matA.ldlt().solve(vecB);
  299. VERIFY_IS_APPROX(matA * vecX, vecB);
  300. }
  301. // LDLT is not guaranteed to work for indefinite matrices, but happens to work fine if matrix is diagonal.
  302. // This test checks that LDLT reports correctly that matrix is indefinite.
  303. // See http://forum.kde.org/viewtopic.php?f=74&t=106942 and bug 736
  304. template<typename MatrixType> void cholesky_definiteness(const MatrixType& m)
  305. {
  306. eigen_assert(m.rows() == 2 && m.cols() == 2);
  307. MatrixType mat;
  308. LDLT<MatrixType> ldlt(2);
  309. {
  310. mat << 1, 0, 0, -1;
  311. ldlt.compute(mat);
  312. VERIFY(ldlt.info()==Success);
  313. VERIFY(!ldlt.isNegative());
  314. VERIFY(!ldlt.isPositive());
  315. VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
  316. }
  317. {
  318. mat << 1, 2, 2, 1;
  319. ldlt.compute(mat);
  320. VERIFY(ldlt.info()==Success);
  321. VERIFY(!ldlt.isNegative());
  322. VERIFY(!ldlt.isPositive());
  323. VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
  324. }
  325. {
  326. mat << 0, 0, 0, 0;
  327. ldlt.compute(mat);
  328. VERIFY(ldlt.info()==Success);
  329. VERIFY(ldlt.isNegative());
  330. VERIFY(ldlt.isPositive());
  331. VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
  332. }
  333. {
  334. mat << 0, 0, 0, 1;
  335. ldlt.compute(mat);
  336. VERIFY(ldlt.info()==Success);
  337. VERIFY(!ldlt.isNegative());
  338. VERIFY(ldlt.isPositive());
  339. VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
  340. }
  341. {
  342. mat << -1, 0, 0, 0;
  343. ldlt.compute(mat);
  344. VERIFY(ldlt.info()==Success);
  345. VERIFY(ldlt.isNegative());
  346. VERIFY(!ldlt.isPositive());
  347. VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
  348. }
  349. }
  350. template<typename>
  351. void cholesky_faillure_cases()
  352. {
  353. MatrixXd mat;
  354. LDLT<MatrixXd> ldlt;
  355. {
  356. mat.resize(2,2);
  357. mat << 0, 1, 1, 0;
  358. ldlt.compute(mat);
  359. VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
  360. VERIFY(ldlt.info()==NumericalIssue);
  361. }
  362. #if (!EIGEN_ARCH_i386) || defined(EIGEN_VECTORIZE_SSE2)
  363. {
  364. mat.resize(3,3);
  365. mat << -1, -3, 3,
  366. -3, -8.9999999999999999999, 1,
  367. 3, 1, 0;
  368. ldlt.compute(mat);
  369. VERIFY(ldlt.info()==NumericalIssue);
  370. VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
  371. }
  372. #endif
  373. {
  374. mat.resize(3,3);
  375. mat << 1, 2, 3,
  376. 2, 4, 1,
  377. 3, 1, 0;
  378. ldlt.compute(mat);
  379. VERIFY(ldlt.info()==NumericalIssue);
  380. VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
  381. }
  382. {
  383. mat.resize(8,8);
  384. mat << 0.1, 0, -0.1, 0, 0, 0, 1, 0,
  385. 0, 4.24667, 0, 2.00333, 0, 0, 0, 0,
  386. -0.1, 0, 0.2, 0, -0.1, 0, 0, 0,
  387. 0, 2.00333, 0, 8.49333, 0, 2.00333, 0, 0,
  388. 0, 0, -0.1, 0, 0.1, 0, 0, 1,
  389. 0, 0, 0, 2.00333, 0, 4.24667, 0, 0,
  390. 1, 0, 0, 0, 0, 0, 0, 0,
  391. 0, 0, 0, 0, 1, 0, 0, 0;
  392. ldlt.compute(mat);
  393. VERIFY(ldlt.info()==NumericalIssue);
  394. VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
  395. }
  396. // bug 1479
  397. {
  398. mat.resize(4,4);
  399. mat << 1, 2, 0, 1,
  400. 2, 4, 0, 2,
  401. 0, 0, 0, 1,
  402. 1, 2, 1, 1;
  403. ldlt.compute(mat);
  404. VERIFY(ldlt.info()==NumericalIssue);
  405. VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
  406. }
  407. }
  408. template<typename MatrixType> void cholesky_verify_assert()
  409. {
  410. MatrixType tmp;
  411. LLT<MatrixType> llt;
  412. VERIFY_RAISES_ASSERT(llt.matrixL())
  413. VERIFY_RAISES_ASSERT(llt.matrixU())
  414. VERIFY_RAISES_ASSERT(llt.solve(tmp))
  415. VERIFY_RAISES_ASSERT(llt.solveInPlace(&tmp))
  416. LDLT<MatrixType> ldlt;
  417. VERIFY_RAISES_ASSERT(ldlt.matrixL())
  418. VERIFY_RAISES_ASSERT(ldlt.permutationP())
  419. VERIFY_RAISES_ASSERT(ldlt.vectorD())
  420. VERIFY_RAISES_ASSERT(ldlt.isPositive())
  421. VERIFY_RAISES_ASSERT(ldlt.isNegative())
  422. VERIFY_RAISES_ASSERT(ldlt.solve(tmp))
  423. VERIFY_RAISES_ASSERT(ldlt.solveInPlace(&tmp))
  424. }
  425. void test_cholesky()
  426. {
  427. int s = 0;
  428. for(int i = 0; i < g_repeat; i++) {
  429. CALL_SUBTEST_1( cholesky(Matrix<double,1,1>()) );
  430. CALL_SUBTEST_3( cholesky(Matrix2d()) );
  431. CALL_SUBTEST_3( cholesky_bug241(Matrix2d()) );
  432. CALL_SUBTEST_3( cholesky_definiteness(Matrix2d()) );
  433. CALL_SUBTEST_4( cholesky(Matrix3f()) );
  434. CALL_SUBTEST_5( cholesky(Matrix4d()) );
  435. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
  436. CALL_SUBTEST_2( cholesky(MatrixXd(s,s)) );
  437. TEST_SET_BUT_UNUSED_VARIABLE(s)
  438. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
  439. CALL_SUBTEST_6( cholesky_cplx(MatrixXcd(s,s)) );
  440. TEST_SET_BUT_UNUSED_VARIABLE(s)
  441. }
  442. // empty matrix, regression test for Bug 785:
  443. CALL_SUBTEST_2( cholesky(MatrixXd(0,0)) );
  444. // This does not work yet:
  445. // CALL_SUBTEST_2( cholesky(Matrix<double,0,0>()) );
  446. CALL_SUBTEST_4( cholesky_verify_assert<Matrix3f>() );
  447. CALL_SUBTEST_7( cholesky_verify_assert<Matrix3d>() );
  448. CALL_SUBTEST_8( cholesky_verify_assert<MatrixXf>() );
  449. CALL_SUBTEST_2( cholesky_verify_assert<MatrixXd>() );
  450. // Test problem size constructors
  451. CALL_SUBTEST_9( LLT<MatrixXf>(10) );
  452. CALL_SUBTEST_9( LDLT<MatrixXf>(10) );
  453. CALL_SUBTEST_2( cholesky_faillure_cases<void>() );
  454. TEST_SET_BUT_UNUSED_VARIABLE(nb_temporaries)
  455. }