runsteady | R Documentation |
Solves the steady-state condition of ordinary differential equations (ODE) in the form:
dy/dt = f(t,y)
by dynamically running till the summed absolute values of the derivatives become smaller than some predefined tolerance.
The R function runsteady
makes use of the FORTRAN ODE solver DLSODE,
written by Alan C. Hindmarsh and Andrew H. Sherman
The system of ODE's is written as an R function or defined in compiled code that has been dynamically loaded. The user has to specify whether or not the problem is stiff and choose the appropriate solution method (e.g. make choices about the type of the Jacobian).
runsteady(y, time = c(0, Inf), func, parms,
stol = 1e-8, rtol = 1e-6, atol = 1e-6,
jacfunc = NULL, jactype = "fullint", mf = NULL,
verbose = FALSE, tcrit = NULL, hmin = 0, hmax = NULL,
hini = 0, ynames = TRUE, maxord = NULL, bandup = NULL,
banddown = NULL, maxsteps = 100000, dllname = NULL,
initfunc = dllname, initpar = parms, rpar = NULL,
ipar = NULL, nout = 0, outnames = NULL,
forcings = NULL, initforc = NULL, fcontrol = NULL,
lrw = NULL, liw = NULL, times = time, ...)
y |
the initial (state) values for the ODE system. If |
time, times |
The simulation time. This should be a 2-valued vector,
consisting of the initial time and the end time.
The last time value should be large enough to make sure that steady-state
is effectively reached in this period.
The simulation will stop either when |
func |
either an R-function that computes the values of the derivatives in the ODE system (the model definition) at time t, or a character string giving the name of a compiled function in a dynamically loaded shared library. If The return value of The derivatives
should be specified in the same order as the state variables |
parms |
vector or list of parameters used in |
stol |
steady-state tolerance; it is assumed that steady-state is reached if the average of absolute values of the derivatives drops below this number. |
rtol |
relative error tolerance of integrator, either a scalar or an
array as long as |
atol |
absolute error tolerance of integrator, either a scalar or an
array as long as |
jacfunc |
if not If the jacobian is a full matrix, If the jacobian is banded, |
jactype |
the structure of the jacobian,
one of "fullint", "fullusr", "bandusr", "bandint", "sparse" - either full,
banded or sparse and estimated internally or by user; overruled if |
mf |
the "method flag" passed to function lsode - overrules
|
verbose |
if |
tcrit |
if not |
hmin |
an optional minimum value of the integration
stepsize. In special situations this parameter may speed up computations
with the cost of precision. Don't use |
hmax |
an optional maximum value of the integration stepsize. If not
specified, |
hini |
initial step size to be attempted; if 0, the initial step size is determined by the solver. |
ynames |
if |
maxord |
the maximum order to be allowed. |
bandup |
number of non-zero bands above the diagonal, in case the jacobian is banded. |
banddown |
number of non-zero bands below the diagonal, in case the jacobian is banded. |
maxsteps |
maximal number of steps. The simulation will stop either
when |
dllname |
a string giving the name of the shared library (without
extension) that contains all the compiled function or subroutine
definitions refered to in |
initfunc |
if not NULL, the name of the initialisation function (which initialises values of parameters), as provided in ‘dllname’. See package vignette. |
initpar |
only when ‘dllname’ is specified and an initialisation
function |
rpar |
only when ‘dllname’ is specified: a vector with double
precision values passed to the dll-functions whose names are specified
by |
ipar |
only when ‘dllname’ is specified: a vector with integer
values passed to the dll-functions whose names are specified by |
nout |
only used if |
outnames |
only used if ‘dllname’ is specified and
|
forcings |
only used if ‘dllname’ is specified: a vector with the
forcing function values, or a list with the forcing function data sets,
each present as a two-columned matrix, with (time,value); interpolation
outside the interval [min( This feature is here for compatibility with models defined in compiled code
from package deSolve; see deSolve's package vignette |
initforc |
if not |
fcontrol |
A list of control parameters for the forcing functions.
See deSolve's package vignette |
lrw |
Only if jactype = 'sparse', the length of the real work array rwork; due to the
sparsicity, this cannot be readily predicted. If For instance, if you get the error: DLSODES- RWORK length is insufficient to proceed. Length needed is .ge. LENRW (=I1), exceeds LRW (=I2) In above message, I1 = 27627 I2 = 25932 set |
liw |
Only if jactype = 'sparse', the length of the integer work array iwork; due to the
sparsicity, this cannot be readily predicted. If |
... |
additional arguments passed to |
The work is done by the Fortran subroutine dlsode
or dlsodes
(if sparse),
whose documentation should be consulted for details (it is included as
comments in the source file ‘src/lsodes.f’). The implementation is
based on the November, 2003 version of lsode, from Netlib.
Before using runsteady
, the user has to decide whether or not the
problem is stiff.
If the problem is nonstiff, use method flag mf
= 10, which selects
a nonstiff (Adams) method, no Jacobian used..
If the problem is stiff, there are four standard choices
which can be specified with jactype
or mf
.
The options for jactype are
jactype = "fullint" : a full jacobian, calculated internally by
lsode
, corresponds to mf
=22.
jactype = "fullusr" : a full jacobian, specified by user function
jacfunc
, corresponds to mf
=21.
jactype = "bandusr" : a banded jacobian, specified by user function
jacfunc
; the size of the bands specified by bandup
and
banddown
, corresponds to mf
=24.
jactype = "bandint" : a banded jacobian, calculated by lsode
;
the size of the bands specified by bandup
and banddown
,
corresponds to mf
=25.
jactype = "sparse" : the soler lsodes
is used, and the sparse jacobian is
calculated by lsodes
- not possible to specify jacfunc
.
More options are available when specifying mf directly.
The legal values of mf
are 10, 11, 12, 13, 14, 15, 20, 21, 22, 23,
24, 25.
mf
is a positive two-digit integer, mf
= (10*METH + MITER),
where
METH indicates the basic linear multistep method: METH = 1 means the implicit Adams method. METH = 2 means the method based on backward differentiation formulas (BDF-s).
MITER indicates the corrector iteration method: MITER = 0 means
functional iteration (no Jacobian matrix is involved).
MITER = 1 means chord iteration with a user-supplied full (NEQ by NEQ)
Jacobian.
MITER = 2 means chord iteration with an internally generated
(difference quotient) full Jacobian
(using NEQ extra calls to func
per df/dy value).
MITER = 3 means chord iteration with an internally generated diagonal
Jacobian approximation
(using 1 extra call to func
per df/dy evaluation).
MITER = 4 means chord iteration with a user-supplied banded Jacobian.
MITER = 5 means chord iteration with an internally generated banded
Jacobian (using ML+MU+1 extra calls to func
per df/dy evaluation).
If MITER = 1 or 4, the user must supply a subroutine jacfunc
.
Inspection of the example below shows how to specify both a banded and full jacobian.
The input parameters rtol
, and atol
determine the error
control performed by the solver.
See stode
for details.
Models may be defined in compiled C or Fortran code, as well as in
an R-function. See function stode
for details.
The output will have the attributes *istate*, and *rstate*, two vectors with several useful elements.
if verbose
= TRUE, the settings of istate and rstate will be written
to the screen.
the following elements of istate are meaningful:
el 1 : gives the conditions under which the last call to the integrator returned. 2 if lsode was successful, -1 if excess work done, -2 means excess accuracy requested. (Tolerances too small), -3 means illegal input detected. (See printed message.), -4 means repeated error test failures. (Check all input), -5 means repeated convergence failures. (Perhaps bad Jacobian supplied or wrong choice of MF or tolerances.), -6 means error weight became zero during problem. (Solution component i vanished, and atol or atol(i) = 0.)
el 12 : The number of steps taken for the problem so far.
el 13 : The number of evaluations for the problem so far.,
el 14 : The number of Jacobian evaluations and LU decompositions so far.,
el 15 : The method order last used (successfully).,
el 16 : The order to be attempted on the next step.,
el 17 : if el 1 =-4,-5: the largest component in the error vector,
rstate contains the following:
1: The step size in t last used (successfully).
2: The step size to be attempted on the next step.
3: The current value of the independent variable which the solver has actually reached, i.e. the current internal mesh point in t.
4: A tolerance scale factor, greater than 1.0, computed when a request for too much accuracy was detected.
For more information, see the comments in the original code lsode.f
A list containing
y |
a vector with the state variable values from the last iteration
during estimation of steady-state condition of the system of equations.
If |
... |
the number of "global" values returned. |
The output will have the attribute steady
, which returns TRUE
,
if steady-state has been reached, the attribute precis
with the
precision attained at the last iteration estimated as the mean absolute
rate of change (sum(abs(dy))/n), the attribute time
with the
simulation time reached and the attribute steps
with the number of
steps performed.
The output will also have the attributes istate
, and rstate
,
two vectors with several useful elements of the dynamic simulation.
See details.
The first element of istate returns the conditions under which the last call
to the integrator returned. Normal is istate[1] = 2
.
If verbose
= TRUE
, the settings of istate and rstate will
be written to the screen.
Karline Soetaert <karline.soetaert@nioz.nl>
Alan C. Hindmarsh, "ODEPACK, A Systematized Collection of ODE Solvers," in Scientific Computing, R. S. Stepleman, et al., Eds. (North-Holland, Amsterdam, 1983), pp. 55-64.
S. C. Eisenstat, M. C. Gursky, M. H. Schultz, and A. H. Sherman, Yale Sparse Matrix Package: I. The Symmetric Codes, Int. J. Num. Meth. Eng., 18 (1982), pp. 1145-1151.
S. C. Eisenstat, M. C. Gursky, M. H. Schultz, and A. H. Sherman, Yale Sparse Matrix Package: II. The Nonsymmetric Codes, Research Report No. 114, Dept. of Computer Sciences, Yale University, 1977.
steady
, for a general interface to most of the steady-state
solvers
steady.band
, to find the steady-state of ODE models with a
banded Jacobian
steady.1D
, steady.2D
,
steady.3D
steady-state solvers for 1-D, 2-D and 3-D
partial differential equations.
stode
, iterative steady-state solver for ODEs with full
or banded Jacobian.
stodes
, iterative steady-state solver for ODEs with arbitrary
sparse Jacobian.
## =======================================================================
## A simple sediment biogeochemical model
## =======================================================================
model<-function(t, y, pars) {
with (as.list(c(y, pars)),{
Min = r*OM
oxicmin = Min*(O2/(O2+ks))
anoxicmin = Min*(1-O2/(O2+ks))* SO4/(SO4+ks2)
dOM = Flux - oxicmin - anoxicmin
dO2 = -oxicmin -2*rox*HS*(O2/(O2+ks)) + D*(BO2-O2)
dSO4 = -0.5*anoxicmin +rox*HS*(O2/(O2+ks)) + D*(BSO4-SO4)
dHS = 0.5*anoxicmin -rox*HS*(O2/(O2+ks)) + D*(BHS-HS)
list(c(dOM, dO2, dSO4, dHS), SumS = SO4+HS)
})
}
# parameter values
pars <- c(D = 1, Flux = 100, r = 0.1, rox = 1,
ks = 1, ks2 = 1, BO2 = 100, BSO4 = 10000, BHS = 0)
# initial conditions
y <- c(OM = 1, O2 = 1, SO4 = 1, HS = 1)
# direct iteration
print( system.time(
ST <- stode(y = y, func = model, parms = pars, pos = TRUE)
))
print( system.time(
ST2 <- runsteady(y = y, func = model, parms = pars, times = c(0, 1000))
))
print( system.time(
ST3 <- runsteady(y = y, func = model, parms = pars, times = c(0, 1000),
jactype = "sparse")
))
rbind("Newton Raphson" = ST$y, "Runsteady" = ST2$y, "Run sparse" = ST3$y)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.