matrixfree_cg.cpp File Reference
#include <iostream>
#include <Eigen/Core>
#include <Eigen/Dense>
#include <Eigen/IterativeLinearSolvers>
#include <unsupported/Eigen/IterativeSolvers>

Classes

struct  Eigen::internal::traits< MatrixReplacement >
 
class  MatrixReplacement
 
struct  Eigen::internal::generic_product_impl< MatrixReplacement, Rhs, SparseShape, DenseShape, GemvProduct >
 

Namespaces

 Eigen
 Namespace containing all symbols from the Eigen library.
 
 Eigen::internal
 Namespace containing low-level routines from the Eigen library.
 

Functions

int main ()
 

Function Documentation

◆ main()

int main ( )
72  {
73  int n = 10;
74  Eigen::SparseMatrix<double> S = Eigen::MatrixXd::Random(n, n).sparseView(0.5, 1);
75  S = S.transpose() * S;
76 
78  A.attachMyMatrix(S);
79 
80  Eigen::VectorXd b(n), x;
81  b.setRandom();
82 
83  // Solve Ax = b using various iterative solver with matrix-free version:
84  {
86  cg.compute(A);
87  x = cg.solve(b);
88  std::cout << "CG: #iterations: " << cg.iterations() << ", estimated error: " << cg.error() << std::endl;
89  }
90 
91  {
93  bicg.compute(A);
94  x = bicg.solve(b);
95  std::cout << "BiCGSTAB: #iterations: " << bicg.iterations() << ", estimated error: " << bicg.error() << std::endl;
96  }
97 
98  {
100  gmres.compute(A);
101  x = gmres.solve(b);
102  std::cout << "GMRES: #iterations: " << gmres.iterations() << ", estimated error: " << gmres.error() << std::endl;
103  }
104 
105  {
107  gmres.compute(A);
108  x = gmres.solve(b);
109  std::cout << "DGMRES: #iterations: " << gmres.iterations() << ", estimated error: " << gmres.error() << std::endl;
110  }
111 
112  {
114  minres.compute(A);
115  x = minres.solve(b);
116  std::cout << "MINRES: #iterations: " << minres.iterations() << ", estimated error: " << minres.error()
117  << std::endl;
118  }
119 }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Scalar * b
Definition: benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
A bi conjugate gradient stabilized solver for sparse square problems.
Definition: BiCGSTAB.h:153
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition: ConjugateGradient.h:152
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
Definition: DGMRES.h:98
A GMRES solver for sparse square problems.
Definition: GMRES.h:255
RealScalar error() const
Definition: IterativeSolverBase.h:270
Derived & compute(const EigenBase< MatrixDerived > &A)
Definition: IterativeSolverBase.h:210
Index iterations() const
Definition: IterativeSolverBase.h:262
A minimal residual solver for sparse symmetric problems.
Definition: MINRES.h:188
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SparseSolverBase.h:84
Definition: matrixfree_cg.cpp:20
bool gmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, const Index &restart, typename Dest::RealScalar &tol_error)
Definition: GMRES.h:59
EIGEN_DONT_INLINE void minres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
Definition: MINRES.h:32
@ S
Definition: quadtree.h:62
list x
Definition: plotDoE.py:28

References b, Eigen::IterativeSolverBase< Derived >::compute(), Eigen::IterativeSolverBase< Derived >::error(), Eigen::internal::gmres(), Eigen::IterativeSolverBase< Derived >::iterations(), Eigen::internal::minres(), n, oomph::QuadTreeNames::S, Eigen::SparseSolverBase< Derived >::solve(), and plotDoE::x.