qr_colpivoting.cpp 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368
  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 <Eigen/SVD>
  13. #include "solverbase.h"
  14. template <typename MatrixType>
  15. void cod() {
  16. STATIC_CHECK(( internal::is_same<typename CompleteOrthogonalDecomposition<MatrixType>::StorageIndex,int>::value ));
  17. Index rows = internal::random<Index>(2, EIGEN_TEST_MAX_SIZE);
  18. Index cols = internal::random<Index>(2, EIGEN_TEST_MAX_SIZE);
  19. Index cols2 = internal::random<Index>(2, EIGEN_TEST_MAX_SIZE);
  20. Index rank = internal::random<Index>(1, (std::min)(rows, cols) - 1);
  21. typedef typename MatrixType::Scalar Scalar;
  22. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime,
  23. MatrixType::RowsAtCompileTime>
  24. MatrixQType;
  25. MatrixType matrix;
  26. createRandomPIMatrixOfRank(rank, rows, cols, matrix);
  27. CompleteOrthogonalDecomposition<MatrixType> cod(matrix);
  28. VERIFY(rank == cod.rank());
  29. VERIFY(cols - cod.rank() == cod.dimensionOfKernel());
  30. VERIFY(!cod.isInjective());
  31. VERIFY(!cod.isInvertible());
  32. VERIFY(!cod.isSurjective());
  33. MatrixQType q = cod.householderQ();
  34. VERIFY_IS_UNITARY(q);
  35. MatrixType z = cod.matrixZ();
  36. VERIFY_IS_UNITARY(z);
  37. MatrixType t;
  38. t.setZero(rows, cols);
  39. t.topLeftCorner(rank, rank) =
  40. cod.matrixT().topLeftCorner(rank, rank).template triangularView<Upper>();
  41. MatrixType c = q * t * z * cod.colsPermutation().inverse();
  42. VERIFY_IS_APPROX(matrix, c);
  43. check_solverbase<MatrixType, MatrixType>(matrix, cod, rows, cols, cols2);
  44. // Verify that we get the same minimum-norm solution as the SVD.
  45. MatrixType exact_solution = MatrixType::Random(cols, cols2);
  46. MatrixType rhs = matrix * exact_solution;
  47. MatrixType cod_solution = cod.solve(rhs);
  48. JacobiSVD<MatrixType> svd(matrix, ComputeThinU | ComputeThinV);
  49. MatrixType svd_solution = svd.solve(rhs);
  50. VERIFY_IS_APPROX(cod_solution, svd_solution);
  51. MatrixType pinv = cod.pseudoInverse();
  52. VERIFY_IS_APPROX(cod_solution, pinv * rhs);
  53. }
  54. template <typename MatrixType, int Cols2>
  55. void cod_fixedsize() {
  56. enum {
  57. Rows = MatrixType::RowsAtCompileTime,
  58. Cols = MatrixType::ColsAtCompileTime
  59. };
  60. typedef typename MatrixType::Scalar Scalar;
  61. typedef CompleteOrthogonalDecomposition<Matrix<Scalar, Rows, Cols> > COD;
  62. int rank = internal::random<int>(1, (std::min)(int(Rows), int(Cols)) - 1);
  63. Matrix<Scalar, Rows, Cols> matrix;
  64. createRandomPIMatrixOfRank(rank, Rows, Cols, matrix);
  65. COD cod(matrix);
  66. VERIFY(rank == cod.rank());
  67. VERIFY(Cols - cod.rank() == cod.dimensionOfKernel());
  68. VERIFY(cod.isInjective() == (rank == Rows));
  69. VERIFY(cod.isSurjective() == (rank == Cols));
  70. VERIFY(cod.isInvertible() == (cod.isInjective() && cod.isSurjective()));
  71. check_solverbase<Matrix<Scalar, Cols, Cols2>, Matrix<Scalar, Rows, Cols2> >(matrix, cod, Rows, Cols, Cols2);
  72. // Verify that we get the same minimum-norm solution as the SVD.
  73. Matrix<Scalar, Cols, Cols2> exact_solution;
  74. exact_solution.setRandom(Cols, Cols2);
  75. Matrix<Scalar, Rows, Cols2> rhs = matrix * exact_solution;
  76. Matrix<Scalar, Cols, Cols2> cod_solution = cod.solve(rhs);
  77. JacobiSVD<MatrixType> svd(matrix, ComputeFullU | ComputeFullV);
  78. Matrix<Scalar, Cols, Cols2> svd_solution = svd.solve(rhs);
  79. VERIFY_IS_APPROX(cod_solution, svd_solution);
  80. typename Inverse<COD>::PlainObject pinv = cod.pseudoInverse();
  81. VERIFY_IS_APPROX(cod_solution, pinv * rhs);
  82. }
  83. template<typename MatrixType> void qr()
  84. {
  85. using std::sqrt;
  86. STATIC_CHECK(( internal::is_same<typename ColPivHouseholderQR<MatrixType>::StorageIndex,int>::value ));
  87. 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);
  88. Index rank = internal::random<Index>(1, (std::min)(rows, cols)-1);
  89. typedef typename MatrixType::Scalar Scalar;
  90. typedef typename MatrixType::RealScalar RealScalar;
  91. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> MatrixQType;
  92. MatrixType m1;
  93. createRandomPIMatrixOfRank(rank,rows,cols,m1);
  94. ColPivHouseholderQR<MatrixType> qr(m1);
  95. VERIFY_IS_EQUAL(rank, qr.rank());
  96. VERIFY_IS_EQUAL(cols - qr.rank(), qr.dimensionOfKernel());
  97. VERIFY(!qr.isInjective());
  98. VERIFY(!qr.isInvertible());
  99. VERIFY(!qr.isSurjective());
  100. MatrixQType q = qr.householderQ();
  101. VERIFY_IS_UNITARY(q);
  102. MatrixType r = qr.matrixQR().template triangularView<Upper>();
  103. MatrixType c = q * r * qr.colsPermutation().inverse();
  104. VERIFY_IS_APPROX(m1, c);
  105. // Verify that the absolute value of the diagonal elements in R are
  106. // non-increasing until they reach the singularity threshold.
  107. RealScalar threshold =
  108. sqrt(RealScalar(rows)) * numext::abs(r(0, 0)) * NumTraits<Scalar>::epsilon();
  109. for (Index i = 0; i < (std::min)(rows, cols) - 1; ++i) {
  110. RealScalar x = numext::abs(r(i, i));
  111. RealScalar y = numext::abs(r(i + 1, i + 1));
  112. if (x < threshold && y < threshold) continue;
  113. if (!test_isApproxOrLessThan(y, x)) {
  114. for (Index j = 0; j < (std::min)(rows, cols); ++j) {
  115. std::cout << "i = " << j << ", |r_ii| = " << numext::abs(r(j, j)) << std::endl;
  116. }
  117. std::cout << "Failure at i=" << i << ", rank=" << rank
  118. << ", threshold=" << threshold << std::endl;
  119. }
  120. VERIFY_IS_APPROX_OR_LESS_THAN(y, x);
  121. }
  122. check_solverbase<MatrixType, MatrixType>(m1, qr, rows, cols, cols2);
  123. {
  124. MatrixType m2, m3;
  125. Index size = rows;
  126. do {
  127. m1 = MatrixType::Random(size,size);
  128. qr.compute(m1);
  129. } while(!qr.isInvertible());
  130. MatrixType m1_inv = qr.inverse();
  131. m3 = m1 * MatrixType::Random(size,cols2);
  132. m2 = qr.solve(m3);
  133. VERIFY_IS_APPROX(m2, m1_inv*m3);
  134. }
  135. }
  136. template<typename MatrixType, int Cols2> void qr_fixedsize()
  137. {
  138. using std::sqrt;
  139. using std::abs;
  140. enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
  141. typedef typename MatrixType::Scalar Scalar;
  142. typedef typename MatrixType::RealScalar RealScalar;
  143. int rank = internal::random<int>(1, (std::min)(int(Rows), int(Cols))-1);
  144. Matrix<Scalar,Rows,Cols> m1;
  145. createRandomPIMatrixOfRank(rank,Rows,Cols,m1);
  146. ColPivHouseholderQR<Matrix<Scalar,Rows,Cols> > qr(m1);
  147. VERIFY_IS_EQUAL(rank, qr.rank());
  148. VERIFY_IS_EQUAL(Cols - qr.rank(), qr.dimensionOfKernel());
  149. VERIFY_IS_EQUAL(qr.isInjective(), (rank == Rows));
  150. VERIFY_IS_EQUAL(qr.isSurjective(), (rank == Cols));
  151. VERIFY_IS_EQUAL(qr.isInvertible(), (qr.isInjective() && qr.isSurjective()));
  152. Matrix<Scalar,Rows,Cols> r = qr.matrixQR().template triangularView<Upper>();
  153. Matrix<Scalar,Rows,Cols> c = qr.householderQ() * r * qr.colsPermutation().inverse();
  154. VERIFY_IS_APPROX(m1, c);
  155. check_solverbase<Matrix<Scalar,Cols,Cols2>, Matrix<Scalar,Rows,Cols2> >(m1, qr, Rows, Cols, Cols2);
  156. // Verify that the absolute value of the diagonal elements in R are
  157. // non-increasing until they reache the singularity threshold.
  158. RealScalar threshold =
  159. sqrt(RealScalar(Rows)) * (std::abs)(r(0, 0)) * NumTraits<Scalar>::epsilon();
  160. for (Index i = 0; i < (std::min)(int(Rows), int(Cols)) - 1; ++i) {
  161. RealScalar x = numext::abs(r(i, i));
  162. RealScalar y = numext::abs(r(i + 1, i + 1));
  163. if (x < threshold && y < threshold) continue;
  164. if (!test_isApproxOrLessThan(y, x)) {
  165. for (Index j = 0; j < (std::min)(int(Rows), int(Cols)); ++j) {
  166. std::cout << "i = " << j << ", |r_ii| = " << numext::abs(r(j, j)) << std::endl;
  167. }
  168. std::cout << "Failure at i=" << i << ", rank=" << rank
  169. << ", threshold=" << threshold << std::endl;
  170. }
  171. VERIFY_IS_APPROX_OR_LESS_THAN(y, x);
  172. }
  173. }
  174. // This test is meant to verify that pivots are chosen such that
  175. // even for a graded matrix, the diagonal of R falls of roughly
  176. // monotonically until it reaches the threshold for singularity.
  177. // We use the so-called Kahan matrix, which is a famous counter-example
  178. // for rank-revealing QR. See
  179. // http://www.netlib.org/lapack/lawnspdf/lawn176.pdf
  180. // page 3 for more detail.
  181. template<typename MatrixType> void qr_kahan_matrix()
  182. {
  183. using std::sqrt;
  184. using std::abs;
  185. typedef typename MatrixType::Scalar Scalar;
  186. typedef typename MatrixType::RealScalar RealScalar;
  187. Index rows = 300, cols = rows;
  188. MatrixType m1;
  189. m1.setZero(rows,cols);
  190. RealScalar s = std::pow(NumTraits<RealScalar>::epsilon(), 1.0 / rows);
  191. RealScalar c = std::sqrt(1 - s*s);
  192. RealScalar pow_s_i(1.0); // pow(s,i)
  193. for (Index i = 0; i < rows; ++i) {
  194. m1(i, i) = pow_s_i;
  195. m1.row(i).tail(rows - i - 1) = -pow_s_i * c * MatrixType::Ones(1, rows - i - 1);
  196. pow_s_i *= s;
  197. }
  198. m1 = (m1 + m1.transpose()).eval();
  199. ColPivHouseholderQR<MatrixType> qr(m1);
  200. MatrixType r = qr.matrixQR().template triangularView<Upper>();
  201. RealScalar threshold =
  202. std::sqrt(RealScalar(rows)) * numext::abs(r(0, 0)) * NumTraits<Scalar>::epsilon();
  203. for (Index i = 0; i < (std::min)(rows, cols) - 1; ++i) {
  204. RealScalar x = numext::abs(r(i, i));
  205. RealScalar y = numext::abs(r(i + 1, i + 1));
  206. if (x < threshold && y < threshold) continue;
  207. if (!test_isApproxOrLessThan(y, x)) {
  208. for (Index j = 0; j < (std::min)(rows, cols); ++j) {
  209. std::cout << "i = " << j << ", |r_ii| = " << numext::abs(r(j, j)) << std::endl;
  210. }
  211. std::cout << "Failure at i=" << i << ", rank=" << qr.rank()
  212. << ", threshold=" << threshold << std::endl;
  213. }
  214. VERIFY_IS_APPROX_OR_LESS_THAN(y, x);
  215. }
  216. }
  217. template<typename MatrixType> void qr_invertible()
  218. {
  219. using std::log;
  220. using std::abs;
  221. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  222. typedef typename MatrixType::Scalar Scalar;
  223. int size = internal::random<int>(10,50);
  224. MatrixType m1(size, size), m2(size, size), m3(size, size);
  225. m1 = MatrixType::Random(size,size);
  226. if (internal::is_same<RealScalar,float>::value)
  227. {
  228. // let's build a matrix more stable to inverse
  229. MatrixType a = MatrixType::Random(size,size*2);
  230. m1 += a * a.adjoint();
  231. }
  232. ColPivHouseholderQR<MatrixType> qr(m1);
  233. check_solverbase<MatrixType, MatrixType>(m1, qr, size, size, size);
  234. // now construct a matrix with prescribed determinant
  235. m1.setZero();
  236. for(int i = 0; i < size; i++) m1(i,i) = internal::random<Scalar>();
  237. RealScalar absdet = abs(m1.diagonal().prod());
  238. m3 = qr.householderQ(); // get a unitary
  239. m1 = m3 * m1 * m3;
  240. qr.compute(m1);
  241. VERIFY_IS_APPROX(absdet, qr.absDeterminant());
  242. VERIFY_IS_APPROX(log(absdet), qr.logAbsDeterminant());
  243. }
  244. template<typename MatrixType> void qr_verify_assert()
  245. {
  246. MatrixType tmp;
  247. ColPivHouseholderQR<MatrixType> qr;
  248. VERIFY_RAISES_ASSERT(qr.matrixQR())
  249. VERIFY_RAISES_ASSERT(qr.solve(tmp))
  250. VERIFY_RAISES_ASSERT(qr.transpose().solve(tmp))
  251. VERIFY_RAISES_ASSERT(qr.adjoint().solve(tmp))
  252. VERIFY_RAISES_ASSERT(qr.householderQ())
  253. VERIFY_RAISES_ASSERT(qr.dimensionOfKernel())
  254. VERIFY_RAISES_ASSERT(qr.isInjective())
  255. VERIFY_RAISES_ASSERT(qr.isSurjective())
  256. VERIFY_RAISES_ASSERT(qr.isInvertible())
  257. VERIFY_RAISES_ASSERT(qr.inverse())
  258. VERIFY_RAISES_ASSERT(qr.absDeterminant())
  259. VERIFY_RAISES_ASSERT(qr.logAbsDeterminant())
  260. }
  261. template<typename MatrixType> void cod_verify_assert()
  262. {
  263. MatrixType tmp;
  264. CompleteOrthogonalDecomposition<MatrixType> cod;
  265. VERIFY_RAISES_ASSERT(cod.matrixQTZ())
  266. VERIFY_RAISES_ASSERT(cod.solve(tmp))
  267. VERIFY_RAISES_ASSERT(cod.transpose().solve(tmp))
  268. VERIFY_RAISES_ASSERT(cod.adjoint().solve(tmp))
  269. VERIFY_RAISES_ASSERT(cod.householderQ())
  270. VERIFY_RAISES_ASSERT(cod.dimensionOfKernel())
  271. VERIFY_RAISES_ASSERT(cod.isInjective())
  272. VERIFY_RAISES_ASSERT(cod.isSurjective())
  273. VERIFY_RAISES_ASSERT(cod.isInvertible())
  274. VERIFY_RAISES_ASSERT(cod.pseudoInverse())
  275. VERIFY_RAISES_ASSERT(cod.absDeterminant())
  276. VERIFY_RAISES_ASSERT(cod.logAbsDeterminant())
  277. }
  278. EIGEN_DECLARE_TEST(qr_colpivoting)
  279. {
  280. for(int i = 0; i < g_repeat; i++) {
  281. CALL_SUBTEST_1( qr<MatrixXf>() );
  282. CALL_SUBTEST_2( qr<MatrixXd>() );
  283. CALL_SUBTEST_3( qr<MatrixXcd>() );
  284. CALL_SUBTEST_4(( qr_fixedsize<Matrix<float,3,5>, 4 >() ));
  285. CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,6,2>, 3 >() ));
  286. CALL_SUBTEST_5(( qr_fixedsize<Matrix<double,1,1>, 1 >() ));
  287. }
  288. for(int i = 0; i < g_repeat; i++) {
  289. CALL_SUBTEST_1( cod<MatrixXf>() );
  290. CALL_SUBTEST_2( cod<MatrixXd>() );
  291. CALL_SUBTEST_3( cod<MatrixXcd>() );
  292. CALL_SUBTEST_4(( cod_fixedsize<Matrix<float,3,5>, 4 >() ));
  293. CALL_SUBTEST_5(( cod_fixedsize<Matrix<double,6,2>, 3 >() ));
  294. CALL_SUBTEST_5(( cod_fixedsize<Matrix<double,1,1>, 1 >() ));
  295. }
  296. for(int i = 0; i < g_repeat; i++) {
  297. CALL_SUBTEST_1( qr_invertible<MatrixXf>() );
  298. CALL_SUBTEST_2( qr_invertible<MatrixXd>() );
  299. CALL_SUBTEST_6( qr_invertible<MatrixXcf>() );
  300. CALL_SUBTEST_3( qr_invertible<MatrixXcd>() );
  301. }
  302. CALL_SUBTEST_7(qr_verify_assert<Matrix3f>());
  303. CALL_SUBTEST_8(qr_verify_assert<Matrix3d>());
  304. CALL_SUBTEST_1(qr_verify_assert<MatrixXf>());
  305. CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
  306. CALL_SUBTEST_6(qr_verify_assert<MatrixXcf>());
  307. CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
  308. CALL_SUBTEST_7(cod_verify_assert<Matrix3f>());
  309. CALL_SUBTEST_8(cod_verify_assert<Matrix3d>());
  310. CALL_SUBTEST_1(cod_verify_assert<MatrixXf>());
  311. CALL_SUBTEST_2(cod_verify_assert<MatrixXd>());
  312. CALL_SUBTEST_6(cod_verify_assert<MatrixXcf>());
  313. CALL_SUBTEST_3(cod_verify_assert<MatrixXcd>());
  314. // Test problem size constructors
  315. CALL_SUBTEST_9(ColPivHouseholderQR<MatrixXf>(10, 20));
  316. CALL_SUBTEST_1( qr_kahan_matrix<MatrixXf>() );
  317. CALL_SUBTEST_2( qr_kahan_matrix<MatrixXd>() );
  318. }