_covariance.py 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649
  1. from functools import cached_property
  2. from types import GenericAlias
  3. import numpy as np
  4. from scipy import linalg
  5. from scipy.stats import _multivariate
  6. __all__ = ["Covariance"]
  7. class Covariance:
  8. """
  9. Representation of a covariance matrix
  10. Calculations involving covariance matrices (e.g. data whitening,
  11. multivariate normal function evaluation) are often performed more
  12. efficiently using a decomposition of the covariance matrix instead of the
  13. covariance matrix itself. This class allows the user to construct an
  14. object representing a covariance matrix using any of several
  15. decompositions and perform calculations using a common interface.
  16. .. note::
  17. The `Covariance` class cannot be instantiated directly. Instead, use
  18. one of the factory methods (e.g. `Covariance.from_diagonal`).
  19. Examples
  20. --------
  21. The `Covariance` class is used by calling one of its
  22. factory methods to create a `Covariance` object, then pass that
  23. representation of the `Covariance` matrix as a shape parameter of a
  24. multivariate distribution.
  25. For instance, the multivariate normal distribution can accept an array
  26. representing a covariance matrix:
  27. >>> from scipy import stats
  28. >>> import numpy as np
  29. >>> d = [1, 2, 3]
  30. >>> A = np.diag(d) # a diagonal covariance matrix
  31. >>> x = [4, -2, 5] # a point of interest
  32. >>> dist = stats.multivariate_normal(mean=[0, 0, 0], cov=A)
  33. >>> dist.pdf(x)
  34. 4.9595685102808205e-08
  35. but the calculations are performed in a very generic way that does not
  36. take advantage of any special properties of the covariance matrix. Because
  37. our covariance matrix is diagonal, we can use ``Covariance.from_diagonal``
  38. to create an object representing the covariance matrix, and
  39. `multivariate_normal` can use this to compute the probability density
  40. function more efficiently.
  41. >>> cov = stats.Covariance.from_diagonal(d)
  42. >>> dist = stats.multivariate_normal(mean=[0, 0, 0], cov=cov)
  43. >>> dist.pdf(x)
  44. 4.9595685102808205e-08
  45. """
  46. # generic type compatibility with scipy-stubs
  47. __class_getitem__ = classmethod(GenericAlias)
  48. def __init__(self):
  49. message = ("The `Covariance` class cannot be instantiated directly. "
  50. "Please use one of the factory methods "
  51. "(e.g. `Covariance.from_diagonal`).")
  52. raise NotImplementedError(message)
  53. @staticmethod
  54. def from_diagonal(diagonal):
  55. r"""
  56. Return a representation of a covariance matrix from its diagonal.
  57. Parameters
  58. ----------
  59. diagonal : array_like
  60. The diagonal elements of a diagonal matrix.
  61. Notes
  62. -----
  63. Let the diagonal elements of a diagonal covariance matrix :math:`D` be
  64. stored in the vector :math:`d`.
  65. When all elements of :math:`d` are strictly positive, whitening of a
  66. data point :math:`x` is performed by computing
  67. :math:`x \cdot d^{-1/2}`, where the inverse square root can be taken
  68. element-wise.
  69. :math:`\log\det{D}` is calculated as :math:`-2 \sum(\log{d})`,
  70. where the :math:`\log` operation is performed element-wise.
  71. This `Covariance` class supports singular covariance matrices. When
  72. computing ``_log_pdet``, non-positive elements of :math:`d` are
  73. ignored. Whitening is not well defined when the point to be whitened
  74. does not lie in the span of the columns of the covariance matrix. The
  75. convention taken here is to treat the inverse square root of
  76. non-positive elements of :math:`d` as zeros.
  77. Examples
  78. --------
  79. Prepare a symmetric positive definite covariance matrix ``A`` and a
  80. data point ``x``.
  81. >>> import numpy as np
  82. >>> from scipy import stats
  83. >>> rng = np.random.default_rng()
  84. >>> n = 5
  85. >>> A = np.diag(rng.random(n))
  86. >>> x = rng.random(size=n)
  87. Extract the diagonal from ``A`` and create the `Covariance` object.
  88. >>> d = np.diag(A)
  89. >>> cov = stats.Covariance.from_diagonal(d)
  90. Compare the functionality of the `Covariance` object against a
  91. reference implementations.
  92. >>> res = cov.whiten(x)
  93. >>> ref = np.diag(d**-0.5) @ x
  94. >>> np.allclose(res, ref)
  95. True
  96. >>> res = cov.log_pdet
  97. >>> ref = np.linalg.slogdet(A)[-1]
  98. >>> np.allclose(res, ref)
  99. True
  100. """
  101. return CovViaDiagonal(diagonal)
  102. @staticmethod
  103. def from_precision(precision, covariance=None):
  104. r"""
  105. Return a representation of a covariance from its precision matrix.
  106. Parameters
  107. ----------
  108. precision : array_like
  109. The precision matrix; that is, the inverse of a square, symmetric,
  110. positive definite covariance matrix.
  111. covariance : array_like, optional
  112. The square, symmetric, positive definite covariance matrix. If not
  113. provided, this may need to be calculated (e.g. to evaluate the
  114. cumulative distribution function of
  115. `scipy.stats.multivariate_normal`) by inverting `precision`.
  116. Notes
  117. -----
  118. Let the covariance matrix be :math:`A`, its precision matrix be
  119. :math:`P = A^{-1}`, and :math:`L` be the lower Cholesky factor such
  120. that :math:`L L^T = P`.
  121. Whitening of a data point :math:`x` is performed by computing
  122. :math:`x^T L`. :math:`\log\det{A}` is calculated as
  123. :math:`-2tr(\log{L})`, where the :math:`\log` operation is performed
  124. element-wise.
  125. This `Covariance` class does not support singular covariance matrices
  126. because the precision matrix does not exist for a singular covariance
  127. matrix.
  128. Examples
  129. --------
  130. Prepare a symmetric positive definite precision matrix ``P`` and a
  131. data point ``x``. (If the precision matrix is not already available,
  132. consider the other factory methods of the `Covariance` class.)
  133. >>> import numpy as np
  134. >>> from scipy import stats
  135. >>> rng = np.random.default_rng()
  136. >>> n = 5
  137. >>> P = rng.random(size=(n, n))
  138. >>> P = P @ P.T # a precision matrix must be positive definite
  139. >>> x = rng.random(size=n)
  140. Create the `Covariance` object.
  141. >>> cov = stats.Covariance.from_precision(P)
  142. Compare the functionality of the `Covariance` object against
  143. reference implementations.
  144. >>> res = cov.whiten(x)
  145. >>> ref = x @ np.linalg.cholesky(P)
  146. >>> np.allclose(res, ref)
  147. True
  148. >>> res = cov.log_pdet
  149. >>> ref = -np.linalg.slogdet(P)[-1]
  150. >>> np.allclose(res, ref)
  151. True
  152. """
  153. return CovViaPrecision(precision, covariance)
  154. @staticmethod
  155. def from_cholesky(cholesky):
  156. r"""
  157. Representation of a covariance provided via the (lower) Cholesky factor
  158. Parameters
  159. ----------
  160. cholesky : array_like
  161. The lower triangular Cholesky factor of the covariance matrix.
  162. Notes
  163. -----
  164. Let the covariance matrix be :math:`A` and :math:`L` be the lower
  165. Cholesky factor such that :math:`L L^T = A`.
  166. Whitening of a data point :math:`x` is performed by computing
  167. :math:`L^{-1} x`. :math:`\log\det{A}` is calculated as
  168. :math:`2tr(\log{L})`, where the :math:`\log` operation is performed
  169. element-wise.
  170. This `Covariance` class does not support singular covariance matrices
  171. because the Cholesky decomposition does not exist for a singular
  172. covariance matrix.
  173. Examples
  174. --------
  175. Prepare a symmetric positive definite covariance matrix ``A`` and a
  176. data point ``x``.
  177. >>> import numpy as np
  178. >>> from scipy import stats
  179. >>> rng = np.random.default_rng()
  180. >>> n = 5
  181. >>> A = rng.random(size=(n, n))
  182. >>> A = A @ A.T # make the covariance symmetric positive definite
  183. >>> x = rng.random(size=n)
  184. Perform the Cholesky decomposition of ``A`` and create the
  185. `Covariance` object.
  186. >>> L = np.linalg.cholesky(A)
  187. >>> cov = stats.Covariance.from_cholesky(L)
  188. Compare the functionality of the `Covariance` object against
  189. reference implementation.
  190. >>> from scipy.linalg import solve_triangular
  191. >>> res = cov.whiten(x)
  192. >>> ref = solve_triangular(L, x, lower=True)
  193. >>> np.allclose(res, ref)
  194. True
  195. >>> res = cov.log_pdet
  196. >>> ref = np.linalg.slogdet(A)[-1]
  197. >>> np.allclose(res, ref)
  198. True
  199. """
  200. return CovViaCholesky(cholesky)
  201. @staticmethod
  202. def from_eigendecomposition(eigendecomposition):
  203. r"""
  204. Representation of a covariance provided via eigendecomposition
  205. Parameters
  206. ----------
  207. eigendecomposition : sequence
  208. A sequence (nominally a tuple) containing the eigenvalue and
  209. eigenvector arrays as computed by `scipy.linalg.eigh` or
  210. `numpy.linalg.eigh`.
  211. Notes
  212. -----
  213. Let the covariance matrix be :math:`A`, let :math:`V` be matrix of
  214. eigenvectors, and let :math:`W` be the diagonal matrix of eigenvalues
  215. such that `V W V^T = A`.
  216. When all of the eigenvalues are strictly positive, whitening of a
  217. data point :math:`x` is performed by computing
  218. :math:`x^T (V W^{-1/2})`, where the inverse square root can be taken
  219. element-wise.
  220. :math:`\log\det{A}` is calculated as :math:`tr(\log{W})`,
  221. where the :math:`\log` operation is performed element-wise.
  222. This `Covariance` class supports singular covariance matrices. When
  223. computing ``_log_pdet``, non-positive eigenvalues are ignored.
  224. Whitening is not well defined when the point to be whitened
  225. does not lie in the span of the columns of the covariance matrix. The
  226. convention taken here is to treat the inverse square root of
  227. non-positive eigenvalues as zeros.
  228. Examples
  229. --------
  230. Prepare a symmetric positive definite covariance matrix ``A`` and a
  231. data point ``x``.
  232. >>> import numpy as np
  233. >>> from scipy import stats
  234. >>> rng = np.random.default_rng()
  235. >>> n = 5
  236. >>> A = rng.random(size=(n, n))
  237. >>> A = A @ A.T # make the covariance symmetric positive definite
  238. >>> x = rng.random(size=n)
  239. Perform the eigendecomposition of ``A`` and create the `Covariance`
  240. object.
  241. >>> w, v = np.linalg.eigh(A)
  242. >>> cov = stats.Covariance.from_eigendecomposition((w, v))
  243. Compare the functionality of the `Covariance` object against
  244. reference implementations.
  245. >>> res = cov.whiten(x)
  246. >>> ref = x @ (v @ np.diag(w**-0.5))
  247. >>> np.allclose(res, ref)
  248. True
  249. >>> res = cov.log_pdet
  250. >>> ref = np.linalg.slogdet(A)[-1]
  251. >>> np.allclose(res, ref)
  252. True
  253. """
  254. return CovViaEigendecomposition(eigendecomposition)
  255. def whiten(self, x):
  256. """
  257. Perform a whitening transformation on data.
  258. "Whitening" ("white" as in "white noise", in which each frequency has
  259. equal magnitude) transforms a set of random variables into a new set of
  260. random variables with unit-diagonal covariance. When a whitening
  261. transform is applied to a sample of points distributed according to
  262. a multivariate normal distribution with zero mean, the covariance of
  263. the transformed sample is approximately the identity matrix.
  264. Parameters
  265. ----------
  266. x : array_like
  267. An array of points. The last dimension must correspond with the
  268. dimensionality of the space, i.e., the number of columns in the
  269. covariance matrix.
  270. Returns
  271. -------
  272. x_ : array_like
  273. The transformed array of points.
  274. References
  275. ----------
  276. .. [1] "Whitening Transformation". Wikipedia.
  277. https://en.wikipedia.org/wiki/Whitening_transformation
  278. .. [2] Novak, Lukas, and Miroslav Vorechovsky. "Generalization of
  279. coloring linear transformation". Transactions of VSB 18.2
  280. (2018): 31-35. :doi:`10.31490/tces-2018-0013`
  281. Examples
  282. --------
  283. >>> import numpy as np
  284. >>> from scipy import stats
  285. >>> rng = np.random.default_rng()
  286. >>> n = 3
  287. >>> A = rng.random(size=(n, n))
  288. >>> cov_array = A @ A.T # make matrix symmetric positive definite
  289. >>> precision = np.linalg.inv(cov_array)
  290. >>> cov_object = stats.Covariance.from_precision(precision)
  291. >>> x = rng.multivariate_normal(np.zeros(n), cov_array, size=(10000))
  292. >>> x_ = cov_object.whiten(x)
  293. >>> np.cov(x_, rowvar=False) # near-identity covariance
  294. array([[0.97862122, 0.00893147, 0.02430451],
  295. [0.00893147, 0.96719062, 0.02201312],
  296. [0.02430451, 0.02201312, 0.99206881]])
  297. """
  298. return self._whiten(np.asarray(x))
  299. def colorize(self, x):
  300. """
  301. Perform a colorizing transformation on data.
  302. "Colorizing" ("color" as in "colored noise", in which different
  303. frequencies may have different magnitudes) transforms a set of
  304. uncorrelated random variables into a new set of random variables with
  305. the desired covariance. When a coloring transform is applied to a
  306. sample of points distributed according to a multivariate normal
  307. distribution with identity covariance and zero mean, the covariance of
  308. the transformed sample is approximately the covariance matrix used
  309. in the coloring transform.
  310. Parameters
  311. ----------
  312. x : array_like
  313. An array of points. The last dimension must correspond with the
  314. dimensionality of the space, i.e., the number of columns in the
  315. covariance matrix.
  316. Returns
  317. -------
  318. x_ : array_like
  319. The transformed array of points.
  320. References
  321. ----------
  322. .. [1] "Whitening Transformation". Wikipedia.
  323. https://en.wikipedia.org/wiki/Whitening_transformation
  324. .. [2] Novak, Lukas, and Miroslav Vorechovsky. "Generalization of
  325. coloring linear transformation". Transactions of VSB 18.2
  326. (2018): 31-35. :doi:`10.31490/tces-2018-0013`
  327. Examples
  328. --------
  329. >>> import numpy as np
  330. >>> from scipy import stats
  331. >>> rng = np.random.default_rng(1638083107694713882823079058616272161)
  332. >>> n = 3
  333. >>> A = rng.random(size=(n, n))
  334. >>> cov_array = A @ A.T # make matrix symmetric positive definite
  335. >>> cholesky = np.linalg.cholesky(cov_array)
  336. >>> cov_object = stats.Covariance.from_cholesky(cholesky)
  337. >>> x = rng.multivariate_normal(np.zeros(n), np.eye(n), size=(10000))
  338. >>> x_ = cov_object.colorize(x)
  339. >>> cov_data = np.cov(x_, rowvar=False)
  340. >>> np.allclose(cov_data, cov_array, rtol=3e-2)
  341. True
  342. """
  343. return self._colorize(np.asarray(x))
  344. @property
  345. def log_pdet(self):
  346. """
  347. Log of the pseudo-determinant of the covariance matrix
  348. """
  349. return np.array(self._log_pdet, dtype=float)[()]
  350. @property
  351. def rank(self):
  352. """
  353. Rank of the covariance matrix
  354. """
  355. return np.array(self._rank, dtype=int)[()]
  356. @property
  357. def covariance(self):
  358. """
  359. Explicit representation of the covariance matrix
  360. """
  361. return self._covariance
  362. @property
  363. def shape(self):
  364. """
  365. Shape of the covariance array
  366. """
  367. return self._shape
  368. def _validate_matrix(self, A, name):
  369. A = np.atleast_2d(A)
  370. m, n = A.shape[-2:]
  371. if m != n or A.ndim != 2 or not (np.issubdtype(A.dtype, np.integer) or
  372. np.issubdtype(A.dtype, np.floating)):
  373. message = (f"The input `{name}` must be a square, "
  374. "two-dimensional array of real numbers.")
  375. raise ValueError(message)
  376. return A
  377. def _validate_vector(self, A, name):
  378. A = np.atleast_1d(A)
  379. if A.ndim != 1 or not (np.issubdtype(A.dtype, np.integer) or
  380. np.issubdtype(A.dtype, np.floating)):
  381. message = (f"The input `{name}` must be a one-dimensional array "
  382. "of real numbers.")
  383. raise ValueError(message)
  384. return A
  385. class CovViaPrecision(Covariance):
  386. __class_getitem__ = None
  387. def __init__(self, precision, covariance=None):
  388. precision = self._validate_matrix(precision, 'precision')
  389. if covariance is not None:
  390. covariance = self._validate_matrix(covariance, 'covariance')
  391. message = "`precision.shape` must equal `covariance.shape`."
  392. if precision.shape != covariance.shape:
  393. raise ValueError(message)
  394. self._chol_P = np.linalg.cholesky(precision)
  395. self._log_pdet = -2*np.log(np.diag(self._chol_P)).sum(axis=-1)
  396. self._rank = precision.shape[-1] # must be full rank if invertible
  397. self._precision = precision
  398. self._cov_matrix = covariance
  399. self._shape = precision.shape
  400. self._allow_singular = False
  401. def _whiten(self, x):
  402. return x @ self._chol_P
  403. @cached_property
  404. def _covariance(self):
  405. n = self._shape[-1]
  406. return (linalg.cho_solve((self._chol_P, True), np.eye(n))
  407. if self._cov_matrix is None else self._cov_matrix)
  408. def _colorize(self, x):
  409. m = x.T.shape[0]
  410. res = linalg.solve_triangular(self._chol_P.T, x.T.reshape(m, -1), lower=False)
  411. return res.reshape(x.T.shape).T
  412. def _dot_diag(x, d):
  413. # If d were a full diagonal matrix, x @ d would always do what we want.
  414. # Special treatment is needed for n-dimensional `d` in which each row
  415. # includes only the diagonal elements of a covariance matrix.
  416. return x * d if x.ndim < 2 else x * np.expand_dims(d, -2)
  417. class CovViaDiagonal(Covariance):
  418. def __init__(self, diagonal):
  419. diagonal = self._validate_vector(diagonal, 'diagonal')
  420. i_zero = diagonal <= 0
  421. positive_diagonal = np.array(diagonal, dtype=np.float64)
  422. positive_diagonal[i_zero] = 1 # ones don't affect determinant
  423. self._log_pdet = np.sum(np.log(positive_diagonal), axis=-1)
  424. psuedo_reciprocals = 1 / np.sqrt(positive_diagonal)
  425. psuedo_reciprocals[i_zero] = 0
  426. self._sqrt_diagonal = np.sqrt(diagonal)
  427. self._LP = psuedo_reciprocals
  428. self._rank = positive_diagonal.shape[-1] - i_zero.sum(axis=-1)
  429. self._covariance = np.apply_along_axis(np.diag, -1, diagonal)
  430. self._i_zero = i_zero
  431. self._shape = self._covariance.shape
  432. self._allow_singular = True
  433. def _whiten(self, x):
  434. return _dot_diag(x, self._LP)
  435. def _colorize(self, x):
  436. return _dot_diag(x, self._sqrt_diagonal)
  437. def _support_mask(self, x):
  438. """
  439. Check whether x lies in the support of the distribution.
  440. """
  441. return ~np.any(_dot_diag(x, self._i_zero), axis=-1)
  442. class CovViaCholesky(Covariance):
  443. __class_getitem__ = None
  444. def __init__(self, cholesky):
  445. L = self._validate_matrix(cholesky, 'cholesky')
  446. self._factor = L
  447. self._log_pdet = 2*np.log(np.diag(self._factor)).sum(axis=-1)
  448. self._rank = L.shape[-1] # must be full rank for cholesky
  449. self._shape = L.shape
  450. self._allow_singular = False
  451. @cached_property
  452. def _covariance(self):
  453. return self._factor @ self._factor.T
  454. def _whiten(self, x):
  455. m = x.T.shape[0]
  456. res = linalg.solve_triangular(self._factor, x.T.reshape(m, -1), lower=True)
  457. return res.reshape(x.T.shape).T
  458. def _colorize(self, x):
  459. return x @ self._factor.T
  460. class CovViaEigendecomposition(Covariance):
  461. __class_getitem__ = None
  462. def __init__(self, eigendecomposition):
  463. eigenvalues, eigenvectors = eigendecomposition
  464. eigenvalues = self._validate_vector(eigenvalues, 'eigenvalues')
  465. eigenvectors = self._validate_matrix(eigenvectors, 'eigenvectors')
  466. message = ("The shapes of `eigenvalues` and `eigenvectors` "
  467. "must be compatible.")
  468. try:
  469. eigenvalues = np.expand_dims(eigenvalues, -2)
  470. eigenvectors, eigenvalues = np.broadcast_arrays(eigenvectors,
  471. eigenvalues)
  472. eigenvalues = eigenvalues[..., 0, :]
  473. except ValueError:
  474. raise ValueError(message)
  475. i_zero = eigenvalues <= 0
  476. positive_eigenvalues = np.array(eigenvalues, dtype=np.float64)
  477. positive_eigenvalues[i_zero] = 1 # ones don't affect determinant
  478. self._log_pdet = np.sum(np.log(positive_eigenvalues), axis=-1)
  479. psuedo_reciprocals = 1 / np.sqrt(positive_eigenvalues)
  480. psuedo_reciprocals[i_zero] = 0
  481. self._LP = eigenvectors * psuedo_reciprocals
  482. self._LA = eigenvectors * np.sqrt(eigenvalues)
  483. self._rank = positive_eigenvalues.shape[-1] - i_zero.sum(axis=-1)
  484. self._w = eigenvalues
  485. self._v = eigenvectors
  486. self._shape = eigenvectors.shape
  487. self._null_basis = eigenvectors * i_zero
  488. # This is only used for `_support_mask`, not to decide whether
  489. # the covariance is singular or not.
  490. self._eps = _multivariate._eigvalsh_to_eps(eigenvalues) * 10**3
  491. self._allow_singular = True
  492. def _whiten(self, x):
  493. return x @ self._LP
  494. def _colorize(self, x):
  495. return x @ self._LA.T
  496. @cached_property
  497. def _covariance(self):
  498. return (self._v * self._w) @ self._v.T
  499. def _support_mask(self, x):
  500. """
  501. Check whether x lies in the support of the distribution.
  502. """
  503. residual = np.linalg.norm(x @ self._null_basis, axis=-1)
  504. in_support = residual < self._eps
  505. return in_support
  506. class CovViaPSD(Covariance):
  507. """
  508. Representation of a covariance provided via an instance of _PSD
  509. """
  510. __class_getitem__ = None
  511. def __init__(self, psd):
  512. self._LP = psd.U
  513. self._log_pdet = psd.log_pdet
  514. self._rank = psd.rank
  515. self._covariance = psd._M
  516. self._shape = psd._M.shape
  517. self._psd = psd
  518. self._allow_singular = False # by default
  519. def _whiten(self, x):
  520. return x @ self._LP
  521. def _support_mask(self, x):
  522. return self._psd._support_mask(x)