sparseqr.cpp 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2012 Desire Nuentsa Wakam <desire.nuentsa_wakam@inria.fr>
  5. // Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
  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. #include "sparse.h"
  10. #include <Eigen/SparseQR>
  11. template<typename MatrixType,typename DenseMat>
  12. int generate_sparse_rectangular_problem(MatrixType& A, DenseMat& dA, int maxRows = 300, int maxCols = 150)
  13. {
  14. eigen_assert(maxRows >= maxCols);
  15. typedef typename MatrixType::Scalar Scalar;
  16. int rows = internal::random<int>(1,maxRows);
  17. int cols = internal::random<int>(1,maxCols);
  18. double density = (std::max)(8./(rows*cols), 0.01);
  19. A.resize(rows,cols);
  20. dA.resize(rows,cols);
  21. initSparse<Scalar>(density, dA, A,ForceNonZeroDiag);
  22. A.makeCompressed();
  23. int nop = internal::random<int>(0, internal::random<double>(0,1) > 0.5 ? cols/2 : 0);
  24. for(int k=0; k<nop; ++k)
  25. {
  26. int j0 = internal::random<int>(0,cols-1);
  27. int j1 = internal::random<int>(0,cols-1);
  28. Scalar s = internal::random<Scalar>();
  29. A.col(j0) = s * A.col(j1);
  30. dA.col(j0) = s * dA.col(j1);
  31. }
  32. // if(rows<cols) {
  33. // A.conservativeResize(cols,cols);
  34. // dA.conservativeResize(cols,cols);
  35. // dA.bottomRows(cols-rows).setZero();
  36. // }
  37. return rows;
  38. }
  39. template<typename Scalar> void test_sparseqr_scalar()
  40. {
  41. typedef typename NumTraits<Scalar>::Real RealScalar;
  42. typedef SparseMatrix<Scalar,ColMajor> MatrixType;
  43. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMat;
  44. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  45. MatrixType A;
  46. DenseMat dA;
  47. DenseVector refX,x,b;
  48. SparseQR<MatrixType, COLAMDOrdering<int> > solver;
  49. generate_sparse_rectangular_problem(A,dA);
  50. b = dA * DenseVector::Random(A.cols());
  51. solver.compute(A);
  52. // Q should be MxM
  53. VERIFY_IS_EQUAL(solver.matrixQ().rows(), A.rows());
  54. VERIFY_IS_EQUAL(solver.matrixQ().cols(), A.rows());
  55. // R should be MxN
  56. VERIFY_IS_EQUAL(solver.matrixR().rows(), A.rows());
  57. VERIFY_IS_EQUAL(solver.matrixR().cols(), A.cols());
  58. // Q and R can be multiplied
  59. DenseMat recoveredA = solver.matrixQ()
  60. * DenseMat(solver.matrixR().template triangularView<Upper>())
  61. * solver.colsPermutation().transpose();
  62. VERIFY_IS_EQUAL(recoveredA.rows(), A.rows());
  63. VERIFY_IS_EQUAL(recoveredA.cols(), A.cols());
  64. // and in the full rank case the original matrix is recovered
  65. if (solver.rank() == A.cols())
  66. {
  67. VERIFY_IS_APPROX(A, recoveredA);
  68. }
  69. if(internal::random<float>(0,1)>0.5f)
  70. solver.factorize(A); // this checks that calling analyzePattern is not needed if the pattern do not change.
  71. if (solver.info() != Success)
  72. {
  73. std::cerr << "sparse QR factorization failed\n";
  74. exit(0);
  75. return;
  76. }
  77. x = solver.solve(b);
  78. if (solver.info() != Success)
  79. {
  80. std::cerr << "sparse QR factorization failed\n";
  81. exit(0);
  82. return;
  83. }
  84. // Compare with a dense QR solver
  85. ColPivHouseholderQR<DenseMat> dqr(dA);
  86. refX = dqr.solve(b);
  87. bool rank_deficient = A.cols()>A.rows() || dqr.rank()<A.cols();
  88. if(rank_deficient)
  89. {
  90. // rank deficient problem -> we might have to increase the threshold
  91. // to get a correct solution.
  92. RealScalar th = RealScalar(20)*dA.colwise().norm().maxCoeff()*(A.rows()+A.cols()) * NumTraits<RealScalar>::epsilon();
  93. for(Index k=0; (k<16) && !test_isApprox(A*x,b); ++k)
  94. {
  95. th *= RealScalar(10);
  96. solver.setPivotThreshold(th);
  97. solver.compute(A);
  98. x = solver.solve(b);
  99. }
  100. }
  101. VERIFY_IS_APPROX(A * x, b);
  102. // For rank deficient problem, the estimated rank might
  103. // be slightly off, so let's only raise a warning in such cases.
  104. if(rank_deficient) ++g_test_level;
  105. VERIFY_IS_EQUAL(solver.rank(), dqr.rank());
  106. if(rank_deficient) --g_test_level;
  107. if(solver.rank()==A.cols()) // full rank
  108. VERIFY_IS_APPROX(x, refX);
  109. // else
  110. // VERIFY((dA * refX - b).norm() * 2 > (A * x - b).norm() );
  111. // Compute explicitly the matrix Q
  112. MatrixType Q, QtQ, idM;
  113. Q = solver.matrixQ();
  114. //Check ||Q' * Q - I ||
  115. QtQ = Q * Q.adjoint();
  116. idM.resize(Q.rows(), Q.rows()); idM.setIdentity();
  117. VERIFY(idM.isApprox(QtQ));
  118. // Q to dense
  119. DenseMat dQ;
  120. dQ = solver.matrixQ();
  121. VERIFY_IS_APPROX(Q, dQ);
  122. }
  123. EIGEN_DECLARE_TEST(sparseqr)
  124. {
  125. for(int i=0; i<g_repeat; ++i)
  126. {
  127. CALL_SUBTEST_1(test_sparseqr_scalar<double>());
  128. CALL_SUBTEST_2(test_sparseqr_scalar<std::complex<double> >());
  129. }
  130. }