scipy.sparse.linalg.gmres¶
- 
scipy.sparse.linalg.gmres(A, b, x0=None, tol=1e-05, restart=None, maxiter=None, M=None, callback=None, restrt=None, atol=None, callback_type=None)¶
- Use Generalized Minimal RESidual iteration to solve - Ax = b.- Parameters
- A{sparse matrix, dense matrix, LinearOperator}
- The real or complex N-by-N matrix of the linear system. Alternatively, - Acan be a linear operator which can produce- Axusing, e.g.,- scipy.sparse.linalg.LinearOperator.
- b{array, matrix}
- Right hand side of the linear system. Has shape (N,) or (N,1). 
 
- Returns
- x{array, matrix}
- The converged solution. 
- infoint
- Provides convergence information:
- 0 : successful exit 
- >0 : convergence to tolerance not achieved, number of iterations 
- <0 : illegal input or breakdown 
 
 
 
- Other Parameters
- x0{array, matrix}
- Starting guess for the solution (a vector of zeros by default). 
- tol, atolfloat, optional
- Tolerances for convergence, - norm(residual) <= max(tol*norm(b), atol). The default for- atolis- 'legacy', which emulates a different legacy behavior.- Warning - The default value for atol will be changed in a future release. For future compatibility, specify atol explicitly. 
- restartint, optional
- Number of iterations between restarts. Larger values increase iteration cost, but may be necessary for convergence. Default is 20. 
- maxiterint, optional
- Maximum number of iterations (restart cycles). Iteration will stop after maxiter steps even if the specified tolerance has not been achieved. 
- M{sparse matrix, dense matrix, 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. 
- 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 restart
- pr_norm: relative (preconditioned) residual norm (float), called on every inner iteration
- legacy(default): same as- pr_norm, but also changes the meaning of ‘maxiter’ to count inner iterations instead of restart cycles.
 
 
- restrtint, optional
- DEPRECATED - use restart instead. 
 
 - 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 - >>> 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) >>> print(exitCode) # 0 indicates successful convergence 0 >>> np.allclose(A.dot(x), b) True 
