radau.py 19 KB

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