scipy.stats.matrix_t#
- scipy.stats.matrix_t = <scipy.stats._multivariate.matrix_t_gen object>[source]#
A matrix t-random variable.
The mean keyword specifies the mean. The row_spread keyword specifies the row-wise spread matrix. The col_spread keyword specifies the column-wise spread matrix.
- Parameters:
- meanarray_like, optional
Mean of the distribution (default: None)
- row_spreadarray_like, optional
Row-wise 2nd order raw central moment matrix (default:
1
)- col_spreadarray_like, optional
Column-wise 2nd order raw central moment matrix (default:
1
)- dfscalar, optional
Degrees of freedom (default:
1
)- seed{None, int, np.random.RandomState, np.random.Generator}, optional
Used for drawing random variates. If seed is None, the RandomState singleton is used. If seed is an int, a new
RandomState
instance is used, seeded with seed. If seed is already aRandomState
orGenerator
instance, then that object is used. Default is None.
Methods
pdf(x, mean=None, row_spread=None, col_spread=None)
Probability density function.
logpdf(x, mean=None, row_spread=None, col_spread=None)
Log of the probability density function.
rvs(mean=None, row_spread=1, col_spread=1, df=1, size=1, random_state=None)
Draw random samples.
Notes
If mean is set to None then a matrix of zeros is used for the mean. The dimensions of this matrix are inferred from the shape of row_spread and col_spread, if these are provided, or set to
1
if ambiguous.row_spread and col_spread can be two-dimensional array_likes specifying the spread matrices directly. Alternatively, a one-dimensional array will be be interpreted as the entries of a diagonal matrix, and a scalar or zero-dimensional array will be interpreted as this value times the identity matrix.
The spread matrices specified by row_spread and col_spread must be (symmetric) positive definite. If the samples in X have shape (m,n) then row_spread must have shape (m,m) and col_spread must have shape (n,n). Spread matrices must be full rank.
The probability density function for
matrix_t
is\[f(X \vert \mathrm{M}, \Sigma, \Omega, \mathrm{df}) = \frac{ \Gamma_n \left( \frac{\mathrm{df} + m + n - 1}{2} \right) \left( \det \left( I_n + (X - \mathrm{M})^T \Sigma^{-1} (X - \mathrm{M}) \Omega^{-1} \right) \right)^{ -\frac{\mathrm{df} + m + n - 1}{2} } }{ \Gamma_n \left( \frac{\mathrm{df} + n - 1}{2} \right) \pi^{mn / 2} \left( \det \Sigma \right)^{n/2} \left( \det \Omega \right)^{m/2} }\]or, alternatively,
\[f(X \vert \mathrm{M}, \Sigma, \Omega, \mathrm{df}) = \frac{ \Gamma_m \left( \frac{\mathrm{df} + m + n - 1}{2} \right) \left( \det \left( I_m + \Sigma^{-1} (X - \mathrm{M}) \Omega^{-1} (X - \mathrm{M})^T \right) \right)^{ -\frac{\mathrm{df} + m + n - 1}{2} } }{ \Gamma_m \left( \frac{\mathrm{df} + n - 1}{2} \right) \pi^{mn / 2} \left( \det \Sigma \right)^{n/2} \left( \det \Omega \right)^{m/2} }\]where \(\mathrm{M}\) is the mean, \(\Sigma\) is the row-wise spread matrix, \(\Omega\) is the column-wise matrix, \(\mathrm{df}\) is the degrees of freedom, and \(\Gamma_n\) is the multivariate gamma function.
These equivalent formulations come from the identity
\[\det\left( I_m + A B \right) = \det\left( I_n + B A \right)\]for \(m \times n\) arrays \(A\) and \(B^T\) and the fact that \(\gamma_n(\mathrm{df} + m) / \gamma_n(\mathrm{df})\) is equal to \(\gamma_m(\mathrm{df} + n) / \gamma_m(\mathrm{df})\), where
\[\gamma_m(\mathrm{df}) = 2^{m(m-1)/2} \Gamma_m\left( (\mathrm{df} + m - 1) / 2 \right)\]denotes a normalized multivariate gamma function.
When \(\mathrm{df} = 1\) this distribution is known as the matrix variate Cauchy.
Added in version 1.17.0.
References
[1]Gupta, A.K., & Nagar, D.K. (2000). Matrix Variate Distributions (1st ed.). Chapman and Hall/CRC.
Examples
>>> import numpy as np >>> from scipy.stats import matrix_t >>> M = np.arange(6).reshape(3,2) >>> M array([[0, 1], [2, 3], [4, 5]]) >>> Sigma = np.diag([1,2,3]) >>> Sigma array([[1, 0, 0], [0, 2, 0], [0, 0, 3]]) >>> Omega = 0.3*np.identity(2) >>> Omega array([[ 0.3, 0. ], [ 0. , 0.3]]) >>> X = M + 0.1 >>> X array([[ 0.1, 1.1], [ 2.1, 3.1], [ 4.1, 5.1]]) >>> df = 3 >>> matrix_t.pdf(X, mean=M, row_spread=Sigma, col_spread=Omega, df=df) 0.9972880280135796
Alternatively, the object may be called (as a function) to fix the mean and spread parameters, returning a “frozen” matrix t random variable:
>>> rv = matrix_t(mean=None, row_spread=1, col_spread=1, df=1) >>> # Frozen object with the same methods but holding the given >>> # mean and spreads and degrees of freedom fixed.