scipy.sparse.linalg.tfqmr#

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

Use Transpose-Free Quasi-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 produce Ax using, e.g., scipy.sparse.linalg.LinearOperator.

b{ndarray}

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

x0{ndarray}

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 rtol=1e-5, the default for atol is rtol.

Warning

The default value for atol will be changed to 0.0 in SciPy 1.14.0.

maxiterint, optional

Maximum number of iterations. Iteration will stop after maxiter steps even if the specified tolerance has not been achieved. Default is min(10000, ndofs * 10), where ndofs = A.shape[0].

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.

callbackfunction, optional

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

showbool, optional

Specify show = True to show the convergence, show = False is to close the output of the convergence. Default is False.

tolfloat, optional, deprecated

Deprecated since version 1.12.0: tfqmr keyword argument tol is deprecated in favor of rtol 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

  • <0 : illegal input or breakdown

Notes

The Transpose-Free QMR algorithm is derived from the CGS algorithm. However, unlike CGS, the convergence curves for the TFQMR method is smoothed by computing a quasi minimization of the residual norm. The implementation supports left preconditioner, and the “residual norm” to compute in convergence criterion is actually an upper bound on the actual residual norm ||b - Axk||.

References

[1]

R. W. Freund, A Transpose-Free Quasi-Minimal Residual Algorithm for Non-Hermitian Linear Systems, SIAM J. Sci. Comput., 14(2), 470-482, 1993.

[2]

Y. Saad, Iterative Methods for Sparse Linear Systems, 2nd edition, SIAM, Philadelphia, 2003.

[3]

C. T. Kelley, Iterative Methods for Linear and Nonlinear Equations, number 16 in Frontiers in Applied Mathematics, SIAM, Philadelphia, 1995.

Examples

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