spsolve_triangular — SciPy v1.15.3 Manual (original) (raw)

scipy.sparse.linalg.

scipy.sparse.linalg.spsolve_triangular(A, b, lower=True, overwrite_A=False, overwrite_b=False, unit_diagonal=False)[source]#

Solve the equation A x = b for x, assuming A is a triangular matrix.

Parameters:

A(M, M) sparse array or matrix

A sparse square triangular matrix. Should be in CSR or CSC format.

b(M,) or (M, N) array_like

Right-hand side matrix in A x = b

lowerbool, optional

Whether A is a lower or upper triangular matrix. Default is lower triangular matrix.

overwrite_Abool, optional

Allow changing A. Enabling gives a performance gain. Default is False.

overwrite_bbool, optional

Allow overwriting data in b. Enabling gives a performance gain. Default is False. If overwrite_b is True, it should be ensured that_b_ has an appropriate dtype to be able to store the result.

unit_diagonalbool, optional

If True, diagonal elements of a are assumed to be 1.

Added in version 1.4.0.

Returns:

x(M,) or (M, N) ndarray

Solution to the system A x = b. Shape of return matches shape of b.

Raises:

LinAlgError

If A is singular or not triangular.

ValueError

If shape of A or shape of b do not match the requirements.

Notes

Added in version 0.19.0.

Examples

import numpy as np from scipy.sparse import csc_array from scipy.sparse.linalg import spsolve_triangular A = csc_array([[3, 0, 0], [1, -1, 0], [2, 0, 1]], dtype=float) B = np.array([[2, 0], [-1, 0], [2, 0]], dtype=float) x = spsolve_triangular(A, B) np.allclose(A.dot(x), B) True