runsteady {rootSolve}R Documentation

Dynamically runs a system of ordinary differential equations (ODE) to steady-state

Description

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).

Usage

runsteady(y, times=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, ...)

Arguments

y the initial (state) values for the ODE system. If y has a name attribute, the names will be used to label the output matrix.
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 times[2] has been reached or when maxsteps have been performed.
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 func is an R-function, it must be defined as: yprime = func(t, y, parms,...). t is the current time point in the integration, y is the current estimate of the variables in the ODE system. If the initial values y has a names attribute, the names will be available inside func. parms is a vector or list of parameters; ... (optional) are any other arguments passed to the function.
The return value of func should be a list, whose first element is a vector containing the derivatives of y with respect to time, and whose next elements are global values that are required at each point in times.
parms vector or list of parameters used in func or jacfunc.
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 y. See details.
atol absolute error tolerance of integrator, either a scalar or an array as long as y. See details.
jacfunc if not NULL, an R function that computes the jacobian of the system of differential equations dydot(i)/dy(j), or a string giving the name of a function or subroutine in ‘dllname’ that computes the jacobian (see Details below for more about this option). In some circumstances, supplying jacfunc can speed up the computations, if the system is stiff. The R calling sequence for jacfunc is identical to that of func.
If the jacobian is a full matrix, jacfunc should return a matrix dydot/dy, where the ith row contains the derivative of dy_i/dt with respect to y_j, or a vector containing the matrix elements by columns (the way R and Fortran store matrices).
If the jacobian is banded, jacfunc should return a matrix containing only the nonzero bands of the jacobian, rotated row-wise. See first example of lsode.
jactype the structure of the jacobian, one of "fullint", "fullusr", "bandusr" or "bandint" - either full or banded and estimated internally or by user; overruled if mf is not NULL
mf the "method flag" passed to function lsode - overrules jactype - provides more options than jactype - see details
verbose if TRUE: full output to the screen, e.g. will output the settings of vectors *istate* and *rstate* - see details
tcrit if not NULL, then lsode cannot integrate past tcrit. The Fortran routine lsode overshoots its targets (times points in the vector times), and interpolates values for the desired time points. If there is a time beyond which integration should not proceed (perhaps because of a singularity), that should be provided in tcrit.
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 hmin if you don't know why!
hmax an optional maximum value of the integration stepsize. If not specified, hmax is set to the largest difference in times, to avoid that the simulation possibly ignores short-term events. If 0, no maximal size is specified
hini initial step size to be attempted; if 0, the initial step size is determined by the solver
ynames if FALSE: names of state variables are not passed to function func ; this may speed up the simulation
maxord the maximum order to be allowed. NULL uses the default, i.e. order 12 if implicit Adams method (meth=1), order 5 if BDF method (meth=2). Reduce maxord to save storage space
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 maxsteps have been performed or when times[2] has been reached.
dllname a string giving the name of the shared library (without extension) that contains all the compiled function or subroutine definitions refered to in func and jacfunc. See help of stode.
initfunc if not NULL, the name of the initialisation function (which initialises values of parameters), as provided in ‘dllname’. See help of stode.
initpar only when ‘dllname’ is specified and an initialisation function initfunc is in the dll: the parameters passed to the initialiser, to initialise the common blocks (fortran) or global variables (C, C++)
rpar only when ‘dllname’ is specified: a vector with double precision values passed to the dll-functions whose names are specified by func and jacfunc
ipar only when ‘dllname’ is specified: a vector with integer values passed to the dll-functions whose names are specified by func and jacfunc
nout only used if dllname is specified and the model is defined in compiled code: the number of output variables calculated in the compiled function func, present in the shared library. Note: it is not automatically checked whether this is indeed the number of output variables calculed in the dll - you have to perform this check in the code - See package vignette.
outnames only used if ‘dllname’ is specified and nout > 0: the names of output variables calculated in the compiled function func, present in the shared library
... additional arguments passed to func and jacfunc allowing this to be a generic function

Details

The work is done by the Fortran subroutine dlsode, whose documentation should be consulted for details (it is included as comments in the source file ‘src/lsode.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

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

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:

rstate contains the following:

For more information, see the comments in the original code lsode.f

Value

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 y has a names attribute, it will be used to label the output values.
... 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

Author(s)

Karline Soetaert <k.soetaert@nioo.knaw.nl>

References

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.

See Also

stode, for steady-state estimation using the Newton-Raphson method, when the jacobian matrix is banded or full.

stodes, for steady-state estimation using the Newton-Raphson method, when the jacobian matrix is sparse.

steady.band, for steady-state estimation, when the jacobian matrix is banded, and where the state variables need NOT to be rearranged

steady.1D, for steady-state estimation, when the jacobian matrix is banded, and where the state variables need to be rearranged

Examples


# 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))))

rbind("Newton Raphson"=ST$y,"Runsteady"=ST2$y)

[Package rootSolve version 1.3 Index]