scipy.integrate.Radau¶
-
class
scipy.integrate.
Radau
(fun, t0, y0, t_bound, max_step=inf, rtol=0.001, atol=1e-06, jac=None, jac_sparsity=None, vectorized=False, **extraneous)[source]¶ Implicit Runge-Kutta method of Radau IIA family of order 5.
Implementation follows [R61]. The error is controlled for a 3rd order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
Parameters: fun : callable
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, k), thenfun
must return array_like with shape (n, k), i.e. each column corresponds to a single column iny
. The choice between the two options is determined by vectorized argument (see below). The vectorized implementation allows faster approximation of the Jacobian by finite differences.t0 : float
Initial time.
y0 : array_like, shape (n,)
Initial state.
t_bound : float
Boundary time — the integration won’t continue beyond it. It also determines the direction of the integration.
max_step : float, optional
Maximum allowed step size. Default is np.inf, i.e. the step is not bounded and determined solely by the solver.
rtol, atol : float 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 then 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.jac : {None, array_like, sparse_matrix, callable}, optional
Jacobian matrix of the right-hand side of the system with respect to y, required only by ‘Radau’ and ‘BDF’ methods. The Jacobian matrix has shape (n, n) and its element (i, j) is equal to
d f_i / d y_j
. There are 3 ways to define the Jacobian:- If array_like or sparse_matrix, then the Jacobian is assumed to be constant.
- If callable, then the Jacobian is assumed to depend on both
t and y, and will be called as
jac(t, y)
as necessary. The return value might be a sparse matrix. - If None (default), then the Jacobian will be approximated by finite differences.
It is generally recommended to provide the Jacobian rather than relying on a finite difference approximation.
jac_sparsity : {None, array_like, sparse matrix}, optional
Defines a sparsity structure of the Jacobian matrix for a finite difference approximation, its shape must be (n, n). If the Jacobian has only few non-zero elements in each row, providing the sparsity structure will greatly speed up the computations [R62]. A zero entry means that a corresponding element in the Jacobian is identically zero. If None (default), the Jacobian is assumed to be dense.
vectorized : bool, optional
Whether fun is implemented in a vectorized fashion. Default is False.
References
[R61] (1, 2) E. Hairer, G. Wanner, “Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems”, Sec. IV.8. [R62] (1, 2) A. Curtis, M. J. D. Powell, and J. Reid, “On the estimation of sparse Jacobian matrices”, Journal of the Institute of Mathematics and its Applications, 13, pp. 117-120, 1974. Attributes
n (int) Number of equations. status (string) Current status of the solver: ‘running’, ‘finished’ or ‘failed’. t_bound (float) Boundary time. direction (float) Integration direction: +1 or -1. t (float) Current time. y (ndarray) Current state. t_old (float) Previous time. None if no steps were made yet. step_size (float) Size of the last successful step. None if no steps were made yet. nfev (int) Number of the system’s rhs evaluations. njev (int) Number of the Jacobian evaluations. nlu (int) Number of LU decompositions. Methods
dense_output
()Compute a local interpolant over the last successful step. step
()Perform one integration step.