scipy.linalg.

cholesky#

scipy.linalg.cholesky(a, lower=False, overwrite_a=False, check_finite=True)[source]#

Compute the Cholesky decomposition of a matrix.

Returns the Cholesky decomposition, \(A = L L^*\) or \(A = U^* U\) of a Hermitian positive-definite matrix A.

The documentation is written assuming array arguments are of specified “core” shapes. However, array argument(s) of this function may have additional “batch” dimensions prepended to the core shape. In this case, the array is treated as a batch of lower-dimensional slices; see Batched Linear Operations for details.

Parameters:
a(M, M) array_like

Matrix to be decomposed

lowerbool, optional

Whether to compute the upper- or lower-triangular Cholesky factorization. During decomposition, only the selected half of the matrix is referenced. Default is upper-triangular.

overwrite_abool, optional

Whether to overwrite data in a (may improve performance).

check_finitebool, optional

Whether to check that the entire input matrix contains only finite numbers. Disabling may give a performance gain, but may result in problems (crashes, non-termination) if the inputs do contain infinities or NaNs.

Returns:
c(M, M) ndarray

Upper- or lower-triangular Cholesky factor of a.

Raises:
LinAlgErrorif decomposition fails.

Notes

During the finiteness check (if selected), the entire matrix a is checked. During decomposition, a is assumed to be symmetric or Hermitian (as applicable), and only the half selected by option lower is referenced. Consequently, if a is asymmetric/non-Hermitian, cholesky may still succeed if the symmetric/Hermitian matrix represented by the selected half is positive definite, yet it may fail if an element in the other half is non-finite.

Examples

>>> import numpy as np
>>> from scipy.linalg import cholesky
>>> a = np.array([[1,-2j],[2j,5]])
>>> L = cholesky(a, lower=True)
>>> L
array([[ 1.+0.j,  0.+0.j],
       [ 0.+2.j,  1.+0.j]])
>>> L @ L.T.conj()
array([[ 1.+0.j,  0.-2.j],
       [ 0.+2.j,  5.+0.j]])