_decomp_cholesky.py 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358
  1. """Cholesky decomposition functions."""
  2. from numpy import asarray_chkfinite, asarray, atleast_2d
  3. # Local imports
  4. from ._misc import LinAlgError, _datacopied
  5. from .lapack import get_lapack_funcs
  6. __all__ = ['cholesky', 'cho_factor', 'cho_solve', 'cholesky_banded',
  7. 'cho_solve_banded']
  8. def _cholesky(a, lower=False, overwrite_a=False, clean=True,
  9. check_finite=True):
  10. """Common code for cholesky() and cho_factor()."""
  11. a1 = asarray_chkfinite(a) if check_finite else asarray(a)
  12. a1 = atleast_2d(a1)
  13. # Dimension check
  14. if a1.ndim != 2:
  15. raise ValueError('Input array needs to be 2D but received '
  16. 'a {}d-array.'.format(a1.ndim))
  17. # Squareness check
  18. if a1.shape[0] != a1.shape[1]:
  19. raise ValueError('Input array is expected to be square but has '
  20. 'the shape: {}.'.format(a1.shape))
  21. # Quick return for square empty array
  22. if a1.size == 0:
  23. return a1.copy(), lower
  24. overwrite_a = overwrite_a or _datacopied(a1, a)
  25. potrf, = get_lapack_funcs(('potrf',), (a1,))
  26. c, info = potrf(a1, lower=lower, overwrite_a=overwrite_a, clean=clean)
  27. if info > 0:
  28. raise LinAlgError("%d-th leading minor of the array is not positive "
  29. "definite" % info)
  30. if info < 0:
  31. raise ValueError('LAPACK reported an illegal value in {}-th argument'
  32. 'on entry to "POTRF".'.format(-info))
  33. return c, lower
  34. def cholesky(a, lower=False, overwrite_a=False, check_finite=True):
  35. """
  36. Compute the Cholesky decomposition of a matrix.
  37. Returns the Cholesky decomposition, :math:`A = L L^*` or
  38. :math:`A = U^* U` of a Hermitian positive-definite matrix A.
  39. Parameters
  40. ----------
  41. a : (M, M) array_like
  42. Matrix to be decomposed
  43. lower : bool, optional
  44. Whether to compute the upper- or lower-triangular Cholesky
  45. factorization. Default is upper-triangular.
  46. overwrite_a : bool, optional
  47. Whether to overwrite data in `a` (may improve performance).
  48. check_finite : bool, optional
  49. Whether to check that the input matrix contains only finite numbers.
  50. Disabling may give a performance gain, but may result in problems
  51. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  52. Returns
  53. -------
  54. c : (M, M) ndarray
  55. Upper- or lower-triangular Cholesky factor of `a`.
  56. Raises
  57. ------
  58. LinAlgError : if decomposition fails.
  59. Examples
  60. --------
  61. >>> import numpy as np
  62. >>> from scipy.linalg import cholesky
  63. >>> a = np.array([[1,-2j],[2j,5]])
  64. >>> L = cholesky(a, lower=True)
  65. >>> L
  66. array([[ 1.+0.j, 0.+0.j],
  67. [ 0.+2.j, 1.+0.j]])
  68. >>> L @ L.T.conj()
  69. array([[ 1.+0.j, 0.-2.j],
  70. [ 0.+2.j, 5.+0.j]])
  71. """
  72. c, lower = _cholesky(a, lower=lower, overwrite_a=overwrite_a, clean=True,
  73. check_finite=check_finite)
  74. return c
  75. def cho_factor(a, lower=False, overwrite_a=False, check_finite=True):
  76. """
  77. Compute the Cholesky decomposition of a matrix, to use in cho_solve
  78. Returns a matrix containing the Cholesky decomposition,
  79. ``A = L L*`` or ``A = U* U`` of a Hermitian positive-definite matrix `a`.
  80. The return value can be directly used as the first parameter to cho_solve.
  81. .. warning::
  82. The returned matrix also contains random data in the entries not
  83. used by the Cholesky decomposition. If you need to zero these
  84. entries, use the function `cholesky` instead.
  85. Parameters
  86. ----------
  87. a : (M, M) array_like
  88. Matrix to be decomposed
  89. lower : bool, optional
  90. Whether to compute the upper or lower triangular Cholesky factorization
  91. (Default: upper-triangular)
  92. overwrite_a : bool, optional
  93. Whether to overwrite data in a (may improve performance)
  94. check_finite : bool, optional
  95. Whether to check that the input matrix contains only finite numbers.
  96. Disabling may give a performance gain, but may result in problems
  97. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  98. Returns
  99. -------
  100. c : (M, M) ndarray
  101. Matrix whose upper or lower triangle contains the Cholesky factor
  102. of `a`. Other parts of the matrix contain random data.
  103. lower : bool
  104. Flag indicating whether the factor is in the lower or upper triangle
  105. Raises
  106. ------
  107. LinAlgError
  108. Raised if decomposition fails.
  109. See Also
  110. --------
  111. cho_solve : Solve a linear set equations using the Cholesky factorization
  112. of a matrix.
  113. Examples
  114. --------
  115. >>> import numpy as np
  116. >>> from scipy.linalg import cho_factor
  117. >>> A = np.array([[9, 3, 1, 5], [3, 7, 5, 1], [1, 5, 9, 2], [5, 1, 2, 6]])
  118. >>> c, low = cho_factor(A)
  119. >>> c
  120. array([[3. , 1. , 0.33333333, 1.66666667],
  121. [3. , 2.44948974, 1.90515869, -0.27216553],
  122. [1. , 5. , 2.29330749, 0.8559528 ],
  123. [5. , 1. , 2. , 1.55418563]])
  124. >>> np.allclose(np.triu(c).T @ np. triu(c) - A, np.zeros((4, 4)))
  125. True
  126. """
  127. c, lower = _cholesky(a, lower=lower, overwrite_a=overwrite_a, clean=False,
  128. check_finite=check_finite)
  129. return c, lower
  130. def cho_solve(c_and_lower, b, overwrite_b=False, check_finite=True):
  131. """Solve the linear equations A x = b, given the Cholesky factorization of A.
  132. Parameters
  133. ----------
  134. (c, lower) : tuple, (array, bool)
  135. Cholesky factorization of a, as given by cho_factor
  136. b : array
  137. Right-hand side
  138. overwrite_b : bool, optional
  139. Whether to overwrite data in b (may improve performance)
  140. check_finite : bool, optional
  141. Whether to check that the input matrices contain only finite numbers.
  142. Disabling may give a performance gain, but may result in problems
  143. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  144. Returns
  145. -------
  146. x : array
  147. The solution to the system A x = b
  148. See Also
  149. --------
  150. cho_factor : Cholesky factorization of a matrix
  151. Examples
  152. --------
  153. >>> import numpy as np
  154. >>> from scipy.linalg import cho_factor, cho_solve
  155. >>> A = np.array([[9, 3, 1, 5], [3, 7, 5, 1], [1, 5, 9, 2], [5, 1, 2, 6]])
  156. >>> c, low = cho_factor(A)
  157. >>> x = cho_solve((c, low), [1, 1, 1, 1])
  158. >>> np.allclose(A @ x - [1, 1, 1, 1], np.zeros(4))
  159. True
  160. """
  161. (c, lower) = c_and_lower
  162. if check_finite:
  163. b1 = asarray_chkfinite(b)
  164. c = asarray_chkfinite(c)
  165. else:
  166. b1 = asarray(b)
  167. c = asarray(c)
  168. if c.ndim != 2 or c.shape[0] != c.shape[1]:
  169. raise ValueError("The factored matrix c is not square.")
  170. if c.shape[1] != b1.shape[0]:
  171. raise ValueError("incompatible dimensions ({} and {})"
  172. .format(c.shape, b1.shape))
  173. overwrite_b = overwrite_b or _datacopied(b1, b)
  174. potrs, = get_lapack_funcs(('potrs',), (c, b1))
  175. x, info = potrs(c, b1, lower=lower, overwrite_b=overwrite_b)
  176. if info != 0:
  177. raise ValueError('illegal value in %dth argument of internal potrs'
  178. % -info)
  179. return x
  180. def cholesky_banded(ab, overwrite_ab=False, lower=False, check_finite=True):
  181. """
  182. Cholesky decompose a banded Hermitian positive-definite matrix
  183. The matrix a is stored in ab either in lower-diagonal or upper-
  184. diagonal ordered form::
  185. ab[u + i - j, j] == a[i,j] (if upper form; i <= j)
  186. ab[ i - j, j] == a[i,j] (if lower form; i >= j)
  187. Example of ab (shape of a is (6,6), u=2)::
  188. upper form:
  189. * * a02 a13 a24 a35
  190. * a01 a12 a23 a34 a45
  191. a00 a11 a22 a33 a44 a55
  192. lower form:
  193. a00 a11 a22 a33 a44 a55
  194. a10 a21 a32 a43 a54 *
  195. a20 a31 a42 a53 * *
  196. Parameters
  197. ----------
  198. ab : (u + 1, M) array_like
  199. Banded matrix
  200. overwrite_ab : bool, optional
  201. Discard data in ab (may enhance performance)
  202. lower : bool, optional
  203. Is the matrix in the lower form. (Default is upper form)
  204. check_finite : bool, optional
  205. Whether to check that the input matrix contains only finite numbers.
  206. Disabling may give a performance gain, but may result in problems
  207. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  208. Returns
  209. -------
  210. c : (u + 1, M) ndarray
  211. Cholesky factorization of a, in the same banded format as ab
  212. See Also
  213. --------
  214. cho_solve_banded :
  215. Solve a linear set equations, given the Cholesky factorization
  216. of a banded Hermitian.
  217. Examples
  218. --------
  219. >>> import numpy as np
  220. >>> from scipy.linalg import cholesky_banded
  221. >>> from numpy import allclose, zeros, diag
  222. >>> Ab = np.array([[0, 0, 1j, 2, 3j], [0, -1, -2, 3, 4], [9, 8, 7, 6, 9]])
  223. >>> A = np.diag(Ab[0,2:], k=2) + np.diag(Ab[1,1:], k=1)
  224. >>> A = A + A.conj().T + np.diag(Ab[2, :])
  225. >>> c = cholesky_banded(Ab)
  226. >>> C = np.diag(c[0, 2:], k=2) + np.diag(c[1, 1:], k=1) + np.diag(c[2, :])
  227. >>> np.allclose(C.conj().T @ C - A, np.zeros((5, 5)))
  228. True
  229. """
  230. if check_finite:
  231. ab = asarray_chkfinite(ab)
  232. else:
  233. ab = asarray(ab)
  234. pbtrf, = get_lapack_funcs(('pbtrf',), (ab,))
  235. c, info = pbtrf(ab, lower=lower, overwrite_ab=overwrite_ab)
  236. if info > 0:
  237. raise LinAlgError("%d-th leading minor not positive definite" % info)
  238. if info < 0:
  239. raise ValueError('illegal value in %d-th argument of internal pbtrf'
  240. % -info)
  241. return c
  242. def cho_solve_banded(cb_and_lower, b, overwrite_b=False, check_finite=True):
  243. """
  244. Solve the linear equations ``A x = b``, given the Cholesky factorization of
  245. the banded Hermitian ``A``.
  246. Parameters
  247. ----------
  248. (cb, lower) : tuple, (ndarray, bool)
  249. `cb` is the Cholesky factorization of A, as given by cholesky_banded.
  250. `lower` must be the same value that was given to cholesky_banded.
  251. b : array_like
  252. Right-hand side
  253. overwrite_b : bool, optional
  254. If True, the function will overwrite the values in `b`.
  255. check_finite : bool, optional
  256. Whether to check that the input matrices contain only finite numbers.
  257. Disabling may give a performance gain, but may result in problems
  258. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  259. Returns
  260. -------
  261. x : array
  262. The solution to the system A x = b
  263. See Also
  264. --------
  265. cholesky_banded : Cholesky factorization of a banded matrix
  266. Notes
  267. -----
  268. .. versionadded:: 0.8.0
  269. Examples
  270. --------
  271. >>> import numpy as np
  272. >>> from scipy.linalg import cholesky_banded, cho_solve_banded
  273. >>> Ab = np.array([[0, 0, 1j, 2, 3j], [0, -1, -2, 3, 4], [9, 8, 7, 6, 9]])
  274. >>> A = np.diag(Ab[0,2:], k=2) + np.diag(Ab[1,1:], k=1)
  275. >>> A = A + A.conj().T + np.diag(Ab[2, :])
  276. >>> c = cholesky_banded(Ab)
  277. >>> x = cho_solve_banded((c, False), np.ones(5))
  278. >>> np.allclose(A @ x - np.ones(5), np.zeros(5))
  279. True
  280. """
  281. (cb, lower) = cb_and_lower
  282. if check_finite:
  283. cb = asarray_chkfinite(cb)
  284. b = asarray_chkfinite(b)
  285. else:
  286. cb = asarray(cb)
  287. b = asarray(b)
  288. # Validate shapes.
  289. if cb.shape[-1] != b.shape[0]:
  290. raise ValueError("shapes of cb and b are not compatible.")
  291. pbtrs, = get_lapack_funcs(('pbtrs',), (cb, b))
  292. x, info = pbtrs(cb, b, lower=lower, overwrite_b=overwrite_b)
  293. if info > 0:
  294. raise LinAlgError("%dth leading minor not positive definite" % info)
  295. if info < 0:
  296. raise ValueError('illegal value in %dth argument of internal pbtrs'
  297. % -info)
  298. return x