scipy.sparse.linalg.gmres#
- scipy.sparse.linalg.gmres(A, b, x0=None, *, tol=<object object>, restart=None, maxiter=None, M=None, callback=None, restrt=<object object>, atol=0.0, callback_type=None, rtol=1e-05)[source]#
Use Generalized Minimal RESidual iteration to solve
Ax = b
.- Parameters:
- A{sparse matrix, ndarray, LinearOperator}
The real or complex N-by-N matrix of the linear system. Alternatively,
A
can be a linear operator which can produceAx
using, e.g.,scipy.sparse.linalg.LinearOperator
.- bndarray
Right hand side of the linear system. Has shape (N,) or (N,1).
- x0ndarray
Starting guess for the solution (a vector of zeros by default).
- atol, rtolfloat
Parameters for the convergence test. For convergence,
norm(b - A @ x) <= max(rtol*norm(b), atol)
should be satisfied. The default isatol=0.
andrtol=1e-5
.- restartint, optional
Number of iterations between restarts. Larger values increase iteration cost, but may be necessary for convergence. If omitted,
min(20, n)
is used.- maxiterint, optional
Maximum number of iterations (restart cycles). Iteration will stop after maxiter steps even if the specified tolerance has not been achieved. See callback_type.
- M{sparse matrix, ndarray, LinearOperator}
Inverse of the preconditioner of A. M should approximate the inverse of A and be easy to solve for (see Notes). Effective preconditioning dramatically improves the rate of convergence, which implies that fewer iterations are needed to reach a given error tolerance. By default, no preconditioner is used. In this implementation, left preconditioning is used, and the preconditioned residual is minimized. However, the final convergence is tested with respect to the
b - A @ x
residual.- callbackfunction
User-supplied function to call after each iteration. It is called as callback(args), where args are selected by callback_type.
- callback_type{‘x’, ‘pr_norm’, ‘legacy’}, optional
- Callback function argument requested:
x
: current iterate (ndarray), called on every restartpr_norm
: relative (preconditioned) residual norm (float), called on every inner iterationlegacy
(default): same aspr_norm
, but also changes the meaning of maxiter to count inner iterations instead of restart cycles.
This keyword has no effect if callback is not set.
- restrtint, optional, deprecated
Deprecated since version 0.11.0:
gmres
keyword argumentrestrt
is deprecated in favor ofrestart
and will be removed in SciPy 1.14.0.- tolfloat, optional, deprecated
Deprecated since version 1.12.0:
gmres
keyword argumenttol
is deprecated in favor ofrtol
and will be removed in SciPy 1.14.0
- Returns:
- xndarray
The converged solution.
- infoint
- Provides convergence information:
0 : successful exit >0 : convergence to tolerance not achieved, number of iterations
See also
Notes
A preconditioner, P, is chosen such that P is close to A but easy to solve for. The preconditioner parameter required by this routine is
M = P^-1
. The inverse should preferably not be calculated explicitly. Rather, use the following template to produce M:# Construct a linear operator that computes P^-1 @ x. import scipy.sparse.linalg as spla M_x = lambda x: spla.spsolve(P, x) M = spla.LinearOperator((n, n), M_x)
Examples
>>> import numpy as np >>> from scipy.sparse import csc_matrix >>> from scipy.sparse.linalg import gmres >>> A = csc_matrix([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float) >>> b = np.array([2, 4, -1], dtype=float) >>> x, exitCode = gmres(A, b, atol=1e-5) >>> print(exitCode) # 0 indicates successful convergence 0 >>> np.allclose(A.dot(x), b) True