pdmpModel-class: Class pdmpModel

pdmpModel-classR Documentation

Class pdmpModel

Description

An S4 class to represent a piecewise deterministic markov process (PDMP).
This class provides a possibility to simulate piecewise deterministic markov processes. These processes are described in [Zei09] and [Ben+15]. There are some restrictions compared to the original definition given in [Dav84]: In our case there are no borders allowed and the number of continous variables should be independent of the state (of the discrete variable) of the process.

Usage

pdmpModel(
  obj = NULL,
  descr = character(0),
  dynfunc = function(t, x, parms) 0 * x,
  jumpfunc = function(t, x, parms, jtype) x,
  ratefunc = function(t, x, parms) c(0),
  times = c(from = 0, to = 10, by = 1),
  init = c(0, 0),
  parms = list(),
  discStates = list(0),
  initfunc = NULL,
  out = NULL,
  solver = "lsodar"
)

Arguments

obj

pdmpModel object that is being built.

descr

a string containing a short description of the model. This parameter is optional and only used in plot methods.

dynfunc

a function(time, variables, parms) that returns a vector with odes for every variable. The order and length of the vector should be the same as in slot "init", discrete variables should have 0 as entry.

jumpfunc

a function(t, x, parms, jtype) that returns the next discrete state the process will jump to. This state depends on parameter jtype. The number of possible jtypes is determined by function ratefunc. The value for jtype will be chosen randomly during simulation, depending ot the rates given in ratefunc.

ratefunc

a function(t, x, parms) that returns a vector with transition rates from the actual state to another state. Only non zero rates are given. The length of the returned vector determines the number of different jumptypes.

times

vector of time steps or vector with three named values "from", "to", "by" specifying the simulation time steps. The from-to-by can be edited with fixParms.

init

initial state of the simulation. This is a named vector giving the names of all variables and their start value.

parms

a vector or list with constant model parameters.

discStates

a list. For every discrete variable, discStates contains a vector with all its possible state values. This entry should have the same name as the discrete variable.

initfunc

this parameter can hold an optional function which has a pdmpModel as only parameter and returnes a (modified) pdmp. This function is called automatically when a new object is created by new or when it is reinitialized by initialize(obj) or before starting a simulation with sim(obj, initialize = TRUE).

out

NULL or an object of class deSolve. If a simulation is done with method sim, the result will be stored in this parameter.

solver

a function or a character string specifying the numerical algorithm used, e.g. "lsodar" from package deSolve. The solver should include root-finding.

Slots

descr

a string containing a short description of the model. This slot is optional and only used in plot methods.

parms

a list with constant model parameters.

times

vector of time steps or vector with three named values "from", "to", "by" specifying the simulation time steps. The from-to-by can be edited with fixParms.

init

initial state of the simulation. This is a named vector giving the names of all variables and their start value.

discStates

a list. For every discrete variable, discStates contains a vector with all its possible state values. This entry should have the same name as the discrete variable.

dynfunc

a function(time, variables, parms) that returns a vector with odes for every variable. The order and length of the vector should be the same as in slot "init", discrete variables should have 0 as entry.

ratefunc

a function(t, x, parms) that returns a vector with transition rates from the actual state to another state. Only non zero rates are given. The length of the returned vector determines the number of different jumptypes.

jumpfunc

a function(t, x, parms, jtype) that returns the next discrete state the process will jump to. This state depends on parameter jtype. The number of possible jtypes is determined by function ratefunc. The value for jtype will be chosen randomly during simulation, depending ot the rates given in ratefunc.

solver

a function or a character string specifying the numerical algorithm used, e.g. "lsodar" from package deSolve. The solver should include root-finding.

initfunc

this slot can hold an optional function which has a pdmpModel as only parameter and returnes an object of class pdmpModel. This function is called automatically when a new object is created by new or when it is reinitialized by initialize(obj) or before starting a simulation with sim(obj, initialize = TRUE).

out

NULL or an object of class deSolve. If a simulation is done with method sim, the result will be stored in this slot.

References

\,[Dav84] Davis, M. H. (1984). Piecewise-deterministic Markov processes: A general class of
non-diffusion stochastic models. Journal of the Royal Statistical Society. Series B
(Methodological), 353-388.
[Zei09] S. Zeiser. Classical and Hybrid Modeling of Gene Regulatory Networks. 2009.
[Ben+15]\,\,\,\, Benaïm, M., Le Borgne, S., Malrieu, F., & Zitt, P. A. (2015). Qualitative properties
of certain piecewise deterministic Markov processes. In Annales de l'Institut Henri
Poincaré, Probabilités et Statistiques (Vol. 51, No. 3, pp. 1040-1075). Institut
Henri Poincaré.

See Also

See simplePdmp and toggleSwitch for two examples that have a detailed documentation explaining every slot. Class pdmpModel provides a method sim for simulation, pdmp-accessors{accessor functions} (with names identical to the slot names) to get or set model parameters, time steps, initial values, the vectorfields, the transition rates and the solver. See multSim and multSimCsv to perform multiple simulations for a pdmpModel.

Examples


# a pdmp with different jumptypes
simplePdmp <- pdmpModel(
  descr = "a PDMP with 2 jumptypes",
  init = c(f = 0, d = 0),
  times = c(from = 0, to = 10, by = 0.01),
  discStates = list(d = -1:1),
  dynfunc = function(t, x, parms) c(x["d"], 0),
  ratefunc = function(t, x, parms) c(1+x["d"], 1-x["d"]),
  jumpfunc = function(t, x, parms, jtype){
    c(0, switch(jtype, x["d"]-1, x["d"]+1))
  }
)
simplePdmp <- sim(simplePdmp, seed = 10) # simulate
plot(simplePdmp)

# the same pdmp defined with two discrete variables
discPdmp <- pdmpModel(
  descr = "a PDMP with 2 discrete variables and 2 jumptypes",
  init = c(f = 0, d1 = 0, d2 = 0),
  discStates = list(d1 = 0:1, d2 = 0:1),
  times = c(from = 0, to = 10, by = 0.01),
  dynfunc = function(t, x, parms) c(x["d1"]-x["d2"], 0, 0),
  ratefunc = function(t, x, parms) c(1,1),
  jumpfunc = function(t, x, parms, jtype){
    c(0, switch(jtype, c((x["d1"]-1)^2, x["d2"]),
                c(x["d1"], (x["d2"]-1)^2)))
  }
)
discPdmp <- sim(discPdmp, seed = 10) # simulate
print(head(out(discPdmp)))

# see ?simplePdmp for an explanation of this example
# and ?toggleSwitch for a more sophisticated example


CharlotteJana/pdmpsim documentation built on Oct. 21, 2024, 4:54 p.m.