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.

106 lines
3.0 KiB

  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 SparseMatrix<Scalar,ColMajor> MatrixType;
  42. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMat;
  43. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  44. MatrixType A;
  45. DenseMat dA;
  46. DenseVector refX,x,b;
  47. SparseQR<MatrixType, COLAMDOrdering<int> > solver;
  48. generate_sparse_rectangular_problem(A,dA);
  49. b = dA * DenseVector::Random(A.cols());
  50. solver.compute(A);
  51. if(internal::random<float>(0,1)>0.5)
  52. solver.factorize(A); // this checks that calling analyzePattern is not needed if the pattern do not change.
  53. if (solver.info() != Success)
  54. {
  55. std::cerr << "sparse QR factorization failed\n";
  56. exit(0);
  57. return;
  58. }
  59. x = solver.solve(b);
  60. if (solver.info() != Success)
  61. {
  62. std::cerr << "sparse QR factorization failed\n";
  63. exit(0);
  64. return;
  65. }
  66. VERIFY_IS_APPROX(A * x, b);
  67. //Compare with a dense QR solver
  68. ColPivHouseholderQR<DenseMat> dqr(dA);
  69. refX = dqr.solve(b);
  70. VERIFY_IS_EQUAL(dqr.rank(), solver.rank());
  71. if(solver.rank()==A.cols()) // full rank
  72. VERIFY_IS_APPROX(x, refX);
  73. // else
  74. // VERIFY((dA * refX - b).norm() * 2 > (A * x - b).norm() );
  75. // Compute explicitly the matrix Q
  76. MatrixType Q, QtQ, idM;
  77. Q = solver.matrixQ();
  78. //Check ||Q' * Q - I ||
  79. QtQ = Q * Q.adjoint();
  80. idM.resize(Q.rows(), Q.rows()); idM.setIdentity();
  81. VERIFY(idM.isApprox(QtQ));
  82. // Q to dense
  83. DenseMat dQ;
  84. dQ = solver.matrixQ();
  85. VERIFY_IS_APPROX(Q, dQ);
  86. }
  87. void test_sparseqr()
  88. {
  89. for(int i=0; i<g_repeat; ++i)
  90. {
  91. CALL_SUBTEST_1(test_sparseqr_scalar<double>());
  92. CALL_SUBTEST_2(test_sparseqr_scalar<std::complex<double> >());
  93. }
  94. }