lsoda.py 7.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188
  1. import numpy as np
  2. from scipy.integrate import ode
  3. from .common import validate_tol, validate_first_step, warn_extraneous
  4. from .base import OdeSolver, DenseOutput
  5. class LSODA(OdeSolver):
  6. """Adams/BDF method with automatic stiffness detection and switching.
  7. This is a wrapper to the Fortran solver from ODEPACK [1]_. It switches
  8. automatically between the nonstiff Adams method and the stiff BDF method.
  9. The method was originally detailed in [2]_.
  10. Parameters
  11. ----------
  12. fun : callable
  13. Right-hand side of the system. The calling signature is ``fun(t, y)``.
  14. Here ``t`` is a scalar, and there are two options for the ndarray ``y``:
  15. It can either have shape (n,); then ``fun`` must return array_like with
  16. shape (n,). Alternatively it can have shape (n, k); then ``fun``
  17. must return an array_like with shape (n, k), i.e. each column
  18. corresponds to a single column in ``y``. The choice between the two
  19. options is determined by `vectorized` argument (see below). The
  20. vectorized implementation allows a faster approximation of the Jacobian
  21. by finite differences (required for this solver).
  22. t0 : float
  23. Initial time.
  24. y0 : array_like, shape (n,)
  25. Initial state.
  26. t_bound : float
  27. Boundary time - the integration won't continue beyond it. It also
  28. determines the direction of the integration.
  29. first_step : float or None, optional
  30. Initial step size. Default is ``None`` which means that the algorithm
  31. should choose.
  32. min_step : float, optional
  33. Minimum allowed step size. Default is 0.0, i.e. the step size is not
  34. bounded and determined solely by the solver.
  35. max_step : float, optional
  36. Maximum allowed step size. Default is np.inf, i.e. the step size is not
  37. bounded and determined solely by the solver.
  38. rtol, atol : float and array_like, optional
  39. Relative and absolute tolerances. The solver keeps the local error
  40. estimates less than ``atol + rtol * abs(y)``. Here `rtol` controls a
  41. relative accuracy (number of correct digits). But if a component of `y`
  42. is approximately below `atol`, the error only needs to fall within
  43. the same `atol` threshold, and the number of correct digits is not
  44. guaranteed. If components of y have different scales, it might be
  45. beneficial to set different `atol` values for different components by
  46. passing array_like with shape (n,) for `atol`. Default values are
  47. 1e-3 for `rtol` and 1e-6 for `atol`.
  48. jac : None or callable, optional
  49. Jacobian matrix of the right-hand side of the system with respect to
  50. ``y``. The Jacobian matrix has shape (n, n) and its element (i, j) is
  51. equal to ``d f_i / d y_j``. The function will be called as
  52. ``jac(t, y)``. If None (default), the Jacobian will be
  53. approximated by finite differences. It is generally recommended to
  54. provide the Jacobian rather than relying on a finite-difference
  55. approximation.
  56. lband, uband : int or None
  57. Parameters defining the bandwidth of the Jacobian,
  58. i.e., ``jac[i, j] != 0 only for i - lband <= j <= i + uband``. Setting
  59. these requires your jac routine to return the Jacobian in the packed format:
  60. the returned array must have ``n`` columns and ``uband + lband + 1``
  61. rows in which Jacobian diagonals are written. Specifically
  62. ``jac_packed[uband + i - j , j] = jac[i, j]``. The same format is used
  63. in `scipy.linalg.solve_banded` (check for an illustration).
  64. These parameters can be also used with ``jac=None`` to reduce the
  65. number of Jacobian elements estimated by finite differences.
  66. vectorized : bool, optional
  67. Whether `fun` is implemented in a vectorized fashion. A vectorized
  68. implementation offers no advantages for this solver. Default is False.
  69. Attributes
  70. ----------
  71. n : int
  72. Number of equations.
  73. status : string
  74. Current status of the solver: 'running', 'finished' or 'failed'.
  75. t_bound : float
  76. Boundary time.
  77. direction : float
  78. Integration direction: +1 or -1.
  79. t : float
  80. Current time.
  81. y : ndarray
  82. Current state.
  83. t_old : float
  84. Previous time. None if no steps were made yet.
  85. nfev : int
  86. Number of evaluations of the right-hand side.
  87. njev : int
  88. Number of evaluations of the Jacobian.
  89. References
  90. ----------
  91. .. [1] A. C. Hindmarsh, "ODEPACK, A Systematized Collection of ODE
  92. Solvers," IMACS Transactions on Scientific Computation, Vol 1.,
  93. pp. 55-64, 1983.
  94. .. [2] L. Petzold, "Automatic selection of methods for solving stiff and
  95. nonstiff systems of ordinary differential equations", SIAM Journal
  96. on Scientific and Statistical Computing, Vol. 4, No. 1, pp. 136-148,
  97. 1983.
  98. """
  99. def __init__(self, fun, t0, y0, t_bound, first_step=None, min_step=0.0,
  100. max_step=np.inf, rtol=1e-3, atol=1e-6, jac=None, lband=None,
  101. uband=None, vectorized=False, **extraneous):
  102. warn_extraneous(extraneous)
  103. super(LSODA, self).__init__(fun, t0, y0, t_bound, vectorized)
  104. if first_step is None:
  105. first_step = 0 # LSODA value for automatic selection.
  106. else:
  107. first_step = validate_first_step(first_step, t0, t_bound)
  108. first_step *= self.direction
  109. if max_step == np.inf:
  110. max_step = 0 # LSODA value for infinity.
  111. elif max_step <= 0:
  112. raise ValueError("`max_step` must be positive.")
  113. if min_step < 0:
  114. raise ValueError("`min_step` must be nonnegative.")
  115. rtol, atol = validate_tol(rtol, atol, self.n)
  116. solver = ode(self.fun, jac)
  117. solver.set_integrator('lsoda', rtol=rtol, atol=atol, max_step=max_step,
  118. min_step=min_step, first_step=first_step,
  119. lband=lband, uband=uband)
  120. solver.set_initial_value(y0, t0)
  121. # Inject t_bound into rwork array as needed for itask=5.
  122. solver._integrator.rwork[0] = self.t_bound
  123. solver._integrator.call_args[4] = solver._integrator.rwork
  124. self._lsoda_solver = solver
  125. def _step_impl(self):
  126. solver = self._lsoda_solver
  127. integrator = solver._integrator
  128. # From lsoda.step and lsoda.integrate itask=5 means take a single
  129. # step and do not go past t_bound.
  130. itask = integrator.call_args[2]
  131. integrator.call_args[2] = 5
  132. solver._y, solver.t = integrator.run(
  133. solver.f, solver.jac or (lambda: None), solver._y, solver.t,
  134. self.t_bound, solver.f_params, solver.jac_params)
  135. integrator.call_args[2] = itask
  136. if solver.successful():
  137. self.t = solver.t
  138. self.y = solver._y
  139. # From LSODA Fortran source njev is equal to nlu.
  140. self.njev = integrator.iwork[12]
  141. self.nlu = integrator.iwork[12]
  142. return True, None
  143. else:
  144. return False, 'Unexpected istate in LSODA.'
  145. def _dense_output_impl(self):
  146. iwork = self._lsoda_solver._integrator.iwork
  147. rwork = self._lsoda_solver._integrator.rwork
  148. order = iwork[14]
  149. h = rwork[11]
  150. yh = np.reshape(rwork[20:20 + (order + 1) * self.n],
  151. (self.n, order + 1), order='F').copy()
  152. return LsodaDenseOutput(self.t_old, self.t, h, order, yh)
  153. class LsodaDenseOutput(DenseOutput):
  154. def __init__(self, t_old, t, h, order, yh):
  155. super(LsodaDenseOutput, self).__init__(t_old, t)
  156. self.h = h
  157. self.yh = yh
  158. self.p = np.arange(order + 1)
  159. def _call_impl(self, t):
  160. if t.ndim == 0:
  161. x = ((t - self.t) / self.h) ** self.p
  162. else:
  163. x = ((t - self.t) / self.h) ** self.p[:, None]
  164. return np.dot(self.yh, x)