solve_ivp

xitorch.integrate.solve_ivp(fcn: Union[Callable[[], torch.Tensor], Callable[[], Sequence[torch.Tensor]]], ts: torch.Tensor, y0: torch.Tensor, params: Sequence[Any] = [], bck_options: Mapping[str, Any] = {}, method: Optional[Union[str, Callable]] = None, **fwd_options) → Union[torch.Tensor, Sequence[torch.Tensor]][source]

Solve the initial value problem (IVP) or also commonly known as ordinary differential equations (ODE), where given the initial value \(\mathbf{y_0}\), it then solves

\[\mathbf{y}(t) = \mathbf{y_0} + \int_{t_0}^{t} \mathbf{f}(t', \mathbf{y}, \theta)\ \mathrm{d}t'\]
Parameters
  • fcn (callable) – The function that represents dy/dt. The function takes an input of a single time t and tensor y with shape (*ny) and produce \(\mathrm{d}\mathbf{y}/\mathrm{d}t\) with shape (*ny). The output of the function must be a tensor with shape (*ny) or a list of tensors.

  • ts (torch.tensor) – The time points where the value of y will be returned. It must be monotonically increasing or decreasing. It is a tensor with shape (nt,).

  • y0 (torch.tensor) – The initial value of y, i.e. y(t[0]) == y0. It is a tensor with shape (*ny) or a list of tensors.

  • params (list) – Sequence of other parameters required in the function.

  • bck_options (dict) – Options for the backward solve_ivp method. If not specified, it will take the same options as fwd_options.

  • method (str or callable or None) – Initial value problem solver. If None, it will choose "rk45".

  • **fwd_options – Method-specific option (see method section below).

Returns

The values of y for each time step in ts. It is a tensor with shape (nt,*ny) or a list of tensors

Return type

torch.tensor or a list of tensors

method="rk45"
solve_ivp(..., method="rk45", *, atol=1e-08, rtol=1e-05)

Perform the adaptive Runge-Kutta steps with order 4 and 5.

Keyword Arguments
  • atol (float) – The absolute error tolerance in deciding the steps

  • rtol (float) – The relative error tolerance in deciding the steps

method="rk23"
solve_ivp(..., method="rk23", *, atol=1e-08, rtol=1e-05)

Perform the adaptive Runge-Kutta steps with order 2 and 3.

Keyword Arguments
  • atol (float) – The absolute error tolerance in deciding the steps

  • rtol (float) – The relative error tolerance in deciding the steps

method="rk4"
solve_ivp(..., method="rk4")

Perform the Runge-Kutta steps of order 4 with a fixed step size.