logdet#

detkit.logdet(A, sym_pos=False, overwrite_A=False, use_scipy=True)#

Computes the logdet of a matrix.

The logdet function is defined by

\[\mathrm{logdet}(\mathbf{A}) := \log_{e} |\mathrm{det}(\mathbf{A})|.\]
Parameters:
A(n, n) array_like

Square matrix. The matrix type can be float32, float64, or float128. If a matrix of the type int32 or int64 is given, the type is cast to float64.

sym_posbool, default=False

If True, the matrix A is assumed to be symmetric and positive-definite (SPD). The computation can be twice as fast as when the matrix is not SPD. This function does not verify whether A is symmetric or positive-definite.

overwrite_Abool, default=False

If True, the input matrix A will be overwritten during the computation. It uses less memory and could potentially be slightly faster.

use_scipybool, default=False

If True, it uses scipy functions which are the wrappers around Fortran routines in BLAS and LAPACK. If False, it uses a C++ library developed in this package.

Returns:
logdetfloat

logdet of A. If A is singular, returns -numpy.inf.

signint

Sign of the determinant of A and can be +1 for positive or -1 for negative determinant. If A is singular, returns 0.

Raises:
RuntimeError

Error raised when sym_pos=True and matrix A is not symmetric positive-definite.

See also

loggdet

Log of determinant terms used in Gaussian process regression.

logpdet

Log of pseudo-determinant of the precision matrix in Gaussian process regression.

Notes

The function logdet is computed using the following algorithms:

  • When sym_pos=False, the logdet function is computed using the PLU decomposition of A.

  • When sym_pos=True, the logdet function is computed using the Cholesky decomposition of A.

This python function is a wrapper to a C++ implementation.

Examples

>>> import numpy
>>> from detkit import logdet

>>> # Generate a random matrix
>>> n = 1000
>>> rng = numpy.random.RandomState(0)
>>> A = rng.rand(n, n)

>>> # Compute logdet of matrix
>>> logdet(A)
(1710.9576831500378, -1)

>>> # Compute logdet of a symmetric and positive-definite matrix
>>> B = A.T @ A
>>> logdet(B, sym_pos=True)
(3421.9153663693114, 1)

>>> # Compute logdet of a singular matrix
>>> A[:, 0] = 0
>>> logdet(A)
(-inf, 0)