qr_fullpivoting.cpp 5.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159
  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. // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #include "main.h"
  11. #include <Eigen/QR>
  12. #include "solverbase.h"
  13. template<typename MatrixType> void qr()
  14. {
  15. STATIC_CHECK(( internal::is_same<typename FullPivHouseholderQR<MatrixType>::StorageIndex,int>::value ));
  16. static const int Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime;
  17. Index max_size = EIGEN_TEST_MAX_SIZE;
  18. Index min_size = numext::maxi(1,EIGEN_TEST_MAX_SIZE/10);
  19. Index rows = Rows == Dynamic ? internal::random<Index>(min_size,max_size) : Rows,
  20. cols = Cols == Dynamic ? internal::random<Index>(min_size,max_size) : Cols,
  21. cols2 = Cols == Dynamic ? internal::random<Index>(min_size,max_size) : Cols,
  22. rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
  23. typedef typename MatrixType::Scalar Scalar;
  24. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
  25. MatrixType m1;
  26. createRandomPIMatrixOfRank(rank,rows,cols,m1);
  27. FullPivHouseholderQR<MatrixType> qr(m1);
  28. VERIFY_IS_EQUAL(rank, qr.rank());
  29. VERIFY_IS_EQUAL(cols - qr.rank(), qr.dimensionOfKernel());
  30. VERIFY(!qr.isInjective());
  31. VERIFY(!qr.isInvertible());
  32. VERIFY(!qr.isSurjective());
  33. MatrixType r = qr.matrixQR();
  34. MatrixQType q = qr.matrixQ();
  35. VERIFY_IS_UNITARY(q);
  36. // FIXME need better way to construct trapezoid
  37. for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
  38. MatrixType c = qr.matrixQ() * r * qr.colsPermutation().inverse();
  39. VERIFY_IS_APPROX(m1, c);
  40. // stress the ReturnByValue mechanism
  41. MatrixType tmp;
  42. VERIFY_IS_APPROX(tmp.noalias() = qr.matrixQ() * r, (qr.matrixQ() * r).eval());
  43. check_solverbase<MatrixType, MatrixType>(m1, qr, rows, cols, cols2);
  44. {
  45. MatrixType m2, m3;
  46. Index size = rows;
  47. do {
  48. m1 = MatrixType::Random(size,size);
  49. qr.compute(m1);
  50. } while(!qr.isInvertible());
  51. MatrixType m1_inv = qr.inverse();
  52. m3 = m1 * MatrixType::Random(size,cols2);
  53. m2 = qr.solve(m3);
  54. VERIFY_IS_APPROX(m2, m1_inv*m3);
  55. }
  56. }
  57. template<typename MatrixType> void qr_invertible()
  58. {
  59. using std::log;
  60. using std::abs;
  61. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  62. typedef typename MatrixType::Scalar Scalar;
  63. Index max_size = numext::mini(50,EIGEN_TEST_MAX_SIZE);
  64. Index min_size = numext::maxi(1,EIGEN_TEST_MAX_SIZE/10);
  65. Index size = internal::random<Index>(min_size,max_size);
  66. MatrixType m1(size, size), m2(size, size), m3(size, size);
  67. m1 = MatrixType::Random(size,size);
  68. if (internal::is_same<RealScalar,float>::value)
  69. {
  70. // let's build a matrix more stable to inverse
  71. MatrixType a = MatrixType::Random(size,size*2);
  72. m1 += a * a.adjoint();
  73. }
  74. FullPivHouseholderQR<MatrixType> qr(m1);
  75. VERIFY(qr.isInjective());
  76. VERIFY(qr.isInvertible());
  77. VERIFY(qr.isSurjective());
  78. check_solverbase<MatrixType, MatrixType>(m1, qr, size, size, size);
  79. // now construct a matrix with prescribed determinant
  80. m1.setZero();
  81. for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>();
  82. RealScalar absdet = abs(m1.diagonal().prod());
  83. m3 = qr.matrixQ(); // get a unitary
  84. m1 = m3 * m1 * m3;
  85. qr.compute(m1);
  86. VERIFY_IS_APPROX(absdet, qr.absDeterminant());
  87. VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant());
  88. }
  89. template<typename MatrixType> void qr_verify_assert()
  90. {
  91. MatrixType tmp;
  92. FullPivHouseholderQR<MatrixType> qr;
  93. VERIFY_RAISES_ASSERT(qr.matrixQR())
  94. VERIFY_RAISES_ASSERT(qr.solve(tmp))
  95. VERIFY_RAISES_ASSERT(qr.transpose().solve(tmp))
  96. VERIFY_RAISES_ASSERT(qr.adjoint().solve(tmp))
  97. VERIFY_RAISES_ASSERT(qr.matrixQ())
  98. VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
  99. VERIFY_RAISES_ASSERT(qr.isInjective())
  100. VERIFY_RAISES_ASSERT(qr.isSurjective())
  101. VERIFY_RAISES_ASSERT(qr.isInvertible())
  102. VERIFY_RAISES_ASSERT(qr.inverse())
  103. VERIFY_RAISES_ASSERT(qr.absDeterminant())
  104. VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
  105. }
  106. EIGEN_DECLARE_TEST(qr_fullpivoting)
  107. {
  108. for(int i = 0; i < 1; i++) {
  109. CALL_SUBTEST_5( qr<Matrix3f>() );
  110. CALL_SUBTEST_6( qr<Matrix3d>() );
  111. CALL_SUBTEST_8( qr<Matrix2f>() );
  112. CALL_SUBTEST_1( qr<MatrixXf>() );
  113. CALL_SUBTEST_2( qr<MatrixXd>() );
  114. CALL_SUBTEST_3( qr<MatrixXcd>() );
  115. }
  116. for(int i = 0; i < g_repeat; i++) {
  117. CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
  118. CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
  119. CALL_SUBTEST_4( qr_invertible<MatrixXcf>() );
  120. CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
  121. }
  122. CALL_SUBTEST_5(qr_verify_assert<Matrix3f>());
  123. CALL_SUBTEST_6(qr_verify_assert<Matrix3d>());
  124. CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
  125. CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
  126. CALL_SUBTEST_4(qr_verify_assert<MatrixXcf>());
  127. CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
  128. // Test problem size constructors
  129. CALL_SUBTEST_7(FullPivHouseholderQR<MatrixXf>(10, 20));
  130. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(10,20)));
  131. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(Matrix<float,10,20>::Random())));
  132. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(20,10)));
  133. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(Matrix<float,20,10>::Random())));
  134. }