scipy.sparse.linalg.minres#

scipy.sparse.linalg.minres(A, b, x0=None, *, shift=0.0, tol=<object object>, maxiter=None, M=None, callback=None, show=False, check=False, rtol=1e-05)[source]#

Use MINimum RESidual iteration to solve Ax=b

MINRES minimizes norm(Ax - b) for a real symmetric matrix A. Unlike the Conjugate Gradient method, A can be indefinite or singular.

If shift != 0 then the method solves (A - shift*I)x = b

Parameters:
A{sparse matrix, ndarray, LinearOperator}

The real symmetric N-by-N matrix of the linear system Alternatively, A can be a linear operator which can produce Ax using, e.g., scipy.sparse.linalg.LinearOperator.

bndarray

Right hand side of the linear system. Has shape (N,) or (N,1).

Returns:
xndarray

The converged solution.

infointeger
Provides convergence information:

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

Other Parameters:
x0ndarray

Starting guess for the solution.

shiftfloat

Value to apply to the system (A - shift * I)x = b. Default is 0.

rtolfloat

Tolerance to achieve. The algorithm terminates when the relative residual is below rtol.

maxiterinteger

Maximum number of iterations. Iteration will stop after maxiter steps even if the specified tolerance has not been achieved.

M{sparse matrix, ndarray, LinearOperator}

Preconditioner for A. The preconditioner should approximate the inverse of A. Effective preconditioning dramatically improves the rate of convergence, which implies that fewer iterations are needed to reach a given error tolerance.

callbackfunction

User-supplied function to call after each iteration. It is called as callback(xk), where xk is the current solution vector.

showbool

If True, print out a summary and metrics related to the solution during iterations. Default is False.

checkbool

If True, run additional input validation to check that A and M (if specified) are symmetric. Default is False.

tolfloat, optional, deprecated

Deprecated since version 1.12.0: minres keyword argument tol is deprecated in favor of rtol and will be removed in SciPy 1.14.0.

References

Solution of sparse indefinite systems of linear equations,

C. C. Paige and M. A. Saunders (1975), SIAM J. Numer. Anal. 12(4), pp. 617-629. https://web.stanford.edu/group/SOL/software/minres/

This file is a translation of the following MATLAB implementation:

https://web.stanford.edu/group/SOL/software/minres/minres-matlab.zip

Examples

>>> import numpy as np
>>> from scipy.sparse import csc_matrix
>>> from scipy.sparse.linalg import minres
>>> A = csc_matrix([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float)
>>> A = A + A.T
>>> b = np.array([2, 4, -1], dtype=float)
>>> x, exitCode = minres(A, b)
>>> print(exitCode)            # 0 indicates successful convergence
0
>>> np.allclose(A.dot(x), b)
True