Next: , Previous: Singular Value Decomposition, Up: Linear Algebra


14.5 Cholesky Decomposition

A symmetric, positive definite square matrix A has a Cholesky decomposition into a product of a lower triangular matrix L and its transpose L^T,

     A = L L^T

This is sometimes referred to as taking the square-root of a matrix. The Cholesky decomposition can only be carried out when all the eigenvalues of the matrix are positive. This decomposition can be used to convert the linear system A x = b into a pair of triangular systems (L y = b, L^T x = y), which can be solved by forward and back-substitution.

— Function: int gsl_linalg_cholesky_decomp (gsl_matrix * A)
— Function: int gsl_linalg_complex_cholesky_decomp (gsl_matrix_complex * A)

These functions factorize the symmetric, positive-definite square matrix A into the Cholesky decomposition A = L L^T (or A = L L^H for the complex case). On input, the values from the diagonal and lower-triangular part of the matrix A are used (the upper triangular part is ignored). On output the diagonal and lower triangular part of the input matrix A contain the matrix L, while the upper triangular part of the input matrix is overwritten with L^T (the diagonal terms being identical for both L and L^T). If the matrix is not positive-definite then the decomposition will fail, returning the error code GSL_EDOM.

When testing whether a matrix is positive-definite, disable the error handler first to avoid triggering an error.

— Function: int gsl_linalg_cholesky_solve (const gsl_matrix * cholesky, const gsl_vector * b, gsl_vector * x)
— Function: int gsl_linalg_complex_cholesky_solve (const gsl_matrix_complex * cholesky, const gsl_vector_complex * b, gsl_vector_complex * x)

These functions solve the system A x = b using the Cholesky decomposition of A held in the matrix cholesky which must have been previously computed by gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp.

— Function: int gsl_linalg_cholesky_svx (const gsl_matrix * cholesky, gsl_vector * x)
— Function: int gsl_linalg_complex_cholesky_svx (const gsl_matrix_complex * cholesky, gsl_vector_complex * x)

These functions solve the system A x = b in-place using the Cholesky decomposition of A held in the matrix cholesky which must have been previously computed by gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp. On input x should contain the right-hand side b, which is replaced by the solution on output.

— Function: int gsl_linalg_cholesky_invert (gsl_matrix * cholesky)
— Function: int gsl_linalg_complex_cholesky_invert (gsl_matrix_complex * cholesky)

These functions compute the inverse of a matrix from its Cholesky decomposition cholesky, which must have been previously computed by gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp. On output, the inverse is stored in-place in cholesky.