SciPy

scipy.integrate.RK45

class scipy.integrate.RK45(fun, t0, y0, t_bound, max_step=inf, rtol=0.001, atol=1e-06, vectorized=False, first_step=None, **extraneous)[source]

Explicit Runge-Kutta method of order 5(4).

This uses the Dormand-Prince pair of formulas [1]. The error is controlled assuming accuracy of the fourth-order method accuracy, but steps are taken using the fifth-order accurate formula (local extrapolation is done). A quartic interpolation polynomial is used for the dense output [2].

Can be applied in the complex domain.

Parameters
funcallable

Right-hand side of the system. The calling signature is fun(t, y). Here t is a scalar, and there are two options for the ndarray y: It can either have shape (n,); then fun must return array_like with shape (n,). Alternatively it can have shape (n, k); then fun must 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).

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 None which means that the algorithm should choose.

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.

vectorizedbool, optional

Whether fun is implemented in a vectorized fashion. Default is False.

References

1

J. R. Dormand, P. J. Prince, “A family of embedded Runge-Kutta formulae”, Journal of Computational and Applied Mathematics, Vol. 6, No. 1, pp. 19-26, 1980.

2

L. W. Shampine, “Some Practical Runge-Kutta Formulas”, Mathematics of Computation,, Vol. 46, No. 173, pp. 135-150, 1986.

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 evaluations of the system’s right-hand side.

njevint

Number of evaluations of the Jacobian. Is always 0 for this solver as it does not use the Jacobian.

nluint

Number of LU decompositions. Is always 0 for this solver.

Methods

dense_output(self)

Compute a local interpolant over the last successful step.

step(self)

Perform one integration step.

Previous topic

scipy.integrate.RK23.step

Next topic

scipy.integrate.RK45.dense_output