scipy.integrate.LSODA¶
-
class
scipy.integrate.
LSODA
(fun, t0, y0, t_bound, first_step=None, min_step=0.0, max_step=inf, rtol=0.001, atol=1e-06, jac=None, lband=None, uband=None, vectorized=False, **extraneous)[source]¶ Adams/BDF method with automatic stiffness detection and switching.
This is a wrapper to the Fortran solver from ODEPACK [R56]. It switches automatically between the nonstiff Adams method and the stiff BDF method. The method was originally detailed in [R57].
Parameters: fun : callable
Right-hand side of the system. The calling signature is
fun(t, y)
. Heret
is a scalar and there are two options for ndarrayy
. It can either have shape (n,), thenfun
must return array_like with shape (n,). Or alternatively it can have shape (n, k), thenfun
must return array_like with shape (n, k), i.e. each column corresponds to a single column iny
. The choice between the two options is determined by vectorized argument (see below). The vectorized implementation allows faster approximation of the Jacobian by finite differences.t0 : float
Initial time.
y0 : array_like, shape (n,)
Initial state.
t_bound : float
Boundary time — the integration won’t continue beyond it. It also determines the direction of the integration.
first_step : float or None, optional
Initial step size. Default is
None
which means that the algorithm should choose.min_step : float, optional
Minimum allowed step size. Default is 0.0, i.e. the step is not bounded and determined solely by the solver.
max_step : float, optional
Maximum allowed step size. Default is
np.inf
, i.e. the step is not bounded and determined solely by the solver.rtol, atol : float and array_like, optional
Relative and absolute tolerances. The solver keeps the local error estimates less than
atol + rtol * abs(y)
. Here rtol controls a relative accuracy (number of correct digits). But if a component of y is approximately below atol then the error only needs to fall within the same atol threshold, and the number of correct digits is not guaranteed. If components of y have different scales, it might be beneficial to set different atol values for different components by passing array_like with shape (n,) for atol. Default values are 1e-3 for rtol and 1e-6 for atol.jac : None or callable, optional
Jacobian matrix of the right-hand side of the system with respect to
y
. The Jacobian matrix has shape (n, n) and its element (i, j) is equal tod f_i / d y_j
. The function will be called asjac(t, y)
. If None (default), then the Jacobian will be approximated by finite differences. It is generally recommended to provide the Jacobian rather than relying on a finite difference approximation.lband, uband : int or None, optional
Jacobian band width:
jac[i, j] != 0 only for i - lband <= j <= i + uband
. Setting these requires your jac routine to return the Jacobian in the packed format: the returned array must haven
columns anduband + lband + 1
rows in which Jacobian diagonals are written. Specificallyjac_packed[uband + i - j , j] = jac[i, j]
. The same format is used inscipy.linalg.solve_banded
(check for an illustration). These parameters can be also used withjac=None
to reduce the number of Jacobian elements estimated by finite differences.vectorized : bool, optional
Whether fun is implemented in a vectorized fashion. A vectorized implementation offers no advantages for this solver. Default is False.
References
[R56] (1, 2) A. C. Hindmarsh, “ODEPACK, A Systematized Collection of ODE Solvers,” IMACS Transactions on Scientific Computation, Vol 1., pp. 55-64, 1983. [R57] (1, 2) L. Petzold, “Automatic selection of methods for solving stiff and nonstiff systems of ordinary differential equations”, SIAM Journal on Scientific and Statistical Computing, Vol. 4, No. 1, pp. 136-148, 1983. Attributes
n (int) Number of equations. status (string) Current status of the solver: ‘running’, ‘finished’ or ‘failed’. t_bound (float) Boundary time. direction (float) Integration direction: +1 or -1. t (float) Current time. y (ndarray) Current state. t_old (float) Previous time. None if no steps were made yet. nfev (int) Number of the system’s rhs evaluations. njev (int) Number of the Jacobian evaluations. Methods
dense_output
()Compute a local interpolant over the last successful step. step
()Perform one integration step.