nlmixr (version 1.0.0-7)

foceiControl: Control Options for FOCEi

Description

Control Options for FOCEi

Usage

foceiControl(sigdig = 4, epsilon = NULL, maxInnerIterations = 1000,
  maxOuterIterations = 5000, n1qn1nsim = NULL, method = c("liblsoda",
  "lsoda", "dop853"), transitAbs = NULL, atol = NULL, rtol = NULL,
  maxstepsOde = 5000L, hmin = 0L, hmax = NULL, hini = 0,
  maxordn = 12L, maxords = 5L, cores, covsInterpolation = c("locf",
  "linear", "nocb", "midpoint"), printInner = 0L, print = 1L,
  printNcol = floor((getOption("width") - 23)/12), scaleTo = 1,
  scaleObjective = 1, derivEps = rep(20 * sqrt(.Machine$double.eps),
  2), derivMethod = c("switch", "forward", "central"),
  derivSwitchTol = NULL, covDerivMethod = c("central", "forward"),
  covMethod = c("r,s", "r", "s", ""), hessEps = 0.001,
  covDerivEps = c(0, 0.001), lbfgsLmm = 7L, lbfgsPgtol = 0,
  lbfgsFactr = NULL, eigen = TRUE, addPosthoc = TRUE,
  diagXform = c("sqrt", "log", "identity"), sumProd = FALSE,
  optExpression = TRUE, ci = 0.95, useColor = crayon::has_color(),
  boundTol = NULL, calcTables = TRUE, noAbort = TRUE,
  interaction = TRUE, cholSEtol = (.Machine$double.eps)^(1/3),
  cholAccept = 0.001, resetEtaP = 0.1, diagOmegaBoundUpper = 5,
  diagOmegaBoundLower = 100, cholSEOpt = FALSE, cholSECov = FALSE,
  fo = FALSE, covTryHarder = FALSE, outerOpt = c("bobyqa",
  "L-BFGS-B", "lbfgsb3", "nlminb", "mma", "lbfgsbLG", "slsqp"),
  innerOpt = c("n1qn1", "BFGS"), rhobeg = 0.2, rhoend = NULL,
  npt = NULL, rel.tol = NULL, x.tol = NULL, eval.max = 4000,
  iter.max = 2000, abstol = NULL, reltol = NULL,
  resetHessianAndEta = FALSE, stateTrim = Inf, ..., stiff)

Arguments

sigdig

Optimization significant digits. This controls:

  • The tolerance of the inner and outer optimization is 10^-sigdig

  • The tolerance of the ODE solvers is 10^(-sigdig-1)

  • The tolerance of the boundary check is 5 * 10 ^ (-sigdig + 1)

  • The significant figures that some tables are rounded to.

epsilon

Precision of estimate for n1qn1 optimization.

maxInnerIterations

Number of iterations for n1qn1 optimization.

maxOuterIterations

Maximum number of L-BFGS-B optimization for outer problem.

n1qn1nsim

Number of function evaluations for n1qn1 optimization.

method

The method for solving ODEs. Currently this supports:

  • "liblsoda" thread safe lsoda. This supports parallel thread-based solving, and ignores user Jacobian specification.

  • "lsoda" -- LSODA solver. Does not support parallel thread-based solving, but allows user Jacobian specification.

  • "dop853" -- DOP853 solver. Does not support parallel thread-based solving nor user Jacobain specification

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 (1e-6 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.

maxstepsOde

Maximum number of steps for ODE solver.

hmin

The minimum absolute step size allowed. The default value is 0.

hmax

The maximum absolute step size allowed. The default checks for the maximum difference in times in your sampling and events, and uses this value. The value 0 is equivalent to infinite maximum absolute step size.

hini

The step size to be attempted on the first step. The default value is determined by the solver (when hini = 0)

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.

cores

Number of cores used in parallel ODE solving. This defaults to the number or system cores determined by rxCores for methods that support parallel solving (ie thread-safe methods like "liblsoda").

covsInterpolation

specifies the interpolation method for time-varying covariates. When solving ODEs it often samples times outside the sampling time specified in events. When this happens, the time varying covariates are interpolated. Currently this can be:

  • "linear" interpolation (the default), which interpolates the covariate by solving the line between the observed covariates and extrapolating the new covariate value.

  • "constant" -- Last observation carried forward.

  • "NOCB" -- Next Observation Carried Backward. This is the same method that NONMEM uses.

  • "midpoint" Last observation carried forward to midpoint; Next observation carried backward to midpoint.

printInner

Integer representing when the inner step is printed. By default this is 0 or do not print. 1 is print every function evaluation, 5 is print every 5 evaluations.

print

Integer representing when the outer step is printed. When this is 0 or do not print the iterations. 1 is print every function evaluation (default), 5 is print every 5 evaluations.

printNcol

Number of columns to printout before wrapping parameter estimates/gradient

scaleTo

Scale the initial parameter estimate to this value. By default this is 1. When zero or below, no scaling is performed.

scaleObjective

Scale the initial objective function to this value. By default this is 1. When scaleObjective is greater than zero, this scaling is performed by:

scaledObj = currentObj / \|initialObj\| * scaleObjective

Therefore, if the initial objective function is negative, the initial scaled objective function would be negative as well. When scaleObjective is less than zero, no scaling is performed.

derivEps

Central/Forward difference tolerances, which is a vector of relative difference and absolute difference. The central/forward difference step size h is calculated as:

h = abs(x)*derivEps[1] + derivEps[2]

derivMethod

indicates the method for calculating derivatives of the outer problem. Currently supports "switch", "central" and "forward" difference methods. Switch starts with forward differences. This will switch to central differences when abs(delta(OFV)) <= derivSwitchTol and switch back to forward differences when abs(delta(OFV)) > derivSwitchTol.

derivSwitchTol

The tolerance to switch forward to central differences.

covDerivMethod

indicates the method for calculating the derivatives while calculating the covariance components (Hessian and S).

covMethod

Method for calculating covariance. In this discussion, R is the Hessian matrix of the objective function. The S matrix is the sum of individual gradient cross-product (evaluated at the individual empirical Bayes estimates).

  • "r,s" Uses the sandwich matrix to calculate the covariance, that is: solve(R) %*% S %*% solve(R)

  • "r" Uses the Hessian matrix to calculate the covariance as 2 %*% solve(R)

  • "s" Uses the crossproduct matrix to calculate the covariance as 4 %*% solve(S)

  • "" Does not calculate the covariance step.

hessEps

is a double value representing the epsilon for the Hessian calculation.

covDerivEps

Central/Forward difference tolerances while calculating the covariance matrices. This is a numeric vector of relative difference and absolute difference. The central/forward difference step size h is calculated as:

h = abs(x)*derivEps[1] + derivEps[2]

lbfgsLmm

An integer giving the number of BFGS updates retained in the "L-BFGS-B" method, It defaults to 7.

lbfgsPgtol

is a double precision variable.

On entry pgtol >= 0 is specified by the user. The iteration will stop when:

max(\| proj g_i \| i = 1, ..., n) <= lbfgsPgtol

where pg_i is the ith component of the projected gradient.

On exit pgtol is unchanged. This defaults to zero, when the check is suppressed.

lbfgsFactr

Controls the convergence of the "L-BFGS-B" method. Convergence occurs when the reduction in the objective is within this factor of the machine tolerance. Default is 1e10, which gives a tolerance of about 2e-6, approximately 4 sigdigs. You can check your exact tolerance by multiplying this value by .Machine$double.eps

eigen

A boolean indicating if eigenvectors are calculated to include a condition number calculation.

addPosthoc

Boolean indicating if posthoc parameters are added to the table output.

diagXform

This is the transformation used on the diagonal of the chol(solve(omega)). This matrix and values are the parameters estimated in FOCEi. The possibilities are:

  • sqrt Estimates the sqrt of the diagonal elements of chol(solve(omega)). This is the default method.

  • log Estimates the log of the diagonal elements of chol(solve(omega))

  • identity Estimates the diagonal elements without any transformations

sumProd

Is a boolean indicating if the model should change multiplication to high precision multiplication and sums to high precision sums using the PreciseSums package. By default this is FALSE.

optExpression

Optimize the RxODE expression to speed up calculation. By default this is turned on.

ci

Confidence level for some tables. By default this is 0.95 or 95% confidence.

useColor

Boolean indicating if focei can use ASCII color codes

boundTol

Tolerance for boundary issues.

calcTables

This boolean is to determine if the foceiFit will calculate tables. By default this is TRUE

noAbort

Boolean to indicate if you should abort the FOCEi evaluation if it runs into troubles. (default TRUE)

interaction

Boolean indicate FOCEi should be used (TRUE) instead of FOCE (FALSE)

cholSEtol

tolerance for Generalized Cholesky Decomposition. Defaults to suggested (.Machine$double.eps)^(1/3)

cholAccept

Tolerance to accept a Generalized Cholesky Decomposition for a R or S matrix.

resetEtaP

represents the p-value for reseting the individual ETA to 0 during optimization (instead of the saved value). The two test statistics used in the z-test are either chol(omega^-1) indicates the ETAs never reset. A p-value of 1 indicates the ETAs always reset.

diagOmegaBoundUpper

This represents the upper bound of the diagonal omega matrix. The upper bound is given by diag(omega)*diagOmegaBoundUpper. If diagOmegaBoundUpper is 1, there is no upper bound on Omega.

diagOmegaBoundLower

This represents the lower bound of the diagonal omega matrix. The lower bound is given by diag(omega)/diagOmegaBoundUpper. If diagOmegaBoundLower is 1, there is no lower bound on Omega.

cholSEOpt

Boolean indicating if the generalized Cholesky should be used while optimizing.

cholSECov

Boolean indicating if the generalized Cholesky should be used while calculating the Covariance Matrix.

fo

is a boolean indicating if this is a fo approximation routine.

covTryHarder

If the R matrix is non-positive definite and cannot be corrected to be non-positive definite try estimating the Hessian on the unscaled parameter space.

outerOpt

optimization method for the outer problem

innerOpt

optimization method for the inner problem (not implemented yet.)

rhobeg

Beginning change in parameters for bobyqa algorithm (trust region). By default this is 0.2 or 20 parameters when the parameters are scaled to 1. rhobeg and rhoend must be set to the initial and final values of a trust region radius, so both must be positive with 0 < rhoend < rhobeg. Typically rhobeg should be about one tenth of the greatest expected change to a variable. Note also that smallest difference abs(upper-lower) should be greater than or equal to rhobeg*2. If this is not the case then rhobeg will be adjusted.

rhoend

The smallest value of the trust region radius that is allowed. If not defined, then 10^(-sigdig-1) will be used.

npt

The number of points used to approximate the objective function via a quadratic approximation for bobyqa. The value of npt must be in the interval [n+2,(n+1)(n+2)/2] where n is the number of parameters in par. Choices that exceed 2*n+1 are not recommended. If not defined, it will be set to 2*n + 1

rel.tol

Relative tolerance before nlminb stops.

x.tol

X tolerance

eval.max

Number of maximum evaluations of the objective function

iter.max

Maximum number of iterations allowed.

abstol

Absolute tolerance for nlmixr

reltol

tolerance for nlmixr

resetHessianAndEta

is a boolean representing if the individual Hessian is reset when ETAs are reset using the option resetEtaP.

stateTrim

Trim state amounts/concentrations to this value.

...

Ignored parameters

stiff

a logical (TRUE by default) indicating whether the ODE system is stiff or not.

For stiff ODE sytems (stiff = TRUE), RxODE uses the LSODA (Livermore Solver for Ordinary Differential Equations) Fortran package, which implements an automatic method switching for stiff and non-stiff problems along the integration interval, authored by Hindmarsh and Petzold (2003).

For non-stiff systems (stiff = FALSE), RxODE uses DOP853, an explicit Runge-Kutta method of order 8(5, 3) of Dormand and Prince as implemented in C by Hairer and Wanner (1993).

Details

Note this uses the R's L-BFGS-B in optim for the outer problem and the BFGS n1qn1 with that allows restoring the prior individual Hessian (for faster optimization speed).

By default FOCEi scales the outer problem parameters to 1.0 for the initial parameter estimates and scales the objective function to 1.0, as suggested by the NAG library and https://www.scipy-lectures.org/advanced/mathematical_optimization/.

However the inner problem is not scaled. Since most eta estimates start near zero, scaling for these parameters do not make sense.

This process of scaling can fix some ill conditioning for the unscaled problem. The covariance step is performed on the unscaled problem, so the condition number of that matrix may not be reflective of the scaled problem's condition-number.

See Also

optim

n1qn1

rxSolve