scipy.sparse.linalg.cg#

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

Use Conjugate Gradient iteration to solve Ax = b.

Parameters:
A{sparse matrix, ndarray, LinearOperator}

The real or complex N-by-N matrix of the linear system. A must represent a hermitian, positive definite matrix. 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).

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 atol=0. and rtol=1e-5.

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.

tolfloat, optional, deprecated

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

Returns:
xndarray

The converged solution.

infointeger
Provides convergence information:

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

Examples

>>> import numpy as np
>>> from scipy.sparse import csc_matrix
>>> from scipy.sparse.linalg import cg
>>> P = np.array([[4, 0, 1, 0],
...               [0, 5, 0, 0],
...               [1, 0, 3, 2],
...               [0, 0, 2, 4]])
>>> A = csc_matrix(P)
>>> b = np.array([-1, -0.5, -1, 2])
>>> x, exit_code = cg(A, b, atol=1e-5)
>>> print(exit_code)    # 0 indicates successful convergence
0
>>> np.allclose(A.dot(x), b)
True