scipy.integrate.OdeSolver¶
-
class
scipy.integrate.
OdeSolver
(fun, t0, y0, t_bound, vectorized, support_complex=False)[source]¶ Base class for ODE solvers.
In order to implement a new solver you need to follow the guidelines:
A constructor must accept parameters presented in the base class (listed below) along with any other parameters specific to a solver.
A constructor must accept arbitrary extraneous arguments
**extraneous
, but warn that these arguments are irrelevant using common.warn_extraneous function. Do not pass these arguments to the base class.A solver must implement a private method _step_impl(self) which propagates a solver one step further. It must return tuple
(success, message)
, wheresuccess
is a boolean indicating whether a step was successful, andmessage
is a string containing description of a failure if a step failed or None otherwise.A solver must implement a private method _dense_output_impl(self) which returns a
DenseOutput
object covering the last successful step.A solver must have attributes listed below in Attributes section. Note that
t_old
andstep_size
are updated automatically.Use fun(self, t, y) method for the system rhs evaluation, this way the number of function evaluations (nfev) will be tracked automatically.
For convenience a base class provides fun_single(self, t, y) and fun_vectorized(self, t, y) for evaluating the rhs in non-vectorized and vectorized fashions respectively (regardless of how fun from the constructor is implemented). These calls don’t increment nfev.
If a solver uses a Jacobian matrix and LU decompositions, it should track the number of Jacobian evaluations (njev) and the number of LU decompositions (nlu).
By convention the function evaluations used to compute a finite difference approximation of the Jacobian should not be counted in nfev, thus use fun_single(self, t, y) or fun_vectorized(self, t, y) when computing a finite difference approximation of the Jacobian.
- Parameters
- funcallable
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, n_points), thenfun
must return array_like with shape (n, n_points) (each column corresponds to a single column iny
). The choice between the two options is determined by vectorized argument (see below).- 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.
- vectorizedbool
Whether fun is implemented in a vectorized fashion.
- support_complexbool, optional
Whether integration in a complex domain should be supported. Generally determined by a derived solver class capabilities. Default is False.
- 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.
- step_sizefloat
Size of the last successful step. None if no steps were made yet.
- nfevint
Number of the system’s rhs evaluations.
- njevint
Number of the Jacobian evaluations.
- nluint
Number of LU decompositions.
Methods
dense_output
(self)Compute a local interpolant over the last successful step.
step
(self)Perform one integration step.