R/rxsolve.R
rxSolve.Rd
This uses RxODE family of objects, file, or model specification to
solve a ODE system. There are many options for a solved RxODE
model, the first are the required object
, and events
with the
some-times optional params
and inits
.
rxSolve(
object,
params = NULL,
events = NULL,
inits = NULL,
scale = NULL,
method = c("liblsoda", "lsoda", "dop853", "indLin"),
transitAbs = NULL,
atol = 1e-08,
rtol = 1e-06,
maxsteps = 70000L,
hmin = 0,
hmax = NA_real_,
hmaxSd = 0,
hini = 0,
maxordn = 12L,
maxords = 5L,
...,
cores,
covsInterpolation = c("locf", "linear", "nocb", "midpoint"),
addCov = FALSE,
matrix = FALSE,
sigma = NULL,
sigmaDf = NULL,
sigmaLower = -Inf,
sigmaUpper = Inf,
nCoresRV = 1L,
sigmaIsChol = FALSE,
sigmaSeparation = c("auto", "lkj", "separation"),
sigmaXform = c("identity", "variance", "log", "nlmixrSqrt", "nlmixrLog",
"nlmixrIdentity"),
nDisplayProgress = 10000L,
amountUnits = NA_character_,
timeUnits = "hours",
stiff,
theta = NULL,
thetaLower = -Inf,
thetaUpper = Inf,
eta = NULL,
addDosing = FALSE,
stateTrim = Inf,
updateObject = FALSE,
omega = NULL,
omegaDf = NULL,
omegaIsChol = FALSE,
omegaSeparation = c("auto", "lkj", "separation"),
omegaXform = c("variance", "identity", "log", "nlmixrSqrt", "nlmixrLog",
"nlmixrIdentity"),
omegaLower = -Inf,
omegaUpper = Inf,
nSub = 1L,
thetaMat = NULL,
thetaDf = NULL,
thetaIsChol = FALSE,
nStud = 1L,
dfSub = 0,
dfObs = 0,
returnType = c("rxSolve", "matrix", "data.frame", "data.frame.TBS", "data.table",
"tbl", "tibble"),
seed = NULL,
nsim = NULL,
minSS = 10L,
maxSS = 1000L,
infSSstep = 12,
strictSS = TRUE,
istateReset = TRUE,
subsetNonmem = TRUE,
maxAtolRtolFactor = 0.1,
from = NULL,
to = NULL,
by = NULL,
length.out = NULL,
iCov = NULL,
keep = NULL,
indLinPhiTol = 1e-07,
indLinPhiM = 0L,
indLinMatExpType = c("expokit", "Al-Mohy", "arma"),
indLinMatExpOrder = 6L,
drop = NULL,
idFactor = TRUE,
mxhnil = 0,
hmxi = 0,
warnIdSort = TRUE,
warnDrop = TRUE,
ssAtol = 1e-08,
ssRtol = 1e-06,
safeZero = TRUE,
sumType = c("pairwise", "fsum", "kahan", "neumaier", "c"),
prodType = c("long double", "double", "logify"),
sensType = c("advan", "autodiff", "forward", "central"),
linDiff = c(tlag = 1.5e-05, f = 1.5e-05, rate = 1.5e-05, dur = 1.5e-05, tlag2 =
1.5e-05, f2 = 1.5e-05, rate2 = 1.5e-05, dur2 = 1.5e-05),
linDiffCentral = c(tlag = TRUE, f = TRUE, rate = TRUE, dur = TRUE, tlag2 = TRUE, f2 =
TRUE, rate2 = TRUE, dur2 = TRUE),
resample = NULL,
resampleID = TRUE,
maxwhile = 1e+05
)
# S3 method for default
rxSolve(
object,
params = NULL,
events = NULL,
inits = NULL,
...,
theta = NULL,
eta = NULL
)
# S3 method for rxSolve
update(object, ...)
# S3 method for RxODE
predict(object, ...)
# S3 method for rxSolve
predict(object, ...)
# S3 method for rxEt
predict(object, ...)
# S3 method for rxParams
predict(object, ...)
# S3 method for RxODE
simulate(object, nsim = 1L, seed = NULL, ...)
# S3 method for rxSolve
simulate(object, nsim = 1L, seed = NULL, ...)
# S3 method for rxParams
simulate(object, nsim = 1L, seed = NULL, ...)
# S3 method for rxSolve
solve(a, b, ...)
# S3 method for RxODE
solve(a, b, ...)
# S3 method for rxParams
solve(a, b, ...)
# S3 method for rxEt
solve(a, b, ...)
rxControl(..., params = NULL, events = NULL, inits = NULL)
object | is a either a RxODE family of objects, or a file-name with a RxODE model specification, or a string with a RxODE model specification. |
---|---|
params | a numeric named vector with values for every parameter in the ODE system; the names must correspond to the parameter identifiers used in the ODE specification; |
events | an |
inits | a vector of initial values of the state variables (e.g., amounts in each compartment), and the order in this vector must be the same as the state variables (e.g., PK/PD compartments); |
scale | a numeric named vector with scaling for ode
parameters of the system. The names must correspond to the
parameter identifiers in the ODE specification. Each of the
ODE variables will be divided by the scaling factor. For
example |
method | The method for solving ODEs. Currently this supports:
|
transitAbs | boolean indicating if this is a transit compartment absorption |
atol | a numeric absolute tolerance (1e-8 by default) used by the ODE solver to determine if a good solution has been achieved; This is also used in the solved linear model to check if prior doses do not add anything to the solution. |
rtol | a numeric relative tolerance ( |
maxsteps | maximum number of (internally defined) steps allowed during one call to the solver. (5000 by default) |
hmin | The minimum absolute step size allowed. The default value is 0. |
hmax | The maximum absolute step size allowed. When
|
hmaxSd | The number of standard deviations of the time difference to add to hmax. The default is 0 |
hini | The step size to be attempted on the first step. The
default value is determined by the solver (when |
maxordn | The maximum order to be allowed for the nonstiff (Adams) method. The default is 12. It can be between 1 and 12. |
maxords | The maximum order to be allowed for the stiff (BDF) method. The default value is 5. This can be between 1 and 5. |
... | Other arguments including scaling factors for each compartment. This includes S# = numeric will scale a compartment # by a dividing the compartment amount by the scale factor, like NONMEM. |
cores | Number of cores used in parallel ODE solving. This
is equivalent to calling |
covsInterpolation | specifies the interpolation method for
time-varying covariates. When solving ODEs it often samples
times outside the sampling time specified in
|
addCov | A boolean indicating if covariates should be added to the output matrix or data frame. By default this is disabled. |
matrix | A boolean indicating if a matrix should be returned instead of the RxODE's solved object. |
sigma | Named sigma covariance or Cholesky decomposition of a covariance matrix. The names of the columns indicate parameters that are simulated. These are simulated for every observation in the solved system. |
sigmaDf | Degrees of freedom of the sigma t-distribution. By
default it is equivalent to |
sigmaLower | Lower bounds for simulated unexplained variability (by default -Inf) |
sigmaUpper | Upper bounds for simulated unexplained variability (by default Inf) |
nCoresRV | Number of cores used for the simulation of the sigma variables. By default this is 1. To reproduce the results you need to run on the same platform with the same number of cores. This is the reason this is set to be one, regardless of what the number of cores are used in threaded ODE solving. |
sigmaIsChol | Boolean indicating if the sigma is in the Cholesky decomposition instead of a symmetric covariance |
sigmaSeparation | separation strategy for sigma; Tells the type of separation strategy when
simulating covariance with parameter uncertainty with standard
deviations modeled in the
|
sigmaXform | When taking
|
nDisplayProgress | An integer indicating the minimum number of c-based solves before a progress bar is shown. By default this is 10,000. |
amountUnits | This supplies the dose units of a data frame supplied instead of an event table. This is for importing the data as an RxODE event table. |
timeUnits | This supplies the time units of a data frame supplied instead of an event table. This is for importing the data as an RxODE event table. |
stiff | a logical (
|
theta | A vector of parameters that will be named |
thetaLower | Lower bounds for simulated population parameter
variability (by default |
thetaUpper | Upper bounds for simulated population unexplained
variability (by default |
eta | A vector of parameters that will be named |
addDosing | Boolean indicating if the solve should add RxODE
EVID and related columns. This will also include dosing
information and estimates at the doses. Be default, RxODE
only includes estimates at the observations. (default
|
stateTrim | When amounts/concentrations in one of the states
are above this value, trim them to be this value. By default
Inf. Also trims to -stateTrim for large negative
amounts/concentrations. If you want to trim between a range
say |
updateObject | This is an internally used flag to update the
RxODE solved object (when supplying an RxODE solved object) as
well as returning a new object. You probably should not
modify it's |
omega | Estimate of Covariance matrix. When omega is a list, assume it is a block matrix and convert it to a full matrix for simulations. |
omegaDf | The degrees of freedom of a t-distribution for
simulation. By default this is |
omegaIsChol | Indicates if the |
omegaSeparation | Omega separation strategy Tells the type of separation strategy when
simulating covariance with parameter uncertainty with standard
deviations modeled in the
|
omegaXform | When taking
|
omegaLower | Lower bounds for simulated ETAs (by default -Inf) |
omegaUpper | Upper bounds for simulated ETAs (by default Inf) |
nSub | Number between subject variabilities ( |
thetaMat | Named theta matrix. |
thetaDf | The degrees of freedom of a t-distribution for
simulation. By default this is |
thetaIsChol | Indicates if the |
nStud | Number virtual studies to characterize uncertainty in estimated parameters. |
dfSub | Degrees of freedom to sample the between subject variability matrix from the inverse Wishart distribution (scaled) or scaled inverse chi squared distribution. |
dfObs | Degrees of freedom to sample the unexplained variability matrix from the inverse Wishart distribution (scaled) or scaled inverse chi squared distribution. |
returnType | This tells what type of object is returned. The currently supported types are:
|
seed | an object specifying if and how the random number generator should be initialized |
nsim | represents the number of simulations. For RxODE, if
you supply single subject event tables (created with
|
minSS | Minimum number of iterations for a steady-state dose |
maxSS | Maximum number of iterations for a steady-state dose |
infSSstep | Step size for determining if a constant infusion has reached steady state. By default this is large value, 420. |
strictSS | Boolean indicating if a strict steady-state is
required. If a strict steady-state is ( |
istateReset | When |
subsetNonmem | subset to NONMEM compatible EVIDs only. By
default |
maxAtolRtolFactor | The maximum |
from | When there is no observations in the event table, start observations at this value. By default this is zero. |
to | When there is no observations in the event table, end observations at this value. By default this is 24 + maximum dose time. |
by | When there are no observations in the event table, this
is the amount to increment for the observations between |
length.out | The number of observations to create if there isn't any observations in the event table. By default this is 200. |
iCov | A data frame of individual non-time varying covariates
to combine with the |
keep | Columns to keep from either the input dataset or the
|
indLinPhiTol | the requested accuracy tolerance on exponential matrix. |
indLinPhiM | the maximum size for the Krylov basis |
indLinMatExpType | This is them matrix exponential type that is use for RxODE. Currently the following are supported:
|
indLinMatExpOrder | an integer, the order of approximation to
be used, for the |
drop | Columns to drop from the output |
idFactor | This boolean indicates if original ID values should be maintained. This changes the default sequentially ordered ID to a factor with the original ID values in the original dataset. By default this is enabled. |
mxhnil | maximum number of messages printed (per problem)
warning that |
hmxi | inverse of the maximum absolute value of |
warnIdSort | Warn if the ID is not present and RxODE assumes the order of the parameters/iCov are the same as the order of the parameters in the input dataset. |
warnDrop | Warn if column(s) were supposed to be dropped, but were not present. |
ssAtol | Steady state atol convergence factor. Can be a vector based on each state. |
ssRtol | Steady state rtol convergence factor. Can be a vector based on each state. |
safeZero | Use safe zero divide and log routines. By default this is turned on but you may turn it off if you wish. |
sumType | Sum type to use for
|
prodType | Product to use for
|
sensType | Sensitivity type for
|
linDiff | This gives the linear difference amount for all the types of linear compartment model parameters where sensitivities are not calculated. The named components of this numeric vector are:
|
linDiffCentral | This gives the which parameters use central
differences for the linear compartment model parameters. The
are the same components as |
resample | A character vector of model variables to resample
from the input dataset; This sampling is done with replacement.
When |
resampleID | boolean representing if the resampling should be
done on an individual basis |
maxwhile | represents the maximum times a while loop is evaluated before exiting. By default this is 100000 |
a | when using |
b | when using |
An “rxSolve” solve object that stores the solved
value in a special data.frame or other type as determined by
returnType
. By default this has 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). It also stores information about the call to allow
dynamic updating of the solved object.
The operations for the object are similar to a data-frame, but
expand the $
and [[""]]
access operators and assignment
operators to resolve based on different parameter values, initial
conditions, solver parameters, or events (by updating the time
variable).
You can call the eventTable()
methods on the solved object to
update the event table and resolve the system of equations.
The rest of the document focus on the different ODE solving methods, followed by the core solving method's options, RxODE event handling options, RxODE's numerical stability options, RxODE's output options, and finally internal RxODE options or compatibility options.
"New Scaling and Squaring Algorithm for the Matrix Exponential", by Awad H. Al-Mohy and Nicholas J. Higham, August 2009
Roger B. Sidje (1998). EXPOKIT: Software package for computing matrix exponentials. ACM - Transactions on Mathematical Software 24(1), 130-156.
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).
Matthew Fidler, Melissa Hallow and Wenping Wang