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.

203 lines
6.6 KiB

  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> void lu_non_invertible()
  13. {
  14. typedef typename MatrixType::Index Index;
  15. typedef typename MatrixType::RealScalar RealScalar;
  16. /* this test covers the following files:
  17. LU.h
  18. */
  19. Index rows, cols, cols2;
  20. if(MatrixType::RowsAtCompileTime==Dynamic)
  21. {
  22. rows = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
  23. }
  24. else
  25. {
  26. rows = MatrixType::RowsAtCompileTime;
  27. }
  28. if(MatrixType::ColsAtCompileTime==Dynamic)
  29. {
  30. cols = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
  31. cols2 = internal::random<int>(2,EIGEN_TEST_MAX_SIZE);
  32. }
  33. else
  34. {
  35. cols2 = cols = MatrixType::ColsAtCompileTime;
  36. }
  37. enum {
  38. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  39. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  40. };
  41. typedef typename internal::kernel_retval_base<FullPivLU<MatrixType> >::ReturnType KernelMatrixType;
  42. typedef typename internal::image_retval_base<FullPivLU<MatrixType> >::ReturnType ImageMatrixType;
  43. typedef Matrix<typename MatrixType::Scalar, ColsAtCompileTime, ColsAtCompileTime>
  44. CMatrixType;
  45. typedef Matrix<typename MatrixType::Scalar, RowsAtCompileTime, RowsAtCompileTime>
  46. RMatrixType;
  47. Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
  48. // The image of the zero matrix should consist of a single (zero) column vector
  49. VERIFY((MatrixType::Zero(rows,cols).fullPivLu().image(MatrixType::Zero(rows,cols)).cols() == 1));
  50. MatrixType m1(rows, cols), m3(rows, cols2);
  51. CMatrixType m2(cols, cols2);
  52. createRandomPIMatrixOfRank(rank, rows, cols, m1);
  53. FullPivLU<MatrixType> lu;
  54. // The special value 0.01 below works well in tests. Keep in mind that we're only computing the rank
  55. // of singular values are either 0 or 1.
  56. // So it's not clear at all that the epsilon should play any role there.
  57. lu.setThreshold(RealScalar(0.01));
  58. lu.compute(m1);
  59. MatrixType u(rows,cols);
  60. u = lu.matrixLU().template triangularView<Upper>();
  61. RMatrixType l = RMatrixType::Identity(rows,rows);
  62. l.block(0,0,rows,(std::min)(rows,cols)).template triangularView<StrictlyLower>()
  63. = lu.matrixLU().block(0,0,rows,(std::min)(rows,cols));
  64. VERIFY_IS_APPROX(lu.permutationP() * m1 * lu.permutationQ(), l*u);
  65. KernelMatrixType m1kernel = lu.kernel();
  66. ImageMatrixType m1image = lu.image(m1);
  67. VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
  68. VERIFY(rank == lu.rank());
  69. VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
  70. VERIFY(!lu.isInjective());
  71. VERIFY(!lu.isInvertible());
  72. VERIFY(!lu.isSurjective());
  73. VERIFY((m1 * m1kernel).isMuchSmallerThan(m1));
  74. VERIFY(m1image.fullPivLu().rank() == rank);
  75. VERIFY_IS_APPROX(m1 * m1.adjoint() * m1image, m1image);
  76. m2 = CMatrixType::Random(cols,cols2);
  77. m3 = m1*m2;
  78. m2 = CMatrixType::Random(cols,cols2);
  79. // test that the code, which does resize(), may be applied to an xpr
  80. m2.block(0,0,m2.rows(),m2.cols()) = lu.solve(m3);
  81. VERIFY_IS_APPROX(m3, m1*m2);
  82. }
  83. template<typename MatrixType> void lu_invertible()
  84. {
  85. /* this test covers the following files:
  86. LU.h
  87. */
  88. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  89. int size = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
  90. MatrixType m1(size, size), m2(size, size), m3(size, size);
  91. FullPivLU<MatrixType> lu;
  92. lu.setThreshold(RealScalar(0.01));
  93. do {
  94. m1 = MatrixType::Random(size,size);
  95. lu.compute(m1);
  96. } while(!lu.isInvertible());
  97. VERIFY_IS_APPROX(m1, lu.reconstructedMatrix());
  98. VERIFY(0 == lu.dimensionOfKernel());
  99. VERIFY(lu.kernel().cols() == 1); // the kernel() should consist of a single (zero) column vector
  100. VERIFY(size == lu.rank());
  101. VERIFY(lu.isInjective());
  102. VERIFY(lu.isSurjective());
  103. VERIFY(lu.isInvertible());
  104. VERIFY(lu.image(m1).fullPivLu().isInvertible());
  105. m3 = MatrixType::Random(size,size);
  106. m2 = lu.solve(m3);
  107. VERIFY_IS_APPROX(m3, m1*m2);
  108. VERIFY_IS_APPROX(m2, lu.inverse()*m3);
  109. }
  110. template<typename MatrixType> void lu_partial_piv()
  111. {
  112. /* this test covers the following files:
  113. PartialPivLU.h
  114. */
  115. typedef typename MatrixType::Index Index;
  116. Index rows = internal::random<Index>(1,4);
  117. Index cols = rows;
  118. MatrixType m1(cols, rows);
  119. m1.setRandom();
  120. PartialPivLU<MatrixType> plu(m1);
  121. VERIFY_IS_APPROX(m1, plu.reconstructedMatrix());
  122. }
  123. template<typename MatrixType> void lu_verify_assert()
  124. {
  125. MatrixType tmp;
  126. FullPivLU<MatrixType> lu;
  127. VERIFY_RAISES_ASSERT(lu.matrixLU())
  128. VERIFY_RAISES_ASSERT(lu.permutationP())
  129. VERIFY_RAISES_ASSERT(lu.permutationQ())
  130. VERIFY_RAISES_ASSERT(lu.kernel())
  131. VERIFY_RAISES_ASSERT(lu.image(tmp))
  132. VERIFY_RAISES_ASSERT(lu.solve(tmp))
  133. VERIFY_RAISES_ASSERT(lu.determinant())
  134. VERIFY_RAISES_ASSERT(lu.rank())
  135. VERIFY_RAISES_ASSERT(lu.dimensionOfKernel())
  136. VERIFY_RAISES_ASSERT(lu.isInjective())
  137. VERIFY_RAISES_ASSERT(lu.isSurjective())
  138. VERIFY_RAISES_ASSERT(lu.isInvertible())
  139. VERIFY_RAISES_ASSERT(lu.inverse())
  140. PartialPivLU<MatrixType> plu;
  141. VERIFY_RAISES_ASSERT(plu.matrixLU())
  142. VERIFY_RAISES_ASSERT(plu.permutationP())
  143. VERIFY_RAISES_ASSERT(plu.solve(tmp))
  144. VERIFY_RAISES_ASSERT(plu.determinant())
  145. VERIFY_RAISES_ASSERT(plu.inverse())
  146. }
  147. void test_lu()
  148. {
  149. for(int i = 0; i < g_repeat; i++) {
  150. CALL_SUBTEST_1( lu_non_invertible<Matrix3f>() );
  151. CALL_SUBTEST_1( lu_verify_assert<Matrix3f>() );
  152. CALL_SUBTEST_2( (lu_non_invertible<Matrix<double, 4, 6> >()) );
  153. CALL_SUBTEST_2( (lu_verify_assert<Matrix<double, 4, 6> >()) );
  154. CALL_SUBTEST_3( lu_non_invertible<MatrixXf>() );
  155. CALL_SUBTEST_3( lu_invertible<MatrixXf>() );
  156. CALL_SUBTEST_3( lu_verify_assert<MatrixXf>() );
  157. CALL_SUBTEST_4( lu_non_invertible<MatrixXd>() );
  158. CALL_SUBTEST_4( lu_invertible<MatrixXd>() );
  159. CALL_SUBTEST_4( lu_partial_piv<MatrixXd>() );
  160. CALL_SUBTEST_4( lu_verify_assert<MatrixXd>() );
  161. CALL_SUBTEST_5( lu_non_invertible<MatrixXcf>() );
  162. CALL_SUBTEST_5( lu_invertible<MatrixXcf>() );
  163. CALL_SUBTEST_5( lu_verify_assert<MatrixXcf>() );
  164. CALL_SUBTEST_6( lu_non_invertible<MatrixXcd>() );
  165. CALL_SUBTEST_6( lu_invertible<MatrixXcd>() );
  166. CALL_SUBTEST_6( lu_partial_piv<MatrixXcd>() );
  167. CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
  168. CALL_SUBTEST_7(( lu_non_invertible<Matrix<float,Dynamic,16> >() ));
  169. // Test problem size constructors
  170. CALL_SUBTEST_9( PartialPivLU<MatrixXf>(10) );
  171. CALL_SUBTEST_9( FullPivLU<MatrixXf>(10, 20); );
  172. }
  173. }