|
template<typename Rhs , typename Dest > |
void | _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const |
|
Index | get_restart () |
|
| GMRES () |
|
template<typename MatrixDerived > |
| GMRES (const EigenBase< MatrixDerived > &A) |
|
void | set_restart (const Index restart) |
|
| ~GMRES () |
|
void | _solve_impl (const Rhs &b, Dest &x) const |
|
void | _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const |
|
internal::enable_if< Rhs::ColsAtCompileTime!=1 &&DestDerived::ColsAtCompileTime!=1 >::type | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const |
|
internal::enable_if< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 >::type | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const |
|
GMRES< _MatrixType, _Preconditioner > & | analyzePattern (const EigenBase< MatrixDerived > &A) |
|
EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
|
GMRES< _MatrixType, _Preconditioner > & | compute (const EigenBase< MatrixDerived > &A) |
|
RealScalar | error () const |
|
GMRES< _MatrixType, _Preconditioner > & | factorize (const EigenBase< MatrixDerived > &A) |
|
ComputationInfo | info () const |
|
Index | iterations () const |
|
| IterativeSolverBase () |
|
| IterativeSolverBase (const EigenBase< MatrixDerived > &A) |
|
Index | maxIterations () const |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () const |
|
EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
|
GMRES< _MatrixType, _Preconditioner > & | setMaxIterations (Index maxIters) |
|
GMRES< _MatrixType, _Preconditioner > & | setTolerance (const RealScalar &tolerance) |
|
const SolveWithGuess< GMRES< _MatrixType, _Preconditioner >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
|
RealScalar | tolerance () const |
|
| ~IterativeSolverBase () |
|
void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
|
GMRES< _MatrixType, _Preconditioner > & | derived () |
|
const GMRES< _MatrixType, _Preconditioner > & | derived () const |
|
const Solve< GMRES< _MatrixType, _Preconditioner >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
|
const Solve< GMRES< _MatrixType, _Preconditioner >, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
|
| SparseSolverBase () |
|
| ~SparseSolverBase () |
|
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;
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
Definition at line 221 of file GMRES.h.
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.
Definition at line 304 of file GMRES.h.