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.

141 lines
4.6 KiB

  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. template<typename MatrixType> void qr()
  13. {
  14. typedef typename MatrixType::Index Index;
  15. Index rows = internal::random<Index>(20,200), cols = internal::random<int>(20,200), cols2 = internal::random<int>(20,200);
  16. Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
  17. typedef typename MatrixType::Scalar Scalar;
  18. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
  19. MatrixType m1;
  20. createRandomPIMatrixOfRank(rank,rows,cols,m1);
  21. FullPivHouseholderQR<MatrixType> qr(m1);
  22. VERIFY_IS_EQUAL(rank, qr.rank());
  23. VERIFY_IS_EQUAL(cols - qr.rank(), qr.dimensionOfKernel());
  24. VERIFY(!qr.isInjective());
  25. VERIFY(!qr.isInvertible());
  26. VERIFY(!qr.isSurjective());
  27. MatrixType r = qr.matrixQR();
  28. MatrixQType q = qr.matrixQ();
  29. VERIFY_IS_UNITARY(q);
  30. // FIXME need better way to construct trapezoid
  31. for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
  32. MatrixType c = qr.matrixQ() * r * qr.colsPermutation().inverse();
  33. VERIFY_IS_APPROX(m1, c);
  34. // stress the ReturnByValue mechanism
  35. MatrixType tmp;
  36. VERIFY_IS_APPROX(tmp.noalias() = qr.matrixQ() * r, (qr.matrixQ() * r).eval());
  37. MatrixType m2 = MatrixType::Random(cols,cols2);
  38. MatrixType m3 = m1*m2;
  39. m2 = MatrixType::Random(cols,cols2);
  40. m2 = qr.solve(m3);
  41. VERIFY_IS_APPROX(m3, m1*m2);
  42. }
  43. template<typename MatrixType> void qr_invertible()
  44. {
  45. using std::log;
  46. using std::abs;
  47. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  48. typedef typename MatrixType::Scalar Scalar;
  49. int size = internal::random<int>(10,50);
  50. MatrixType m1(size, size), m2(size, size), m3(size, size);
  51. m1 = MatrixType::Random(size,size);
  52. if (internal::is_same<RealScalar,float>::value)
  53. {
  54. // let's build a matrix more stable to inverse
  55. MatrixType a = MatrixType::Random(size,size*2);
  56. m1 += a * a.adjoint();
  57. }
  58. FullPivHouseholderQR<MatrixType> qr(m1);
  59. VERIFY(qr.isInjective());
  60. VERIFY(qr.isInvertible());
  61. VERIFY(qr.isSurjective());
  62. m3 = MatrixType::Random(size,size);
  63. m2 = qr.solve(m3);
  64. VERIFY_IS_APPROX(m3, m1*m2);
  65. // now construct a matrix with prescribed determinant
  66. m1.setZero();
  67. for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>();
  68. RealScalar absdet = abs(m1.diagonal().prod());
  69. m3 = qr.matrixQ(); // get a unitary
  70. m1 = m3 * m1 * m3;
  71. qr.compute(m1);
  72. VERIFY_IS_APPROX(absdet, qr.absDeterminant());
  73. VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant());
  74. }
  75. template<typename MatrixType> void qr_verify_assert()
  76. {
  77. MatrixType tmp;
  78. FullPivHouseholderQR<MatrixType> qr;
  79. VERIFY_RAISES_ASSERT(qr.matrixQR())
  80. VERIFY_RAISES_ASSERT(qr.solve(tmp))
  81. VERIFY_RAISES_ASSERT(qr.matrixQ())
  82. VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
  83. VERIFY_RAISES_ASSERT(qr.isInjective())
  84. VERIFY_RAISES_ASSERT(qr.isSurjective())
  85. VERIFY_RAISES_ASSERT(qr.isInvertible())
  86. VERIFY_RAISES_ASSERT(qr.inverse())
  87. VERIFY_RAISES_ASSERT(qr.absDeterminant())
  88. VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
  89. }
  90. void test_qr_fullpivoting()
  91. {
  92. for(int i = 0; i < 1; i++) {
  93. // FIXME : very weird bug here
  94. // CALL_SUBTEST(qr(Matrix2f()) );
  95. CALL_SUBTEST_1( qr<MatrixXf>() );
  96. CALL_SUBTEST_2( qr<MatrixXd>() );
  97. CALL_SUBTEST_3( qr<MatrixXcd>() );
  98. }
  99. for(int i = 0; i < g_repeat; i++) {
  100. CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
  101. CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
  102. CALL_SUBTEST_4( qr_invertible<MatrixXcf>() );
  103. CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
  104. }
  105. CALL_SUBTEST_5(qr_verify_assert<Matrix3f>());
  106. CALL_SUBTEST_6(qr_verify_assert<Matrix3d>());
  107. CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
  108. CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
  109. CALL_SUBTEST_4(qr_verify_assert<MatrixXcf>());
  110. CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
  111. // Test problem size constructors
  112. CALL_SUBTEST_7(FullPivHouseholderQR<MatrixXf>(10, 20));
  113. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(10,20)));
  114. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,10,20> >(Matrix<float,10,20>::Random())));
  115. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(20,10)));
  116. CALL_SUBTEST_7((FullPivHouseholderQR<Matrix<float,20,10> >(Matrix<float,20,10>::Random())));
  117. }