| RxODE {RxODE} | R Documentation |
Create a dynamic ODE-based model object suitably for translation into fast C code
RxODE( model, modName = basename(wd), wd = getwd(), filename = NULL, extraC = NULL, debug = FALSE, calcJac = NULL, calcSens = NULL, collapseModel = FALSE, package = NULL, ... )
model |
This is the ODE model specification. It can be:
(see also the |
modName |
a string to be used as the model name. This string
is used for naming various aspects of the computations,
including generating C symbol names, dynamic libraries,
etc. Therefore, it is necessary that |
wd |
character string with a working directory where to
create a subdirectory according to |
filename |
A file name or connection object where the
ODE-based model specification resides. Only one of |
extraC |
Extra c code to include in the model. This can be
useful to specify functions in the model. These C functions
should usually take |
debug |
is a boolean indicating if the executable should be compiled with verbose debugging information turned on. |
calcJac |
boolean indicating if RxODE will calculate the Jacobain according to the specified ODEs. |
calcSens |
boolean indicating if RxODE will calculate the sensitivities according to the specified ODEs. |
collapseModel |
boolean indicating if RxODE will remove all LHS variables when calculating sensitivities. |
package |
Package name for pre-compiled binaries. |
... |
ignored arguments. The “Rx” in the name The ODE-based model specification may be coded inside a character
string or in a text file, see Section RxODE Syntax below for
coding details. An internal For evaluating |
An object (closure) of class “RxODE” (see Chambers and Temple Lang (2001))
consisting of the following list of strings and functions:
modName |
the name of the model (a copy of the input argument). |
model |
a character string holding the source model specification. |
get.modelVars |
a function that returns a list with 3 character
vectors, |
solve |
this function solves (integrates) the ODE. This
is done by passing the code to
For stiff ODE systems ( For non-stiff systems (
The output of “solve” is a matrix with as many rows as there are sampled time points and as many columns as system variables (as defined by the ODEs and additional assignments in the RxODE model code). |
isValid |
a function that (naively) checks for model validity, namely that the C object code reflects the latest model specification. |
version |
a string with the version of the |
dynLoad |
a function with one |
dynUnload |
a function with no argument that unloads the model object code. |
delete |
removes all created model files, including C and DLL files.
The model object is no longer valid and should be removed, e.g.,
|
run |
deprecated, use |
parse |
deprecated. |
compile |
deprecated. |
get.index |
deprecated. |
getObj |
internal (not user callable) function. |
An RxODE model specification consists of one or more
statements terminated by semi-colons, ‘;’, and
optional comments (comments are delimited by # and an
end-of-line marker). NB: Comments are not allowed inside
statements.
A block of statements is a set of statements delimited by curly
braces, ‘{ ... }’. Statements can be either
assignments or conditional if statements. Assignment
statements can be: (1) “simple” assignments, where the left
hand is an identifier (i.e., variable), (2) special
“time-derivative” assignments, where the left hand specifies
the change of that variable with respect to time e.g.,
d/dt(depot), or (3) special “jacobian” assignments,
where the left hand specifies the change of of the ODE with respect
to one of the parameters, e.g. df(depot)/dy(kel). The
“jacobian” assignments are not required, and are only useful
for very stiff differential systems.
Expressions in assignment and ‘if’ statements can be
numeric or logical (no character expressions are currently
supported). Numeric expressions can include the following numeric
operators (‘+’, ‘-’, ‘*’,
‘/’, ‘^’), and those mathematical
functions defined in the C or the R math libraries (e.g.,
fabs, exp, log, sin). (Notice that the
modulo operator ‘%’ is currently not supported.)
Identifiers in an RxODE model specification can refer to:
state variables in the dynamic system (e.g., compartments in a pharmacokinetics/pharamcodynamics model);
implied input variable, t (time),
podo (oral dose, for absorption models), and
tlast (last time point);
model parameters, (ka rate of absorption, CL
clearance, etc.);
pi, for the constant pi.
others, as created by assignments as part of the model specification.
Identifiers consists of case-sensitive alphanumeric characters, plus the underscore ‘_’ character. NB: the dot ‘.’ character is not a valid character identifier.
The values of these variables at pre-specified time points are
saved as part of the fitted/integrated/solved model (see
eventTable, in particular its member function
add.sampling that defines a set of time points at which to
capture a snapshot of the system via the values of these variables).
The ODE specification mini-language is parsed with the help of the
open source tool dparser, Plevyak (2015).
Melissa Hallow, Wenping Wang and Matthew Fidler
Chamber, J. M. and Temple Lang, D. (2001) Object Oriented Programming in R. R News, Vol. 1, No. 3, September 2001. https://cran.r-project.org/doc/Rnews/Rnews_2001-3.pdf.
Hindmarsh, A. C. ODEPACK, A Systematized Collection of ODE Solvers. Scientific Computing, R. S. Stepleman et al. (Eds.), North-Holland, Amsterdam, 1983, pp. 55-64.
Petzold, L. R. Automatic Selection of Methods for Solving Stiff and Nonstiff Systems of Ordinary Differential Equations. Siam J. Sci. Stat. Comput. 4 (1983), pp. 136-148.
Hairer, E., Norsett, S. P., and Wanner, G. Solving ordinary differential equations I, nonstiff problems. 2nd edition, Springer Series in Computational Mathematics, Springer-Verlag (1993).
Plevyak, J.
dparser, http://dparser.sourceforge.net. Web. 12 Oct. 2015.
eventTable, et, add.sampling, add.dosing
# Step 1 - Create a model specification
ode <- "
# A 4-compartment model, 3 PK and a PD (effect) compartment
# (notice state variable names 'depot', 'centr', 'peri', 'eff')
C2 = centr/V2;
C3 = peri/V3;
d/dt(depot) =-KA*depot;
d/dt(centr) = KA*depot - CL*C2 - Q*C2 + Q*C3;
d/dt(peri) = Q*C2 - Q*C3;
d/dt(eff) = Kin - Kout*(1-C2/(EC50+C2))*eff;
"
m1 <- RxODE(model = ode)
print(m1)
# Step 2 - Create the model input as an EventTable,
# including dosing and observation (sampling) events
# QD (once daily) dosing for 5 days.
qd <- eventTable(amount.units = "ug", time.units = "hours")
qd$add.dosing(dose = 10000, nbr.doses = 5, dosing.interval = 24)
# Sample the system hourly during the first day, every 8 hours
# then after
qd$add.sampling(0:24)
qd$add.sampling(seq(from = 24+8, to = 5*24, by = 8))
# Step 3 - set starting parameter estimates and initial
# values of the state
theta <-
c(KA = .291, CL = 18.6,
V2 = 40.2, Q = 10.5, V3 = 297.0,
Kin = 1.0, Kout = 1.0, EC50 = 200.0)
# init state variable
inits <- c(0, 0, 0, 1);
# Step 4 - Fit the model to the data
qd.cp <- m1$solve(theta, events = qd, inits)
head(qd.cp)
# This returns a matrix. Note that you can also
# solve using name initial values. For example:
inits <- c(eff = 1);
qd.cp <- solve(m1, theta, events = qd, inits);
print(qd.cp)
plot(qd.cp)