Learn R Programming

deTestSet (version 1.1.7.4)

mebdfi: Solver for Differential Algebraic Equations (DAE) up to index 3

Description

Solves either:

  • a system of ordinary differential equations (ODE) of the form $$y' = f(t,y,...)$$ or

  • a system of differential algebraic equations (DAE) of the form $$F(t,y,y') = 0$$ or

  • a system of linearly implicit DAES in the form $$M y' = f(t,y)$$

using the Modified Extended Backward Differentiation formulas for stiff fully implicit inital value problems

These formulas increase the absolute stability regions of the classical BDFs.

The orders of the implemented formulae range from 1 to 8.

The R function mebdfi provides an interface to the Fortran DAE solver of the same name, written by T.J. Abdulla and J.R. Cash.

The system of DE's is written as an R function or can be defined in compiled code that has been dynamically loaded.

Usage

mebdfi(y, times, func = NULL, parms, dy = NULL, res = NULL,
  nind=c(length(y),0,0), rtol = 1e-6, atol = 1e-6, jacfunc = NULL,
  jacres = NULL, jactype = "fullint", mass = NULL, verbose = FALSE,
  tcrit = NULL,  hini = 0, ynames = TRUE, maxord = 7, bandup = NULL,
  banddown = NULL, maxsteps = 5000, dllname = NULL,
  initfunc = dllname, initpar = parms, rpar = NULL,
  ipar = NULL, nout = 0, outnames = NULL,
  forcings=NULL, initforc = NULL, fcontrol=NULL,  ...)

Value

A matrix of class deSolve with up to as many rows as elements in

times and as many columns as elements in y plus the number of "global" values returned in the next elements of the return from func or

res, plus an additional column (the first) for the time value. There will be one row for each element in times unless the Fortran routine `mebdfi' returns with an unrecoverable error. If

y has a names attribute, it will be used to label the columns of the output value.

Arguments

y

the initial (state) values for the DE system. If y has a name attribute, the names will be used to label the output matrix.

times

time sequence for which output is wanted; the first value of times must be the initial time; if only one step is to be taken; set times = NULL.

func

cannot be used if the model is a DAE system. If an ODE system, func should be an R-function that computes the values of the derivatives in the ODE system (the model definition) at time t.

func must be defined as: func <- function(t, y, parms,...).
t is the current time point in the integration, y is the current estimate of the variables in the ODE or DAE system. If the initial values y has a names attribute, the names will be available inside func, unless ynames is FALSE. 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. The derivatives should be specified in the same order as the specification of the state variables y.

Note that it is not possible to define func as a compiled function in a dynamically loaded shared library. Use res instead.

parms

vector or list of parameters used in func, jacfunc, or res

dy

the initial derivatives of the state variables of the DE system. Ignored if an ODE.

res

if a DAE system: either an R-function that computes the residual function F(t,y,y') of the DAE system (the model defininition) at time t, or a character string giving the name of a compiled function in a dynamically loaded shared library.

If res is a user-supplied R-function, it must be defined as: res <- function(t, y, dy, parms, ...).

Here t is the current time point in the integration, y is the current estimate of the variables in the DAE system, dy are the corresponding derivatives. If the initial y or dy have a names attribute, the names will be available inside res, unless ynames is FALSE. parms is a vector of parameters.

The return value of res should be a list, whose first element is a vector containing the residuals of the DAE system, i.e. delta = F(t,y,y'), and whose next elements contain output variables that are required at each point in times.

If res is a string, then dllname must give the name of the shared library (without extension) which must be loaded before mebdfi() is called (see package vignette "compiledCode" for more information).

nind

if a DAE system: a three-valued vector with the number of variables of index 1, 2, 3 respectively. The equations must be defined such that the index 1 variables precede the index 2 variables which in turn precede the index 3 variables. The sum of the variables of different index should equal N, the total number of variables.

rtol

relative error tolerance, either a scalar or a vector, one value for each y,

atol

absolute error tolerance, either a scalar or a vector, one value for each y.

jacfunc

if not NULL, an R function that computes the Jacobian of the system of differential equations. Only used in case the system is an ODE (y' = f(t,y)), specified by func. 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.

jacres

jacres and not jacfunc should be used if the system is specified by the residual function F(t,y,y'), i.e. jacres is used in conjunction with res.

If jacres is an R-function, the calling sequence for jacres is identical to that of res, but with extra parameter cj. Thus it should be defined as: jacres <- function(t, y, dy, parms, cj, ...). Here t is the current time point in the integration, y is the current estimate of the variables in the ODE system, \(y'\) are the corresponding derivatives and cj is a scalar, which is normally proportional to the inverse of the stepsize. If the initial y or dy have a names attribute, the names will be available inside jacres, unless ynames is FALSE. parms is a vector of parameters (which may have a names attribute).

If the Jacobian is a full matrix, jacres should return the matrix dG/dy + cj*dG/dyprime, where the ith row is the sum of the derivatives of \(G_i\) with respect to \(y_j\) and the scaled derivatives of \(G_i\) with respect to \(dy_j\).

If the Jacobian is banded, jacres should return only the nonzero bands of the Jacobian, rotated rowwise. See details for the calling sequence when jacres is a string.

jactype

the structure of the Jacobian, one of "fullint", "fullusr", "bandusr" or "bandint" - either full or banded and estimated internally or by the user.

mass

the mass matrix. If not NULL, the problem is a linearly implicit DAE and defined as \(mass dy/dt = f(t,y)\). The mass-matrix should be of dimension n*n where n is the number of y-values.

If mass=NULL then the model is either an ODE or a DAE, specified with res

verbose

if TRUE: full output to the screen, e.g. will print the diagnostiscs of the integration - see details.

tcrit

the Fortran routine mebdfi 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.

hini

initial step size to be attempted; if 0, the initial step size is set to 1e-6, but it may be better to set it equal to rtol. The solver is quite sensitive to values of hini; sometimes if it fails, it helps to decrease/increase hini

ynames

logical, if FALSE names of state variables are not passed to function func; this may speed up the simulation especially for large models.

maxord

the maximum order to be allowed, an integer between 2 and 7. The default is maxord = 7, but values of 4-5 may be better for difficult problems; hihger order methods are more efficient but less stable.

bandup

number of non-zero bands above the diagonal, in case the Jacobian is banded (and jactype one of "bandint","bandusr")

banddown

number of non-zero bands below the diagonal, in case the Jacobian is banded (and jactype one of "bandint","bandusr")

maxsteps

maximal number of steps per output interval taken by the solver.

dllname

a string giving the name of the shared library (without extension) that contains all the compiled function or subroutine definitions referred to in res and jacres. See package vignette "compiledCode".

initfunc

if not NULL, the name of the initialisation function (which initialises values of parameters), as provided in dllname. See package vignette "compiledCode".

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 res and jacres.

ipar

only when dllname is specified: a vector with integer values passed to the dll-functions whose names are specified by res and jacres.

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 res, 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 "compiledCode".

outnames

only used if dllname is specified and nout > 0: the names of output variables calculated in the compiled function res, present in the shared library. These names will be used to label the output matrix.

forcings

only used if dllname is specified: a list with the forcing function data sets, each present as a two-columned matrix, with (time,value); interpolation outside the interval [min(times), max(times)] is done by taking the value at the closest data extreme.

See package vignette "compiledCode".

initforc

if not NULL, the name of the forcing function initialisation function, as provided in dllname. It MUST be present if forcings has been given a value. See package vignette "compiledCode".

fcontrol

A list of control parameters for the forcing functions. vignette compiledCode from package deSolve.

...

additional arguments passed to func, jacfunc, res and jacres, allowing this to be a generic function.

Author

Karline Soetaert <karline.soetaert@nioz.nl>

Jeff Cash

Details

The mebdfi solver uses modified extended backward differentiation formulas of orders one through eight (specified with maxord) to solve either:

  • an ODE system of the form $$y' = f(t,y,...)$$ for y = Y, or

  • a DAE system of the form $$F(t,y,y') = 0$$ for y = Y and y' = YPRIME.

The recommended value of maxord is eight, unless it is believed that there are severe stability problems in which case maxord = 4 or 5 should be tried instead.

ODEs are specified in func, DAEs are specified in res.

If a DAE system, Values for Y and YPRIME at the initial time must be given as input. Ideally,these values should be consistent, that is, if T, Y, YPRIME are the given initial values, they should satisfy F(T,Y,YPRIME) = 0.

The form of the Jacobian can be specified by jactype. This is one of:

jactype = "fullint":

a full Jacobian, calculated internally by mebdfi, the default,

jactype = "fullusr":

a full Jacobian, specified by user function jacfunc or jacres,

jactype = "bandusr":

a banded Jacobian, specified by user function jacfunc or jacres; the size of the bands specified by bandup and banddown,

jactype = "bandint":

a banded Jacobian, calculated by mebdfi; the size of the bands specified by bandup and banddown.

If jactype = "fullusr" or "bandusr" then the user must supply a subroutine jacfunc.

If jactype = "fullusr" or "bandusr" then the user must supply a subroutine jacfunc or jacres.

The input parameters rtol, and atol determine the error control performed by the solver. If the request for precision exceeds the capabilities of the machine, mebdfi will return an error code.

res and jacres may be defined in compiled C or Fortran code, as well as in an R-function. See deSolve's vignette "compiledCode" for details. Examples in Fortran are in the dynload subdirectory of the deSolve package directory.

The diagnostics of the integration can be printed to screen by calling diagnostics. If verbose = TRUE, the diagnostics will written to the screen at the end of the integration.

See vignette("deSolve") for an explanation of each element in the vectors containing the diagnostic properties and how to directly access them.

References

J. R. Cash, The integration of stiff initial value problems in O.D.E.S using modified extended backward differentiation formulae, Comp. and Maths. with applics., 9, 645-657, (1983).

J.R. Cash and S. Considine, an MEBDF code for stiff initial value problems, ACM Trans Math Software, 142-158, (1992).

J.R. Cash, Stable recursions with applications to the numerical solution of stiff systems, Academic Press,(1979).

See Also

  • gamd and bimd two other DAE solvers,

  • daspk another DAE solver from package deSolve,

diagnostics to print diagnostic messages.

Examples

Run this code
## =======================================================================
## Coupled chemical reactions including an equilibrium
## modeled as (1) an ODE and (2) as a DAE
##
## The model describes three chemical species A,B,D:
## subjected to equilibrium reaction D <- > A + B
## D is produced at a constant rate, prod
## B is consumed at 1s-t order rate, r
## Chemical problem formulation 1: ODE
## =======================================================================

## Dissociation constant
K <- 1 

## parameters
pars <- c(
        ka   = 1e6,     # forward rate
        r    = 1,
        prod = 0.1)


Fun_ODE <- function (t, y, pars)
{
  with (as.list(c(y, pars)), {
    ra  <- ka*D        # forward rate
    rb  <- ka/K *A*B   # backward rate

    ## rates of changes
    dD  <- -ra + rb + prod
    dA  <-  ra - rb
    dB  <-  ra - rb - r*B
    return(list(dy = c(dA, dB, dD),
                CONC = A+B+D))
  })
}

## =======================================================================
## Chemical problem formulation 2: DAE
## 1. get rid of the fast reactions ra and rb by taking
## linear combinations   : dD+dA = prod (res1) and
##                         dB-dA = -r*B (res2)
## 2. In addition, the equilibrium condition (eq) reads:
## as ra = rb : ka*D = ka/K*A*B = >      K*D = A*B
## =======================================================================

Res_DAE <- function (t, y, yprime, pars)
{
  with (as.list(c(y, yprime, pars)), {
    ## residuals of lumped rates of changes
    res1 <- -dD - dA + prod
    res2 <- -dB + dA - r*B
    
    ## and the equilibrium equation
    eq   <- K*D - A*B

    return(list(c(res1, res2, eq),
                CONC = A+B+D))
  })
}

times <- seq(0, 100, by = 1)

## Initial conc; D is in equilibrium with A,B
y     <- c(A = 2, B = 3, D = 2*3/K)

## ODE model solved with mebdfi
ODE <- as.data.frame(mebdfi(y = y, times = times, func = Fun_ODE,
                     parms = pars, atol = 1e-8, rtol = 1e-8))

## Initial rate of change
dy  <- c(dA = 0, dB = 0, dD = 0) 
## DAE model solved with mebdfi
DAE <- as.data.frame(mebdfi(y = y, dy = dy, times = times,  
         res = Res_DAE, parms = pars, atol = 1e-8, rtol = 1e-8))

         
## =======================================================================
## Chemical problem formulation 3: Mass * Func
## Based on the DAE formulation
## =======================================================================

Mass_FUN <- function (t, y, pars)
{
  with (as.list(c(y, pars)), {

    ## as above, but without the 
    f1 <- prod
    f2 <- - r*B
    
    ## and the equilibrium equation
    f3   <- K*D - A*B

    return(list(c(f1, f2, f3),
                CONC = A+B+D))
  })
}
Mass <- matrix(nr=3, nc=3, byrow = TRUE, 
  data=c(1,  0, 1,         # dA + 0 + dB
        -1,  1, 0,         # -dA + dB +0
         0,  0, 0))        # algebraic
         
times <- seq(0, 100, by = 2)

## Initial conc; D is in equilibrium with A,B
y     <- c(A = 2, B = 3, D = 2*3/K)

## ODE model solved with daspk
ODE <- as.data.frame(daspk(y = y, times = times, func = Fun_ODE,
                     parms = pars, atol = 1e-10, rtol = 1e-10))

## Initial rate of change
dy  <- c(dA = 0, dB = 0, dD = 0) 

## DAE model solved with daspk
DAE <- as.data.frame(daspk(y = y, dy = dy, times = times,
         res = Res_DAE, parms = pars, atol = 1e-10, rtol = 1e-10))

MASS<- mebdfi(y = y, times = times, func = Mass_FUN, 
              parms = pars, mass = Mass)

## ================
## plotting output
## ================
opa <- par(mfrow = c(2, 2))

for (i in 2:5)
{
plot(ODE$time, ODE[, i], xlab = "time",
     ylab = "conc", main = names(ODE)[i], type = "l")
points(DAE$time, DAE[,i], col = "red")
}
legend("bottomright",lty = c(1,NA),pch = c(NA,1),
       col = c("black","red"),legend = c("ODE","DAE"))      
       
# difference between both implementations:
max(abs(ODE-DAE))

par(mfrow = opa)

## =============================================================================
##
## Example 3: higher index DAE
##
## Car axis problem, index 3 DAE, 8 differential, 2 algebraic equations
## from
## F. Mazzia and C. Magherini. Test Set for Initial Value Problem Solvers,
## release 2.4. Department
## of Mathematics, University of Bari and INdAM, Research Unit of Bari,
## February 2008.
## Available at http://www.dm.uniba.it/~testset.
## =============================================================================

# car returns the residuals of the implicit DAE
car <- function(t, y, dy, pars){
  with(as.list(c(pars, y)), {
      f <- rep(0, 10)

      yb  <- r*sin(w*t)
      xb  <- sqrt(L*L - yb*yb)
      Ll  <- sqrt(xl^2 + yl^2)
      Lr  <- sqrt((xr-xb)^2 + (yr-yb)^2)

      f[1:4] <- y[5:8]
      k <- M*eps*eps/2

      f[5]  <- (L0-Ll)*xl/Ll + lam1*xb+2*lam2*(xl-xr)
      f[6]  <- (L0-Ll)*yl/Ll + lam1*yb+2*lam2*(yl-yr)-k*g
      f[7]  <- (L0-Lr)*(xr-xb)/Lr - 2*lam2*(xl-xr)
      f[8]  <- (L0-Lr)*(yr-yb)/Lr - 2*lam2*(yl-yr)-k*g

      f[9]  <- xb*xl+yb*yl
      f[10] <- (xl-xr)^2+(yl-yr)^2-L*L

      delt       <- dy-f
      delt[5:8]  <- k*dy[5:8]-f[5:8]
      delt[9:10] <- -f[9:10]

      list(delt=delt,f=f)
  })
}


# parameters
pars <- c(eps = 1e-2, M = 10, L = 1, L0 = 0.5,
          r   = 0.1,  w = 10, g = 1)

# initial conditions: state variables
yini <-  with (as.list(pars),
   c(xl = 0, yl = L0, xr = L, yr = L0, xla = -L0/L,
     yla = 0, xra = -L0/L, yra = 0, lam1 = 0, lam2 = 0)
              )

# initial conditions: derivates
dyini <- rep(0, 10)
FF    <- car(0, yini, dyini, pars)
dyini[1:4] <- yini[5:8]
dyini[5:8] <- 2/pars["M"]/(pars["eps"])^2*FF$f[5:8]

# check consistency of initial condition: delt should be = 0.
car(0, yini, dyini, pars)

# running the model
times <- seq(0, 3, by = 0.01)
nind  <- c(4, 4, 2)   # index 1, 2 and 3 variables
out   <- mebdfi(y = yini, dy = dyini, times, res = car, parms = pars,
                nind = nind, rtol = 1e-5, atol = 1e-5)

plot(out, which = 1:4, type = "l", lwd=2)

mtext(outer = TRUE, side = 3, line = -0.5, cex = 1.5, "car axis")

Run the code above in your browser using DataLab