dpnp.linalg.cholesky

dpnp.linalg.cholesky(a, upper=False)[source]

Cholesky decomposition.

Return the lower or upper Cholesky decomposition, L * L.H or U.H * U, of the square matrix a, where L is lower-triangular, U is upper-triangular, and .H is the conjugate transpose operator (which is the ordinary transpose if a is real-valued). a must be Hermitian (symmetric if real-valued) and positive-definite. No checking is performed to verify whether a is Hermitian or not. In addition, only the lower or upper-triangular and diagonal elements of a are used. Only L or U is actually returned.

For full documentation refer to numpy.linalg.cholesky.

Parameters:
  • a ((..., M, M) {dpnp.ndarray, usm_ndarray}) -- Hermitian (symmetric if all elements are real), positive-definite input matrix.

  • upper ({bool}, optional) -- If True, the result must be the upper-triangular Cholesky factor. If False, the result must be the lower-triangular Cholesky factor. Default: False.

Returns:

L -- Lower or upper-triangular Cholesky factor of a.

Return type:

(..., M, M) dpnp.ndarray

Examples

>>> import dpnp as np
>>> A = np.array([[1.0, 2.0],[2.0, 5.0]])
>>> A
array([[1., 2.],
       [2., 5.]])
>>> L = np.linalg.cholesky(A)
>>> L
array([[1., 0.],
       [2., 1.]])
>>> np.dot(L, L.T.conj()) # verify that L * L.H = A
array([[1., 2.],
       [2., 5.]])

The upper-triangular Cholesky factor can also be obtained:

>>> np.linalg.cholesky(A, upper=True)
array([[ 1.+0.j, -0.-2.j],
       [ 0.+0.j,  1.+0.j]]