common.py 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451
  1. from itertools import groupby
  2. from warnings import warn
  3. import numpy as np
  4. from scipy.sparse import find, coo_matrix
  5. EPS = np.finfo(float).eps
  6. def validate_first_step(first_step, t0, t_bound):
  7. """Assert that first_step is valid and return it."""
  8. if first_step <= 0:
  9. raise ValueError("`first_step` must be positive.")
  10. if first_step > np.abs(t_bound - t0):
  11. raise ValueError("`first_step` exceeds bounds.")
  12. return first_step
  13. def validate_max_step(max_step):
  14. """Assert that max_Step is valid and return it."""
  15. if max_step <= 0:
  16. raise ValueError("`max_step` must be positive.")
  17. return max_step
  18. def warn_extraneous(extraneous):
  19. """Display a warning for extraneous keyword arguments.
  20. The initializer of each solver class is expected to collect keyword
  21. arguments that it doesn't understand and warn about them. This function
  22. prints a warning for each key in the supplied dictionary.
  23. Parameters
  24. ----------
  25. extraneous : dict
  26. Extraneous keyword arguments
  27. """
  28. if extraneous:
  29. warn("The following arguments have no effect for a chosen solver: "
  30. f"{', '.join(f'`{x}`' for x in extraneous)}.",
  31. stacklevel=3)
  32. def validate_tol(rtol, atol, n):
  33. """Validate tolerance values."""
  34. if np.any(rtol < 100 * EPS):
  35. warn("At least one element of `rtol` is too small. "
  36. f"Setting `rtol = np.maximum(rtol, {100 * EPS})`.",
  37. stacklevel=3)
  38. rtol = np.maximum(rtol, 100 * EPS)
  39. atol = np.asarray(atol)
  40. if atol.ndim > 0 and atol.shape != (n,):
  41. raise ValueError("`atol` has wrong shape.")
  42. if np.any(atol < 0):
  43. raise ValueError("`atol` must be positive.")
  44. return rtol, atol
  45. def norm(x):
  46. """Compute RMS norm."""
  47. return np.linalg.norm(x) / x.size ** 0.5
  48. def select_initial_step(fun, t0, y0, t_bound,
  49. max_step, f0, direction, order, rtol, atol):
  50. """Empirically select a good initial step.
  51. The algorithm is described in [1]_.
  52. Parameters
  53. ----------
  54. fun : callable
  55. Right-hand side of the system.
  56. t0 : float
  57. Initial value of the independent variable.
  58. y0 : ndarray, shape (n,)
  59. Initial value of the dependent variable.
  60. t_bound : float
  61. End-point of integration interval; used to ensure that t0+step<=tbound
  62. and that fun is only evaluated in the interval [t0,tbound]
  63. max_step : float
  64. Maximum allowable step size.
  65. f0 : ndarray, shape (n,)
  66. Initial value of the derivative, i.e., ``fun(t0, y0)``.
  67. direction : float
  68. Integration direction.
  69. order : float
  70. Error estimator order. It means that the error controlled by the
  71. algorithm is proportional to ``step_size ** (order + 1)`.
  72. rtol : float
  73. Desired relative tolerance.
  74. atol : float
  75. Desired absolute tolerance.
  76. Returns
  77. -------
  78. h_abs : float
  79. Absolute value of the suggested initial step.
  80. References
  81. ----------
  82. .. [1] E. Hairer, S. P. Norsett G. Wanner, "Solving Ordinary Differential
  83. Equations I: Nonstiff Problems", Sec. II.4.
  84. """
  85. if y0.size == 0:
  86. return np.inf
  87. interval_length = abs(t_bound - t0)
  88. if interval_length == 0.0:
  89. return 0.0
  90. scale = atol + np.abs(y0) * rtol
  91. d0 = norm(y0 / scale)
  92. d1 = norm(f0 / scale)
  93. if d0 < 1e-5 or d1 < 1e-5:
  94. h0 = 1e-6
  95. else:
  96. h0 = 0.01 * d0 / d1
  97. # Check t0+h0*direction doesn't take us beyond t_bound
  98. h0 = min(h0, interval_length)
  99. y1 = y0 + h0 * direction * f0
  100. f1 = fun(t0 + h0 * direction, y1)
  101. d2 = norm((f1 - f0) / scale) / h0
  102. if d1 <= 1e-15 and d2 <= 1e-15:
  103. h1 = max(1e-6, h0 * 1e-3)
  104. else:
  105. h1 = (0.01 / max(d1, d2)) ** (1 / (order + 1))
  106. return min(100 * h0, h1, interval_length, max_step)
  107. class OdeSolution:
  108. """Continuous ODE solution.
  109. It is organized as a collection of `DenseOutput` objects which represent
  110. local interpolants. It provides an algorithm to select a right interpolant
  111. for each given point.
  112. The interpolants cover the range between `t_min` and `t_max` (see
  113. Attributes below). Evaluation outside this interval is not forbidden, but
  114. the accuracy is not guaranteed.
  115. When evaluating at a breakpoint (one of the values in `ts`) a segment with
  116. the lower index is selected.
  117. Parameters
  118. ----------
  119. ts : array_like, shape (n_segments + 1,)
  120. Time instants between which local interpolants are defined. Must
  121. be strictly increasing or decreasing (zero segment with two points is
  122. also allowed).
  123. interpolants : list of DenseOutput with n_segments elements
  124. Local interpolants. An i-th interpolant is assumed to be defined
  125. between ``ts[i]`` and ``ts[i + 1]``.
  126. alt_segment : boolean
  127. Requests the alternative interpolant segment selection scheme. At each
  128. solver integration point, two interpolant segments are available. The
  129. default (False) and alternative (True) behaviours select the segment
  130. for which the requested time corresponded to ``t`` and ``t_old``,
  131. respectively. This functionality is only relevant for testing the
  132. interpolants' accuracy: different integrators use different
  133. construction strategies.
  134. Attributes
  135. ----------
  136. t_min, t_max : float
  137. Time range of the interpolation.
  138. """
  139. def __init__(self, ts, interpolants, alt_segment=False):
  140. ts = np.asarray(ts)
  141. d = np.diff(ts)
  142. # The first case covers integration on zero segment.
  143. if not ((ts.size == 2 and ts[0] == ts[-1])
  144. or np.all(d > 0) or np.all(d < 0)):
  145. raise ValueError("`ts` must be strictly increasing or decreasing.")
  146. self.n_segments = len(interpolants)
  147. if ts.shape != (self.n_segments + 1,):
  148. raise ValueError("Numbers of time stamps and interpolants "
  149. "don't match.")
  150. self.ts = ts
  151. self.interpolants = interpolants
  152. if ts[-1] >= ts[0]:
  153. self.t_min = ts[0]
  154. self.t_max = ts[-1]
  155. self.ascending = True
  156. self.side = "right" if alt_segment else "left"
  157. self.ts_sorted = ts
  158. else:
  159. self.t_min = ts[-1]
  160. self.t_max = ts[0]
  161. self.ascending = False
  162. self.side = "left" if alt_segment else "right"
  163. self.ts_sorted = ts[::-1]
  164. def _call_single(self, t):
  165. # Here we preserve a certain symmetry that when t is in self.ts,
  166. # if alt_segment=False, then we prioritize a segment with a lower
  167. # index.
  168. ind = np.searchsorted(self.ts_sorted, t, side=self.side)
  169. segment = min(max(ind - 1, 0), self.n_segments - 1)
  170. if not self.ascending:
  171. segment = self.n_segments - 1 - segment
  172. return self.interpolants[segment](t)
  173. def __call__(self, t):
  174. """Evaluate the solution.
  175. Parameters
  176. ----------
  177. t : float or array_like with shape (n_points,)
  178. Points to evaluate at.
  179. Returns
  180. -------
  181. y : ndarray, shape (n_states,) or (n_states, n_points)
  182. Computed values. Shape depends on whether `t` is a scalar or a
  183. 1-D array.
  184. """
  185. t = np.asarray(t)
  186. if t.ndim == 0:
  187. return self._call_single(t)
  188. order = np.argsort(t)
  189. reverse = np.empty_like(order)
  190. reverse[order] = np.arange(order.shape[0])
  191. t_sorted = t[order]
  192. # See comment in self._call_single.
  193. segments = np.searchsorted(self.ts_sorted, t_sorted, side=self.side)
  194. segments -= 1
  195. segments[segments < 0] = 0
  196. segments[segments > self.n_segments - 1] = self.n_segments - 1
  197. if not self.ascending:
  198. segments = self.n_segments - 1 - segments
  199. ys = []
  200. group_start = 0
  201. for segment, group in groupby(segments):
  202. group_end = group_start + len(list(group))
  203. y = self.interpolants[segment](t_sorted[group_start:group_end])
  204. ys.append(y)
  205. group_start = group_end
  206. ys = np.hstack(ys)
  207. ys = ys[:, reverse]
  208. return ys
  209. NUM_JAC_DIFF_REJECT = EPS ** 0.875
  210. NUM_JAC_DIFF_SMALL = EPS ** 0.75
  211. NUM_JAC_DIFF_BIG = EPS ** 0.25
  212. NUM_JAC_MIN_FACTOR = 1e3 * EPS
  213. NUM_JAC_FACTOR_INCREASE = 10
  214. NUM_JAC_FACTOR_DECREASE = 0.1
  215. def num_jac(fun, t, y, f, threshold, factor, sparsity=None):
  216. """Finite differences Jacobian approximation tailored for ODE solvers.
  217. This function computes finite difference approximation to the Jacobian
  218. matrix of `fun` with respect to `y` using forward differences.
  219. The Jacobian matrix has shape (n, n) and its element (i, j) is equal to
  220. ``d f_i / d y_j``.
  221. A special feature of this function is the ability to correct the step
  222. size from iteration to iteration. The main idea is to keep the finite
  223. difference significantly separated from its round-off error which
  224. approximately equals ``EPS * np.abs(f)``. It reduces a possibility of a
  225. huge error and assures that the estimated derivative are reasonably close
  226. to the true values (i.e., the finite difference approximation is at least
  227. qualitatively reflects the structure of the true Jacobian).
  228. Parameters
  229. ----------
  230. fun : callable
  231. Right-hand side of the system implemented in a vectorized fashion.
  232. t : float
  233. Current time.
  234. y : ndarray, shape (n,)
  235. Current state.
  236. f : ndarray, shape (n,)
  237. Value of the right hand side at (t, y).
  238. threshold : float
  239. Threshold for `y` value used for computing the step size as
  240. ``factor * np.maximum(np.abs(y), threshold)``. Typically, the value of
  241. absolute tolerance (atol) for a solver should be passed as `threshold`.
  242. factor : ndarray with shape (n,) or None
  243. Factor to use for computing the step size. Pass None for the very
  244. evaluation, then use the value returned from this function.
  245. sparsity : tuple (structure, groups) or None
  246. Sparsity structure of the Jacobian, `structure` must be csc_matrix.
  247. Returns
  248. -------
  249. J : ndarray or csc_matrix, shape (n, n)
  250. Jacobian matrix.
  251. factor : ndarray, shape (n,)
  252. Suggested `factor` for the next evaluation.
  253. """
  254. y = np.asarray(y)
  255. n = y.shape[0]
  256. if n == 0:
  257. return np.empty((0, 0)), factor
  258. if factor is None:
  259. factor = np.full(n, EPS ** 0.5)
  260. else:
  261. factor = factor.copy()
  262. # Direct the step as ODE dictates, hoping that such a step won't lead to
  263. # a problematic region. For complex ODEs it makes sense to use the real
  264. # part of f as we use steps along real axis.
  265. f_sign = 2 * (np.real(f) >= 0).astype(float) - 1
  266. y_scale = f_sign * np.maximum(threshold, np.abs(y))
  267. h = (y + factor * y_scale) - y
  268. # Make sure that the step is not 0 to start with. Not likely it will be
  269. # executed often.
  270. for i in np.nonzero(h == 0)[0]:
  271. while h[i] == 0:
  272. factor[i] *= 10
  273. h[i] = (y[i] + factor[i] * y_scale[i]) - y[i]
  274. if sparsity is None:
  275. return _dense_num_jac(fun, t, y, f, h, factor, y_scale)
  276. else:
  277. structure, groups = sparsity
  278. return _sparse_num_jac(fun, t, y, f, h, factor, y_scale,
  279. structure, groups)
  280. def _dense_num_jac(fun, t, y, f, h, factor, y_scale):
  281. n = y.shape[0]
  282. h_vecs = np.diag(h)
  283. f_new = fun(t, y[:, None] + h_vecs)
  284. diff = f_new - f[:, None]
  285. max_ind = np.argmax(np.abs(diff), axis=0)
  286. r = np.arange(n)
  287. max_diff = np.abs(diff[max_ind, r])
  288. scale = np.maximum(np.abs(f[max_ind]), np.abs(f_new[max_ind, r]))
  289. diff_too_small = max_diff < NUM_JAC_DIFF_REJECT * scale
  290. if np.any(diff_too_small):
  291. ind, = np.nonzero(diff_too_small)
  292. new_factor = NUM_JAC_FACTOR_INCREASE * factor[ind]
  293. h_new = (y[ind] + new_factor * y_scale[ind]) - y[ind]
  294. h_vecs[ind, ind] = h_new
  295. f_new = fun(t, y[:, None] + h_vecs[:, ind])
  296. diff_new = f_new - f[:, None]
  297. max_ind = np.argmax(np.abs(diff_new), axis=0)
  298. r = np.arange(ind.shape[0])
  299. max_diff_new = np.abs(diff_new[max_ind, r])
  300. scale_new = np.maximum(np.abs(f[max_ind]), np.abs(f_new[max_ind, r]))
  301. update = max_diff[ind] * scale_new < max_diff_new * scale[ind]
  302. if np.any(update):
  303. update, = np.nonzero(update)
  304. update_ind = ind[update]
  305. factor[update_ind] = new_factor[update]
  306. h[update_ind] = h_new[update]
  307. diff[:, update_ind] = diff_new[:, update]
  308. scale[update_ind] = scale_new[update]
  309. max_diff[update_ind] = max_diff_new[update]
  310. diff /= h
  311. factor[max_diff < NUM_JAC_DIFF_SMALL * scale] *= NUM_JAC_FACTOR_INCREASE
  312. factor[max_diff > NUM_JAC_DIFF_BIG * scale] *= NUM_JAC_FACTOR_DECREASE
  313. factor = np.maximum(factor, NUM_JAC_MIN_FACTOR)
  314. return diff, factor
  315. def _sparse_num_jac(fun, t, y, f, h, factor, y_scale, structure, groups):
  316. n = y.shape[0]
  317. n_groups = np.max(groups) + 1
  318. h_vecs = np.empty((n_groups, n))
  319. for group in range(n_groups):
  320. e = np.equal(group, groups)
  321. h_vecs[group] = h * e
  322. h_vecs = h_vecs.T
  323. f_new = fun(t, y[:, None] + h_vecs)
  324. df = f_new - f[:, None]
  325. i, j, _ = find(structure)
  326. diff = coo_matrix((df[i, groups[j]], (i, j)), shape=(n, n)).tocsc()
  327. max_ind = np.array(abs(diff).argmax(axis=0)).ravel()
  328. r = np.arange(n)
  329. max_diff = np.asarray(np.abs(diff[max_ind, r])).ravel()
  330. scale = np.maximum(np.abs(f[max_ind]),
  331. np.abs(f_new[max_ind, groups[r]]))
  332. diff_too_small = max_diff < NUM_JAC_DIFF_REJECT * scale
  333. if np.any(diff_too_small):
  334. ind, = np.nonzero(diff_too_small)
  335. new_factor = NUM_JAC_FACTOR_INCREASE * factor[ind]
  336. h_new = (y[ind] + new_factor * y_scale[ind]) - y[ind]
  337. h_new_all = np.zeros(n)
  338. h_new_all[ind] = h_new
  339. groups_unique = np.unique(groups[ind])
  340. groups_map = np.empty(n_groups, dtype=int)
  341. h_vecs = np.empty((groups_unique.shape[0], n))
  342. for k, group in enumerate(groups_unique):
  343. e = np.equal(group, groups)
  344. h_vecs[k] = h_new_all * e
  345. groups_map[group] = k
  346. h_vecs = h_vecs.T
  347. f_new = fun(t, y[:, None] + h_vecs)
  348. df = f_new - f[:, None]
  349. i, j, _ = find(structure[:, ind])
  350. diff_new = coo_matrix((df[i, groups_map[groups[ind[j]]]],
  351. (i, j)), shape=(n, ind.shape[0])).tocsc()
  352. max_ind_new = np.array(abs(diff_new).argmax(axis=0)).ravel()
  353. r = np.arange(ind.shape[0])
  354. max_diff_new = np.asarray(np.abs(diff_new[max_ind_new, r])).ravel()
  355. scale_new = np.maximum(
  356. np.abs(f[max_ind_new]),
  357. np.abs(f_new[max_ind_new, groups_map[groups[ind]]]))
  358. update = max_diff[ind] * scale_new < max_diff_new * scale[ind]
  359. if np.any(update):
  360. update, = np.nonzero(update)
  361. update_ind = ind[update]
  362. factor[update_ind] = new_factor[update]
  363. h[update_ind] = h_new[update]
  364. diff[:, update_ind] = diff_new[:, update]
  365. scale[update_ind] = scale_new[update]
  366. max_diff[update_ind] = max_diff_new[update]
  367. diff.data /= np.repeat(h, np.diff(diff.indptr))
  368. factor[max_diff < NUM_JAC_DIFF_SMALL * scale] *= NUM_JAC_FACTOR_INCREASE
  369. factor[max_diff > NUM_JAC_DIFF_BIG * scale] *= NUM_JAC_FACTOR_DECREASE
  370. factor = np.maximum(factor, NUM_JAC_MIN_FACTOR)
  371. return diff, factor