Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1"""Cholesky decomposition functions.""" 

2 

3from numpy import asarray_chkfinite, asarray, atleast_2d 

4 

5# Local imports 

6from .misc import LinAlgError, _datacopied 

7from .lapack import get_lapack_funcs 

8 

9__all__ = ['cholesky', 'cho_factor', 'cho_solve', 'cholesky_banded', 

10 'cho_solve_banded'] 

11 

12 

13def _cholesky(a, lower=False, overwrite_a=False, clean=True, 

14 check_finite=True): 

15 """Common code for cholesky() and cho_factor().""" 

16 

17 a1 = asarray_chkfinite(a) if check_finite else asarray(a) 

18 a1 = atleast_2d(a1) 

19 

20 # Dimension check 

21 if a1.ndim != 2: 

22 raise ValueError('Input array needs to be 2D but received ' 

23 'a {}d-array.'.format(a1.ndim)) 

24 # Squareness check 

25 if a1.shape[0] != a1.shape[1]: 

26 raise ValueError('Input array is expected to be square but has ' 

27 'the shape: {}.'.format(a1.shape)) 

28 

29 # Quick return for square empty array 

30 if a1.size == 0: 

31 return a1.copy(), lower 

32 

33 overwrite_a = overwrite_a or _datacopied(a1, a) 

34 potrf, = get_lapack_funcs(('potrf',), (a1,)) 

35 c, info = potrf(a1, lower=lower, overwrite_a=overwrite_a, clean=clean) 

36 if info > 0: 

37 raise LinAlgError("%d-th leading minor of the array is not positive " 

38 "definite" % info) 

39 if info < 0: 

40 raise ValueError('LAPACK reported an illegal value in {}-th argument' 

41 'on entry to "POTRF".'.format(-info)) 

42 return c, lower 

43 

44 

45def cholesky(a, lower=False, overwrite_a=False, check_finite=True): 

46 """ 

47 Compute the Cholesky decomposition of a matrix. 

48 

49 Returns the Cholesky decomposition, :math:`A = L L^*` or 

50 :math:`A = U^* U` of a Hermitian positive-definite matrix A. 

51 

52 Parameters 

53 ---------- 

54 a : (M, M) array_like 

55 Matrix to be decomposed 

56 lower : bool, optional 

57 Whether to compute the upper- or lower-triangular Cholesky 

58 factorization. Default is upper-triangular. 

59 overwrite_a : bool, optional 

60 Whether to overwrite data in `a` (may improve performance). 

61 check_finite : bool, optional 

62 Whether to check that the input matrix contains only finite numbers. 

63 Disabling may give a performance gain, but may result in problems 

64 (crashes, non-termination) if the inputs do contain infinities or NaNs. 

65 

66 Returns 

67 ------- 

68 c : (M, M) ndarray 

69 Upper- or lower-triangular Cholesky factor of `a`. 

70 

71 Raises 

72 ------ 

73 LinAlgError : if decomposition fails. 

74 

75 Examples 

76 -------- 

77 >>> from scipy.linalg import cholesky 

78 >>> a = np.array([[1,-2j],[2j,5]]) 

79 >>> L = cholesky(a, lower=True) 

80 >>> L 

81 array([[ 1.+0.j, 0.+0.j], 

82 [ 0.+2.j, 1.+0.j]]) 

83 >>> L @ L.T.conj() 

84 array([[ 1.+0.j, 0.-2.j], 

85 [ 0.+2.j, 5.+0.j]]) 

86 

87 """ 

88 c, lower = _cholesky(a, lower=lower, overwrite_a=overwrite_a, clean=True, 

89 check_finite=check_finite) 

90 return c 

91 

92 

93def cho_factor(a, lower=False, overwrite_a=False, check_finite=True): 

94 """ 

95 Compute the Cholesky decomposition of a matrix, to use in cho_solve 

96 

97 Returns a matrix containing the Cholesky decomposition, 

98 ``A = L L*`` or ``A = U* U`` of a Hermitian positive-definite matrix `a`. 

99 The return value can be directly used as the first parameter to cho_solve. 

100 

101 .. warning:: 

102 The returned matrix also contains random data in the entries not 

103 used by the Cholesky decomposition. If you need to zero these 

104 entries, use the function `cholesky` instead. 

105 

106 Parameters 

107 ---------- 

108 a : (M, M) array_like 

109 Matrix to be decomposed 

110 lower : bool, optional 

111 Whether to compute the upper or lower triangular Cholesky factorization 

112 (Default: upper-triangular) 

113 overwrite_a : bool, optional 

114 Whether to overwrite data in a (may improve performance) 

115 check_finite : bool, optional 

116 Whether to check that the input matrix contains only finite numbers. 

117 Disabling may give a performance gain, but may result in problems 

118 (crashes, non-termination) if the inputs do contain infinities or NaNs. 

119 

120 Returns 

121 ------- 

122 c : (M, M) ndarray 

123 Matrix whose upper or lower triangle contains the Cholesky factor 

124 of `a`. Other parts of the matrix contain random data. 

125 lower : bool 

126 Flag indicating whether the factor is in the lower or upper triangle 

127 

128 Raises 

129 ------ 

130 LinAlgError 

131 Raised if decomposition fails. 

132 

133 See also 

134 -------- 

135 cho_solve : Solve a linear set equations using the Cholesky factorization 

136 of a matrix. 

137 

138 Examples 

139 -------- 

140 >>> from scipy.linalg import cho_factor 

141 >>> A = np.array([[9, 3, 1, 5], [3, 7, 5, 1], [1, 5, 9, 2], [5, 1, 2, 6]]) 

142 >>> c, low = cho_factor(A) 

143 >>> c 

144 array([[3. , 1. , 0.33333333, 1.66666667], 

145 [3. , 2.44948974, 1.90515869, -0.27216553], 

146 [1. , 5. , 2.29330749, 0.8559528 ], 

147 [5. , 1. , 2. , 1.55418563]]) 

148 >>> np.allclose(np.triu(c).T @ np. triu(c) - A, np.zeros((4, 4))) 

149 True 

150 

151 """ 

152 c, lower = _cholesky(a, lower=lower, overwrite_a=overwrite_a, clean=False, 

153 check_finite=check_finite) 

154 return c, lower 

155 

156 

157def cho_solve(c_and_lower, b, overwrite_b=False, check_finite=True): 

158 """Solve the linear equations A x = b, given the Cholesky factorization of A. 

159 

160 Parameters 

161 ---------- 

162 (c, lower) : tuple, (array, bool) 

163 Cholesky factorization of a, as given by cho_factor 

164 b : array 

165 Right-hand side 

166 overwrite_b : bool, optional 

167 Whether to overwrite data in b (may improve performance) 

168 check_finite : bool, optional 

169 Whether to check that the input matrices contain only finite numbers. 

170 Disabling may give a performance gain, but may result in problems 

171 (crashes, non-termination) if the inputs do contain infinities or NaNs. 

172 

173 Returns 

174 ------- 

175 x : array 

176 The solution to the system A x = b 

177 

178 See also 

179 -------- 

180 cho_factor : Cholesky factorization of a matrix 

181 

182 Examples 

183 -------- 

184 >>> from scipy.linalg import cho_factor, cho_solve 

185 >>> A = np.array([[9, 3, 1, 5], [3, 7, 5, 1], [1, 5, 9, 2], [5, 1, 2, 6]]) 

186 >>> c, low = cho_factor(A) 

187 >>> x = cho_solve((c, low), [1, 1, 1, 1]) 

188 >>> np.allclose(A @ x - [1, 1, 1, 1], np.zeros(4)) 

189 True 

190 

191 """ 

192 (c, lower) = c_and_lower 

193 if check_finite: 

194 b1 = asarray_chkfinite(b) 

195 c = asarray_chkfinite(c) 

196 else: 

197 b1 = asarray(b) 

198 c = asarray(c) 

199 if c.ndim != 2 or c.shape[0] != c.shape[1]: 

200 raise ValueError("The factored matrix c is not square.") 

201 if c.shape[1] != b1.shape[0]: 

202 raise ValueError("incompatible dimensions.") 

203 

204 overwrite_b = overwrite_b or _datacopied(b1, b) 

205 

206 potrs, = get_lapack_funcs(('potrs',), (c, b1)) 

207 x, info = potrs(c, b1, lower=lower, overwrite_b=overwrite_b) 

208 if info != 0: 

209 raise ValueError('illegal value in %dth argument of internal potrs' 

210 % -info) 

211 return x 

212 

213 

214def cholesky_banded(ab, overwrite_ab=False, lower=False, check_finite=True): 

215 """ 

216 Cholesky decompose a banded Hermitian positive-definite matrix 

217 

218 The matrix a is stored in ab either in lower-diagonal or upper- 

219 diagonal ordered form:: 

220 

221 ab[u + i - j, j] == a[i,j] (if upper form; i <= j) 

222 ab[ i - j, j] == a[i,j] (if lower form; i >= j) 

223 

224 Example of ab (shape of a is (6,6), u=2):: 

225 

226 upper form: 

227 * * a02 a13 a24 a35 

228 * a01 a12 a23 a34 a45 

229 a00 a11 a22 a33 a44 a55 

230 

231 lower form: 

232 a00 a11 a22 a33 a44 a55 

233 a10 a21 a32 a43 a54 * 

234 a20 a31 a42 a53 * * 

235 

236 Parameters 

237 ---------- 

238 ab : (u + 1, M) array_like 

239 Banded matrix 

240 overwrite_ab : bool, optional 

241 Discard data in ab (may enhance performance) 

242 lower : bool, optional 

243 Is the matrix in the lower form. (Default is upper form) 

244 check_finite : bool, optional 

245 Whether to check that the input matrix contains only finite numbers. 

246 Disabling may give a performance gain, but may result in problems 

247 (crashes, non-termination) if the inputs do contain infinities or NaNs. 

248 

249 Returns 

250 ------- 

251 c : (u + 1, M) ndarray 

252 Cholesky factorization of a, in the same banded format as ab 

253 

254 See also 

255 -------- 

256 cho_solve_banded : Solve a linear set equations, given the Cholesky factorization 

257 of a banded hermitian. 

258 

259 Examples 

260 -------- 

261 >>> from scipy.linalg import cholesky_banded 

262 >>> from numpy import allclose, zeros, diag 

263 >>> Ab = np.array([[0, 0, 1j, 2, 3j], [0, -1, -2, 3, 4], [9, 8, 7, 6, 9]]) 

264 >>> A = np.diag(Ab[0,2:], k=2) + np.diag(Ab[1,1:], k=1) 

265 >>> A = A + A.conj().T + np.diag(Ab[2, :]) 

266 >>> c = cholesky_banded(Ab) 

267 >>> C = np.diag(c[0, 2:], k=2) + np.diag(c[1, 1:], k=1) + np.diag(c[2, :]) 

268 >>> np.allclose(C.conj().T @ C - A, np.zeros((5, 5))) 

269 True 

270 

271 """ 

272 if check_finite: 

273 ab = asarray_chkfinite(ab) 

274 else: 

275 ab = asarray(ab) 

276 

277 pbtrf, = get_lapack_funcs(('pbtrf',), (ab,)) 

278 c, info = pbtrf(ab, lower=lower, overwrite_ab=overwrite_ab) 

279 if info > 0: 

280 raise LinAlgError("%d-th leading minor not positive definite" % info) 

281 if info < 0: 

282 raise ValueError('illegal value in %d-th argument of internal pbtrf' 

283 % -info) 

284 return c 

285 

286 

287def cho_solve_banded(cb_and_lower, b, overwrite_b=False, check_finite=True): 

288 """ 

289 Solve the linear equations ``A x = b``, given the Cholesky factorization of 

290 the banded hermitian ``A``. 

291 

292 Parameters 

293 ---------- 

294 (cb, lower) : tuple, (ndarray, bool) 

295 `cb` is the Cholesky factorization of A, as given by cholesky_banded. 

296 `lower` must be the same value that was given to cholesky_banded. 

297 b : array_like 

298 Right-hand side 

299 overwrite_b : bool, optional 

300 If True, the function will overwrite the values in `b`. 

301 check_finite : bool, optional 

302 Whether to check that the input matrices contain only finite numbers. 

303 Disabling may give a performance gain, but may result in problems 

304 (crashes, non-termination) if the inputs do contain infinities or NaNs. 

305 

306 Returns 

307 ------- 

308 x : array 

309 The solution to the system A x = b 

310 

311 See also 

312 -------- 

313 cholesky_banded : Cholesky factorization of a banded matrix 

314 

315 Notes 

316 ----- 

317 

318 .. versionadded:: 0.8.0 

319 

320 Examples 

321 -------- 

322 >>> from scipy.linalg import cholesky_banded, cho_solve_banded 

323 >>> Ab = np.array([[0, 0, 1j, 2, 3j], [0, -1, -2, 3, 4], [9, 8, 7, 6, 9]]) 

324 >>> A = np.diag(Ab[0,2:], k=2) + np.diag(Ab[1,1:], k=1) 

325 >>> A = A + A.conj().T + np.diag(Ab[2, :]) 

326 >>> c = cholesky_banded(Ab) 

327 >>> x = cho_solve_banded((c, False), np.ones(5)) 

328 >>> np.allclose(A @ x - np.ones(5), np.zeros(5)) 

329 True 

330 

331 """ 

332 (cb, lower) = cb_and_lower 

333 if check_finite: 

334 cb = asarray_chkfinite(cb) 

335 b = asarray_chkfinite(b) 

336 else: 

337 cb = asarray(cb) 

338 b = asarray(b) 

339 

340 # Validate shapes. 

341 if cb.shape[-1] != b.shape[0]: 

342 raise ValueError("shapes of cb and b are not compatible.") 

343 

344 pbtrs, = get_lapack_funcs(('pbtrs',), (cb, b)) 

345 x, info = pbtrs(cb, b, lower=lower, overwrite_b=overwrite_b) 

346 if info > 0: 

347 raise LinAlgError("%dth leading minor not positive definite" % info) 

348 if info < 0: 

349 raise ValueError('illegal value in %dth argument of internal pbtrs' 

350 % -info) 

351 return x