scipy.stats.Covariance#

class scipy.stats.Covariance[source]#

Representation of a covariance matrix

Calculations involving covariance matrices (e.g. data whitening, multivariate normal function evaluation) are often performed more efficiently using a decomposition of the covariance matrix instead of the covariance matrix itself. This class allows the user to construct an object representing a covariance matrix using any of several decompositions and perform calculations using a common interface.

Note

The Covariance class cannot be instantiated directly. Instead, use one of the factory methods (e.g. Covariance.from_diagonal).

Examples

The Covariance class is is used by calling one of its factory methods to create a Covariance object, then pass that representation of the Covariance matrix as a shape parameter of a multivariate distribution.

For instance, the multivariate normal distribution can accept an array representing a covariance matrix:

>>> from scipy import stats
>>> import numpy as np
>>> d = [1, 2, 3]
>>> A = np.diag(d)  # a diagonal covariance matrix
>>> x = [4, -2, 5]  # a point of interest
>>> dist = stats.multivariate_normal(mean=[0, 0, 0], cov=A)
>>> dist.pdf(x)
4.9595685102808205e-08

but the calculations are performed in a very generic way that does not take advantage of any special properties of the covariance matrix. Because our covariance matrix is diagonal, we can use Covariance.from_diagonal to create an object representing the covariance matrix, and multivariate_normal can use this to compute the probability density function more efficiently.

>>> cov = stats.Covariance.from_diagonal(d)
>>> dist = stats.multivariate_normal(mean=[0, 0, 0], cov=cov)
>>> dist.pdf(x)
4.9595685102808205e-08
Attributes:
covariance

Explicit representation of the covariance matrix

log_pdet

Log of the pseudo-determinant of the covariance matrix

rank

Rank of the covariance matrix

shape

Shape of the covariance array

Methods

colorize(x)

Perform a colorizing transformation on data.

from_cholesky(cholesky)

Representation of a covariance provided via the (lower) Cholesky factor

from_diagonal(diagonal)

Return a representation of a covariance matrix from its diagonal.

from_eigendecomposition(eigendecomposition)

Representation of a covariance provided via eigendecomposition

from_precision(precision[, covariance])

Return a representation of a covariance from its precision matrix.

whiten(x)

Perform a whitening transformation on data.