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, first_step=None, **extraneous)[source]#
Implicit Runge-Kutta method of Radau IIA family of order 5.
The implementation follows [1]. The error is controlled with a third-order accurate embedded formula. A cubic polynomial which satisfies the collocation conditions is used for the dense output.
- Parameters:
- funcallable
Right-hand side of the system: the time derivative of the state
y
at timet
. The calling signature isfun(t, y)
, wheret
is a scalar andy
is an ndarray withlen(y) = len(y0)
.fun
must return an array of the same shape asy
. See vectorized for more information.- 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)
. HHere rtol controls a relative accuracy (number of correct digits), while atol controls absolute accuracy (number of correct decimal places). To achieve the desired rtol, set atol to be smaller than the smallest value that can be expected fromrtol * abs(y)
so that rtol dominates the allowable error. If atol is larger thanrtol * abs(y)
the number of correct digits is not guaranteed. Conversely, to achieve the desired atol set rtol such thatrtol * abs(y)
is always smaller than atol. 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 by this method. The Jacobian matrix has shape (n, n) and its element (i, j) is equal to
d f_i / d y_j
. There are three ways to define the Jacobian:If array_like or sparse_matrix, the Jacobian is assumed to be constant.
If callable, the Jacobian is assumed to depend on both t and y; it will be called as
jac(t, y)
as necessary. For the ‘Radau’ and ‘BDF’ methods, the return value might be a sparse matrix.If None (default), 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). This argument is ignored if jac is not None. If the Jacobian has only few non-zero elements in each row, providing the sparsity structure will greatly speed up the computations [2]. A zero entry means that a corresponding element in the Jacobian is always zero. If None (default), the Jacobian is assumed to be dense.
- vectorizedbool, optional
Whether fun can be called in a vectorized fashion. Default is False.
If
vectorized
is False, fun will always be called withy
of shape(n,)
, wheren = len(y0)
.If
vectorized
is True, fun may be called withy
of shape(n, k)
, wherek
is an integer. In this case, fun must behave such thatfun(t, y)[:, i] == fun(t, y[:, i])
(i.e. each column of the returned array is the time derivative of the state corresponding with a column ofy
).Setting
vectorized=True
allows for faster finite difference approximation of the Jacobian by this method, but may result in slower execution overall in some circumstances (e.g. smalllen(y0)
).
References
[1]E. Hairer, G. Wanner, “Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems”, Sec. IV.8.
[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:
- 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 evaluations of the right-hand side.
- njevint
Number of evaluations of the Jacobian.
- nluint
Number of LU decompositions.
Methods
Compute a local interpolant over the last successful step.
step
()Perform one integration step.