covmats.CovViaCholesky#

class covmats.CovViaCholesky(*args, **kwargs)[source]#

Representation of a covariance via a Cholesky factorization.

Parameters:

cholesky (array_like) – The lower triangular Cholesky factor of the covariance matrix.

Notes

Let the covariance matrix be \(A\) and \(L\) be the lower Cholesky factor such that \(L L^T = A\). Whitening of a data point \(x\) is performed by computing \(L^{-1} x\). \(\log\det{A}\) is calculated as \(2tr(\log{L})\), where the \(\log\) operation is performed element-wise.

This Covariance class does not support singular covariance matrices because the Cholesky decomposition does not exist for a singular covariance matrix.

Examples

Prepare a symmetric positive definite covariance matrix A and a data point x.

>>> import numpy as np
>>> import covmats
>>> rng = np.random.default_rng()
>>> n = 5
>>> A = rng.random(size=(n, n))
>>> A = A @ A.T  # make the covariance symmetric positive definite
>>> x = rng.random(size=n)

Perform the Cholesky decomposition of A and create the Covariance object.

>>> L = np.linalg.cholesky(A)
>>> cov = covmats.CovViaCholesky(L)

Compare the functionality of the Covariance object against reference implementation.

>>> from scipy.linalg import solve_triangular
>>> res = cov.whiten(x)
>>> ref = solve_triangular(L, x, lower=True)
>>> np.allclose(res, ref)
True
>>> res = cov.log_pdet
>>> ref = np.linalg.slogdet(A)[-1]
>>> np.allclose(res, ref)
True
__init__(L: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str], covariance: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str] | None = None) None[source]#

Initialize the instance.

Parameters:
  • L (NDArrayFloat) – Lower triangle of the covariance matrix Cholesky factorization.

  • covariance (Optional[NDArrayFloat], optional) – Dense covariance matrix, by default None

Properties

H

Hermitian adjoint.

L

Lower triangle of the covariance matrix Cholesky factorization.

T

Transpose this linear operator.

covariance

Explicit dense representation of the covariance matrix.

log_pdet

Log of the pseudo-determinant of the covariance matrix.

n_pts

Number of points in the domain (n).

ndim

precision

Explicit dense representation of the precision matrix.

rank

Rank of the covariance matrix.

shape

Shape of the covariance matrix (n, n).

subspace_size

Subspace size of the covariance matrix.

Methods

adjoint

Hermitian adjoint.

colorize

Perform a colorizing transformation on data.

dot

Matrix-matrix or matrix-vector multiplication.

from_cholesky

Representation of a covariance provided via choleksy factorization.

from_diagonal

Representation of a covariance provided via diagonal.

from_eigendecomposition

Representation of a covariance provided via eigendecomposition.

from_precision

Return a representation of a covariance from its precision matrix.

get_diagonal

Return the diagonal entries of the matrix (variances).

get_trace

Return the trace of the covariance matrix (sum of diagonal elements).

matmat

Matrix-matrix multiplication.

matvec

Matrix-vector multiplication.

rmatmat

Adjoint matrix-matrix multiplication.

rmatvec

Adjoint matrix-vector multiplication.

sample_mvnormal

Draw samples from the multivariate normal N(0, A).

solve

Solve Ax = b, with A, the current covariance matrix instance.

todense

Explicit dense representation of the covariance matrix with shape (n, n).

transpose

Transpose this linear operator.

whiten

Perform a whitening transformation on data.