|
| GMRES () |
|
template<typename MatrixDerived > |
| GMRES (const EigenBase< MatrixDerived > &A) |
|
| ~GMRES () |
|
Index | get_restart () |
|
void | set_restart (const Index restart) |
|
template<typename Rhs , typename Dest > |
void | _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const |
|
template<typename Rhs , typename Dest > |
void | _solve_impl (const Rhs &b, Dest &x) const |
|
| IterativeSolverBase () |
|
| IterativeSolverBase (const EigenBase< MatrixDerived > &A) |
|
| IterativeSolverBase (IterativeSolverBase &&)=default |
|
| ~IterativeSolverBase () |
|
GMRES< MatrixType_, Preconditioner_ > & | analyzePattern (const EigenBase< MatrixDerived > &A) |
|
GMRES< MatrixType_, Preconditioner_ > & | factorize (const EigenBase< MatrixDerived > &A) |
|
GMRES< MatrixType_, Preconditioner_ > & | compute (const EigenBase< MatrixDerived > &A) |
|
EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
|
EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
|
RealScalar | tolerance () const |
|
GMRES< MatrixType_, Preconditioner_ > & | setTolerance (const RealScalar &tolerance) |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () const |
|
Index | maxIterations () const |
|
GMRES< MatrixType_, Preconditioner_ > & | setMaxIterations (Index maxIters) |
|
Index | iterations () const |
|
RealScalar | error () const |
|
const SolveWithGuess< GMRES< MatrixType_, Preconditioner_ >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
|
ComputationInfo | info () const |
|
void | _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const |
|
std::enable_if_t< Rhs::ColsAtCompileTime !=1 &&DestDerived::ColsAtCompileTime !=1 > | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const |
|
std::enable_if_t< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 > | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const |
|
void | _solve_impl (const Rhs &b, Dest &x) const |
|
GMRES< MatrixType_, Preconditioner_ > & | derived () |
|
const GMRES< MatrixType_, Preconditioner_ > & | derived () const |
|
| SparseSolverBase () |
|
| SparseSolverBase (SparseSolverBase &&other) |
|
| ~SparseSolverBase () |
|
Derived & | derived () |
|
const Derived & | derived () const |
|
template<typename Rhs > |
const Solve< Derived, Rhs > | solve (const MatrixBase< Rhs > &b) const |
|
template<typename Rhs > |
const Solve< Derived, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
|
template<typename Rhs , typename Dest > |
void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
|
template<typename MatrixType_, typename Preconditioner_>
class Eigen::GMRES< MatrixType_, Preconditioner_ >
A GMRES solver for sparse square problems.
This class allows to solve for A.x = b sparse linear problems using a generalized minimal residual method. The vectors x and b can be either dense or sparse.
- Template Parameters
-
MatrixType_ | the type of the sparse matrix A, can be a dense or a sparse matrix. |
Preconditioner_ | the type of the preconditioner. Default is DiagonalPreconditioner |
The maximal number of iterations and tolerance value can be controlled via the setMaxIterations() and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations and NumTraits<Scalar>::epsilon() for the tolerance.
This class can be used as the direct solver classes. Here is a typical usage example:
SparseMatrix<double>
A(
n,
n);
GMRES<SparseMatrix<double> >
solver(
A);
std::cout <<
"#iterations: " <<
solver.iterations() << std::endl;
std::cout <<
"estimated error: " <<
solver.error() << std::endl;
BiCGSTAB< SparseMatrix< double > > solver
Definition: BiCGSTAB_simple.cpp:5
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
list x
Definition: plotDoE.py:28
By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method.
GMRES can also be used in a matrix-free context, see the following example .
- See also
- class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
template<typename MatrixType_ , typename Preconditioner_ >
template<typename MatrixDerived >
Initialize the solver with matrix A for further Ax=b
solving.
This constructor is a shortcut for the default constructor followed by a call to compute().
- Warning
- this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if A is changed this class becomes invalid. Call compute() to update it with the new matrix A, or modify a copy of A.
template<typename MatrixType_ , typename Preconditioner_ >
template<typename Rhs , typename Dest >
void Eigen::GMRES< MatrixType_, Preconditioner_ >::_solve_vector_with_guess_impl |
( |
const Rhs & |
b, |
|
|
Dest & |
x |
|
) |
| const |
|
inline |
ComputationInfo m_info
Definition: IterativeSolverBase.h:389
RealScalar m_error
Definition: IterativeSolverBase.h:387
Index m_iterations
Definition: IterativeSolverBase.h:388
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:374
Index maxIterations() const
Definition: IterativeSolverBase.h:251
Preconditioner m_preconditioner
Definition: IterativeSolverBase.h:382
RealScalar m_tolerance
Definition: IterativeSolverBase.h:385
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
Eigen::DenseIndex ret
Definition: level1_cplx_impl.h:43
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
References b, Eigen::internal::gmres(), Eigen::GMRES< MatrixType_, Preconditioner_ >::m_error, Eigen::GMRES< MatrixType_, Preconditioner_ >::m_info, Eigen::GMRES< MatrixType_, Preconditioner_ >::m_iterations, Eigen::IterativeSolverBase< Derived >::m_preconditioner, Eigen::GMRES< MatrixType_, Preconditioner_ >::m_restart, Eigen::IterativeSolverBase< Derived >::m_tolerance, Eigen::GMRES< MatrixType_, Preconditioner_ >::matrix(), Eigen::IterativeSolverBase< Derived >::maxIterations(), Eigen::NoConvergence, Eigen::NumericalIssue, ret, Eigen::Success, and plotDoE::x.