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.

407 lines
13 KiB

  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 STORMEIGEN_NO_ASSERTION_CHECKING
  10. #define STORMEIGEN_NO_ASSERTION_CHECKING
  11. #endif
  12. #define TEST_ENABLE_TEMPORARY_TRACKING
  13. #include "main.h"
  14. #include <StormEigen/Cholesky>
  15. #include <StormEigen/QR>
  16. template<typename MatrixType,template <typename,int> class CholType> void test_chol_update(const MatrixType& symm)
  17. {
  18. typedef typename MatrixType::Scalar Scalar;
  19. typedef typename MatrixType::RealScalar RealScalar;
  20. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  21. MatrixType symmLo = symm.template triangularView<Lower>();
  22. MatrixType symmUp = symm.template triangularView<Upper>();
  23. MatrixType symmCpy = symm;
  24. CholType<MatrixType,Lower> chollo(symmLo);
  25. CholType<MatrixType,Upper> cholup(symmUp);
  26. for (int k=0; k<10; ++k)
  27. {
  28. VectorType vec = VectorType::Random(symm.rows());
  29. RealScalar sigma = internal::random<RealScalar>();
  30. symmCpy += sigma * vec * vec.adjoint();
  31. // we are doing some downdates, so it might be the case that the matrix is not SPD anymore
  32. CholType<MatrixType,Lower> chol(symmCpy);
  33. if(chol.info()!=Success)
  34. break;
  35. chollo.rankUpdate(vec, sigma);
  36. VERIFY_IS_APPROX(symmCpy, chollo.reconstructedMatrix());
  37. cholup.rankUpdate(vec, sigma);
  38. VERIFY_IS_APPROX(symmCpy, cholup.reconstructedMatrix());
  39. }
  40. }
  41. template<typename MatrixType> void cholesky(const MatrixType& m)
  42. {
  43. typedef typename MatrixType::Index Index;
  44. /* this test covers the following files:
  45. LLT.h LDLT.h
  46. */
  47. Index rows = m.rows();
  48. Index cols = m.cols();
  49. typedef typename MatrixType::Scalar Scalar;
  50. typedef typename NumTraits<Scalar>::Real RealScalar;
  51. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
  52. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  53. MatrixType a0 = MatrixType::Random(rows,cols);
  54. VectorType vecB = VectorType::Random(rows), vecX(rows);
  55. MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
  56. SquareMatrixType symm = a0 * a0.adjoint();
  57. // let's make sure the matrix is not singular or near singular
  58. for (int k=0; k<3; ++k)
  59. {
  60. MatrixType a1 = MatrixType::Random(rows,cols);
  61. symm += a1 * a1.adjoint();
  62. }
  63. {
  64. SquareMatrixType symmUp = symm.template triangularView<Upper>();
  65. SquareMatrixType symmLo = symm.template triangularView<Lower>();
  66. LLT<SquareMatrixType,Lower> chollo(symmLo);
  67. VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
  68. vecX = chollo.solve(vecB);
  69. VERIFY_IS_APPROX(symm * vecX, vecB);
  70. matX = chollo.solve(matB);
  71. VERIFY_IS_APPROX(symm * matX, matB);
  72. // test the upper mode
  73. LLT<SquareMatrixType,Upper> cholup(symmUp);
  74. VERIFY_IS_APPROX(symm, cholup.reconstructedMatrix());
  75. vecX = cholup.solve(vecB);
  76. VERIFY_IS_APPROX(symm * vecX, vecB);
  77. matX = cholup.solve(matB);
  78. VERIFY_IS_APPROX(symm * matX, matB);
  79. MatrixType neg = -symmLo;
  80. chollo.compute(neg);
  81. VERIFY(chollo.info()==NumericalIssue);
  82. VERIFY_IS_APPROX(MatrixType(chollo.matrixL().transpose().conjugate()), MatrixType(chollo.matrixU()));
  83. VERIFY_IS_APPROX(MatrixType(chollo.matrixU().transpose().conjugate()), MatrixType(chollo.matrixL()));
  84. VERIFY_IS_APPROX(MatrixType(cholup.matrixL().transpose().conjugate()), MatrixType(cholup.matrixU()));
  85. VERIFY_IS_APPROX(MatrixType(cholup.matrixU().transpose().conjugate()), MatrixType(cholup.matrixL()));
  86. // test some special use cases of SelfCwiseBinaryOp:
  87. MatrixType m1 = MatrixType::Random(rows,cols), m2(rows,cols);
  88. m2 = m1;
  89. m2 += symmLo.template selfadjointView<Lower>().llt().solve(matB);
  90. VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
  91. m2 = m1;
  92. m2 -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
  93. VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
  94. m2 = m1;
  95. m2.noalias() += symmLo.template selfadjointView<Lower>().llt().solve(matB);
  96. VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
  97. m2 = m1;
  98. m2.noalias() -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
  99. VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
  100. }
  101. // LDLT
  102. {
  103. int sign = internal::random<int>()%2 ? 1 : -1;
  104. if(sign == -1)
  105. {
  106. symm = -symm; // test a negative matrix
  107. }
  108. SquareMatrixType symmUp = symm.template triangularView<Upper>();
  109. SquareMatrixType symmLo = symm.template triangularView<Lower>();
  110. LDLT<SquareMatrixType,Lower> ldltlo(symmLo);
  111. VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
  112. vecX = ldltlo.solve(vecB);
  113. VERIFY_IS_APPROX(symm * vecX, vecB);
  114. matX = ldltlo.solve(matB);
  115. VERIFY_IS_APPROX(symm * matX, matB);
  116. LDLT<SquareMatrixType,Upper> ldltup(symmUp);
  117. VERIFY_IS_APPROX(symm, ldltup.reconstructedMatrix());
  118. vecX = ldltup.solve(vecB);
  119. VERIFY_IS_APPROX(symm * vecX, vecB);
  120. matX = ldltup.solve(matB);
  121. VERIFY_IS_APPROX(symm * matX, matB);
  122. VERIFY_IS_APPROX(MatrixType(ldltlo.matrixL().transpose().conjugate()), MatrixType(ldltlo.matrixU()));
  123. VERIFY_IS_APPROX(MatrixType(ldltlo.matrixU().transpose().conjugate()), MatrixType(ldltlo.matrixL()));
  124. VERIFY_IS_APPROX(MatrixType(ldltup.matrixL().transpose().conjugate()), MatrixType(ldltup.matrixU()));
  125. VERIFY_IS_APPROX(MatrixType(ldltup.matrixU().transpose().conjugate()), MatrixType(ldltup.matrixL()));
  126. if(MatrixType::RowsAtCompileTime==Dynamic)
  127. {
  128. // note : each inplace permutation requires a small temporary vector (mask)
  129. // check inplace solve
  130. matX = matB;
  131. VERIFY_EVALUATION_COUNT(matX = ldltlo.solve(matX), 0);
  132. VERIFY_IS_APPROX(matX, ldltlo.solve(matB).eval());
  133. matX = matB;
  134. VERIFY_EVALUATION_COUNT(matX = ldltup.solve(matX), 0);
  135. VERIFY_IS_APPROX(matX, ldltup.solve(matB).eval());
  136. }
  137. // restore
  138. if(sign == -1)
  139. symm = -symm;
  140. // check matrices coming from linear constraints with Lagrange multipliers
  141. if(rows>=3)
  142. {
  143. SquareMatrixType A = symm;
  144. Index c = internal::random<Index>(0,rows-2);
  145. A.bottomRightCorner(c,c).setZero();
  146. // Make sure a solution exists:
  147. vecX.setRandom();
  148. vecB = A * vecX;
  149. vecX.setZero();
  150. ldltlo.compute(A);
  151. VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
  152. vecX = ldltlo.solve(vecB);
  153. VERIFY_IS_APPROX(A * vecX, vecB);
  154. }
  155. // check non-full rank matrices
  156. if(rows>=3)
  157. {
  158. Index r = internal::random<Index>(1,rows-1);
  159. Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,r);
  160. SquareMatrixType A = a * a.adjoint();
  161. // Make sure a solution exists:
  162. vecX.setRandom();
  163. vecB = A * vecX;
  164. vecX.setZero();
  165. ldltlo.compute(A);
  166. VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
  167. vecX = ldltlo.solve(vecB);
  168. VERIFY_IS_APPROX(A * vecX, vecB);
  169. }
  170. // check matrices with a wide spectrum
  171. if(rows>=3)
  172. {
  173. RealScalar s = (std::min)(16,std::numeric_limits<RealScalar>::max_exponent10/8);
  174. Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,rows);
  175. Matrix<RealScalar,Dynamic,1> d = Matrix<RealScalar,Dynamic,1>::Random(rows);
  176. for(Index k=0; k<rows; ++k)
  177. d(k) = d(k)*std::pow(RealScalar(10),internal::random<RealScalar>(-s,s));
  178. SquareMatrixType A = a * d.asDiagonal() * a.adjoint();
  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. if(ldltlo.vectorD().real().cwiseAbs().minCoeff()>RealScalar(0))
  187. {
  188. VERIFY_IS_APPROX(A * vecX,vecB);
  189. }
  190. else
  191. {
  192. RealScalar large_tol = std::sqrt(test_precision<RealScalar>());
  193. VERIFY((A * vecX).isApprox(vecB, large_tol));
  194. ++g_test_level;
  195. VERIFY_IS_APPROX(A * vecX,vecB);
  196. --g_test_level;
  197. }
  198. }
  199. }
  200. // update/downdate
  201. CALL_SUBTEST(( test_chol_update<SquareMatrixType,LLT>(symm) ));
  202. CALL_SUBTEST(( test_chol_update<SquareMatrixType,LDLT>(symm) ));
  203. }
  204. template<typename MatrixType> void cholesky_cplx(const MatrixType& m)
  205. {
  206. // classic test
  207. cholesky(m);
  208. // test mixing real/scalar types
  209. typedef typename MatrixType::Index Index;
  210. Index rows = m.rows();
  211. Index cols = m.cols();
  212. typedef typename MatrixType::Scalar Scalar;
  213. typedef typename NumTraits<Scalar>::Real RealScalar;
  214. typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RealMatrixType;
  215. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  216. RealMatrixType a0 = RealMatrixType::Random(rows,cols);
  217. VectorType vecB = VectorType::Random(rows), vecX(rows);
  218. MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
  219. RealMatrixType symm = a0 * a0.adjoint();
  220. // let's make sure the matrix is not singular or near singular
  221. for (int k=0; k<3; ++k)
  222. {
  223. RealMatrixType a1 = RealMatrixType::Random(rows,cols);
  224. symm += a1 * a1.adjoint();
  225. }
  226. {
  227. RealMatrixType symmLo = symm.template triangularView<Lower>();
  228. LLT<RealMatrixType,Lower> chollo(symmLo);
  229. VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
  230. vecX = chollo.solve(vecB);
  231. VERIFY_IS_APPROX(symm * vecX, vecB);
  232. // matX = chollo.solve(matB);
  233. // VERIFY_IS_APPROX(symm * matX, matB);
  234. }
  235. // LDLT
  236. {
  237. int sign = internal::random<int>()%2 ? 1 : -1;
  238. if(sign == -1)
  239. {
  240. symm = -symm; // test a negative matrix
  241. }
  242. RealMatrixType symmLo = symm.template triangularView<Lower>();
  243. LDLT<RealMatrixType,Lower> ldltlo(symmLo);
  244. VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
  245. vecX = ldltlo.solve(vecB);
  246. VERIFY_IS_APPROX(symm * vecX, vecB);
  247. // matX = ldltlo.solve(matB);
  248. // VERIFY_IS_APPROX(symm * matX, matB);
  249. }
  250. }
  251. // regression test for bug 241
  252. template<typename MatrixType> void cholesky_bug241(const MatrixType& m)
  253. {
  254. eigen_assert(m.rows() == 2 && m.cols() == 2);
  255. typedef typename MatrixType::Scalar Scalar;
  256. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  257. MatrixType matA;
  258. matA << 1, 1, 1, 1;
  259. VectorType vecB;
  260. vecB << 1, 1;
  261. VectorType vecX = matA.ldlt().solve(vecB);
  262. VERIFY_IS_APPROX(matA * vecX, vecB);
  263. }
  264. // LDLT is not guaranteed to work for indefinite matrices, but happens to work fine if matrix is diagonal.
  265. // This test checks that LDLT reports correctly that matrix is indefinite.
  266. // See http://forum.kde.org/viewtopic.php?f=74&t=106942 and bug 736
  267. template<typename MatrixType> void cholesky_definiteness(const MatrixType& m)
  268. {
  269. eigen_assert(m.rows() == 2 && m.cols() == 2);
  270. MatrixType mat;
  271. LDLT<MatrixType> ldlt(2);
  272. {
  273. mat << 1, 0, 0, -1;
  274. ldlt.compute(mat);
  275. VERIFY(!ldlt.isNegative());
  276. VERIFY(!ldlt.isPositive());
  277. }
  278. {
  279. mat << 1, 2, 2, 1;
  280. ldlt.compute(mat);
  281. VERIFY(!ldlt.isNegative());
  282. VERIFY(!ldlt.isPositive());
  283. }
  284. {
  285. mat << 0, 0, 0, 0;
  286. ldlt.compute(mat);
  287. VERIFY(ldlt.isNegative());
  288. VERIFY(ldlt.isPositive());
  289. }
  290. {
  291. mat << 0, 0, 0, 1;
  292. ldlt.compute(mat);
  293. VERIFY(!ldlt.isNegative());
  294. VERIFY(ldlt.isPositive());
  295. }
  296. {
  297. mat << -1, 0, 0, 0;
  298. ldlt.compute(mat);
  299. VERIFY(ldlt.isNegative());
  300. VERIFY(!ldlt.isPositive());
  301. }
  302. }
  303. template<typename MatrixType> void cholesky_verify_assert()
  304. {
  305. MatrixType tmp;
  306. LLT<MatrixType> llt;
  307. VERIFY_RAISES_ASSERT(llt.matrixL())
  308. VERIFY_RAISES_ASSERT(llt.matrixU())
  309. VERIFY_RAISES_ASSERT(llt.solve(tmp))
  310. VERIFY_RAISES_ASSERT(llt.solveInPlace(&tmp))
  311. LDLT<MatrixType> ldlt;
  312. VERIFY_RAISES_ASSERT(ldlt.matrixL())
  313. VERIFY_RAISES_ASSERT(ldlt.permutationP())
  314. VERIFY_RAISES_ASSERT(ldlt.vectorD())
  315. VERIFY_RAISES_ASSERT(ldlt.isPositive())
  316. VERIFY_RAISES_ASSERT(ldlt.isNegative())
  317. VERIFY_RAISES_ASSERT(ldlt.solve(tmp))
  318. VERIFY_RAISES_ASSERT(ldlt.solveInPlace(&tmp))
  319. }
  320. void test_cholesky()
  321. {
  322. int s = 0;
  323. for(int i = 0; i < g_repeat; i++) {
  324. CALL_SUBTEST_1( cholesky(Matrix<double,1,1>()) );
  325. CALL_SUBTEST_3( cholesky(Matrix2d()) );
  326. CALL_SUBTEST_3( cholesky_bug241(Matrix2d()) );
  327. CALL_SUBTEST_3( cholesky_definiteness(Matrix2d()) );
  328. CALL_SUBTEST_4( cholesky(Matrix3f()) );
  329. CALL_SUBTEST_5( cholesky(Matrix4d()) );
  330. s = internal::random<int>(1,STORMEIGEN_TEST_MAX_SIZE);
  331. CALL_SUBTEST_2( cholesky(MatrixXd(s,s)) );
  332. TEST_SET_BUT_UNUSED_VARIABLE(s)
  333. s = internal::random<int>(1,STORMEIGEN_TEST_MAX_SIZE/2);
  334. CALL_SUBTEST_6( cholesky_cplx(MatrixXcd(s,s)) );
  335. TEST_SET_BUT_UNUSED_VARIABLE(s)
  336. }
  337. CALL_SUBTEST_4( cholesky_verify_assert<Matrix3f>() );
  338. CALL_SUBTEST_7( cholesky_verify_assert<Matrix3d>() );
  339. CALL_SUBTEST_8( cholesky_verify_assert<MatrixXf>() );
  340. CALL_SUBTEST_2( cholesky_verify_assert<MatrixXd>() );
  341. // Test problem size constructors
  342. CALL_SUBTEST_9( LLT<MatrixXf>(10) );
  343. CALL_SUBTEST_9( LDLT<MatrixXf>(10) );
  344. TEST_SET_BUT_UNUSED_VARIABLE(nb_temporaries)
  345. }