cholesky.cpp 18 KB

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