legate_sparse.linalg.gmres#

legate_sparse.linalg.gmres(
A,
b,
x0=None,
tol=None,
restart=None,
maxiter=None,
M=None,
callback=None,
restrt=None,
atol=0.0,
callback_type=None,
rtol=1e-05,
)#

Solve a linear system using the Generalized Minimal Residual method.

Parameters:
  • A (sparse matrix or LinearOperator) – The coefficient matrix of the linear system.

  • b (cupynumeric.ndarray) – Right-hand side of the linear system with shape (n,) or (n, 1).

  • x0 (cupynumeric.ndarray, optional) – Starting guess for the solution. If None, uses zero vector.

  • tol (float, optional) – Legacy tolerance parameter. If provided, overrides rtol.

  • restart (int, optional) – Number of iterations between restarts. Larger values increase iteration cost but may be necessary for convergence. Default is 20.

  • maxiter (int, optional) – Maximum number of iterations. If None, uses 10 * n.

  • M (sparse matrix or LinearOperator, optional) – Preconditioner for A. The preconditioner should approximate the inverse of A. If None, uses identity.

  • callback (callable, optional) – User-specified function called on every restart.

  • restrt (int, optional) – Deprecated alias for restart parameter.

  • atol (float, optional) – Absolute tolerance for convergence. Default is 0.0.

  • callback_type (str, optional) – Type of callback argument: ‘x’ for current solution vector, ‘pr_norm’ for relative preconditioned residual norm. Default is ‘pr_norm’.

  • rtol (float, optional) – Relative tolerance for convergence. Default is 1e-5.

Returns:

(x, info) where x is the converged solution and info provides convergence information.

Return type:

tuple

Raises:
  • AssertionError – If b is not 1D or A is not square.

  • ValueError – If callback_type is not ‘x’ or ‘pr_norm’.

Notes

This implementation is adapted from CuPy’s GMRES solver. The method uses Arnoldi iteration to build a Krylov subspace and solves the least squares problem in that subspace.

For convergence, the residual must satisfy: norm(b - A @ x) <= max(rtol * norm(b), atol)

The restart parameter controls the trade-off between memory usage and convergence rate. Larger restart values may improve convergence but require more memory.

References

M. Wang, H. Klie, M. Parashar and H. Sudan, “Solving Sparse Linear Systems on NVIDIA Tesla GPUs”, ICCS 2009 (2009).

Examples

>>> import cupynumeric as np
>>> from legate_sparse import csr_array, linalg
>>> A = csr_array([[4, 1, 0], [1, 3, 1], [0, 1, 2]])
>>> b = np.array([1, 2, 3])
>>> x, info = linalg.gmres(A, b, restart=10)
>>> print(f"Solution: {x}, Info: {info}")