lu_factor#
- scipy.linalg.lu_factor(a, overwrite_a=False, check_finite=True)[source]#
- Compute pivoted LU decomposition of a matrix. - The decomposition is: - A = P L U - where P is a permutation matrix, L lower triangular with unit diagonal elements, and U upper triangular. - 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, N) array_like
- Matrix to decompose 
- overwrite_abool, optional
- Whether to overwrite data in A (may increase performance) 
- check_finitebool, optional
- Whether to check that the 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:
- lu(M, N) ndarray
- Matrix containing U in its upper triangle, and L in its lower triangle. The unit diagonal elements of L are not stored. 
- piv(K,) ndarray
- Pivot indices representing the permutation matrix P: row i of matrix was interchanged with row piv[i]. Of shape - (K,), with- K = min(M, N).
 
 - See also - Notes - This is a wrapper to the - *GETRFroutines from LAPACK. Unlike- lu, it outputs the L and U factors into a single array and returns pivot indices instead of a permutation matrix.- While the underlying - *GETRFroutines return 1-based pivot indices, the- pivarray returned by- lu_factorcontains 0-based indices.- Examples - >>> import numpy as np >>> from scipy.linalg import lu_factor >>> A = np.array([[2, 5, 8, 7], [5, 2, 2, 8], [7, 5, 6, 6], [5, 4, 4, 8]]) >>> lu, piv = lu_factor(A) >>> piv array([2, 2, 3, 3], dtype=int32) - Convert LAPACK’s - pivarray to NumPy index and test the permutation- >>> def pivot_to_permutation(piv): ... perm = np.arange(len(piv)) ... for i in range(len(piv)): ... perm[i], perm[piv[i]] = perm[piv[i]], perm[i] ... return perm ... >>> p_inv = pivot_to_permutation(piv) >>> p_inv array([2, 0, 3, 1]) >>> L, U = np.tril(lu, k=-1) + np.eye(4), np.triu(lu) >>> np.allclose(A[p_inv] - L @ U, np.zeros((4, 4))) True - The P matrix in P L U is defined by the inverse permutation and can be recovered using argsort: - >>> p = np.argsort(p_inv) >>> p array([1, 3, 0, 2]) >>> np.allclose(A - L[p] @ U, np.zeros((4, 4))) True - or alternatively: - >>> P = np.eye(4)[p] >>> np.allclose(A - P @ L @ U, np.zeros((4, 4))) True