cho_factor#
- scipy.linalg.cho_factor(a, lower=False, overwrite_a=False, check_finite=True)[source]#
Compute the Cholesky decomposition of a matrix, to use in cho_solve
Returns a matrix containing the Cholesky decomposition,
A = L L*
orA = U* U
of a Hermitian positive-definite matrix a. The return value can be directly used as the first parameter to cho_solve.Warning
The returned matrix also contains random data in the entries not used by the Cholesky decomposition. If you need to zero these entries, use the function
cholesky
instead.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: 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
Matrix whose upper or lower triangle contains the Cholesky factor of a. Other parts of the matrix contain random data.
- lowerbool
Flag indicating whether the factor is in the lower or upper triangle
- Raises:
- LinAlgError
Raised if decomposition fails.
See also
cho_solve
Solve a linear set equations using the Cholesky factorization of a matrix.
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 cho_factor >>> A = np.array([[9, 3, 1, 5], [3, 7, 5, 1], [1, 5, 9, 2], [5, 1, 2, 6]]) >>> c, low = cho_factor(A) >>> c array([[3. , 1. , 0.33333333, 1.66666667], [3. , 2.44948974, 1.90515869, -0.27216553], [1. , 5. , 2.29330749, 0.8559528 ], [5. , 1. , 2. , 1.55418563]]) >>> np.allclose(np.triu(c).T @ np. triu(c) - A, np.zeros((4, 4))) True