_odepack_py.py 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271
  1. # Author: Travis Oliphant
  2. __all__ = ['odeint', 'ODEintWarning']
  3. import numpy as np
  4. from . import _odepack
  5. from copy import copy
  6. import warnings
  7. from scipy._lib._array_api import xp_capabilities
  8. class ODEintWarning(Warning):
  9. """Warning raised during the execution of `odeint`."""
  10. pass
  11. _msgs = {2: "Integration successful.",
  12. 1: "Nothing was done; the integration time was 0.",
  13. -1: "Excess work done on this call (perhaps wrong Dfun type).",
  14. -2: "Excess accuracy requested (tolerances too small).",
  15. -3: "Illegal input detected (internal error).",
  16. -4: "Repeated error test failures (internal error).",
  17. -5: "Repeated convergence failures (perhaps bad Jacobian or tolerances).",
  18. -6: "Error weight became zero during problem.",
  19. -7: "Internal workspace insufficient to finish (internal error).",
  20. -8: "Run terminated (internal error)."
  21. }
  22. @xp_capabilities(out_of_scope=True)
  23. def odeint(func, y0, t, args=(), Dfun=None, col_deriv=0, full_output=0,
  24. ml=None, mu=None, rtol=None, atol=None, tcrit=None, h0=0.0,
  25. hmax=0.0, hmin=0.0, ixpr=0, mxstep=0, mxhnil=0, mxordn=12,
  26. mxords=5, printmessg=0, tfirst=False):
  27. """
  28. Integrate a system of ordinary differential equations.
  29. .. note:: For new code, use `scipy.integrate.solve_ivp` to solve a
  30. differential equation.
  31. Solve a system of ordinary differential equations using lsoda from the
  32. FORTRAN library odepack.
  33. Solves the initial value problem for stiff or non-stiff systems
  34. of first order ode-s::
  35. dy/dt = func(y, t, ...) [or func(t, y, ...)]
  36. where y can be a vector.
  37. .. note:: By default, the required order of the first two arguments of
  38. `func` are in the opposite order of the arguments in the system
  39. definition function used by the `scipy.integrate.ode` class and
  40. the function `scipy.integrate.solve_ivp`. To use a function with
  41. the signature ``func(t, y, ...)``, the argument `tfirst` must be
  42. set to ``True``.
  43. Parameters
  44. ----------
  45. func : callable(y, t, ...) or callable(t, y, ...)
  46. Computes the derivative of y at t.
  47. If the signature is ``callable(t, y, ...)``, then the argument
  48. `tfirst` must be set ``True``.
  49. `func` must not modify the data in `y`, as it is a
  50. view of the data used internally by the ODE solver.
  51. y0 : array
  52. Initial condition on y (can be a vector).
  53. t : array
  54. A sequence of time points for which to solve for y. The initial
  55. value point should be the first element of this sequence.
  56. This sequence must be monotonically increasing or monotonically
  57. decreasing; repeated values are allowed.
  58. args : tuple, optional
  59. Extra arguments to pass to function.
  60. Dfun : callable(y, t, ...) or callable(t, y, ...)
  61. Gradient (Jacobian) of `func`.
  62. If the signature is ``callable(t, y, ...)``, then the argument
  63. `tfirst` must be set ``True``.
  64. `Dfun` must not modify the data in `y`, as it is a
  65. view of the data used internally by the ODE solver.
  66. col_deriv : bool, optional
  67. True if `Dfun` defines derivatives down columns (faster),
  68. otherwise `Dfun` should define derivatives across rows.
  69. full_output : bool, optional
  70. True if to return a dictionary of optional outputs as the second output
  71. printmessg : bool, optional
  72. Whether to print the convergence message
  73. tfirst : bool, optional
  74. If True, the first two arguments of `func` (and `Dfun`, if given)
  75. must ``t, y`` instead of the default ``y, t``.
  76. .. versionadded:: 1.1.0
  77. Returns
  78. -------
  79. y : array, shape (len(t), len(y0))
  80. Array containing the value of y for each desired time in t,
  81. with the initial value `y0` in the first row.
  82. infodict : dict, only returned if full_output == True
  83. Dictionary containing additional output information
  84. ======= ============================================================
  85. key meaning
  86. ======= ============================================================
  87. 'hu' vector of step sizes successfully used for each time step
  88. 'tcur' vector with the value of t reached for each time step
  89. (will always be at least as large as the input times)
  90. 'tolsf' vector of tolerance scale factors, greater than 1.0,
  91. computed when a request for too much accuracy was detected
  92. 'tsw' value of t at the time of the last method switch
  93. (given for each time step)
  94. 'nst' cumulative number of time steps
  95. 'nfe' cumulative number of function evaluations for each time step
  96. 'nje' cumulative number of jacobian evaluations for each time step
  97. 'nqu' a vector of method orders for each successful step
  98. 'imxer' index of the component of largest magnitude in the
  99. weighted local error vector (e / ewt) on an error return, -1
  100. otherwise
  101. 'lenrw' the length of the double work array required
  102. 'leniw' the length of integer work array required
  103. 'mused' a vector of method indicators for each successful time step:
  104. 1: adams (nonstiff), 2: bdf (stiff)
  105. ======= ============================================================
  106. Other Parameters
  107. ----------------
  108. ml, mu : int, optional
  109. If either of these are not None or non-negative, then the
  110. Jacobian is assumed to be banded. These give the number of
  111. lower and upper non-zero diagonals in this banded matrix.
  112. For the banded case, `Dfun` should return a matrix whose
  113. rows contain the non-zero bands (starting with the lowest diagonal).
  114. Thus, the return matrix `jac` from `Dfun` should have shape
  115. ``(ml + mu + 1, len(y0))`` when ``ml >=0`` or ``mu >=0``.
  116. The data in `jac` must be stored such that ``jac[i - j + mu, j]``
  117. holds the derivative of the ``i``\\ th equation with respect to the
  118. ``j``\\ th state variable. If `col_deriv` is True, the transpose of
  119. this `jac` must be returned.
  120. rtol, atol : float, optional
  121. The input parameters `rtol` and `atol` determine the error
  122. control performed by the solver. The solver will control the
  123. vector, e, of estimated local errors in y, according to an
  124. inequality of the form ``max-norm of (e / ewt) <= 1``,
  125. where ewt is a vector of positive error weights computed as
  126. ``ewt = rtol * abs(y) + atol``.
  127. rtol and atol can be either vectors the same length as y or scalars.
  128. Defaults to 1.49012e-8.
  129. tcrit : ndarray, optional
  130. Vector of critical points (e.g., singularities) where integration
  131. care should be taken.
  132. h0 : float, (0: solver-determined), optional
  133. The step size to be attempted on the first step.
  134. hmax : float, (0: solver-determined), optional
  135. The maximum absolute step size allowed.
  136. hmin : float, (0: solver-determined), optional
  137. The minimum absolute step size allowed.
  138. ixpr : bool, optional
  139. Whether to generate extra printing at method switches.
  140. mxstep : int, (0: solver-determined), optional
  141. Maximum number of (internally defined) steps allowed for each
  142. integration point in t.
  143. mxhnil : int, (0: solver-determined), optional
  144. Maximum number of messages printed.
  145. mxordn : int, (0: solver-determined), optional
  146. Maximum order to be allowed for the non-stiff (Adams) method.
  147. mxords : int, (0: solver-determined), optional
  148. Maximum order to be allowed for the stiff (BDF) method.
  149. See Also
  150. --------
  151. solve_ivp : solve an initial value problem for a system of ODEs
  152. ode : a more object-oriented integrator based on VODE
  153. quad : for finding the area under a curve
  154. Examples
  155. --------
  156. The second order differential equation for the angle `theta` of a
  157. pendulum acted on by gravity with friction can be written::
  158. theta''(t) + b*theta'(t) + c*sin(theta(t)) = 0
  159. where `b` and `c` are positive constants, and a prime (') denotes a
  160. derivative. To solve this equation with `odeint`, we must first convert
  161. it to a system of first order equations. By defining the angular
  162. velocity ``omega(t) = theta'(t)``, we obtain the system::
  163. theta'(t) = omega(t)
  164. omega'(t) = -b*omega(t) - c*sin(theta(t))
  165. Let `y` be the vector [`theta`, `omega`]. We implement this system
  166. in Python as:
  167. >>> import numpy as np
  168. >>> def pend(y, t, b, c):
  169. ... theta, omega = y
  170. ... dydt = [omega, -b*omega - c*np.sin(theta)]
  171. ... return dydt
  172. ...
  173. We assume the constants are `b` = 0.25 and `c` = 5.0:
  174. >>> b = 0.25
  175. >>> c = 5.0
  176. For initial conditions, we assume the pendulum is nearly vertical
  177. with `theta(0)` = `pi` - 0.1, and is initially at rest, so
  178. `omega(0)` = 0. Then the vector of initial conditions is
  179. >>> y0 = [np.pi - 0.1, 0.0]
  180. We will generate a solution at 101 evenly spaced samples in the interval
  181. 0 <= `t` <= 10. So our array of times is:
  182. >>> t = np.linspace(0, 10, 101)
  183. Call `odeint` to generate the solution. To pass the parameters
  184. `b` and `c` to `pend`, we give them to `odeint` using the `args`
  185. argument.
  186. >>> from scipy.integrate import odeint
  187. >>> sol = odeint(pend, y0, t, args=(b, c))
  188. The solution is an array with shape (101, 2). The first column
  189. is `theta(t)`, and the second is `omega(t)`. The following code
  190. plots both components.
  191. >>> import matplotlib.pyplot as plt
  192. >>> plt.plot(t, sol[:, 0], 'b', label='theta(t)')
  193. >>> plt.plot(t, sol[:, 1], 'g', label='omega(t)')
  194. >>> plt.legend(loc='best')
  195. >>> plt.xlabel('t')
  196. >>> plt.grid()
  197. >>> plt.show()
  198. """
  199. if ml is None:
  200. ml = -1 # changed to zero inside function call
  201. if mu is None:
  202. mu = -1 # changed to zero inside function call
  203. dt = np.diff(t)
  204. if not ((dt >= 0).all() or (dt <= 0).all()):
  205. raise ValueError("The values in t must be monotonically increasing "
  206. "or monotonically decreasing; repeated values are "
  207. "allowed.")
  208. t = copy(t)
  209. y0 = copy(y0)
  210. output = _odepack.odeint(func, y0, t, args, Dfun, col_deriv, ml, mu,
  211. full_output, rtol, atol, tcrit, h0, hmax, hmin,
  212. ixpr, mxstep, mxhnil, mxordn, mxords,
  213. int(bool(tfirst)))
  214. if output[-1] < 0:
  215. warning_msg = (f"{_msgs[output[-1]]} Run with full_output = 1 to "
  216. f"get quantitative information.")
  217. warnings.warn(warning_msg, ODEintWarning, stacklevel=2)
  218. elif printmessg:
  219. warning_msg = _msgs[output[-1]]
  220. warnings.warn(warning_msg, ODEintWarning, stacklevel=2)
  221. if full_output:
  222. output[1]['message'] = _msgs[output[-1]]
  223. output = output[:-1]
  224. if len(output) == 1:
  225. return output[0]
  226. else:
  227. return output