_covariance.py 22 KB

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