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.

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