matrixfree_cg.cpp 4.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129
  1. #include <iostream>
  2. #include <Eigen/Core>
  3. #include <Eigen/Dense>
  4. #include <Eigen/IterativeLinearSolvers>
  5. #include <unsupported/Eigen/IterativeSolvers>
  6. class MatrixReplacement;
  7. using Eigen::SparseMatrix;
  8. namespace Eigen {
  9. namespace internal {
  10. // MatrixReplacement looks-like a SparseMatrix, so let's inherits its traits:
  11. template<>
  12. struct traits<MatrixReplacement> : public Eigen::internal::traits<Eigen::SparseMatrix<double> >
  13. {};
  14. }
  15. }
  16. // Example of a matrix-free wrapper from a user type to Eigen's compatible type
  17. // For the sake of simplicity, this example simply wrap a Eigen::SparseMatrix.
  18. class MatrixReplacement : public Eigen::EigenBase<MatrixReplacement> {
  19. public:
  20. // Required typedefs, constants, and method:
  21. typedef double Scalar;
  22. typedef double RealScalar;
  23. typedef int StorageIndex;
  24. enum {
  25. ColsAtCompileTime = Eigen::Dynamic,
  26. MaxColsAtCompileTime = Eigen::Dynamic,
  27. IsRowMajor = false
  28. };
  29. Index rows() const { return mp_mat->rows(); }
  30. Index cols() const { return mp_mat->cols(); }
  31. template<typename Rhs>
  32. Eigen::Product<MatrixReplacement,Rhs,Eigen::AliasFreeProduct> operator*(const Eigen::MatrixBase<Rhs>& x) const {
  33. return Eigen::Product<MatrixReplacement,Rhs,Eigen::AliasFreeProduct>(*this, x.derived());
  34. }
  35. // Custom API:
  36. MatrixReplacement() : mp_mat(0) {}
  37. void attachMyMatrix(const SparseMatrix<double> &mat) {
  38. mp_mat = &mat;
  39. }
  40. const SparseMatrix<double> my_matrix() const { return *mp_mat; }
  41. private:
  42. const SparseMatrix<double> *mp_mat;
  43. };
  44. // Implementation of MatrixReplacement * Eigen::DenseVector though a specialization of internal::generic_product_impl:
  45. namespace Eigen {
  46. namespace internal {
  47. template<typename Rhs>
  48. struct generic_product_impl<MatrixReplacement, Rhs, SparseShape, DenseShape, GemvProduct> // GEMV stands for matrix-vector
  49. : generic_product_impl_base<MatrixReplacement,Rhs,generic_product_impl<MatrixReplacement,Rhs> >
  50. {
  51. typedef typename Product<MatrixReplacement,Rhs>::Scalar Scalar;
  52. template<typename Dest>
  53. static void scaleAndAddTo(Dest& dst, const MatrixReplacement& lhs, const Rhs& rhs, const Scalar& alpha)
  54. {
  55. // This method should implement "dst += alpha * lhs * rhs" inplace,
  56. // however, for iterative solvers, alpha is always equal to 1, so let's not bother about it.
  57. assert(alpha==Scalar(1) && "scaling is not implemented");
  58. EIGEN_ONLY_USED_FOR_DEBUG(alpha);
  59. // Here we could simply call dst.noalias() += lhs.my_matrix() * rhs,
  60. // but let's do something fancier (and less efficient):
  61. for(Index i=0; i<lhs.cols(); ++i)
  62. dst += rhs(i) * lhs.my_matrix().col(i);
  63. }
  64. };
  65. }
  66. }
  67. int main()
  68. {
  69. int n = 10;
  70. Eigen::SparseMatrix<double> S = Eigen::MatrixXd::Random(n,n).sparseView(0.5,1);
  71. S = S.transpose()*S;
  72. MatrixReplacement A;
  73. A.attachMyMatrix(S);
  74. Eigen::VectorXd b(n), x;
  75. b.setRandom();
  76. // Solve Ax = b using various iterative solver with matrix-free version:
  77. {
  78. Eigen::ConjugateGradient<MatrixReplacement, Eigen::Lower|Eigen::Upper, Eigen::IdentityPreconditioner> cg;
  79. cg.compute(A);
  80. x = cg.solve(b);
  81. std::cout << "CG: #iterations: " << cg.iterations() << ", estimated error: " << cg.error() << std::endl;
  82. }
  83. {
  84. Eigen::BiCGSTAB<MatrixReplacement, Eigen::IdentityPreconditioner> bicg;
  85. bicg.compute(A);
  86. x = bicg.solve(b);
  87. std::cout << "BiCGSTAB: #iterations: " << bicg.iterations() << ", estimated error: " << bicg.error() << std::endl;
  88. }
  89. {
  90. Eigen::GMRES<MatrixReplacement, Eigen::IdentityPreconditioner> gmres;
  91. gmres.compute(A);
  92. x = gmres.solve(b);
  93. std::cout << "GMRES: #iterations: " << gmres.iterations() << ", estimated error: " << gmres.error() << std::endl;
  94. }
  95. {
  96. Eigen::DGMRES<MatrixReplacement, Eigen::IdentityPreconditioner> gmres;
  97. gmres.compute(A);
  98. x = gmres.solve(b);
  99. std::cout << "DGMRES: #iterations: " << gmres.iterations() << ", estimated error: " << gmres.error() << std::endl;
  100. }
  101. {
  102. Eigen::MINRES<MatrixReplacement, Eigen::Lower|Eigen::Upper, Eigen::IdentityPreconditioner> minres;
  103. minres.compute(A);
  104. x = minres.solve(b);
  105. std::cout << "MINRES: #iterations: " << minres.iterations() << ", estimated error: " << minres.error() << std::endl;
  106. }
  107. }