scipy.sparse.linalg.spsolve#

scipy.sparse.linalg.spsolve(A, b, permc_spec=None, use_umfpack=True)[source]#

Solve the sparse linear system Ax=b, where b may be a vector or a matrix.

Parameters:
Andarray or sparse matrix

The square matrix A will be converted into CSC or CSR form

bndarray or sparse matrix

The matrix or vector representing the right hand side of the equation. If a vector, b.shape must be (n,) or (n, 1).

permc_specstr, optional

How to permute the columns of the matrix for sparsity preservation. (default: ‘COLAMD’)

  • NATURAL: natural ordering.

  • MMD_ATA: minimum degree ordering on the structure of A^T A.

  • MMD_AT_PLUS_A: minimum degree ordering on the structure of A^T+A.

  • COLAMD: approximate minimum degree column ordering [1], [2].

use_umfpackbool, optional

if True (default) then use UMFPACK for the solution [3], [4], [5], [6] . This is only referenced if b is a vector and scikits.umfpack is installed.

Returns:
xndarray or sparse matrix

the solution of the sparse linear equation. If b is a vector, then x is a vector of size A.shape[1] If b is a matrix, then x is a matrix of size (A.shape[1], b.shape[1])

Notes

For solving the matrix expression AX = B, this solver assumes the resulting matrix X is sparse, as is often the case for very sparse inputs. If the resulting X is dense, the construction of this sparse result will be relatively expensive. In that case, consider converting A to a dense matrix and using scipy.linalg.solve or its variants.

References

[1]

T. A. Davis, J. R. Gilbert, S. Larimore, E. Ng, Algorithm 836: COLAMD, an approximate column minimum degree ordering algorithm, ACM Trans. on Mathematical Software, 30(3), 2004, pp. 377–380. DOI:10.1145/1024074.1024080

[2]

T. A. Davis, J. R. Gilbert, S. Larimore, E. Ng, A column approximate minimum degree ordering algorithm, ACM Trans. on Mathematical Software, 30(3), 2004, pp. 353–376. DOI:10.1145/1024074.1024079

[3]

T. A. Davis, Algorithm 832: UMFPACK - an unsymmetric-pattern multifrontal method with a column pre-ordering strategy, ACM Trans. on Mathematical Software, 30(2), 2004, pp. 196–199. https://dl.acm.org/doi/abs/10.1145/992200.992206

[4]

T. A. Davis, A column pre-ordering strategy for the unsymmetric-pattern multifrontal method, ACM Trans. on Mathematical Software, 30(2), 2004, pp. 165–195. https://dl.acm.org/doi/abs/10.1145/992200.992205

[5]

T. A. Davis and I. S. Duff, A combined unifrontal/multifrontal method for unsymmetric sparse matrices, ACM Trans. on Mathematical Software, 25(1), 1999, pp. 1–19. https://doi.org/10.1145/305658.287640

[6]

T. A. Davis and I. S. Duff, An unsymmetric-pattern multifrontal method for sparse LU factorization, SIAM J. Matrix Analysis and Computations, 18(1), 1997, pp. 140–158. https://doi.org/10.1137/S0895479894246905T.

Examples

>>> import numpy as np
>>> from scipy.sparse import csc_matrix
>>> from scipy.sparse.linalg import spsolve
>>> A = csc_matrix([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float)
>>> B = csc_matrix([[2, 0], [-1, 0], [2, 0]], dtype=float)
>>> x = spsolve(A, B)
>>> np.allclose(A.dot(x).toarray(), B.toarray())
True