radau.py 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563
  1. from __future__ import division, print_function, absolute_import
  2. import numpy as np
  3. from scipy.linalg import lu_factor, lu_solve
  4. from scipy.sparse import csc_matrix, issparse, eye
  5. from scipy.sparse.linalg import splu
  6. from scipy.optimize._numdiff import group_columns
  7. from .common import (validate_max_step, validate_tol, select_initial_step,
  8. norm, num_jac, EPS, warn_extraneous,
  9. validate_first_step)
  10. from .base import OdeSolver, DenseOutput
  11. S6 = 6 ** 0.5
  12. # Butcher tableau. A is not used directly, see below.
  13. C = np.array([(4 - S6) / 10, (4 + S6) / 10, 1])
  14. E = np.array([-13 - 7 * S6, -13 + 7 * S6, -1]) / 3
  15. # Eigendecomposition of A is done: A = T L T**-1. There is 1 real eigenvalue
  16. # and a complex conjugate pair. They are written below.
  17. MU_REAL = 3 + 3 ** (2 / 3) - 3 ** (1 / 3)
  18. MU_COMPLEX = (3 + 0.5 * (3 ** (1 / 3) - 3 ** (2 / 3))
  19. - 0.5j * (3 ** (5 / 6) + 3 ** (7 / 6)))
  20. # These are transformation matrices.
  21. T = np.array([
  22. [0.09443876248897524, -0.14125529502095421, 0.03002919410514742],
  23. [0.25021312296533332, 0.20412935229379994, -0.38294211275726192],
  24. [1, 1, 0]])
  25. TI = np.array([
  26. [4.17871859155190428, 0.32768282076106237, 0.52337644549944951],
  27. [-4.17871859155190428, -0.32768282076106237, 0.47662355450055044],
  28. [0.50287263494578682, -2.57192694985560522, 0.59603920482822492]])
  29. # These linear combinations are used in the algorithm.
  30. TI_REAL = TI[0]
  31. TI_COMPLEX = TI[1] + 1j * TI[2]
  32. # Interpolator coefficients.
  33. P = np.array([
  34. [13/3 + 7*S6/3, -23/3 - 22*S6/3, 10/3 + 5 * S6],
  35. [13/3 - 7*S6/3, -23/3 + 22*S6/3, 10/3 - 5 * S6],
  36. [1/3, -8/3, 10/3]])
  37. NEWTON_MAXITER = 6 # Maximum number of Newton iterations.
  38. MIN_FACTOR = 0.2 # Minimum allowed decrease in a step size.
  39. MAX_FACTOR = 10 # Maximum allowed increase in a step size.
  40. def solve_collocation_system(fun, t, y, h, Z0, scale, tol,
  41. LU_real, LU_complex, solve_lu):
  42. """Solve the collocation system.
  43. Parameters
  44. ----------
  45. fun : callable
  46. Right-hand side of the system.
  47. t : float
  48. Current time.
  49. y : ndarray, shape (n,)
  50. Current state.
  51. h : float
  52. Step to try.
  53. Z0 : ndarray, shape (3, n)
  54. Initial guess for the solution. It determines new values of `y` at
  55. ``t + h * C`` as ``y + Z0``, where ``C`` is the Radau method constants.
  56. scale : float
  57. Problem tolerance scale, i.e. ``rtol * abs(y) + atol``.
  58. tol : float
  59. Tolerance to which solve the system. This value is compared with
  60. the normalized by `scale` error.
  61. LU_real, LU_complex
  62. LU decompositions of the system Jacobians.
  63. solve_lu : callable
  64. Callable which solves a linear system given a LU decomposition. The
  65. signature is ``solve_lu(LU, b)``.
  66. Returns
  67. -------
  68. converged : bool
  69. Whether iterations converged.
  70. n_iter : int
  71. Number of completed iterations.
  72. Z : ndarray, shape (3, n)
  73. Found solution.
  74. rate : float
  75. The rate of convergence.
  76. """
  77. n = y.shape[0]
  78. M_real = MU_REAL / h
  79. M_complex = MU_COMPLEX / h
  80. W = TI.dot(Z0)
  81. Z = Z0
  82. F = np.empty((3, n))
  83. ch = h * C
  84. dW_norm_old = None
  85. dW = np.empty_like(W)
  86. converged = False
  87. for k in range(NEWTON_MAXITER):
  88. for i in range(3):
  89. F[i] = fun(t + ch[i], y + Z[i])
  90. if not np.all(np.isfinite(F)):
  91. break
  92. f_real = F.T.dot(TI_REAL) - M_real * W[0]
  93. f_complex = F.T.dot(TI_COMPLEX) - M_complex * (W[1] + 1j * W[2])
  94. dW_real = solve_lu(LU_real, f_real)
  95. dW_complex = solve_lu(LU_complex, f_complex)
  96. dW[0] = dW_real
  97. dW[1] = dW_complex.real
  98. dW[2] = dW_complex.imag
  99. dW_norm = norm(dW / scale)
  100. if dW_norm_old is not None:
  101. rate = dW_norm / dW_norm_old
  102. else:
  103. rate = None
  104. if (rate is not None and (rate >= 1 or
  105. rate ** (NEWTON_MAXITER - k) / (1 - rate) * dW_norm > tol)):
  106. break
  107. W += dW
  108. Z = T.dot(W)
  109. if (dW_norm == 0 or
  110. rate is not None and rate / (1 - rate) * dW_norm < tol):
  111. converged = True
  112. break
  113. dW_norm_old = dW_norm
  114. return converged, k + 1, Z, rate
  115. def predict_factor(h_abs, h_abs_old, error_norm, error_norm_old):
  116. """Predict by which factor to increase/decrease the step size.
  117. The algorithm is described in [1]_.
  118. Parameters
  119. ----------
  120. h_abs, h_abs_old : float
  121. Current and previous values of the step size, `h_abs_old` can be None
  122. (see Notes).
  123. error_norm, error_norm_old : float
  124. Current and previous values of the error norm, `error_norm_old` can
  125. be None (see Notes).
  126. Returns
  127. -------
  128. factor : float
  129. Predicted factor.
  130. Notes
  131. -----
  132. If `h_abs_old` and `error_norm_old` are both not None then a two-step
  133. algorithm is used, otherwise a one-step algorithm is used.
  134. References
  135. ----------
  136. .. [1] E. Hairer, S. P. Norsett G. Wanner, "Solving Ordinary Differential
  137. Equations II: Stiff and Differential-Algebraic Problems", Sec. IV.8.
  138. """
  139. if error_norm_old is None or h_abs_old is None or error_norm == 0:
  140. multiplier = 1
  141. else:
  142. multiplier = h_abs / h_abs_old * (error_norm_old / error_norm) ** 0.25
  143. with np.errstate(divide='ignore'):
  144. factor = min(1, multiplier) * error_norm ** -0.25
  145. return factor
  146. class Radau(OdeSolver):
  147. """Implicit Runge-Kutta method of Radau IIA family of order 5.
  148. The implementation follows [1]_. The error is controlled with a
  149. third-order accurate embedded formula. A cubic polynomial which satisfies
  150. the collocation conditions is used for the dense output.
  151. Parameters
  152. ----------
  153. fun : callable
  154. Right-hand side of the system. The calling signature is ``fun(t, y)``.
  155. Here ``t`` is a scalar, and there are two options for the ndarray ``y``:
  156. It can either have shape (n,); then ``fun`` must return array_like with
  157. shape (n,). Alternatively it can have shape (n, k); then ``fun``
  158. must return an array_like with shape (n, k), i.e. each column
  159. corresponds to a single column in ``y``. The choice between the two
  160. options is determined by `vectorized` argument (see below). The
  161. vectorized implementation allows a faster approximation of the Jacobian
  162. by finite differences (required for this solver).
  163. t0 : float
  164. Initial time.
  165. y0 : array_like, shape (n,)
  166. Initial state.
  167. t_bound : float
  168. Boundary time - the integration won't continue beyond it. It also
  169. determines the direction of the integration.
  170. first_step : float or None, optional
  171. Initial step size. Default is ``None`` which means that the algorithm
  172. should choose.
  173. max_step : float, optional
  174. Maximum allowed step size. Default is np.inf, i.e. the step size is not
  175. bounded and determined solely by the solver.
  176. rtol, atol : float and array_like, optional
  177. Relative and absolute tolerances. The solver keeps the local error
  178. estimates less than ``atol + rtol * abs(y)``. Here `rtol` controls a
  179. relative accuracy (number of correct digits). But if a component of `y`
  180. is approximately below `atol`, the error only needs to fall within
  181. the same `atol` threshold, and the number of correct digits is not
  182. guaranteed. If components of y have different scales, it might be
  183. beneficial to set different `atol` values for different components by
  184. passing array_like with shape (n,) for `atol`. Default values are
  185. 1e-3 for `rtol` and 1e-6 for `atol`.
  186. jac : {None, array_like, sparse_matrix, callable}, optional
  187. Jacobian matrix of the right-hand side of the system with respect to
  188. y, required by this method. The Jacobian matrix has shape (n, n) and
  189. its element (i, j) is equal to ``d f_i / d y_j``.
  190. There are three ways to define the Jacobian:
  191. * If array_like or sparse_matrix, the Jacobian is assumed to
  192. be constant.
  193. * If callable, the Jacobian is assumed to depend on both
  194. t and y; it will be called as ``jac(t, y)`` as necessary.
  195. For the 'Radau' and 'BDF' methods, the return value might be a
  196. sparse matrix.
  197. * If None (default), the Jacobian will be approximated by
  198. finite differences.
  199. It is generally recommended to provide the Jacobian rather than
  200. relying on a finite-difference approximation.
  201. jac_sparsity : {None, array_like, sparse matrix}, optional
  202. Defines a sparsity structure of the Jacobian matrix for a
  203. finite-difference approximation. Its shape must be (n, n). This argument
  204. is ignored if `jac` is not `None`. If the Jacobian has only few non-zero
  205. elements in *each* row, providing the sparsity structure will greatly
  206. speed up the computations [2]_. A zero entry means that a corresponding
  207. element in the Jacobian is always zero. If None (default), the Jacobian
  208. is assumed to be dense.
  209. vectorized : bool, optional
  210. Whether `fun` is implemented in a vectorized fashion. Default is False.
  211. Attributes
  212. ----------
  213. n : int
  214. Number of equations.
  215. status : string
  216. Current status of the solver: 'running', 'finished' or 'failed'.
  217. t_bound : float
  218. Boundary time.
  219. direction : float
  220. Integration direction: +1 or -1.
  221. t : float
  222. Current time.
  223. y : ndarray
  224. Current state.
  225. t_old : float
  226. Previous time. None if no steps were made yet.
  227. step_size : float
  228. Size of the last successful step. None if no steps were made yet.
  229. nfev : int
  230. Number of evaluations of the right-hand side.
  231. njev : int
  232. Number of evaluations of the Jacobian.
  233. nlu : int
  234. Number of LU decompositions.
  235. References
  236. ----------
  237. .. [1] E. Hairer, G. Wanner, "Solving Ordinary Differential Equations II:
  238. Stiff and Differential-Algebraic Problems", Sec. IV.8.
  239. .. [2] A. Curtis, M. J. D. Powell, and J. Reid, "On the estimation of
  240. sparse Jacobian matrices", Journal of the Institute of Mathematics
  241. and its Applications, 13, pp. 117-120, 1974.
  242. """
  243. def __init__(self, fun, t0, y0, t_bound, max_step=np.inf,
  244. rtol=1e-3, atol=1e-6, jac=None, jac_sparsity=None,
  245. vectorized=False, first_step=None, **extraneous):
  246. warn_extraneous(extraneous)
  247. super(Radau, self).__init__(fun, t0, y0, t_bound, vectorized)
  248. self.y_old = None
  249. self.max_step = validate_max_step(max_step)
  250. self.rtol, self.atol = validate_tol(rtol, atol, self.n)
  251. self.f = self.fun(self.t, self.y)
  252. # Select initial step assuming the same order which is used to control
  253. # the error.
  254. if first_step is None:
  255. self.h_abs = select_initial_step(
  256. self.fun, self.t, self.y, self.f, self.direction,
  257. 3, self.rtol, self.atol)
  258. else:
  259. self.h_abs = validate_first_step(first_step, t0, t_bound)
  260. self.h_abs_old = None
  261. self.error_norm_old = None
  262. self.newton_tol = max(10 * EPS / rtol, min(0.03, rtol ** 0.5))
  263. self.sol = None
  264. self.jac_factor = None
  265. self.jac, self.J = self._validate_jac(jac, jac_sparsity)
  266. if issparse(self.J):
  267. def lu(A):
  268. self.nlu += 1
  269. return splu(A)
  270. def solve_lu(LU, b):
  271. return LU.solve(b)
  272. I = eye(self.n, format='csc')
  273. else:
  274. def lu(A):
  275. self.nlu += 1
  276. return lu_factor(A, overwrite_a=True)
  277. def solve_lu(LU, b):
  278. return lu_solve(LU, b, overwrite_b=True)
  279. I = np.identity(self.n)
  280. self.lu = lu
  281. self.solve_lu = solve_lu
  282. self.I = I
  283. self.current_jac = True
  284. self.LU_real = None
  285. self.LU_complex = None
  286. self.Z = None
  287. def _validate_jac(self, jac, sparsity):
  288. t0 = self.t
  289. y0 = self.y
  290. if jac is None:
  291. if sparsity is not None:
  292. if issparse(sparsity):
  293. sparsity = csc_matrix(sparsity)
  294. groups = group_columns(sparsity)
  295. sparsity = (sparsity, groups)
  296. def jac_wrapped(t, y, f):
  297. self.njev += 1
  298. J, self.jac_factor = num_jac(self.fun_vectorized, t, y, f,
  299. self.atol, self.jac_factor,
  300. sparsity)
  301. return J
  302. J = jac_wrapped(t0, y0, self.f)
  303. elif callable(jac):
  304. J = jac(t0, y0)
  305. self.njev = 1
  306. if issparse(J):
  307. J = csc_matrix(J)
  308. def jac_wrapped(t, y, _=None):
  309. self.njev += 1
  310. return csc_matrix(jac(t, y), dtype=float)
  311. else:
  312. J = np.asarray(J, dtype=float)
  313. def jac_wrapped(t, y, _=None):
  314. self.njev += 1
  315. return np.asarray(jac(t, y), dtype=float)
  316. if J.shape != (self.n, self.n):
  317. raise ValueError("`jac` is expected to have shape {}, but "
  318. "actually has {}."
  319. .format((self.n, self.n), J.shape))
  320. else:
  321. if issparse(jac):
  322. J = csc_matrix(jac)
  323. else:
  324. J = np.asarray(jac, dtype=float)
  325. if J.shape != (self.n, self.n):
  326. raise ValueError("`jac` is expected to have shape {}, but "
  327. "actually has {}."
  328. .format((self.n, self.n), J.shape))
  329. jac_wrapped = None
  330. return jac_wrapped, J
  331. def _step_impl(self):
  332. t = self.t
  333. y = self.y
  334. f = self.f
  335. max_step = self.max_step
  336. atol = self.atol
  337. rtol = self.rtol
  338. min_step = 10 * np.abs(np.nextafter(t, self.direction * np.inf) - t)
  339. if self.h_abs > max_step:
  340. h_abs = max_step
  341. h_abs_old = None
  342. error_norm_old = None
  343. elif self.h_abs < min_step:
  344. h_abs = min_step
  345. h_abs_old = None
  346. error_norm_old = None
  347. else:
  348. h_abs = self.h_abs
  349. h_abs_old = self.h_abs_old
  350. error_norm_old = self.error_norm_old
  351. J = self.J
  352. LU_real = self.LU_real
  353. LU_complex = self.LU_complex
  354. current_jac = self.current_jac
  355. jac = self.jac
  356. rejected = False
  357. step_accepted = False
  358. message = None
  359. while not step_accepted:
  360. if h_abs < min_step:
  361. return False, self.TOO_SMALL_STEP
  362. h = h_abs * self.direction
  363. t_new = t + h
  364. if self.direction * (t_new - self.t_bound) > 0:
  365. t_new = self.t_bound
  366. h = t_new - t
  367. h_abs = np.abs(h)
  368. if self.sol is None:
  369. Z0 = np.zeros((3, y.shape[0]))
  370. else:
  371. Z0 = self.sol(t + h * C).T - y
  372. scale = atol + np.abs(y) * rtol
  373. converged = False
  374. while not converged:
  375. if LU_real is None or LU_complex is None:
  376. LU_real = self.lu(MU_REAL / h * self.I - J)
  377. LU_complex = self.lu(MU_COMPLEX / h * self.I - J)
  378. converged, n_iter, Z, rate = solve_collocation_system(
  379. self.fun, t, y, h, Z0, scale, self.newton_tol,
  380. LU_real, LU_complex, self.solve_lu)
  381. if not converged:
  382. if current_jac:
  383. break
  384. J = self.jac(t, y, f)
  385. current_jac = True
  386. LU_real = None
  387. LU_complex = None
  388. if not converged:
  389. h_abs *= 0.5
  390. LU_real = None
  391. LU_complex = None
  392. continue
  393. y_new = y + Z[-1]
  394. ZE = Z.T.dot(E) / h
  395. error = self.solve_lu(LU_real, f + ZE)
  396. scale = atol + np.maximum(np.abs(y), np.abs(y_new)) * rtol
  397. error_norm = norm(error / scale)
  398. safety = 0.9 * (2 * NEWTON_MAXITER + 1) / (2 * NEWTON_MAXITER
  399. + n_iter)
  400. if rejected and error_norm > 1:
  401. error = self.solve_lu(LU_real, self.fun(t, y + error) + ZE)
  402. error_norm = norm(error / scale)
  403. if error_norm > 1:
  404. factor = predict_factor(h_abs, h_abs_old,
  405. error_norm, error_norm_old)
  406. h_abs *= max(MIN_FACTOR, safety * factor)
  407. LU_real = None
  408. LU_complex = None
  409. rejected = True
  410. else:
  411. step_accepted = True
  412. recompute_jac = jac is not None and n_iter > 2 and rate > 1e-3
  413. factor = predict_factor(h_abs, h_abs_old, error_norm, error_norm_old)
  414. factor = min(MAX_FACTOR, safety * factor)
  415. if not recompute_jac and factor < 1.2:
  416. factor = 1
  417. else:
  418. LU_real = None
  419. LU_complex = None
  420. f_new = self.fun(t_new, y_new)
  421. if recompute_jac:
  422. J = jac(t_new, y_new, f_new)
  423. current_jac = True
  424. elif jac is not None:
  425. current_jac = False
  426. self.h_abs_old = self.h_abs
  427. self.error_norm_old = error_norm
  428. self.h_abs = h_abs * factor
  429. self.y_old = y
  430. self.t = t_new
  431. self.y = y_new
  432. self.f = f_new
  433. self.Z = Z
  434. self.LU_real = LU_real
  435. self.LU_complex = LU_complex
  436. self.current_jac = current_jac
  437. self.J = J
  438. self.t_old = t
  439. self.sol = self._compute_dense_output()
  440. return step_accepted, message
  441. def _compute_dense_output(self):
  442. Q = np.dot(self.Z.T, P)
  443. return RadauDenseOutput(self.t_old, self.t, self.y_old, Q)
  444. def _dense_output_impl(self):
  445. return self.sol
  446. class RadauDenseOutput(DenseOutput):
  447. def __init__(self, t_old, t, y_old, Q):
  448. super(RadauDenseOutput, self).__init__(t_old, t)
  449. self.h = t - t_old
  450. self.Q = Q
  451. self.order = Q.shape[1] - 1
  452. self.y_old = y_old
  453. def _call_impl(self, t):
  454. x = (t - self.t_old) / self.h
  455. if t.ndim == 0:
  456. p = np.tile(x, self.order + 1)
  457. p = np.cumprod(p)
  458. else:
  459. p = np.tile(x, (self.order + 1, 1))
  460. p = np.cumprod(p, axis=0)
  461. # Here we don't multiply by h, not a mistake.
  462. y = np.dot(self.Q, p)
  463. if y.ndim == 2:
  464. y += self.y_old[:, None]
  465. else:
  466. y += self.y_old
  467. return y