qmr#
- scipy.sparse.linalg.qmr(A, b, x0=None, *, rtol=1e-05, atol=0.0, maxiter=None, M1=None, M2=None, callback=None)[source]#
Use Quasi-Minimal Residual iteration to solve
Ax = b
.- Parameters:
- A{sparse array, ndarray, LinearOperator}
The real-valued N-by-N matrix of the linear system. Alternatively,
A
can be a linear operator which can produceAx
andA^T x
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.
- atol, rtolfloat, optional
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
.- maxiterinteger
Maximum number of iterations. Iteration will stop after maxiter steps even if the specified tolerance has not been achieved.
- M1{sparse array, ndarray, LinearOperator}
Left preconditioner for A.
- M2{sparse array, ndarray, LinearOperator}
Right preconditioner for A. Used together with the left preconditioner M1. The matrix M1@A@M2 should have better conditioned than A alone.
- callbackfunction
User-supplied function to call after each iteration. It is called as callback(xk), where xk is the current solution vector.
- Returns:
- xndarray
The converged solution.
- infointeger
- Provides convergence information:
0 : successful exit >0 : convergence to tolerance not achieved, number of iterations <0 : parameter breakdown
See also
Examples
>>> import numpy as np >>> from scipy.sparse import csc_array >>> from scipy.sparse.linalg import qmr >>> A = csc_array([[3., 2., 0.], [1., -1., 0.], [0., 5., 1.]]) >>> b = np.array([2., 4., -1.]) >>> x, exitCode = qmr(A, b, atol=1e-5) >>> print(exitCode) # 0 indicates successful convergence 0 >>> np.allclose(A.dot(x), b) True