scipy.sparse.linalg.

qmr#

scipy.sparse.linalg.qmr(A, b, x0=None, *, rtol=1e-05, atol=0.0, maxiter=None, M1=None, M2=None, callback=None)[source]#

Solve Ax = b with the Quasi-Minimal Residual method.

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 produce Ax and A^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.

rtol, atolfloat, optional

Parameters for the convergence test. For convergence, norm(b - A @ x) <= max(rtol*norm(b), atol) should be satisfied. The default is atol=0. and rtol=1e-5.

maxiterint

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.

infoint
Provides convergence information:

0 : successful exit >0 : convergence to tolerance not achieved, number of iterations <0 : parameter breakdown

See also

LinearOperator

Notes

Array API Standard Support

qmr has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.

Library

CPU

GPU

NumPy

n/a

CuPy

n/a

PyTorch

JAX

Dask

n/a

See Support for the array API standard for more information.

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