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.

150 lines
5.3 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>(2,EIGEN_TEST_MAX_SIZE), cols = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE), cols2 = internal::random<Index>(2,EIGEN_TEST_MAX_SIZE);
  16. Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
  17. typedef typename MatrixType::Scalar Scalar;
  18. typedef typename MatrixType::RealScalar RealScalar;
  19. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
  20. typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
  21. MatrixType m1;
  22. createRandomPIMatrixOfRank(rank,rows,cols,m1);
  23. ColPivHouseholderQR<MatrixType> qr(m1);
  24. VERIFY(rank == qr.rank());
  25. VERIFY(cols - qr.rank() == qr.dimensionOfKernel());
  26. VERIFY(!qr.isInjective());
  27. VERIFY(!qr.isInvertible());
  28. VERIFY(!qr.isSurjective());
  29. MatrixQType q = qr.householderQ();
  30. VERIFY_IS_UNITARY(q);
  31. MatrixType r = qr.matrixQR().template triangularView<Upper>();
  32. MatrixType c = q * r * qr.colsPermutation().inverse();
  33. VERIFY_IS_APPROX(m1, c);
  34. MatrixType m2 = MatrixType::Random(cols,cols2);
  35. MatrixType m3 = m1*m2;
  36. m2 = MatrixType::Random(cols,cols2);
  37. m2 = qr.solve(m3);
  38. VERIFY_IS_APPROX(m3, m1*m2);
  39. }
  40. template<typename MatrixType, int Cols2> void qr_fixedsize()
  41. {
  42. enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
  43. typedef typename MatrixType::Scalar Scalar;
  44. int rank = internal::random<int>(1, (std::min)(int(Rows), int(Cols))-1);
  45. Matrix<Scalar,Rows,Cols> m1;
  46. createRandomPIMatrixOfRank(rank,Rows,Cols,m1);
  47. ColPivHouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
  48. VERIFY(rank == qr.rank());
  49. VERIFY(Cols - qr.rank() == qr.dimensionOfKernel());
  50. VERIFY(qr.isInjective() == (rank == Rows));
  51. VERIFY(qr.isSurjective() == (rank == Cols));
  52. VERIFY(qr.isInvertible() == (qr.isInjective() && qr.isSurjective()));
  53. Matrix<Scalar,Rows,Cols> r = qr.matrixQR().template triangularView<Upper>();
  54. Matrix<Scalar,Rows,Cols> c = qr.householderQ() * r * qr.colsPermutation().inverse();
  55. VERIFY_IS_APPROX(m1, c);
  56. Matrix<Scalar,Cols,Cols2> m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
  57. Matrix<Scalar,Rows,Cols2> m3 = m1*m2;
  58. m2 = Matrix<Scalar,Cols,Cols2>::Random(Cols,Cols2);
  59. m2 = qr.solve(m3);
  60. VERIFY_IS_APPROX(m3, m1*m2);
  61. }
  62. template<typename MatrixType> void qr_invertible()
  63. {
  64. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  65. typedef typename MatrixType::Scalar Scalar;
  66. int size = internal::random<int>(10,50);
  67. MatrixType m1(size, size), m2(size, size), m3(size, size);
  68. m1 = MatrixType::Random(size,size);
  69. if (internal::is_same<RealScalar,float>::value)
  70. {
  71. // let's build a matrix more stable to inverse
  72. MatrixType a = MatrixType::Random(size,size*2);
  73. m1 += a * a.adjoint();
  74. }
  75. ColPivHouseholderQR<MatrixType> qr(m1);
  76. m3 = MatrixType::Random(size,size);
  77. m2 = qr.solve(m3);
  78. //VERIFY_IS_APPROX(m3, m1*m2);
  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 = internal::abs(m1.diagonal().prod());
  83. m3 = qr.householderQ(); // get a unitary
  84. m1 = m3 * m1 * m3;
  85. qr.compute(m1);
  86. VERIFY_IS_APPROX(absdet, qr.absDeterminant());
  87. VERIFY_IS_APPROX(internal::log(absdet), qr.logAbsDeterminant());
  88. }
  89. template<typename MatrixType> void qr_verify_assert()
  90. {
  91. MatrixType tmp;
  92. ColPivHouseholderQR<MatrixType> qr;
  93. VERIFY_RAISES_ASSERT(qr.matrixQR())
  94. VERIFY_RAISES_ASSERT(qr.solve(tmp))
  95. VERIFY_RAISES_ASSERT(qr.householderQ())
  96. VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
  97. VERIFY_RAISES_ASSERT(qr.isInjective())
  98. VERIFY_RAISES_ASSERT(qr.isSurjective())
  99. VERIFY_RAISES_ASSERT(qr.isInvertible())
  100. VERIFY_RAISES_ASSERT(qr.inverse())
  101. VERIFY_RAISES_ASSERT(qr.absDeterminant())
  102. VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
  103. }
  104. void test_qr_colpivoting()
  105. {
  106. for(int i = 0; i < g_repeat; i++) {
  107. CALL_SUBTEST_1( qr<MatrixXf>() );
  108. CALL_SUBTEST_2( qr<MatrixXd>() );
  109. CALL_SUBTEST_3( qr<MatrixXcd>() );
  110. CALL_SUBTEST_4(( qr_fixedsize<Matrix<float,3,5>, 4 >() ));
  111. CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,6,2>, 3 >() ));
  112. CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,1,1>, 1 >() ));
  113. }
  114. for(int i = 0; i < g_repeat; i++) {
  115. CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
  116. CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
  117. CALL_SUBTEST_6( qr_invertible<MatrixXcf>() );
  118. CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
  119. }
  120. CALL_SUBTEST_7(qr_verify_assert<Matrix3f>());
  121. CALL_SUBTEST_8(qr_verify_assert<Matrix3d>());
  122. CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
  123. CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
  124. CALL_SUBTEST_6(qr_verify_assert<MatrixXcf>());
  125. CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
  126. // Test problem size constructors
  127. CALL_SUBTEST_9(ColPivHouseholderQR<MatrixXf>(10, 20));
  128. }