brainpy.integrators.ode.RK2#

class brainpy.integrators.ode.RK2(f, beta=0.6666666666666666, var_type=None, dt=None, name=None, show_code=False, state_delays=None, neutral_delays=None)[source]#

Generic second order Runge-Kutta method for ODEs.

Derivation

In the brainpy.integrators.ode.midpoint(), brainpy.integrators.ode.heun2(), and brainpy.integrators.ode.ralston2(), we have already known first-order Euler method brainpy.integrators.ode.euler() has big estimation error.

Here, we seek to derive a generic second order Runge-Kutta method 1 for the given ODE system with a given initial value,

\[y'(t) = f(t,y(t)), \qquad y(t_0)=y_0,\]

we want to get a generic solution:

\[\begin{align} y_{n+1} &= y_{n} + h \left ( a_1 K_1 + a_2 K_2 \right ) \tag{1} \end{align}\]

where \(a_1\) and \(a_2\) are some weights to be determined, and \(K_1\) and \(K_2\) are derivatives on the form:

\[\begin{align} K_1 & = f(t_n,y_n) \qquad \text{and} \qquad K_2 = f(t_n + p_1 h,y_n + p_2 K_1 h ) \tag{2} \end{align}\]

By substitution of (2) in (1) we get:

\[\begin{align} y_{n+1} &= y_{n} + a_1 h f(t_n,y_n) + a_2 h f(t_n + p_1 h,y_n + p_2 K_1 h) \tag{3} \end{align}\]

Now, we may find a Taylor-expansion of \(f(t_n + p_1 h, y_n + p_2 K_1 h )\)

\[\begin{split}\begin{align} f(t_n + p_1 h, y_n + p_2 K_1 h ) &= f + p_1 h f_t + p_2 K_1 h f_y + \text{h.o.t.} \nonumber \\ & = f + p_1 h f_t + p_2 h f f_y + \text{h.o.t.} \tag{4} \end{align}\end{split}\]

where \(f_t \equiv \frac{\partial f}{\partial t}\) and \(f_y \equiv \frac{\partial f}{\partial y}\).

By substitution of (4) in (3) we eliminate the implicit dependency of \(y_{n+1}\)

\[\begin{split}\begin{align} y_{n+1} &= y_{n} + a_1 h f(t_n,y_n) + a_2 h \left (f + p_1 h f_t + p_2 h f f_y \right ) \nonumber \\ &= y_{n} + (a_1 + a_2) h f + \left (a_2 p_1 f_t + a_2 p_2 f f_y \right) h^2 \tag{5} \end{align}\end{split}\]

In the next, we try to get the second order Taylor expansion of the solution:

\[\begin{align} y(t_n+h) = y_n + h y' + \frac{h^2}{2} y'' + O(h^3) \tag{6} \end{align}\]

where the second order derivative is given by

\[\begin{align} y'' = \frac{d^2 y}{dt^2} = \frac{df}{dt} = \frac{\partial{f}}{\partial{t}} \frac{dt}{dt} + \frac{\partial{f}}{\partial{y}} \frac{dy}{dt} = f_t + f f_y \tag{7} \end{align}\]

Substitution of (7) into (6) yields:

\[\begin{align} y(t_n+h) = y_n + h f + \frac{h^2}{2} \left (f_t + f f_y \right ) + O(h^3) \tag{8} \end{align}\]

Finally, in order to approximate (8) by using (5), we get the generic second order Runge-Kutta method, where

\[\begin{split}\begin{aligned} a_1 + a_2 = 1 \\ a_2 p_1 = \frac{1}{2} \\ a_2 p_2 = \frac{1}{2}. \end{aligned}\end{split}\]

Furthermore, let \(p_1=\beta\), we get

\[\begin{split}\begin{aligned} p_1 = & \beta \\ p_2 = & \beta \\ a_2 = &\frac{1}{2\beta} \\ a_1 = &1 - \frac{1}{2\beta} . \end{aligned}\end{split}\]

Therefore, the corresponding Butcher tableau is:

\[\begin{split}\begin{array}{c|cc} 0 & 0 & 0 \\ \beta & \beta & 0 \\ \hline & 1 - {1 \over 2 * \beta} & {1 \over 2 * \beta} \end{array}\end{split}\]

References

1

Chapra, Steven C., and Raymond P. Canale. Numerical methods for engineers. Vol. 1221. New York: Mcgraw-hill, 2011.

__init__(f, beta=0.6666666666666666, var_type=None, dt=None, name=None, show_code=False, state_delays=None, neutral_delays=None)[source]#

Methods

__init__(f[, beta, var_type, dt, name, ...])

build()

cpu()

Move all variable into the CPU device.

cuda()

Move all variables into the GPU device.

load_state_dict(state_dict[, warn])

Copy parameters and buffers from state_dict into this module and its descendants.

load_states(filename[, verbose])

Load the model states.

nodes([method, level, include_self])

Collect all children nodes.

register_implicit_nodes(*nodes[, node_cls])

register_implicit_vars(*variables, ...)

save_states(filename[, variables])

Save the model states.

set_integral(f)

Set the integral function.

state_dict()

Returns a dictionary containing a whole state of the module.

to(device)

Moves all variables into the given device.

tpu()

Move all variables into the TPU device.

train_vars([method, level, include_self])

The shortcut for retrieving all trainable variables.

tree_flatten()

Flattens the object as a PyTree.

tree_unflatten(aux, dynamic_values)

New in version 2.3.1.

unique_name([name, type_])

Get the unique name for this object.

vars([method, level, include_self, ...])

Collect all variables in this node and the children nodes.

Attributes

A

B

C

arguments

All arguments when calling the numer integrator of the differential equation.

dt

The numerical integration precision.

integral

The integral function.

name

Name of the model.

neutral_delays

neutral delays.

parameters

The parameters defined in the differential equation.

state_delays

State delays.

variables

The variables defined in the differential equation.