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 [R838579b36be5-1]. It switches automatically between the nonstiff Adams method and the stiff BDF method. The method was originally detailed in [R838579b36be5-2]. - Parameters
- funcallable
- Right-hand side of the system. The calling signature is - fun(t, y). Here- tis a scalar, and there are two options for the ndarray- y: It can either have shape (n,); then- funmust return array_like with shape (n,). Alternatively it can have shape (n, k); then- funmust return an array_like with shape (n, k), i.e. each column corresponds to a single column in- y. The choice between the two options is determined by vectorized argument (see below). The vectorized implementation allows a faster approximation of the Jacobian by finite differences (required for this solver).
- t0float
- Initial time. 
- y0array_like, shape (n,)
- Initial state. 
- t_boundfloat
- Boundary time - the integration won’t continue beyond it. It also determines the direction of the integration. 
- first_stepfloat or None, optional
- Initial step size. Default is - Nonewhich means that the algorithm should choose.
- min_stepfloat, optional
- Minimum allowed step size. Default is 0.0, i.e. the step size is not bounded and determined solely by the solver. 
- max_stepfloat, optional
- Maximum allowed step size. Default is np.inf, i.e. the step size is not bounded and determined solely by the solver. 
- rtol, atolfloat 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, 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.
- jacNone 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 to- d f_i / d y_j. The function will be called as- jac(t, y). If None (default), 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, ubandint or None
- Parameters defining the bandwidth of the Jacobian, i.e., - 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 have- ncolumns and- uband + lband + 1rows in which Jacobian diagonals are written. Specifically- jac_packed[uband + i - j , j] = jac[i, j]. The same format is used in- scipy.linalg.solve_banded(check for an illustration). These parameters can be also used with- jac=Noneto reduce the number of Jacobian elements estimated by finite differences.
- vectorizedbool, optional
- Whether fun is implemented in a vectorized fashion. A vectorized implementation offers no advantages for this solver. Default is False. 
 
 - References - R838579b36be5-1
- A. C. Hindmarsh, “ODEPACK, A Systematized Collection of ODE Solvers,” IMACS Transactions on Scientific Computation, Vol 1., pp. 55-64, 1983. 
- R838579b36be5-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
- nint
- Number of equations. 
- statusstring
- Current status of the solver: ‘running’, ‘finished’ or ‘failed’. 
- t_boundfloat
- Boundary time. 
- directionfloat
- Integration direction: +1 or -1. 
- tfloat
- Current time. 
- yndarray
- Current state. 
- t_oldfloat
- Previous time. None if no steps were made yet. 
- nfevint
- Number of evaluations of the right-hand side. 
- njevint
- Number of evaluations of the Jacobian. 
 
 - Methods - dense_output(self)- Compute a local interpolant over the last successful step. - step(self)- Perform one integration step. 
