lu.cpp 9.6 KB

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