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

1import numpy as np 

2from scipy.linalg import lu_factor, lu_solve 

3from scipy.sparse import issparse, csc_matrix, eye 

4from scipy.sparse.linalg import splu 

5from scipy.optimize._numdiff import group_columns 

6from .common import (validate_max_step, validate_tol, select_initial_step, 

7 norm, EPS, num_jac, validate_first_step, 

8 warn_extraneous) 

9from .base import OdeSolver, DenseOutput 

10 

11 

12MAX_ORDER = 5 

13NEWTON_MAXITER = 4 

14MIN_FACTOR = 0.2 

15MAX_FACTOR = 10 

16 

17 

18def compute_R(order, factor): 

19 """Compute the matrix for changing the differences array.""" 

20 I = np.arange(1, order + 1)[:, None] 

21 J = np.arange(1, order + 1) 

22 M = np.zeros((order + 1, order + 1)) 

23 M[1:, 1:] = (I - 1 - factor * J) / I 

24 M[0] = 1 

25 return np.cumprod(M, axis=0) 

26 

27 

28def change_D(D, order, factor): 

29 """Change differences array in-place when step size is changed.""" 

30 R = compute_R(order, factor) 

31 U = compute_R(order, 1) 

32 RU = R.dot(U) 

33 D[:order + 1] = np.dot(RU.T, D[:order + 1]) 

34 

35 

36def solve_bdf_system(fun, t_new, y_predict, c, psi, LU, solve_lu, scale, tol): 

37 """Solve the algebraic system resulting from BDF method.""" 

38 d = 0 

39 y = y_predict.copy() 

40 dy_norm_old = None 

41 converged = False 

42 for k in range(NEWTON_MAXITER): 

43 f = fun(t_new, y) 

44 if not np.all(np.isfinite(f)): 

45 break 

46 

47 dy = solve_lu(LU, c * f - psi - d) 

48 dy_norm = norm(dy / scale) 

49 

50 if dy_norm_old is None: 

51 rate = None 

52 else: 

53 rate = dy_norm / dy_norm_old 

54 

55 if (rate is not None and (rate >= 1 or 

56 rate ** (NEWTON_MAXITER - k) / (1 - rate) * dy_norm > tol)): 

57 break 

58 

59 y += dy 

60 d += dy 

61 

62 if (dy_norm == 0 or 

63 rate is not None and rate / (1 - rate) * dy_norm < tol): 

64 converged = True 

65 break 

66 

67 dy_norm_old = dy_norm 

68 

69 return converged, k + 1, y, d 

70 

71 

72class BDF(OdeSolver): 

73 """Implicit method based on backward-differentiation formulas. 

74 

75 This is a variable order method with the order varying automatically from 

76 1 to 5. The general framework of the BDF algorithm is described in [1]_. 

77 This class implements a quasi-constant step size as explained in [2]_. 

78 The error estimation strategy for the constant-step BDF is derived in [3]_. 

79 An accuracy enhancement using modified formulas (NDF) [2]_ is also implemented. 

80 

81 Can be applied in the complex domain. 

82 

83 Parameters 

84 ---------- 

85 fun : callable 

86 Right-hand side of the system. The calling signature is ``fun(t, y)``. 

87 Here ``t`` is a scalar, and there are two options for the ndarray ``y``: 

88 It can either have shape (n,); then ``fun`` must return array_like with 

89 shape (n,). Alternatively it can have shape (n, k); then ``fun`` 

90 must return an array_like with shape (n, k), i.e. each column 

91 corresponds to a single column in ``y``. The choice between the two 

92 options is determined by `vectorized` argument (see below). The 

93 vectorized implementation allows a faster approximation of the Jacobian 

94 by finite differences (required for this solver). 

95 t0 : float 

96 Initial time. 

97 y0 : array_like, shape (n,) 

98 Initial state. 

99 t_bound : float 

100 Boundary time - the integration won't continue beyond it. It also 

101 determines the direction of the integration. 

102 first_step : float or None, optional 

103 Initial step size. Default is ``None`` which means that the algorithm 

104 should choose. 

105 max_step : float, optional 

106 Maximum allowed step size. Default is np.inf, i.e., the step size is not 

107 bounded and determined solely by the solver. 

108 rtol, atol : float and array_like, optional 

109 Relative and absolute tolerances. The solver keeps the local error 

110 estimates less than ``atol + rtol * abs(y)``. Here `rtol` controls a 

111 relative accuracy (number of correct digits). But if a component of `y` 

112 is approximately below `atol`, the error only needs to fall within 

113 the same `atol` threshold, and the number of correct digits is not 

114 guaranteed. If components of y have different scales, it might be 

115 beneficial to set different `atol` values for different components by 

116 passing array_like with shape (n,) for `atol`. Default values are 

117 1e-3 for `rtol` and 1e-6 for `atol`. 

118 jac : {None, array_like, sparse_matrix, callable}, optional 

119 Jacobian matrix of the right-hand side of the system with respect to y, 

120 required by this method. The Jacobian matrix has shape (n, n) and its 

121 element (i, j) is equal to ``d f_i / d y_j``. 

122 There are three ways to define the Jacobian: 

123 

124 * If array_like or sparse_matrix, the Jacobian is assumed to 

125 be constant. 

126 * If callable, the Jacobian is assumed to depend on both 

127 t and y; it will be called as ``jac(t, y)`` as necessary. 

128 For the 'Radau' and 'BDF' methods, the return value might be a 

129 sparse matrix. 

130 * If None (default), the Jacobian will be approximated by 

131 finite differences. 

132 

133 It is generally recommended to provide the Jacobian rather than 

134 relying on a finite-difference approximation. 

135 jac_sparsity : {None, array_like, sparse matrix}, optional 

136 Defines a sparsity structure of the Jacobian matrix for a 

137 finite-difference approximation. Its shape must be (n, n). This argument 

138 is ignored if `jac` is not `None`. If the Jacobian has only few non-zero 

139 elements in *each* row, providing the sparsity structure will greatly 

140 speed up the computations [4]_. A zero entry means that a corresponding 

141 element in the Jacobian is always zero. If None (default), the Jacobian 

142 is assumed to be dense. 

143 vectorized : bool, optional 

144 Whether `fun` is implemented in a vectorized fashion. Default is False. 

145 

146 Attributes 

147 ---------- 

148 n : int 

149 Number of equations. 

150 status : string 

151 Current status of the solver: 'running', 'finished' or 'failed'. 

152 t_bound : float 

153 Boundary time. 

154 direction : float 

155 Integration direction: +1 or -1. 

156 t : float 

157 Current time. 

158 y : ndarray 

159 Current state. 

160 t_old : float 

161 Previous time. None if no steps were made yet. 

162 step_size : float 

163 Size of the last successful step. None if no steps were made yet. 

164 nfev : int 

165 Number of evaluations of the right-hand side. 

166 njev : int 

167 Number of evaluations of the Jacobian. 

168 nlu : int 

169 Number of LU decompositions. 

170 

171 References 

172 ---------- 

173 .. [1] G. D. Byrne, A. C. Hindmarsh, "A Polyalgorithm for the Numerical 

174 Solution of Ordinary Differential Equations", ACM Transactions on 

175 Mathematical Software, Vol. 1, No. 1, pp. 71-96, March 1975. 

176 .. [2] L. F. Shampine, M. W. Reichelt, "THE MATLAB ODE SUITE", SIAM J. SCI. 

177 COMPUTE., Vol. 18, No. 1, pp. 1-22, January 1997. 

178 .. [3] E. Hairer, G. Wanner, "Solving Ordinary Differential Equations I: 

179 Nonstiff Problems", Sec. III.2. 

180 .. [4] A. Curtis, M. J. D. Powell, and J. Reid, "On the estimation of 

181 sparse Jacobian matrices", Journal of the Institute of Mathematics 

182 and its Applications, 13, pp. 117-120, 1974. 

183 """ 

184 def __init__(self, fun, t0, y0, t_bound, max_step=np.inf, 

185 rtol=1e-3, atol=1e-6, jac=None, jac_sparsity=None, 

186 vectorized=False, first_step=None, **extraneous): 

187 warn_extraneous(extraneous) 

188 super(BDF, self).__init__(fun, t0, y0, t_bound, vectorized, 

189 support_complex=True) 

190 self.max_step = validate_max_step(max_step) 

191 self.rtol, self.atol = validate_tol(rtol, atol, self.n) 

192 f = self.fun(self.t, self.y) 

193 if first_step is None: 

194 self.h_abs = select_initial_step(self.fun, self.t, self.y, f, 

195 self.direction, 1, 

196 self.rtol, self.atol) 

197 else: 

198 self.h_abs = validate_first_step(first_step, t0, t_bound) 

199 self.h_abs_old = None 

200 self.error_norm_old = None 

201 

202 self.newton_tol = max(10 * EPS / rtol, min(0.03, rtol ** 0.5)) 

203 

204 self.jac_factor = None 

205 self.jac, self.J = self._validate_jac(jac, jac_sparsity) 

206 if issparse(self.J): 

207 def lu(A): 

208 self.nlu += 1 

209 return splu(A) 

210 

211 def solve_lu(LU, b): 

212 return LU.solve(b) 

213 

214 I = eye(self.n, format='csc', dtype=self.y.dtype) 

215 else: 

216 def lu(A): 

217 self.nlu += 1 

218 return lu_factor(A, overwrite_a=True) 

219 

220 def solve_lu(LU, b): 

221 return lu_solve(LU, b, overwrite_b=True) 

222 

223 I = np.identity(self.n, dtype=self.y.dtype) 

224 

225 self.lu = lu 

226 self.solve_lu = solve_lu 

227 self.I = I 

228 

229 kappa = np.array([0, -0.1850, -1/9, -0.0823, -0.0415, 0]) 

230 self.gamma = np.hstack((0, np.cumsum(1 / np.arange(1, MAX_ORDER + 1)))) 

231 self.alpha = (1 - kappa) * self.gamma 

232 self.error_const = kappa * self.gamma + 1 / np.arange(1, MAX_ORDER + 2) 

233 

234 D = np.empty((MAX_ORDER + 3, self.n), dtype=self.y.dtype) 

235 D[0] = self.y 

236 D[1] = f * self.h_abs * self.direction 

237 self.D = D 

238 

239 self.order = 1 

240 self.n_equal_steps = 0 

241 self.LU = None 

242 

243 def _validate_jac(self, jac, sparsity): 

244 t0 = self.t 

245 y0 = self.y 

246 

247 if jac is None: 

248 if sparsity is not None: 

249 if issparse(sparsity): 

250 sparsity = csc_matrix(sparsity) 

251 groups = group_columns(sparsity) 

252 sparsity = (sparsity, groups) 

253 

254 def jac_wrapped(t, y): 

255 self.njev += 1 

256 f = self.fun_single(t, y) 

257 J, self.jac_factor = num_jac(self.fun_vectorized, t, y, f, 

258 self.atol, self.jac_factor, 

259 sparsity) 

260 return J 

261 J = jac_wrapped(t0, y0) 

262 elif callable(jac): 

263 J = jac(t0, y0) 

264 self.njev += 1 

265 if issparse(J): 

266 J = csc_matrix(J, dtype=y0.dtype) 

267 

268 def jac_wrapped(t, y): 

269 self.njev += 1 

270 return csc_matrix(jac(t, y), dtype=y0.dtype) 

271 else: 

272 J = np.asarray(J, dtype=y0.dtype) 

273 

274 def jac_wrapped(t, y): 

275 self.njev += 1 

276 return np.asarray(jac(t, y), dtype=y0.dtype) 

277 

278 if J.shape != (self.n, self.n): 

279 raise ValueError("`jac` is expected to have shape {}, but " 

280 "actually has {}." 

281 .format((self.n, self.n), J.shape)) 

282 else: 

283 if issparse(jac): 

284 J = csc_matrix(jac, dtype=y0.dtype) 

285 else: 

286 J = np.asarray(jac, dtype=y0.dtype) 

287 

288 if J.shape != (self.n, self.n): 

289 raise ValueError("`jac` is expected to have shape {}, but " 

290 "actually has {}." 

291 .format((self.n, self.n), J.shape)) 

292 jac_wrapped = None 

293 

294 return jac_wrapped, J 

295 

296 def _step_impl(self): 

297 t = self.t 

298 D = self.D 

299 

300 max_step = self.max_step 

301 min_step = 10 * np.abs(np.nextafter(t, self.direction * np.inf) - t) 

302 if self.h_abs > max_step: 

303 h_abs = max_step 

304 change_D(D, self.order, max_step / self.h_abs) 

305 self.n_equal_steps = 0 

306 elif self.h_abs < min_step: 

307 h_abs = min_step 

308 change_D(D, self.order, min_step / self.h_abs) 

309 self.n_equal_steps = 0 

310 else: 

311 h_abs = self.h_abs 

312 

313 atol = self.atol 

314 rtol = self.rtol 

315 order = self.order 

316 

317 alpha = self.alpha 

318 gamma = self.gamma 

319 error_const = self.error_const 

320 

321 J = self.J 

322 LU = self.LU 

323 current_jac = self.jac is None 

324 

325 step_accepted = False 

326 while not step_accepted: 

327 if h_abs < min_step: 

328 return False, self.TOO_SMALL_STEP 

329 

330 h = h_abs * self.direction 

331 t_new = t + h 

332 

333 if self.direction * (t_new - self.t_bound) > 0: 

334 t_new = self.t_bound 

335 change_D(D, order, np.abs(t_new - t) / h_abs) 

336 self.n_equal_steps = 0 

337 LU = None 

338 

339 h = t_new - t 

340 h_abs = np.abs(h) 

341 

342 y_predict = np.sum(D[:order + 1], axis=0) 

343 

344 scale = atol + rtol * np.abs(y_predict) 

345 psi = np.dot(D[1: order + 1].T, gamma[1: order + 1]) / alpha[order] 

346 

347 converged = False 

348 c = h / alpha[order] 

349 while not converged: 

350 if LU is None: 

351 LU = self.lu(self.I - c * J) 

352 

353 converged, n_iter, y_new, d = solve_bdf_system( 

354 self.fun, t_new, y_predict, c, psi, LU, self.solve_lu, 

355 scale, self.newton_tol) 

356 

357 if not converged: 

358 if current_jac: 

359 break 

360 J = self.jac(t_new, y_predict) 

361 LU = None 

362 current_jac = True 

363 

364 if not converged: 

365 factor = 0.5 

366 h_abs *= factor 

367 change_D(D, order, factor) 

368 self.n_equal_steps = 0 

369 LU = None 

370 continue 

371 

372 safety = 0.9 * (2 * NEWTON_MAXITER + 1) / (2 * NEWTON_MAXITER 

373 + n_iter) 

374 

375 scale = atol + rtol * np.abs(y_new) 

376 error = error_const[order] * d 

377 error_norm = norm(error / scale) 

378 

379 if error_norm > 1: 

380 factor = max(MIN_FACTOR, 

381 safety * error_norm ** (-1 / (order + 1))) 

382 h_abs *= factor 

383 change_D(D, order, factor) 

384 self.n_equal_steps = 0 

385 # As we didn't have problems with convergence, we don't 

386 # reset LU here. 

387 else: 

388 step_accepted = True 

389 

390 self.n_equal_steps += 1 

391 

392 self.t = t_new 

393 self.y = y_new 

394 

395 self.h_abs = h_abs 

396 self.J = J 

397 self.LU = LU 

398 

399 # Update differences. The principal relation here is 

400 # D^{j + 1} y_n = D^{j} y_n - D^{j} y_{n - 1}. Keep in mind that D 

401 # contained difference for previous interpolating polynomial and 

402 # d = D^{k + 1} y_n. Thus this elegant code follows. 

403 D[order + 2] = d - D[order + 1] 

404 D[order + 1] = d 

405 for i in reversed(range(order + 1)): 

406 D[i] += D[i + 1] 

407 

408 if self.n_equal_steps < order + 1: 

409 return True, None 

410 

411 if order > 1: 

412 error_m = error_const[order - 1] * D[order] 

413 error_m_norm = norm(error_m / scale) 

414 else: 

415 error_m_norm = np.inf 

416 

417 if order < MAX_ORDER: 

418 error_p = error_const[order + 1] * D[order + 2] 

419 error_p_norm = norm(error_p / scale) 

420 else: 

421 error_p_norm = np.inf 

422 

423 error_norms = np.array([error_m_norm, error_norm, error_p_norm]) 

424 with np.errstate(divide='ignore'): 

425 factors = error_norms ** (-1 / np.arange(order, order + 3)) 

426 

427 delta_order = np.argmax(factors) - 1 

428 order += delta_order 

429 self.order = order 

430 

431 factor = min(MAX_FACTOR, safety * np.max(factors)) 

432 self.h_abs *= factor 

433 change_D(D, order, factor) 

434 self.n_equal_steps = 0 

435 self.LU = None 

436 

437 return True, None 

438 

439 def _dense_output_impl(self): 

440 return BdfDenseOutput(self.t_old, self.t, self.h_abs * self.direction, 

441 self.order, self.D[:self.order + 1].copy()) 

442 

443 

444class BdfDenseOutput(DenseOutput): 

445 def __init__(self, t_old, t, h, order, D): 

446 super(BdfDenseOutput, self).__init__(t_old, t) 

447 self.order = order 

448 self.t_shift = self.t - h * np.arange(self.order) 

449 self.denom = h * (1 + np.arange(self.order)) 

450 self.D = D 

451 

452 def _call_impl(self, t): 

453 if t.ndim == 0: 

454 x = (t - self.t_shift) / self.denom 

455 p = np.cumprod(x) 

456 else: 

457 x = (t - self.t_shift[:, None]) / self.denom[:, None] 

458 p = np.cumprod(x, axis=0) 

459 

460 y = np.dot(self.D[1:].T, p) 

461 if y.ndim == 1: 

462 y += self.D[0] 

463 else: 

464 y += self.D[0, :, None] 

465 

466 return y