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.

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