rkMethod {deSolve} | R Documentation |
This function returns a list specifying coefficients and properties of ODE solver methods from the Runge-Kutta family.
rkMethod(method = NULL, ...)
method |
a string constant naming one of the pre-defined methods
of the Runge-Kutta family of solvers. The most common methods are
the fixed-step methods "euler" , "rk2" , "rk4" or
the variable step methods "rk23bs" (alias "ode23" ) or
"rk45dp7" (alias "ode45" ).
|
... |
specification of a user-defined solver, see Value and example below. |
This function supplies method
settings for rk
or
ode
. If called without arguments, the names of all
implemented solvers of the Runge-Kutta family are returned.
The following comparison gives an idea how the algorithms of deSolve are related to similar algorithms of other simulation languages:
rkMethod | | | Description |
"euler" | | | Euler's Method |
"rk2" | | | 2nd order Runge-Kutta, fixed time step (Heun's method) |
"rk4" | | | classical 4th order Runge-Kutta, fixed time step |
"rk23" | | | Runge-Kutta, order 2(3), Octave: ode23 |
"rk23bs", "ode23" | | | Bogacki-Shampine, order 2(3), Matlab: ode23 |
"rk34f" | | | Runge-Kutta-Fehlberg, order 3(4) |
"rk45ck" | | | Runge-Kutta Cash-Karp, order 4(5) |
"rk45f" | | | Runge-Kutta-Fehlberg, order 4(5), Octave: ode45, pair=1 |
"rk45e" | | | Runge-Kutta-England, order 4(5) |
"rk45dp6" | | | Dormand-Prince, order 4(5), local order 6 |
"rk45dp7", "ode45" | | | Dormand-Prince 4(5), local order 7 |
| | (also known as dopri5, MATLAB: ode45, Octave: ode45, pair=0) |
Note that this table is based on the Runge-Kutta coefficients only, but the algorithms differ also in their implementation, in their stepsize adaption strategy and interpolation methods. The table reflects the state at time of writing and it is of course possible that implementations change.
A list with the following elements:
ID |
name of the method (character) |
varstep |
boolean value specifying if the method allows for
variable time step (TRUE ) or not (FALSE ).
|
FSAL |
(first same as last) optional boolean value specifying if
the method allows re-use of the last function evaluation
(TRUE ) or not (FALSE or NULL ).
|
A |
coefficient matrix of the method. As link{rk} supports
only explicit methods, this matrix must be lower triangular.
A must be a vector for fixed step methods where only the
subdiagonal values are different from zero.
|
b1 |
weighting coefficients for averaging the function evaluations of method 1. |
b2 |
weighting coefficients for averaging the function evaluations of method 2 (optional, for embedded methods that allow variable time step). |
c |
coefficients for calculating the intermediate time steps. |
d |
optional coefficients for built-in polynomial interpolation
of the outputs from internal steps (dense output), currently only
available for method rk45dp7 (Dormand-Prince).
|
stage |
number of function evaluations needed (corresponds to number of rows in A). |
Qerr |
global error order of the method, important for automatic time-step adjustment. |
nknots |
integer value specifying the order of interpolation polynomials
for methods without dense output (default = 5th order). If nknots < 2
internal
interpolation is switched off. Note that this works only if time steps match
exactly. This may not work as expected for all cases because of
floating point rounding errors (see R FAQ about ``Why doesn't R
think these numbers are equal?'').
|
alpha |
optional tuning parameter for stepsize
adjustment. If alpha is omitted, it is set to
1/Qerr - 0.75 beta. The default value is
1/Qerr (for beta = 0). |
beta |
optional tuning parameter for stepsize adjustment. Typical values are 0 (default) or 0.4/Qerr |
"rk23"
(alias
"ode23"
) for simple problems and "rk45dp7"
(alias
"ode45"
) for rough problems. The default solver is
"rk45dp7"
(alias "ode45"), because of its relatively high
order (4), re-use of the last intermediate steps (FSAL = first
same as last) and built-in polynomial interpolation (dense
output).
"rk23bs"
, that supports also FSAL, may be useful for
slightly stiff systems if demands on precision are relatively low.
"rk45ck"
.
"rk4"
is traditionally used in cases where an
adequate stepsize is known a-priori or if external forcing data
are provided for fixed time steps only and frequent interpolation
of external data needs to be avoided.
"rk45dp7"
(alias "ode45"
) contains an
efficient built-in interpolation scheme (dense output) based on
intermediate function evaluations. Interpolation for all other
methods is done with Neville-Aitken polynomials (5th order by
default) which needs at least 6 internal time steps. Local
interpolation inaccuracies can be minimized if external and
internal time steps are not too different.
Thomas Petzoldt thomas.petzoldt@tu-dresden.de
Bogacki, P. and Shampine L.F. (1989) A 3(2) pair of Runge-Kutta formulas, Appl. Math. Lett. 2, 1–9.
Butcher, J. C. (1987) The numerical analysis of ordinary differential equations, Runge-Kutta and general linear methods, Wiley, Chichester and New York.
Cash, J. R. and Karp A. H., 1990. A variable order Runge-Kutta method for initial value problems with rapidly varying right-hand sides, ACM Transactions on Mathematical Software 16, 201–222.
Dormand, J. R. and Prince, P. J. (1980) A family of embedded Runge-Kutta formulae, J. Comput. Appl. Math. 6(1), 19–26.
Dormand, J. R. and Prince, P. J. (1981) High order embedded Runge-Kutta formulae, J. Comput. Appl. Math. 7(1), 67–75.
Engeln-Muellges, G. and Reutter, F. (1996) Numerik Algorithmen: Entscheidungshilfe zur Auswahl und Nutzung. VDI Verlag, Duesseldorf.
Fehlberg, E. (1967) Klassische Runge-Kutta-Formeln fuenfter and siebenter Ordnung mit Schrittweiten-Kontrolle, Computing (Arch. Elektron. Rechnen) 4, 93–106.
Kutta, W. (1901) Beitrag zur naeherungsweisen Integration totaler Differentialgleichungen, Z. Math. Phys. 46, 435–453.
Octave-Forge - Extra Packages for GNU Octave, Package OdePkg. http://octave.sourceforge.net/doc/odepkg.html
Runge, C. (1895) Ueber die numerische Aufloesung von Differentialgleichungen, Math. Ann. 46, 167–178.
MATLAB (R) is a registed property of The Mathworks Inc. http://www.mathworks.com/
rkMethod() # returns the names of all available methods rkMethod("rk45dp7") # parameters of the Dormand-Prince 5(4) method rkMethod("ode45") # an alias for the same method func <- function(t, x, parms) { with(as.list(c(parms, x)),{ dP <- a * P - b * C * P dC <- b * P * C - c * C res <- c(dP, dC) list(res) }) } times <- seq(0, 200, length = 101) parms <- c(a = 0.1, b = 0.1, c = 0.1) x <- c(P = 2, C = 1) ode(x, times, func, parms, method = rkMethod("rk4")) ode(x, times, func, parms, method = "ode45") o0 <- ode(x, times, func, parms, method = "lsoda") o1 <- ode(x, times, func, parms, method = rkMethod("rk45dp7")) ## disable dense-output interpolation ## and fall back to Neville-Aitken polynomials o2 <- ode(x, times, func, parms, method = rkMethod("rk45dp7", d = NULL)) ## show differences between methods par(mfrow=c(3,1)) matplot(o1[,1], o1[,-1], type="l", xlab="Time", main="State Variables") matplot(o1[,1], o2[,-1], type="p", pch=16, add=TRUE) matplot(o0[,1], o1[,-1]-o0[,-1], type="l", lty="solid", xlab="Time", main="Difference: rk45dp7 - lsoda") matplot(o1[,1], o2[,-1]-o0[,-1], type="l", lty="dashed", xlab="Time", main="difference: Neville-Aitken - Dense Output") abline(h=0, col="grey") ## define and use a new rk method ode(x, times, func, parms, method = rkMethod(ID = "midpoint", varstep = FALSE, A = c(0, 1/2), b1 = c(0, 1), c = c(0, 1/2), stage = 2, Qerr = 1 ) )