ivp.py 31 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757
  1. import inspect
  2. import numpy as np
  3. from .bdf import BDF
  4. from .radau import Radau
  5. from .rk import RK23, RK45, DOP853
  6. from .lsoda import LSODA
  7. from scipy.optimize import OptimizeResult
  8. from .common import EPS, OdeSolution
  9. from .base import OdeSolver
  10. from scipy._lib._array_api import xp_capabilities
  11. METHODS = {'RK23': RK23,
  12. 'RK45': RK45,
  13. 'DOP853': DOP853,
  14. 'Radau': Radau,
  15. 'BDF': BDF,
  16. 'LSODA': LSODA}
  17. MESSAGES = {0: "The solver successfully reached the end of the integration interval.",
  18. 1: "A termination event occurred."}
  19. class OdeResult(OptimizeResult):
  20. pass
  21. def prepare_events(events):
  22. """Standardize event functions and extract attributes."""
  23. if callable(events):
  24. events = (events,)
  25. max_events = np.empty(len(events))
  26. direction = np.empty(len(events))
  27. for i, event in enumerate(events):
  28. terminal = getattr(event, 'terminal', None)
  29. direction[i] = getattr(event, 'direction', 0)
  30. message = ('The `terminal` attribute of each event '
  31. 'must be a boolean or positive integer.')
  32. if terminal is None or terminal == 0:
  33. max_events[i] = np.inf
  34. elif int(terminal) == terminal and terminal > 0:
  35. max_events[i] = terminal
  36. else:
  37. raise ValueError(message)
  38. return events, max_events, direction
  39. def solve_event_equation(event, sol, t_old, t):
  40. """Solve an equation corresponding to an ODE event.
  41. The equation is ``event(t, y(t)) = 0``, here ``y(t)`` is known from an
  42. ODE solver using some sort of interpolation. It is solved by
  43. `scipy.optimize.brentq` with xtol=atol=4*EPS.
  44. Parameters
  45. ----------
  46. event : callable
  47. Function ``event(t, y)``.
  48. sol : callable
  49. Function ``sol(t)`` which evaluates an ODE solution between `t_old`
  50. and `t`.
  51. t_old, t : float
  52. Previous and new values of time. They will be used as a bracketing
  53. interval.
  54. Returns
  55. -------
  56. root : float
  57. Found solution.
  58. """
  59. from scipy.optimize import brentq
  60. return brentq(lambda t: event(t, sol(t)), t_old, t,
  61. xtol=4 * EPS, rtol=4 * EPS)
  62. def handle_events(sol, events, active_events, event_count, max_events,
  63. t_old, t):
  64. """Helper function to handle events.
  65. Parameters
  66. ----------
  67. sol : DenseOutput
  68. Function ``sol(t)`` which evaluates an ODE solution between `t_old`
  69. and `t`.
  70. events : list of callables, length n_events
  71. Event functions with signatures ``event(t, y)``.
  72. active_events : ndarray
  73. Indices of events which occurred.
  74. event_count : ndarray
  75. Current number of occurrences for each event.
  76. max_events : ndarray, shape (n_events,)
  77. Number of occurrences allowed for each event before integration
  78. termination is issued.
  79. t_old, t : float
  80. Previous and new values of time.
  81. Returns
  82. -------
  83. root_indices : ndarray
  84. Indices of events which take zero between `t_old` and `t` and before
  85. a possible termination.
  86. roots : ndarray
  87. Values of t at which events occurred.
  88. terminate : bool
  89. Whether a terminal event occurred.
  90. """
  91. roots = [solve_event_equation(events[event_index], sol, t_old, t)
  92. for event_index in active_events]
  93. roots = np.asarray(roots)
  94. if np.any(event_count[active_events] >= max_events[active_events]):
  95. if t > t_old:
  96. order = np.argsort(roots)
  97. else:
  98. order = np.argsort(-roots)
  99. active_events = active_events[order]
  100. roots = roots[order]
  101. t = np.nonzero(event_count[active_events]
  102. >= max_events[active_events])[0][0]
  103. active_events = active_events[:t + 1]
  104. roots = roots[:t + 1]
  105. terminate = True
  106. else:
  107. terminate = False
  108. return active_events, roots, terminate
  109. def find_active_events(g, g_new, direction):
  110. """Find which event occurred during an integration step.
  111. Parameters
  112. ----------
  113. g, g_new : array_like, shape (n_events,)
  114. Values of event functions at a current and next points.
  115. direction : ndarray, shape (n_events,)
  116. Event "direction" according to the definition in `solve_ivp`.
  117. Returns
  118. -------
  119. active_events : ndarray
  120. Indices of events which occurred during the step.
  121. """
  122. g, g_new = np.asarray(g), np.asarray(g_new)
  123. up = (g <= 0) & (g_new >= 0)
  124. down = (g >= 0) & (g_new <= 0)
  125. either = up | down
  126. mask = (up & (direction > 0) |
  127. down & (direction < 0) |
  128. either & (direction == 0))
  129. return np.nonzero(mask)[0]
  130. @xp_capabilities(np_only=True)
  131. def solve_ivp(fun, t_span, y0, method='RK45', t_eval=None, dense_output=False,
  132. events=None, vectorized=False, args=None, **options):
  133. """Solve an initial value problem for a system of ODEs.
  134. This function numerically integrates a system of ordinary differential
  135. equations given an initial value::
  136. dy / dt = f(t, y)
  137. y(t0) = y0
  138. Here t is a 1-D independent variable (time), y(t) is an
  139. N-D vector-valued function (state), and an N-D
  140. vector-valued function f(t, y) determines the differential equations.
  141. The goal is to find y(t) approximately satisfying the differential
  142. equations, given an initial value y(t0)=y0.
  143. Some of the solvers support integration in the complex domain, but note
  144. that for stiff ODE solvers, the right-hand side must be
  145. complex-differentiable (satisfy Cauchy-Riemann equations [11]_).
  146. To solve a problem in the complex domain, pass y0 with a complex data type.
  147. Another option always available is to rewrite your problem for real and
  148. imaginary parts separately.
  149. Parameters
  150. ----------
  151. fun : callable
  152. Right-hand side of the system: the time derivative of the state ``y``
  153. at time ``t``. The calling signature is ``fun(t, y)``, where ``t`` is a
  154. scalar and ``y`` is an ndarray with ``len(y) = len(y0)``. Additional
  155. arguments need to be passed if ``args`` is used (see documentation of
  156. ``args`` argument). ``fun`` must return an array of the same shape as
  157. ``y``. See `vectorized` for more information.
  158. t_span : 2-member sequence
  159. Interval of integration (t0, tf). The solver starts with t=t0 and
  160. integrates until it reaches t=tf. Both t0 and tf must be floats
  161. or values interpretable by the float conversion function.
  162. y0 : array_like, shape (n,)
  163. Initial state. For problems in the complex domain, pass `y0` with a
  164. complex data type (even if the initial value is purely real).
  165. method : string or `OdeSolver`, optional
  166. Integration method to use:
  167. * 'RK45' (default): Explicit Runge-Kutta method of order 5(4) [1]_.
  168. The error is controlled assuming accuracy of the fourth-order
  169. method, but steps are taken using the fifth-order accurate
  170. formula (local extrapolation is done). A quartic interpolation
  171. polynomial is used for the dense output [2]_. Can be applied in
  172. the complex domain.
  173. * 'RK23': Explicit Runge-Kutta method of order 3(2) [3]_. The error
  174. is controlled assuming accuracy of the second-order method, but
  175. steps are taken using the third-order accurate formula (local
  176. extrapolation is done). A cubic Hermite polynomial is used for the
  177. dense output. Can be applied in the complex domain.
  178. * 'DOP853': Explicit Runge-Kutta method of order 8 [13]_.
  179. Python implementation of the "DOP853" algorithm originally
  180. written in Fortran [14]_. A 7-th order interpolation polynomial
  181. accurate to 7-th order is used for the dense output.
  182. Can be applied in the complex domain.
  183. * 'Radau': Implicit Runge-Kutta method of the Radau IIA family of
  184. order 5 [4]_. The error is controlled with a third-order accurate
  185. embedded formula. A cubic polynomial which satisfies the
  186. collocation conditions is used for the dense output.
  187. * 'BDF': Implicit multi-step variable-order (1 to 5) method based
  188. on a backward differentiation formula for the derivative
  189. approximation [5]_. The implementation follows the one described
  190. in [6]_. A quasi-constant step scheme is used and accuracy is
  191. enhanced using the NDF modification. Can be applied in the
  192. complex domain.
  193. * 'LSODA': Adams/BDF method with automatic stiffness detection and
  194. switching [7]_, [8]_. This is a wrapper of the Fortran solver
  195. from ODEPACK.
  196. Explicit Runge-Kutta methods ('RK23', 'RK45', 'DOP853') should be used
  197. for non-stiff problems and implicit methods ('Radau', 'BDF') for
  198. stiff problems [9]_. Among Runge-Kutta methods, 'DOP853' is recommended
  199. for solving with high precision (low values of `rtol` and `atol`).
  200. If not sure, first try to run 'RK45'. If it makes unusually many
  201. iterations, diverges, or fails, your problem is likely to be stiff and
  202. you should use 'Radau' or 'BDF'. 'LSODA' can also be a good universal
  203. choice, but it might be somewhat less convenient to work with as it
  204. wraps old Fortran code.
  205. You can also pass an arbitrary class derived from `OdeSolver` which
  206. implements the solver.
  207. t_eval : array_like or None, optional
  208. Times at which to store the computed solution, must be sorted and lie
  209. within `t_span`. If None (default), use points selected by the solver.
  210. dense_output : bool, optional
  211. Whether to compute a continuous solution. Default is False.
  212. events : callable, or list of callables, optional
  213. Events to track. If None (default), no events will be tracked.
  214. Each event occurs at the zeros of a continuous function of time and
  215. state. Each function must have the signature ``event(t, y)`` where
  216. additional argument have to be passed if ``args`` is used (see
  217. documentation of ``args`` argument). Each function must return a
  218. float. The solver will find an accurate value of `t` at which
  219. ``event(t, y(t)) = 0`` using a root-finding algorithm. By default,
  220. all zeros will be found. The solver looks for a sign change over
  221. each step, so if multiple zero crossings occur within one step,
  222. events may be missed. Additionally each `event` function might
  223. have the following attributes:
  224. terminal: bool or int, optional
  225. When boolean, whether to terminate integration if this event occurs.
  226. When integral, termination occurs after the specified the number of
  227. occurrences of this event.
  228. Implicitly False if not assigned.
  229. direction: float, optional
  230. Direction of a zero crossing. If `direction` is positive,
  231. `event` will only trigger when going from negative to positive,
  232. and vice versa if `direction` is negative. If 0, then either
  233. direction will trigger event. Implicitly 0 if not assigned.
  234. You can assign attributes like ``event.terminal = True`` to any
  235. function in Python.
  236. vectorized : bool, optional
  237. Whether `fun` can be called in a vectorized fashion. Default is False.
  238. If ``vectorized`` is False, `fun` will always be called with ``y`` of
  239. shape ``(n,)``, where ``n = len(y0)``.
  240. If ``vectorized`` is True, `fun` may be called with ``y`` of shape
  241. ``(n, k)``, where ``k`` is an integer. In this case, `fun` must behave
  242. such that ``fun(t, y)[:, i] == fun(t, y[:, i])`` (i.e. each column of
  243. the returned array is the time derivative of the state corresponding
  244. with a column of ``y``).
  245. Setting ``vectorized=True`` allows for faster finite difference
  246. approximation of the Jacobian by methods 'Radau' and 'BDF', but
  247. will result in slower execution for other methods and for 'Radau' and
  248. 'BDF' in some circumstances (e.g. small ``len(y0)``).
  249. args : tuple, optional
  250. Additional arguments to pass to the user-defined functions. If given,
  251. the additional arguments are passed to all user-defined functions.
  252. So if, for example, `fun` has the signature ``fun(t, y, a, b, c)``,
  253. then `jac` (if given) and any event functions must have the same
  254. signature, and `args` must be a tuple of length 3.
  255. **options
  256. Options passed to a chosen solver. All options available for already
  257. implemented solvers are listed below.
  258. first_step : float or None, optional
  259. Initial step size. Default is `None` which means that the algorithm
  260. should choose.
  261. max_step : float, optional
  262. Maximum allowed step size. Default is np.inf, i.e., the step size is not
  263. bounded and determined solely by the solver.
  264. rtol, atol : float or array_like, optional
  265. Relative and absolute tolerances. The solver keeps the local error
  266. estimates less than ``atol + rtol * abs(y)``. Here `rtol` controls a
  267. relative accuracy (number of correct digits), while `atol` controls
  268. absolute accuracy (number of correct decimal places). To achieve the
  269. desired `rtol`, set `atol` to be smaller than the smallest value that
  270. can be expected from ``rtol * abs(y)`` so that `rtol` dominates the
  271. allowable error. If `atol` is larger than ``rtol * abs(y)`` the
  272. number of correct digits is not guaranteed. Conversely, to achieve the
  273. desired `atol` set `rtol` such that ``rtol * abs(y)`` is always smaller
  274. than `atol`. If components of y have different scales, it might be
  275. beneficial to set different `atol` values for different components by
  276. passing array_like with shape (n,) for `atol`. Default values are
  277. 1e-3 for `rtol` and 1e-6 for `atol`.
  278. jac : array_like, sparse_matrix, callable or None, optional
  279. Jacobian matrix of the right-hand side of the system with respect
  280. to y, required by the 'Radau', 'BDF' and 'LSODA' method. The
  281. Jacobian matrix has shape (n, n) and its element (i, j) is equal to
  282. ``d f_i / d y_j``. There are three ways to define the Jacobian:
  283. * If array_like or sparse_matrix, the Jacobian is assumed to
  284. be constant. Not supported by 'LSODA'.
  285. * If callable, the Jacobian is assumed to depend on both
  286. t and y; it will be called as ``jac(t, y)``, as necessary.
  287. Additional arguments have to be passed if ``args`` is
  288. used (see documentation of ``args`` argument).
  289. For 'Radau' and 'BDF' methods, the return value might be a
  290. sparse matrix.
  291. * If None (default), the Jacobian will be approximated by
  292. finite differences.
  293. It is generally recommended to provide the Jacobian rather than
  294. relying on a finite-difference approximation.
  295. jac_sparsity : array_like, sparse matrix or None, optional
  296. Defines a sparsity structure of the Jacobian matrix for a finite-
  297. difference approximation. Its shape must be (n, n). This argument
  298. is ignored if `jac` is not `None`. If the Jacobian has only few
  299. non-zero elements in *each* row, providing the sparsity structure
  300. will greatly speed up the computations [10]_. A zero entry means that
  301. a corresponding element in the Jacobian is always zero. If None
  302. (default), the Jacobian is assumed to be dense.
  303. Not supported by 'LSODA', see `lband` and `uband` instead.
  304. lband, uband : int or None, optional
  305. Parameters defining the bandwidth of the Jacobian for the 'LSODA'
  306. method, i.e., ``jac[i, j] != 0 only for i - lband <= j <= i + uband``.
  307. Default is None. Setting these requires your jac routine to return the
  308. Jacobian in the packed format: the returned array must have ``n``
  309. columns and ``uband + lband + 1`` rows in which Jacobian diagonals are
  310. written. Specifically ``jac_packed[uband + i - j , j] = jac[i, j]``.
  311. The same format is used in `scipy.linalg.solve_banded` (check for an
  312. illustration). These parameters can be also used with ``jac=None`` to
  313. reduce the number of Jacobian elements estimated by finite differences.
  314. min_step : float, optional
  315. The minimum allowed step size for 'LSODA' method.
  316. By default `min_step` is zero.
  317. Returns
  318. -------
  319. Bunch object with the following fields defined:
  320. t : ndarray, shape (n_points,)
  321. Time points.
  322. y : ndarray, shape (n, n_points)
  323. Values of the solution at `t`.
  324. sol : `OdeSolution` or None
  325. Found solution as `OdeSolution` instance; None if `dense_output` was
  326. set to False.
  327. t_events : list of ndarray or None
  328. Contains for each event type a list of arrays at which an event of
  329. that type event was detected. None if `events` was None.
  330. y_events : list of ndarray or None
  331. For each value of `t_events`, the corresponding value of the solution.
  332. None if `events` was None.
  333. nfev : int
  334. Number of evaluations of the right-hand side.
  335. njev : int
  336. Number of evaluations of the Jacobian.
  337. nlu : int
  338. Number of LU decompositions.
  339. status : int
  340. Reason for algorithm termination:
  341. * -1: Integration step failed.
  342. * 0: The solver successfully reached the end of `tspan`.
  343. * 1: A termination event occurred.
  344. message : string
  345. Human-readable description of the termination reason.
  346. success : bool
  347. True if the solver reached the interval end or a termination event
  348. occurred (``status >= 0``).
  349. References
  350. ----------
  351. .. [1] J. R. Dormand, P. J. Prince, "A family of embedded Runge-Kutta
  352. formulae", Journal of Computational and Applied Mathematics, Vol. 6,
  353. No. 1, pp. 19-26, 1980.
  354. .. [2] L. W. Shampine, "Some Practical Runge-Kutta Formulas", Mathematics
  355. of Computation,, Vol. 46, No. 173, pp. 135-150, 1986.
  356. .. [3] P. Bogacki, L.F. Shampine, "A 3(2) Pair of Runge-Kutta Formulas",
  357. Appl. Math. Lett. Vol. 2, No. 4. pp. 321-325, 1989.
  358. .. [4] E. Hairer, G. Wanner, "Solving Ordinary Differential Equations II:
  359. Stiff and Differential-Algebraic Problems", Sec. IV.8.
  360. .. [5] `Backward Differentiation Formula
  361. <https://en.wikipedia.org/wiki/Backward_differentiation_formula>`_
  362. on Wikipedia.
  363. .. [6] L. F. Shampine, M. W. Reichelt, "THE MATLAB ODE SUITE", SIAM J. SCI.
  364. COMPUTE., Vol. 18, No. 1, pp. 1-22, January 1997.
  365. .. [7] A. C. Hindmarsh, "ODEPACK, A Systematized Collection of ODE
  366. Solvers," IMACS Transactions on Scientific Computation, Vol 1.,
  367. pp. 55-64, 1983.
  368. .. [8] L. Petzold, "Automatic selection of methods for solving stiff and
  369. nonstiff systems of ordinary differential equations", SIAM Journal
  370. on Scientific and Statistical Computing, Vol. 4, No. 1, pp. 136-148,
  371. 1983.
  372. .. [9] `Stiff equation <https://en.wikipedia.org/wiki/Stiff_equation>`_ on
  373. Wikipedia.
  374. .. [10] A. Curtis, M. J. D. Powell, and J. Reid, "On the estimation of
  375. sparse Jacobian matrices", Journal of the Institute of Mathematics
  376. and its Applications, 13, pp. 117-120, 1974.
  377. .. [11] `Cauchy-Riemann equations
  378. <https://en.wikipedia.org/wiki/Cauchy-Riemann_equations>`_ on
  379. Wikipedia.
  380. .. [12] `Lotka-Volterra equations
  381. <https://en.wikipedia.org/wiki/Lotka%E2%80%93Volterra_equations>`_
  382. on Wikipedia.
  383. .. [13] E. Hairer, S. P. Norsett G. Wanner, "Solving Ordinary Differential
  384. Equations I: Nonstiff Problems", Sec. II.
  385. .. [14] `Page with original Fortran code of DOP853
  386. <http://www.unige.ch/~hairer/software.html>`_.
  387. Examples
  388. --------
  389. Basic exponential decay showing automatically chosen time points.
  390. >>> import numpy as np
  391. >>> from scipy.integrate import solve_ivp
  392. >>> def exponential_decay(t, y): return -0.5 * y
  393. >>> sol = solve_ivp(exponential_decay, [0, 10], [2, 4, 8])
  394. >>> print(sol.t)
  395. [ 0. 0.11487653 1.26364188 3.06061781 4.81611105 6.57445806
  396. 8.33328988 10. ]
  397. >>> print(sol.y)
  398. [[2. 1.88836035 1.06327177 0.43319312 0.18017253 0.07483045
  399. 0.03107158 0.01350781]
  400. [4. 3.7767207 2.12654355 0.86638624 0.36034507 0.14966091
  401. 0.06214316 0.02701561]
  402. [8. 7.5534414 4.25308709 1.73277247 0.72069014 0.29932181
  403. 0.12428631 0.05403123]]
  404. Specifying points where the solution is desired.
  405. >>> sol = solve_ivp(exponential_decay, [0, 10], [2, 4, 8],
  406. ... t_eval=[0, 1, 2, 4, 10])
  407. >>> print(sol.t)
  408. [ 0 1 2 4 10]
  409. >>> print(sol.y)
  410. [[2. 1.21305369 0.73534021 0.27066736 0.01350938]
  411. [4. 2.42610739 1.47068043 0.54133472 0.02701876]
  412. [8. 4.85221478 2.94136085 1.08266944 0.05403753]]
  413. Cannon fired upward with terminal event upon impact. The ``terminal`` and
  414. ``direction`` fields of an event are applied by monkey patching a function.
  415. Here ``y[0]`` is position and ``y[1]`` is velocity. The projectile starts
  416. at position 0 with velocity +10. Note that the integration never reaches
  417. t=100 because the event is terminal.
  418. >>> def upward_cannon(t, y): return [y[1], -0.5]
  419. >>> def hit_ground(t, y): return y[0]
  420. >>> hit_ground.terminal = True
  421. >>> hit_ground.direction = -1
  422. >>> sol = solve_ivp(upward_cannon, [0, 100], [0, 10], events=hit_ground)
  423. >>> print(sol.t_events)
  424. [array([40.])]
  425. >>> print(sol.t)
  426. [0.00000000e+00 9.99900010e-05 1.09989001e-03 1.10988901e-02
  427. 1.11088891e-01 1.11098890e+00 1.11099890e+01 4.00000000e+01]
  428. Use `dense_output` and `events` to find position, which is 100, at the apex
  429. of the cannonball's trajectory. Apex is not defined as terminal, so both
  430. apex and hit_ground are found. There is no information at t=20, so the sol
  431. attribute is used to evaluate the solution. The sol attribute is returned
  432. by setting ``dense_output=True``. Alternatively, the `y_events` attribute
  433. can be used to access the solution at the time of the event.
  434. >>> def apex(t, y): return y[1]
  435. >>> sol = solve_ivp(upward_cannon, [0, 100], [0, 10],
  436. ... events=(hit_ground, apex), dense_output=True)
  437. >>> print(sol.t_events)
  438. [array([40.]), array([20.])]
  439. >>> print(sol.t)
  440. [0.00000000e+00 9.99900010e-05 1.09989001e-03 1.10988901e-02
  441. 1.11088891e-01 1.11098890e+00 1.11099890e+01 4.00000000e+01]
  442. >>> print(sol.sol(sol.t_events[1][0]))
  443. [100. 0.]
  444. >>> print(sol.y_events)
  445. [array([[-5.68434189e-14, -1.00000000e+01]]),
  446. array([[1.00000000e+02, 1.77635684e-15]])]
  447. As an example of a system with additional parameters, we'll implement
  448. the Lotka-Volterra equations [12]_.
  449. >>> def lotkavolterra(t, z, a, b, c, d):
  450. ... x, y = z
  451. ... return [a*x - b*x*y, -c*y + d*x*y]
  452. ...
  453. We pass in the parameter values a=1.5, b=1, c=3 and d=1 with the `args`
  454. argument.
  455. >>> sol = solve_ivp(lotkavolterra, [0, 15], [10, 5], args=(1.5, 1, 3, 1),
  456. ... dense_output=True)
  457. Compute a dense solution and plot it.
  458. >>> t = np.linspace(0, 15, 300)
  459. >>> z = sol.sol(t)
  460. >>> import matplotlib.pyplot as plt
  461. >>> plt.plot(t, z.T)
  462. >>> plt.xlabel('t')
  463. >>> plt.legend(['x', 'y'], shadow=True)
  464. >>> plt.title('Lotka-Volterra System')
  465. >>> plt.show()
  466. A couple examples of using solve_ivp to solve the differential
  467. equation ``y' = Ay`` with complex matrix ``A``.
  468. >>> A = np.array([[-0.25 + 0.14j, 0, 0.33 + 0.44j],
  469. ... [0.25 + 0.58j, -0.2 + 0.14j, 0],
  470. ... [0, 0.2 + 0.4j, -0.1 + 0.97j]])
  471. Solving an IVP with ``A`` from above and ``y`` as 3x1 vector:
  472. >>> def deriv_vec(t, y):
  473. ... return A @ y
  474. >>> result = solve_ivp(deriv_vec, [0, 25],
  475. ... np.array([10 + 0j, 20 + 0j, 30 + 0j]),
  476. ... t_eval=np.linspace(0, 25, 101))
  477. >>> print(result.y[:, 0])
  478. [10.+0.j 20.+0.j 30.+0.j]
  479. >>> print(result.y[:, -1])
  480. [18.46291039+45.25653651j 10.01569306+36.23293216j
  481. -4.98662741+80.07360388j]
  482. Solving an IVP with ``A`` from above with ``y`` as 3x3 matrix :
  483. >>> def deriv_mat(t, y):
  484. ... return (A @ y.reshape(3, 3)).flatten()
  485. >>> y0 = np.array([[2 + 0j, 3 + 0j, 4 + 0j],
  486. ... [5 + 0j, 6 + 0j, 7 + 0j],
  487. ... [9 + 0j, 34 + 0j, 78 + 0j]])
  488. >>> result = solve_ivp(deriv_mat, [0, 25], y0.flatten(),
  489. ... t_eval=np.linspace(0, 25, 101))
  490. >>> print(result.y[:, 0].reshape(3, 3))
  491. [[ 2.+0.j 3.+0.j 4.+0.j]
  492. [ 5.+0.j 6.+0.j 7.+0.j]
  493. [ 9.+0.j 34.+0.j 78.+0.j]]
  494. >>> print(result.y[:, -1].reshape(3, 3))
  495. [[ 5.67451179 +12.07938445j 17.2888073 +31.03278837j
  496. 37.83405768 +63.25138759j]
  497. [ 3.39949503 +11.82123994j 21.32530996 +44.88668871j
  498. 53.17531184+103.80400411j]
  499. [ -2.26105874 +22.19277664j -15.1255713 +70.19616341j
  500. -38.34616845+153.29039931j]]
  501. """
  502. if method not in METHODS and not (
  503. inspect.isclass(method) and issubclass(method, OdeSolver)):
  504. raise ValueError(f"`method` must be one of {METHODS} or OdeSolver class.")
  505. t0, tf = map(float, t_span)
  506. if args is not None:
  507. # Wrap the user's fun (and jac, if given) in lambdas to hide the
  508. # additional parameters. Pass in the original fun as a keyword
  509. # argument to keep it in the scope of the lambda.
  510. try:
  511. _ = [*(args)]
  512. except TypeError as exp:
  513. suggestion_tuple = (
  514. "Supplied 'args' cannot be unpacked. Please supply `args`"
  515. f" as a tuple (e.g. `args=({args},)`)"
  516. )
  517. raise TypeError(suggestion_tuple) from exp
  518. def fun(t, x, fun=fun):
  519. return fun(t, x, *args)
  520. jac = options.get('jac')
  521. if callable(jac):
  522. options['jac'] = lambda t, x: jac(t, x, *args)
  523. if t_eval is not None:
  524. t_eval = np.asarray(t_eval)
  525. if t_eval.ndim != 1:
  526. raise ValueError("`t_eval` must be 1-dimensional.")
  527. if np.any(t_eval < min(t0, tf)) or np.any(t_eval > max(t0, tf)):
  528. raise ValueError("Values in `t_eval` are not within `t_span`.")
  529. d = np.diff(t_eval)
  530. if tf > t0 and np.any(d <= 0) or tf < t0 and np.any(d >= 0):
  531. raise ValueError("Values in `t_eval` are not properly sorted.")
  532. if tf > t0:
  533. t_eval_i = 0
  534. else:
  535. # Make order of t_eval decreasing to use np.searchsorted.
  536. t_eval = t_eval[::-1]
  537. # This will be an upper bound for slices.
  538. t_eval_i = t_eval.shape[0]
  539. if method in METHODS:
  540. method = METHODS[method]
  541. solver = method(fun, t0, y0, tf, vectorized=vectorized, **options)
  542. if t_eval is None:
  543. ts = [t0]
  544. ys = [y0]
  545. elif t_eval is not None and dense_output:
  546. ts = []
  547. ti = [t0]
  548. ys = []
  549. else:
  550. ts = []
  551. ys = []
  552. interpolants = []
  553. if events is not None:
  554. events, max_events, event_dir = prepare_events(events)
  555. event_count = np.zeros(len(events))
  556. if args is not None:
  557. # Wrap user functions in lambdas to hide the additional parameters.
  558. # The original event function is passed as a keyword argument to the
  559. # lambda to keep the original function in scope (i.e., avoid the
  560. # late binding closure "gotcha").
  561. events = [lambda t, x, event=event: event(t, x, *args)
  562. for event in events]
  563. g = [event(t0, y0) for event in events]
  564. t_events = [[] for _ in range(len(events))]
  565. y_events = [[] for _ in range(len(events))]
  566. else:
  567. t_events = None
  568. y_events = None
  569. status = None
  570. while status is None:
  571. message = solver.step()
  572. if solver.status == 'finished':
  573. status = 0
  574. elif solver.status == 'failed':
  575. status = -1
  576. break
  577. t_old = solver.t_old
  578. t = solver.t
  579. y = solver.y
  580. if dense_output:
  581. sol = solver.dense_output()
  582. interpolants.append(sol)
  583. else:
  584. sol = None
  585. if events is not None:
  586. g_new = [event(t, y) for event in events]
  587. active_events = find_active_events(g, g_new, event_dir)
  588. if active_events.size > 0:
  589. if sol is None:
  590. sol = solver.dense_output()
  591. event_count[active_events] += 1
  592. root_indices, roots, terminate = handle_events(
  593. sol, events, active_events, event_count, max_events,
  594. t_old, t)
  595. for e, te in zip(root_indices, roots):
  596. t_events[e].append(te)
  597. y_events[e].append(sol(te))
  598. if terminate:
  599. status = 1
  600. t = roots[-1]
  601. y = sol(t)
  602. g = g_new
  603. if t_eval is None:
  604. donot_append = (len(ts) > 1 and
  605. ts[-1] == t and
  606. dense_output)
  607. if not donot_append:
  608. ts.append(t)
  609. ys.append(y)
  610. else:
  611. if len(interpolants) > 0:
  612. interpolants.pop()
  613. else:
  614. # The value in t_eval equal to t will be included.
  615. if solver.direction > 0:
  616. t_eval_i_new = np.searchsorted(t_eval, t, side='right')
  617. t_eval_step = t_eval[t_eval_i:t_eval_i_new]
  618. else:
  619. t_eval_i_new = np.searchsorted(t_eval, t, side='left')
  620. # It has to be done with two slice operations, because
  621. # you can't slice to 0th element inclusive using backward
  622. # slicing.
  623. t_eval_step = t_eval[t_eval_i_new:t_eval_i][::-1]
  624. if t_eval_step.size > 0:
  625. if sol is None:
  626. sol = solver.dense_output()
  627. ts.append(t_eval_step)
  628. ys.append(sol(t_eval_step))
  629. t_eval_i = t_eval_i_new
  630. if t_eval is not None and dense_output:
  631. ti.append(t)
  632. message = MESSAGES.get(status, message)
  633. if t_events is not None:
  634. t_events = [np.asarray(te) for te in t_events]
  635. y_events = [np.asarray(ye) for ye in y_events]
  636. if t_eval is None:
  637. ts = np.array(ts)
  638. ys = np.vstack(ys).T
  639. elif ts:
  640. ts = np.hstack(ts)
  641. ys = np.hstack(ys)
  642. if dense_output:
  643. if t_eval is None:
  644. sol = OdeSolution(
  645. ts, interpolants, alt_segment=True if method in [BDF, LSODA] else False
  646. )
  647. else:
  648. sol = OdeSolution(
  649. ti, interpolants, alt_segment=True if method in [BDF, LSODA] else False
  650. )
  651. else:
  652. sol = None
  653. return OdeResult(t=ts, y=ys, sol=sol, t_events=t_events, y_events=y_events,
  654. nfev=solver.nfev, njev=solver.njev, nlu=solver.nlu,
  655. status=status, message=message, success=status >= 0)