legate_sparse.linalg.cg#

legate_sparse.linalg.cg(
A,
b,
x0=None,
tol=None,
maxiter=None,
M=None,
callback=None,
atol=0.0,
rtol=1e-05,
conv_test_iters=25,
)#

Solve a linear system using the Conjugate Gradient method.

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

  • b (cupynumeric.ndarray) – Right-hand side of the linear system.

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

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

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

  • M (sparse matrix or LinearOperator, optional) – Preconditioner for A. If None, uses identity.

  • callback (callable, optional) – User-specified function called after each iteration.

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

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

  • conv_test_iters (int, optional) – Number of iterations between convergence tests. Default is 25.

Returns:

(x, info) where x is the solution and info is zero if solution is converged else number of iterations

Return type:

tuple

Raises:

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

Notes

This implementation follows SciPy’s CG solver semantics closely. The method uses fused vector operations to avoid unnecessary memory allocations and improve performance.

Convergence is tested every conv_test_iters iterations to avoid the overhead of computing the residual norm in every iteration.

Examples

>>> import cupynumeric as np
>>> from legate_sparse import csr_array, linalg
>>> A = csr_array([[4, 1], [1, 3]])
>>> b = np.array([1, 2])
>>> x, iters = linalg.cg(A, b)
>>> print(f"Solution: {x}, Iterations: {iters}")