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.

310 lines
10 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 typename NumTraits<Scalar>::Real RealScalar;
  58. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
  59. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  60. MatrixType a0 = MatrixType::Random(rows,cols);
  61. VectorType vecB = VectorType::Random(rows), vecX(rows);
  62. MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
  63. SquareMatrixType symm = a0 * a0.adjoint();
  64. // let's make sure the matrix is not singular or near singular
  65. for (int k=0; k<3; ++k)
  66. {
  67. MatrixType a1 = MatrixType::Random(rows,cols);
  68. symm += a1 * a1.adjoint();
  69. }
  70. SquareMatrixType symmUp = symm.template triangularView<Upper>();
  71. SquareMatrixType symmLo = symm.template triangularView<Lower>();
  72. // to test if really Cholesky only uses the upper triangular part, uncomment the following
  73. // FIXME: currently that fails !!
  74. //symm.template part<StrictlyLower>().setZero();
  75. {
  76. LLT<SquareMatrixType,Lower> chollo(symmLo);
  77. VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
  78. vecX = chollo.solve(vecB);
  79. VERIFY_IS_APPROX(symm * vecX, vecB);
  80. matX = chollo.solve(matB);
  81. VERIFY_IS_APPROX(symm * matX, matB);
  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. MatrixType neg = -symmLo;
  90. chollo.compute(neg);
  91. VERIFY(chollo.info()==NumericalIssue);
  92. VERIFY_IS_APPROX(MatrixType(chollo.matrixL().transpose().conjugate()), MatrixType(chollo.matrixU()));
  93. VERIFY_IS_APPROX(MatrixType(chollo.matrixU().transpose().conjugate()), MatrixType(chollo.matrixL()));
  94. VERIFY_IS_APPROX(MatrixType(cholup.matrixL().transpose().conjugate()), MatrixType(cholup.matrixU()));
  95. VERIFY_IS_APPROX(MatrixType(cholup.matrixU().transpose().conjugate()), MatrixType(cholup.matrixL()));
  96. }
  97. // LDLT
  98. {
  99. int sign = internal::random<int>()%2 ? 1 : -1;
  100. if(sign == -1)
  101. {
  102. symm = -symm; // test a negative matrix
  103. }
  104. SquareMatrixType symmUp = symm.template triangularView<Upper>();
  105. SquareMatrixType symmLo = symm.template triangularView<Lower>();
  106. LDLT<SquareMatrixType,Lower> ldltlo(symmLo);
  107. VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
  108. vecX = ldltlo.solve(vecB);
  109. VERIFY_IS_APPROX(symm * vecX, vecB);
  110. matX = ldltlo.solve(matB);
  111. VERIFY_IS_APPROX(symm * matX, matB);
  112. LDLT<SquareMatrixType,Upper> ldltup(symmUp);
  113. VERIFY_IS_APPROX(symm, ldltup.reconstructedMatrix());
  114. vecX = ldltup.solve(vecB);
  115. VERIFY_IS_APPROX(symm * vecX, vecB);
  116. matX = ldltup.solve(matB);
  117. VERIFY_IS_APPROX(symm * matX, matB);
  118. VERIFY_IS_APPROX(MatrixType(ldltlo.matrixL().transpose().conjugate()), MatrixType(ldltlo.matrixU()));
  119. VERIFY_IS_APPROX(MatrixType(ldltlo.matrixU().transpose().conjugate()), MatrixType(ldltlo.matrixL()));
  120. VERIFY_IS_APPROX(MatrixType(ldltup.matrixL().transpose().conjugate()), MatrixType(ldltup.matrixU()));
  121. VERIFY_IS_APPROX(MatrixType(ldltup.matrixU().transpose().conjugate()), MatrixType(ldltup.matrixL()));
  122. if(MatrixType::RowsAtCompileTime==Dynamic)
  123. {
  124. // note : each inplace permutation requires a small temporary vector (mask)
  125. // check inplace solve
  126. matX = matB;
  127. VERIFY_EVALUATION_COUNT(matX = ldltlo.solve(matX), 0);
  128. VERIFY_IS_APPROX(matX, ldltlo.solve(matB).eval());
  129. matX = matB;
  130. VERIFY_EVALUATION_COUNT(matX = ldltup.solve(matX), 0);
  131. VERIFY_IS_APPROX(matX, ldltup.solve(matB).eval());
  132. }
  133. // restore
  134. if(sign == -1)
  135. symm = -symm;
  136. }
  137. // test some special use cases of SelfCwiseBinaryOp:
  138. MatrixType m1 = MatrixType::Random(rows,cols), m2(rows,cols);
  139. m2 = m1;
  140. m2 += symmLo.template selfadjointView<Lower>().llt().solve(matB);
  141. VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
  142. m2 = m1;
  143. m2 -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
  144. VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
  145. m2 = m1;
  146. m2.noalias() += symmLo.template selfadjointView<Lower>().llt().solve(matB);
  147. VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
  148. m2 = m1;
  149. m2.noalias() -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
  150. VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
  151. // update/downdate
  152. CALL_SUBTEST(( test_chol_update<SquareMatrixType,LLT>(symm) ));
  153. CALL_SUBTEST(( test_chol_update<SquareMatrixType,LDLT>(symm) ));
  154. }
  155. template<typename MatrixType> void cholesky_cplx(const MatrixType& m)
  156. {
  157. // classic test
  158. cholesky(m);
  159. // test mixing real/scalar types
  160. typedef typename MatrixType::Index Index;
  161. Index rows = m.rows();
  162. Index cols = m.cols();
  163. typedef typename MatrixType::Scalar Scalar;
  164. typedef typename NumTraits<Scalar>::Real RealScalar;
  165. typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RealMatrixType;
  166. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  167. RealMatrixType a0 = RealMatrixType::Random(rows,cols);
  168. VectorType vecB = VectorType::Random(rows), vecX(rows);
  169. MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
  170. RealMatrixType symm = a0 * a0.adjoint();
  171. // let's make sure the matrix is not singular or near singular
  172. for (int k=0; k<3; ++k)
  173. {
  174. RealMatrixType a1 = RealMatrixType::Random(rows,cols);
  175. symm += a1 * a1.adjoint();
  176. }
  177. {
  178. RealMatrixType symmLo = symm.template triangularView<Lower>();
  179. LLT<RealMatrixType,Lower> chollo(symmLo);
  180. VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
  181. vecX = chollo.solve(vecB);
  182. VERIFY_IS_APPROX(symm * vecX, vecB);
  183. // matX = chollo.solve(matB);
  184. // VERIFY_IS_APPROX(symm * matX, matB);
  185. }
  186. // LDLT
  187. {
  188. int sign = internal::random<int>()%2 ? 1 : -1;
  189. if(sign == -1)
  190. {
  191. symm = -symm; // test a negative matrix
  192. }
  193. RealMatrixType symmLo = symm.template triangularView<Lower>();
  194. LDLT<RealMatrixType,Lower> ldltlo(symmLo);
  195. VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
  196. vecX = ldltlo.solve(vecB);
  197. VERIFY_IS_APPROX(symm * vecX, vecB);
  198. // matX = ldltlo.solve(matB);
  199. // VERIFY_IS_APPROX(symm * matX, matB);
  200. }
  201. }
  202. // regression test for bug 241
  203. template<typename MatrixType> void cholesky_bug241(const MatrixType& m)
  204. {
  205. eigen_assert(m.rows() == 2 && m.cols() == 2);
  206. typedef typename MatrixType::Scalar Scalar;
  207. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  208. MatrixType matA;
  209. matA << 1, 1, 1, 1;
  210. VectorType vecB;
  211. vecB << 1, 1;
  212. VectorType vecX = matA.ldlt().solve(vecB);
  213. VERIFY_IS_APPROX(matA * vecX, vecB);
  214. }
  215. template<typename MatrixType> void cholesky_verify_assert()
  216. {
  217. MatrixType tmp;
  218. LLT<MatrixType> llt;
  219. VERIFY_RAISES_ASSERT(llt.matrixL())
  220. VERIFY_RAISES_ASSERT(llt.matrixU())
  221. VERIFY_RAISES_ASSERT(llt.solve(tmp))
  222. VERIFY_RAISES_ASSERT(llt.solveInPlace(&tmp))
  223. LDLT<MatrixType> ldlt;
  224. VERIFY_RAISES_ASSERT(ldlt.matrixL())
  225. VERIFY_RAISES_ASSERT(ldlt.permutationP())
  226. VERIFY_RAISES_ASSERT(ldlt.vectorD())
  227. VERIFY_RAISES_ASSERT(ldlt.isPositive())
  228. VERIFY_RAISES_ASSERT(ldlt.isNegative())
  229. VERIFY_RAISES_ASSERT(ldlt.solve(tmp))
  230. VERIFY_RAISES_ASSERT(ldlt.solveInPlace(&tmp))
  231. }
  232. void test_cholesky()
  233. {
  234. int s;
  235. for(int i = 0; i < g_repeat; i++) {
  236. CALL_SUBTEST_1( cholesky(Matrix<double,1,1>()) );
  237. CALL_SUBTEST_3( cholesky(Matrix2d()) );
  238. CALL_SUBTEST_3( cholesky_bug241(Matrix2d()) );
  239. CALL_SUBTEST_4( cholesky(Matrix3f()) );
  240. CALL_SUBTEST_5( cholesky(Matrix4d()) );
  241. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
  242. CALL_SUBTEST_2( cholesky(MatrixXd(s,s)) );
  243. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
  244. CALL_SUBTEST_6( cholesky_cplx(MatrixXcd(s,s)) );
  245. }
  246. CALL_SUBTEST_4( cholesky_verify_assert<Matrix3f>() );
  247. CALL_SUBTEST_7( cholesky_verify_assert<Matrix3d>() );
  248. CALL_SUBTEST_8( cholesky_verify_assert<MatrixXf>() );
  249. CALL_SUBTEST_2( cholesky_verify_assert<MatrixXd>() );
  250. // Test problem size constructors
  251. CALL_SUBTEST_9( LLT<MatrixXf>(10) );
  252. CALL_SUBTEST_9( LDLT<MatrixXf>(10) );
  253. EIGEN_UNUSED_VARIABLE(s)
  254. }