lu.cpp 8.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
  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. #include "main.h"
  10. #include <Eigen/LU>
  11. #include "solverbase.h"
  12. using namespace std;
  13. template<typename MatrixType>
  14. typename MatrixType::RealScalar matrix_l1_norm(const MatrixType& m) {
  15. return m.cwiseAbs().colwise().sum().maxCoeff();
  16. }
  17. template<typename MatrixType> void lu_non_invertible()
  18. {
  19. STATIC_CHECK(( internal::is_same<typename FullPivLU<MatrixType>::StorageIndex,int>::value ));
  20. typedef typename MatrixType::RealScalar RealScalar;
  21. /* this test covers the following files:
  22. LU.h
  23. */
  24. Index rows, cols, cols2;
  25. if(MatrixType::RowsAtCompileTime==Dynamic)
  26. {
  27. rows = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
  28. }
  29. else
  30. {
  31. rows = MatrixType::RowsAtCompileTime;
  32. }
  33. if(MatrixType::ColsAtCompileTime==Dynamic)
  34. {
  35. cols = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
  36. cols2 = internal::random<int>(2,EIGEN_TEST_MAX_SIZE);
  37. }
  38. else
  39. {
  40. cols2 = cols = MatrixType::ColsAtCompileTime;
  41. }
  42. enum {
  43. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  44. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  45. };
  46. typedef typename internal::kernel_retval_base<FullPivLU<MatrixType> >::ReturnType KernelMatrixType;
  47. typedef typename internal::image_retval_base<FullPivLU<MatrixType> >::ReturnType ImageMatrixType;
  48. typedef Matrix<typename MatrixType::Scalar, ColsAtCompileTime, ColsAtCompileTime>
  49. CMatrixType;
  50. typedef Matrix<typename MatrixType::Scalar, RowsAtCompileTime, RowsAtCompileTime>
  51. RMatrixType;
  52. Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
  53. // The image of the zero matrix should consist of a single (zero) column vector
  54. VERIFY((MatrixType::Zero(rows,cols).fullPivLu().image(MatrixType::Zero(rows,cols)).cols() == 1));
  55. // The kernel of the zero matrix is the entire space, and thus is an invertible matrix of dimensions cols.
  56. KernelMatrixType kernel = MatrixType::Zero(rows,cols).fullPivLu().kernel();
  57. VERIFY((kernel.fullPivLu().isInvertible()));
  58. MatrixType m1(rows, cols), m3(rows, cols2);
  59. CMatrixType m2(cols, cols2);
  60. createRandomPIMatrixOfRank(rank, rows, cols, m1);
  61. FullPivLU<MatrixType> lu;
  62. // The special value 0.01 below works well in tests. Keep in mind that we're only computing the rank
  63. // of singular values are either 0 or 1.
  64. // So it's not clear at all that the epsilon should play any role there.
  65. lu.setThreshold(RealScalar(0.01));
  66. lu.compute(m1);
  67. MatrixType u(rows,cols);
  68. u = lu.matrixLU().template triangularView<Upper>();
  69. RMatrixType l = RMatrixType::Identity(rows,rows);
  70. l.block(0,0,rows,(std::min)(rows,cols)).template triangularView<StrictlyLower>()
  71. = lu.matrixLU().block(0,0,rows,(std::min)(rows,cols));
  72. VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l*u);
  73. KernelMatrixType m1kernel = lu.kernel();
  74. ImageMatrixType m1image = lu.image(m1);
  75. VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
  76. VERIFY(rank == lu.rank());
  77. VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
  78. VERIFY(!lu.isInjective());
  79. VERIFY(!lu.isInvertible());
  80. VERIFY(!lu.isSurjective());
  81. VERIFY_IS_MUCH_SMALLER_THAN((m1 * m1kernel), m1);
  82. VERIFY(m1image.fullPivLu().rank() == rank);
  83. VERIFY_IS_APPROX(m1 * m1.adjoint() * m1image, m1image);
  84. check_solverbase<CMatrixType, MatrixType>(m1, lu, rows, cols, cols2);
  85. m2 = CMatrixType::Random(cols,cols2);
  86. m3 = m1*m2;
  87. m2 = CMatrixType::Random(cols,cols2);
  88. // test that the code, which does resize(), may be applied to an xpr
  89. m2.block(0,0,m2.rows(),m2.cols()) = lu.solve(m3);
  90. VERIFY_IS_APPROX(m3, m1*m2);
  91. }
  92. template<typename MatrixType> void lu_invertible()
  93. {
  94. /* this test covers the following files:
  95. FullPivLU.h
  96. */
  97. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  98. Index size = MatrixType::RowsAtCompileTime;
  99. if( size==Dynamic)
  100. size = internal::random<Index>(1,EIGEN_TEST_MAX_SIZE);
  101. MatrixType m1(size, size), m2(size, size), m3(size, size);
  102. FullPivLU<MatrixType> lu;
  103. lu.setThreshold(RealScalar(0.01));
  104. do {
  105. m1 = MatrixType::Random(size,size);
  106. lu.compute(m1);
  107. } while(!lu.isInvertible());
  108. VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
  109. VERIFY(0 == lu.dimensionOfKernel());
  110. VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector
  111. VERIFY(size == lu.rank());
  112. VERIFY(lu.isInjective());
  113. VERIFY(lu.isSurjective());
  114. VERIFY(lu.isInvertible());
  115. VERIFY(lu.image(m1).fullPivLu().isInvertible());
  116. check_solverbase<MatrixType, MatrixType>(m1, lu, size, size, size);
  117. MatrixType m1_inverse = lu.inverse();
  118. m3 = MatrixType::Random(size,size);
  119. m2 = lu.solve(m3);
  120. VERIFY_IS_APPROX(m2, m1_inverse*m3);
  121. RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
  122. const RealScalar rcond_est = lu.rcond();
  123. // Verify that the estimated condition number is within a factor of 10 of the
  124. // truth.
  125. VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
  126. // Regression test for Bug 302
  127. MatrixType m4 = MatrixType::Random(size,size);
  128. VERIFY_IS_APPROX(lu.solve(m3*m4), lu.solve(m3)*m4);
  129. }
  130. template<typename MatrixType> void lu_partial_piv(Index size = MatrixType::ColsAtCompileTime)
  131. {
  132. /* this test covers the following files:
  133. PartialPivLU.h
  134. */
  135. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  136. MatrixType m1(size, size), m2(size, size), m3(size, size);
  137. m1.setRandom();
  138. PartialPivLU<MatrixType> plu(m1);
  139. STATIC_CHECK(( internal::is_same<typename PartialPivLU<MatrixType>::StorageIndex,int>::value ));
  140. VERIFY_IS_APPROX(m1, plu.reconstructedMatrix());
  141. check_solverbase<MatrixType, MatrixType>(m1, plu, size, size, size);
  142. MatrixType m1_inverse = plu.inverse();
  143. m3 = MatrixType::Random(size,size);
  144. m2 = plu.solve(m3);
  145. VERIFY_IS_APPROX(m2, m1_inverse*m3);
  146. RealScalar rcond = (RealScalar(1) / matrix_l1_norm(m1)) / matrix_l1_norm(m1_inverse);
  147. const RealScalar rcond_est = plu.rcond();
  148. // Verify that the estimate is within a factor of 10 of the truth.
  149. VERIFY(rcond_est > rcond / 10 && rcond_est < rcond * 10);
  150. }
  151. template<typename MatrixType> void lu_verify_assert()
  152. {
  153. MatrixType tmp;
  154. FullPivLU<MatrixType> lu;
  155. VERIFY_RAISES_ASSERT(lu.matrixLU())
  156. VERIFY_RAISES_ASSERT(lu.permutationP())
  157. VERIFY_RAISES_ASSERT(lu.permutationQ())
  158. VERIFY_RAISES_ASSERT(lu.kernel())
  159. VERIFY_RAISES_ASSERT(lu.image(tmp))
  160. VERIFY_RAISES_ASSERT(lu.solve(tmp))
  161. VERIFY_RAISES_ASSERT(lu.transpose().solve(tmp))
  162. VERIFY_RAISES_ASSERT(lu.adjoint().solve(tmp))
  163. VERIFY_RAISES_ASSERT(lu.determinant())
  164. VERIFY_RAISES_ASSERT(lu.rank())
  165. VERIFY_RAISES_ASSERT(lu.dimensionOfKernel())
  166. VERIFY_RAISES_ASSERT(lu.isInjective())
  167. VERIFY_RAISES_ASSERT(lu.isSurjective())
  168. VERIFY_RAISES_ASSERT(lu.isInvertible())
  169. VERIFY_RAISES_ASSERT(lu.inverse())
  170. PartialPivLU<MatrixType> plu;
  171. VERIFY_RAISES_ASSERT(plu.matrixLU())
  172. VERIFY_RAISES_ASSERT(plu.permutationP())
  173. VERIFY_RAISES_ASSERT(plu.solve(tmp))
  174. VERIFY_RAISES_ASSERT(plu.transpose().solve(tmp))
  175. VERIFY_RAISES_ASSERT(plu.adjoint().solve(tmp))
  176. VERIFY_RAISES_ASSERT(plu.determinant())
  177. VERIFY_RAISES_ASSERT(plu.inverse())
  178. }
  179. EIGEN_DECLARE_TEST(lu)
  180. {
  181. for(int i = 0; i < g_repeat; i++) {
  182. CALL_SUBTEST_1( lu_non_invertible<Matrix3f>() );
  183. CALL_SUBTEST_1( lu_invertible<Matrix3f>() );
  184. CALL_SUBTEST_1( lu_verify_assert<Matrix3f>() );
  185. CALL_SUBTEST_1( lu_partial_piv<Matrix3f>() );
  186. CALL_SUBTEST_2( (lu_non_invertible<Matrix<double, 4, 6> >()) );
  187. CALL_SUBTEST_2( (lu_verify_assert<Matrix<double, 4, 6> >()) );
  188. CALL_SUBTEST_2( lu_partial_piv<Matrix2d>() );
  189. CALL_SUBTEST_2( lu_partial_piv<Matrix4d>() );
  190. CALL_SUBTEST_2( (lu_partial_piv<Matrix<double,6,6> >()) );
  191. CALL_SUBTEST_3( lu_non_invertible<MatrixXf>() );
  192. CALL_SUBTEST_3( lu_invertible<MatrixXf>() );
  193. CALL_SUBTEST_3( lu_verify_assert<MatrixXf>() );
  194. CALL_SUBTEST_4( lu_non_invertible<MatrixXd>() );
  195. CALL_SUBTEST_4( lu_invertible<MatrixXd>() );
  196. CALL_SUBTEST_4( lu_partial_piv<MatrixXd>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)) );
  197. CALL_SUBTEST_4( lu_verify_assert<MatrixXd>() );
  198. CALL_SUBTEST_5( lu_non_invertible<MatrixXcf>() );
  199. CALL_SUBTEST_5( lu_invertible<MatrixXcf>() );
  200. CALL_SUBTEST_5( lu_verify_assert<MatrixXcf>() );
  201. CALL_SUBTEST_6( lu_non_invertible<MatrixXcd>() );
  202. CALL_SUBTEST_6( lu_invertible<MatrixXcd>() );
  203. CALL_SUBTEST_6( lu_partial_piv<MatrixXcd>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)) );
  204. CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
  205. CALL_SUBTEST_7(( lu_non_invertible<Matrix<float,Dynamic,16> >() ));
  206. // Test problem size constructors
  207. CALL_SUBTEST_9( PartialPivLU<MatrixXf>(10) );
  208. CALL_SUBTEST_9( FullPivLU<MatrixXf>(10, 20); );
  209. }
  210. }