scipy.sparse.linalg.

gcrotmk#

scipy.sparse.linalg.gcrotmk(A, b, x0=None, *, rtol=1e-05, atol=0.0, maxiter=1000, M=None, callback=None, m=20, k=None, CU=None, discard_C=False, truncate='oldest')[source]#

Solve a matrix equation using flexible GCROT(m,k) algorithm.

Parameters:
A{sparse array, 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., 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 rtol=1e-5 and atol=0.0.

maxiterint, optional

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

M{sparse array, ndarray, LinearOperator}, optional

Preconditioner for A. The preconditioner should approximate the inverse of A. gcrotmk is a ‘flexible’ algorithm and the preconditioner can vary from iteration to iteration. Effective preconditioning dramatically improves the rate of convergence, which implies that fewer iterations are needed to reach a given error tolerance.

callbackfunction, optional

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

mint, optional

Number of inner FGMRES iterations per each outer iteration. Default: 20

kint, optional

Number of vectors to carry between inner FGMRES iterations. According to [2], good values are around m. Default: m

CUlist of tuples, optional

List of tuples (c, u) which contain the columns of the matrices C and U in the GCROT(m,k) algorithm. For details, see [2]. The list given and vectors contained in it are modified in-place. If not given, start from empty matrices. The c elements in the tuples can be None, in which case the vectors are recomputed via c = A u on start and orthogonalized as described in [3].

discard_Cbool, optional

Discard the C-vectors at the end. Useful if recycling Krylov subspaces for different linear systems.

truncate{‘oldest’, ‘smallest’}, optional

Truncation scheme to use. Drop: oldest vectors, or vectors with smallest singular values using the scheme discussed in [1,2]. See [2] for detailed comparison. Default: ‘oldest’

Returns:
xndarray

The solution found.

infoint

Provides convergence information:

  • 0 : successful exit

  • >0 : convergence to tolerance not achieved, number of iterations

References

[1]

E. de Sturler, ‘’Truncation strategies for optimal Krylov subspace methods’’, SIAM J. Numer. Anal. 36, 864 (1999).

[2] (1,2,3)

J.E. Hicken and D.W. Zingg, ‘’A simplified and flexible variant of GCROT for solving nonsymmetric linear systems’’, SIAM J. Sci. Comput. 32, 172 (2010).

[3]

M.L. Parks, E. de Sturler, G. Mackey, D.D. Johnson, S. Maiti, ‘’Recycling Krylov subspaces for sequences of linear systems’’, SIAM J. Sci. Comput. 28, 1651 (2006).

Examples

>>> import numpy as np
>>> from scipy.sparse import csc_array
>>> from scipy.sparse.linalg import gcrotmk
>>> R = np.random.randn(5, 5)
>>> A = csc_array(R)
>>> b = np.random.randn(5)
>>> x, exit_code = gcrotmk(A, b, atol=1e-5)
>>> print(exit_code)
0
>>> np.allclose(A.dot(x), b)
True