common.py 14 KB

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