aboutsummaryrefslogtreecommitdiff
path: root/R
diff options
context:
space:
mode:
authorJohannes Ranke <jranke@uni-bremen.de>2020-11-06 00:03:29 +0100
committerJohannes Ranke <jranke@uni-bremen.de>2020-11-06 00:03:29 +0100
commitb5b446b718b15ccaae5b197e147fc1358f0f564e (patch)
treea36f32ee664c6925b5afdb812daca41075968152 /R
parent2f24fe0ce70d040e491619d7f2413fc902e433f1 (diff)
Fast analytical solutions for saemix, update.mmkin
Also, use logit transformation for g and for solitary formation fractions, addressing #10.
Diffstat (limited to 'R')
-rw-r--r--R/mkinerrmin.R2
-rw-r--r--R/mkinfit.R18
-rw-r--r--R/mmkin.R32
-rw-r--r--R/saemix.R213
-rw-r--r--R/transform_odeparms.R61
5 files changed, 215 insertions, 111 deletions
diff --git a/R/mkinerrmin.R b/R/mkinerrmin.R
index 1994aed0..f52692ba 100644
--- a/R/mkinerrmin.R
+++ b/R/mkinerrmin.R
@@ -107,7 +107,7 @@ mkinerrmin <- function(fit, alpha = 0.05)
if (obs_var == fit$obs_vars[[1]]) {
special_parms = c("alpha", "log_alpha", "beta", "log_beta",
"k1", "log_k1", "k2", "log_k2",
- "g", "g_ilr", "tb", "log_tb")
+ "g", "g_ilr", "g_qlogis", "tb", "log_tb")
n.optim <- n.optim + length(intersect(special_parms, names(parms.optim)))
}
diff --git a/R/mkinfit.R b/R/mkinfit.R
index 1b1bb73d..7fa1c56e 100644
--- a/R/mkinfit.R
+++ b/R/mkinfit.R
@@ -89,12 +89,11 @@ if(getRversion() >= '2.15.1') utils::globalVariables(c("name", "time", "value"))
#' models and the break point tb of the HS model. If FALSE, zero is used as
#' a lower bound for the rates in the optimisation.
#' @param transform_fractions Boolean specifying if formation fractions
-#' constants should be transformed in the model specification used in the
-#' fitting for better compliance with the assumption of normal distribution
-#' of the estimator. The default (TRUE) is to do transformations. If TRUE,
-#' the g parameter of the DFOP and HS models are also transformed, as they
-#' can also be seen as compositional data. The transformation used for these
-#' transformations is the [ilr()] transformation.
+#' should be transformed in the model specification used in the fitting for
+#' better compliance with the assumption of normal distribution of the
+#' estimator. The default (TRUE) is to do transformations. If TRUE,
+#' the g parameter of the DFOP model is also transformed. Transformations
+#' are described in [transform_odeparms].
#' @param quiet Suppress printing out the current value of the negative
#' log-likelihood after each improvement?
#' @param atol Absolute error tolerance, passed to [deSolve::ode()]. Default
@@ -187,15 +186,14 @@ if(getRversion() >= '2.15.1') utils::globalVariables(c("name", "time", "value"))
#'
#' # Fit the model quietly to the FOCUS example dataset D using defaults
#' fit <- mkinfit(SFO_SFO, FOCUS_D, quiet = TRUE)
-#' # Since mkin 0.9.50.3, we get a warning about non-normality of residuals,
-#' # so we try an alternative error model
+#' plot_sep(fit)
+#' # As lower parent values appear to have lower variance, we try an alternative error model
#' fit.tc <- mkinfit(SFO_SFO, FOCUS_D, quiet = TRUE, error_model = "tc")
#' # This avoids the warning, and the likelihood ratio test confirms it is preferable
#' lrtest(fit.tc, fit)
#' # We can also allow for different variances of parent and metabolite as error model
#' fit.obs <- mkinfit(SFO_SFO, FOCUS_D, quiet = TRUE, error_model = "obs")
-#' # This also avoids the warning about non-normality, but the two-component error model
-#' # has significantly higher likelihood
+#' # The two-component error model has significantly higher likelihood
#' lrtest(fit.obs, fit.tc)
#' parms(fit.tc)
#' endpoints(fit.tc)
diff --git a/R/mmkin.R b/R/mmkin.R
index 6f088de0..f3c07e42 100644
--- a/R/mmkin.R
+++ b/R/mmkin.R
@@ -64,8 +64,9 @@
#'
#' @export mmkin
mmkin <- function(models = c("SFO", "FOMC", "DFOP"), datasets,
- cores = detectCores(), cluster = NULL, ...)
+ cores = parallel::detectCores(), cluster = NULL, ...)
{
+ call <- match.call()
parent_models_available = c("SFO", "FOMC", "DFOP", "HS", "SFORB", "IORE", "logistic")
n.m <- length(models)
n.d <- length(datasets)
@@ -100,12 +101,14 @@ mmkin <- function(models = c("SFO", "FOMC", "DFOP"), datasets,
}
if (is.null(cluster)) {
- results <- parallel::mclapply(as.list(1:n.fits), fit_function, mc.cores = cores)
+ results <- parallel::mclapply(as.list(1:n.fits), fit_function,
+ mc.cores = cores, mc.preschedule = FALSE)
} else {
results <- parallel::parLapply(cluster, as.list(1:n.fits), fit_function)
}
attributes(results) <- attributes(fit_indices)
+ attr(results, "call") <- call
class(results) <- "mmkin"
return(results)
}
@@ -187,3 +190,28 @@ print.mmkin <- function(x, ...) {
}
}
+
+#' @export
+update.mmkin <- function(object, ..., evaluate = TRUE)
+{
+ call <- attr(object, "call")
+
+ update_arguments <- match.call(expand.dots = FALSE)$...
+
+ if (length(update_arguments) > 0) {
+ update_arguments_in_call <- !is.na(match(names(update_arguments), names(call)))
+ }
+
+ for (a in names(update_arguments)[update_arguments_in_call]) {
+ call[[a]] <- update_arguments[[a]]
+ }
+
+ update_arguments_not_in_call <- !update_arguments_in_call
+ if(any(update_arguments_not_in_call)) {
+ call <- c(as.list(call), update_arguments[update_arguments_not_in_call])
+ call <- as.call(call)
+ }
+
+ if(evaluate) eval(call, parent.frame())
+ else call
+}
diff --git a/R/saemix.R b/R/saemix.R
index f8714cf2..8632c1a4 100644
--- a/R/saemix.R
+++ b/R/saemix.R
@@ -20,47 +20,40 @@
#' @importFrom saemix saemixData saemixModel
#' @importFrom stats var
#' @examples
+#' \dontrun{
+#' library(saemix)
#' ds <- lapply(experimental_data_for_UBA_2019[6:10],
#' function(x) subset(x$data[c("name", "time", "value")]))
#' names(ds) <- paste("Dataset", 6:10)
-#' sfo_sfo <- mkinmod(parent = mkinsub("SFO", "A1"),
+#' f_mmkin_parent_p0_fixed <- mmkin("FOMC", ds, cores = 1,
+#' state.ini = c(parent = 100), fixed_initials = "parent", quiet = TRUE)
+#' m_saemix_p0_fixed <- saemix_model(f_mmkin_parent_p0_fixed["FOMC", ])
+#' d_saemix_parent <- saemix_data(f_mmkin_parent_p0_fixed)
+#' saemix_options <- list(seed = 123456, displayProgress = FALSE,
+#' save = FALSE, save.graphs = FALSE, nbiter.saemix = c(200, 80))
+#' f_saemix_p0_fixed <- saemix(m_saemix_p0_fixed, d_saemix_parent, saemix_options)
+#'
+#' f_mmkin_parent <- mmkin(c("SFO", "FOMC", "DFOP"), ds, quiet = TRUE)
+#' m_saemix_sfo <- saemix_model(f_mmkin_parent["SFO", ])
+#' m_saemix_fomc <- saemix_model(f_mmkin_parent["FOMC", ])
+#' m_saemix_dfop <- saemix_model(f_mmkin_parent["DFOP", ])
+#' d_saemix_parent <- saemix_data(f_mmkin_parent["SFO", ])
+#' f_saemix_sfo <- saemix(m_saemix_sfo, d_saemix_parent, saemix_options)
+#' f_saemix_fomc <- saemix(m_saemix_fomc, d_saemix_parent, saemix_options)
+#' f_saemix_dfop <- saemix(m_saemix_dfop, d_saemix_parent, saemix_options)
+#' compare.saemix(list(f_saemix_sfo, f_saemix_fomc, f_saemix_dfop))
+#' f_mmkin_parent_tc <- update(f_mmkin_parent, error_model = "tc")
+#' m_saemix_fomc_tc <- saemix_model(f_mmkin_parent_tc["FOMC", ])
+#' f_saemix_fomc_tc <- saemix(m_saemix_fomc_tc, d_saemix_parent, saemix_options)
+#' compare.saemix(list(f_saemix_fomc, f_saemix_fomc_tc))
+#'
+#' dfop_sfo <- mkinmod(parent = mkinsub("DFOP", "A1"),
#' A1 = mkinsub("SFO"))
-#' \dontrun{
-#' f_mmkin <- mmkin(list("SFO-SFO" = sfo_sfo), ds, quiet = TRUE)
-#' library(saemix)
-#' m_saemix <- saemix_model(f_mmkin, cores = 1)
+#' f_mmkin <- mmkin(list("DFOP-SFO" = dfop_sfo), ds, quiet = TRUE)
+#' m_saemix <- saemix_model(f_mmkin)
#' d_saemix <- saemix_data(f_mmkin)
-#' saemix_options <- list(seed = 123456,
-#' save = FALSE, save.graphs = FALSE, displayProgress = FALSE,
-#' nbiter.saemix = c(200, 80))
#' f_saemix <- saemix(m_saemix, d_saemix, saemix_options)
-#' plot(f_saemix, plot.type = "convergence")
-#' }
-#' # Synthetic data with two-component error
-#' sampling_times = c(0, 1, 3, 7, 14, 28, 60, 90, 120)
-#' dt50_sfo_in <- c(80, 90, 100, 111.111, 125)
-#' k_in <- log(2) / dt50_sfo_in
#'
-#' SFO <- mkinmod(parent = mkinsub("SFO"))
-#'
-#' pred_sfo <- function(k) {
-#' mkinpredict(SFO, c(k_parent = k),
-#' c(parent = 100), sampling_times)
-#' }
-#'
-#' ds_sfo_mean <- lapply(k_in, pred_sfo)
-#' set.seed(123456L)
-#' ds_sfo_syn <- lapply(ds_sfo_mean, function(ds) {
-#' add_err(ds, sdfunc = function(value) sqrt(1^2 + value^2 * 0.07^2),
-#' n = 1)[[1]]
-#' })
-#' \dontrun{
-#' f_mmkin_syn <- mmkin("SFO", ds_sfo_syn, error_model = "tc", quiet = TRUE)
-#' # plot(f_mmkin_syn)
-#' m_saemix_tc <- saemix_model(f_mmkin_syn, cores = 1)
-#' d_saemix_tc <- saemix_data(f_mmkin_syn)
-#' f_saemix_tc <- saemix(m_saemix_tc, d_saemix_tc, saemix_options)
-#' plot(f_saemix_tc, plot.type = "convergence")
#' }
#' @return An [saemix::SaemixModel] object.
#' @export
@@ -68,14 +61,14 @@ saemix_model <- function(object, cores = 1) {
if (nrow(object) > 1) stop("Only row objects allowed")
mkin_model <- object[[1]]$mkinmod
- analytical <- is.function(mkin_model$deg_func)
+ solution_type <- object[[1]]$solution_type
degparms_optim <- mean_degparms(object)
- psi0 <- matrix(degparms_optim, nrow = 1)
- colnames(psi0) <- names(degparms_optim)
-
degparms_fixed <- object[[1]]$bparms.fixed
+ # Transformations are done in the degradation function
+ transform.par = rep(0, length(degparms_optim))
+
odeini_optim_parm_names <- grep('_0$', names(degparms_optim), value = TRUE)
odeini_fixed_parm_names <- grep('_0$', names(degparms_fixed), value = TRUE)
@@ -85,50 +78,114 @@ saemix_model <- function(object, cores = 1) {
odeini_fixed <- degparms_fixed[odeini_fixed_parm_names]
names(odeini_fixed) <- gsub('_0$', '', odeini_fixed_parm_names)
- model_function <- function(psi, id, xidep) {
+ model_function <- FALSE
+
+ if (length(mkin_model$spec) == 1 & mkin_model$use_of_ff == "max") {
+ parent_type <- mkin_model$spec[[1]]$type
+ if (length(odeini_fixed) == 1) {
+ if (parent_type == "SFO") {
+ stop("saemix needs at least two parameters to work on.")
+ }
+ if (parent_type == "FOMC") {
+ model_function <- function(psi, id, xidep) {
+ odeini_fixed / (xidep[, "time"]/exp(psi[id, 2]) + 1)^exp(psi[id, 1])
+ }
+ }
+ if (parent_type == "DFOP") {
+ model_function <- function(psi, id, xidep) {
+ g <- plogis(psi[id, 3])
+ t = xidep[, "time"]
+ odeini_fixed * (g * exp(- exp(psi[id, 1]) * t) +
+ (1 - g) * exp(- exp(psi[id, 2]) * t))
+ }
+ }
+ if (parent_type == "HS") {
+ model_function <- function(psi, id, xidep) {
+ tb <- exp(psi[id, 3])
+ t = xidep[, "time"]
+ k1 = exp(psi[id, 1])
+ odeini_fixed * ifelse(t <= tb,
+ exp(- k1 * t),
+ exp(- k1 * t) * exp(- exp(psi[id, 2]) * (t - tb)))
+ }
+ }
+ } else {
+ if (length(odeparms_fixed) == 0) {
+ if (parent_type == "SFO") {
+ model_function <- function(psi, id, xidep) {
+ psi[id, 1] * exp( - exp(psi[id, 2]) * xidep[, "time"])
+ }
+ }
+ if (parent_type == "FOMC") {
+ model_function <- function(psi, id, xidep) {
+ psi[id, 1] / (xidep[, "time"]/exp(psi[id, 3]) + 1)^exp(psi[id, 2])
+ }
+ }
+ if (parent_type == "DFOP") {
+ model_function <- function(psi, id, xidep) {
+ g <- plogis(psi[id, 4])
+ t = xidep[, "time"]
+ psi[id, 1] * (g * exp(- exp(psi[id, 2]) * t) +
+ (1 - g) * exp(- exp(psi[id, 3]) * t))
+ }
+ }
+ if (parent_type == "HS") {
+ model_function <- function(psi, id, xidep) {
+ tb <- exp(psi[id, 4])
+ t = xidep[, "time"]
+ k1 = exp(psi[id, 2])
+ psi[id, 1] * ifelse(t <= tb,
+ exp(- k1 * t),
+ exp(- k1 * t) * exp(- exp(psi[id, 3]) * (t - tb)))
+ }
+ }
+ }
+ }
+ }
- uid <- unique(id)
+ if (!is.function(model_function)) {
+ model_function <- function(psi, id, xidep) {
- res_list <- parallel::mclapply(uid, function(i) {
- transparms_optim <- psi[i, ]
- names(transparms_optim) <- names(degparms_optim)
+ uid <- unique(id)
- odeini_optim <- transparms_optim[odeini_optim_parm_names]
- names(odeini_optim) <- gsub('_0$', '', odeini_optim_parm_names)
+ res_list <- parallel::mclapply(uid, function(i) {
+ transparms_optim <- psi[i, ]
+ names(transparms_optim) <- names(degparms_optim)
- odeini <- c(odeini_optim, odeini_fixed)[names(mkin_model$diffs)]
+ odeini_optim <- transparms_optim[odeini_optim_parm_names]
+ names(odeini_optim) <- gsub('_0$', '', odeini_optim_parm_names)
- ode_transparms_optim_names <- setdiff(names(transparms_optim), odeini_optim_parm_names)
- odeparms_optim <- backtransform_odeparms(transparms_optim[ode_transparms_optim_names], mkin_model,
- transform_rates = object[[1]]$transform_rates,
- transform_fractions = object[[1]]$transform_fractions)
- odeparms <- c(odeparms_optim, odeparms_fixed)
+ odeini <- c(odeini_optim, odeini_fixed)[names(mkin_model$diffs)]
- xidep_i <- subset(xidep, id == i)
+ ode_transparms_optim_names <- setdiff(names(transparms_optim), odeini_optim_parm_names)
+ odeparms_optim <- backtransform_odeparms(transparms_optim[ode_transparms_optim_names], mkin_model,
+ transform_rates = object[[1]]$transform_rates,
+ transform_fractions = object[[1]]$transform_fractions)
+ odeparms <- c(odeparms_optim, odeparms_fixed)
- if (analytical) {
- out_values <- mkin_model$deg_func(xidep_i, odeini, odeparms)
- } else {
+ xidep_i <- subset(xidep, id == i)
- i_time <- xidep_i$time
- i_name <- xidep_i$name
+ if (solution_type == "analytical") {
+ out_values <- mkin_model$deg_func(xidep_i, odeini, odeparms)
+ } else {
- out_wide <- mkinpredict(mkin_model,
- odeparms = odeparms, odeini = odeini,
- solution_type = object[[1]]$solution_type,
- outtimes = sort(unique(i_time)))
+ i_time <- xidep_i$time
+ i_name <- xidep_i$name
- out_index <- cbind(as.character(i_time), as.character(i_name))
- out_values <- out_wide[out_index]
- }
- return(out_values)
- }, mc.cores = cores)
- res <- unlist(res_list)
- return(res)
- }
+ out_wide <- mkinpredict(mkin_model,
+ odeparms = odeparms, odeini = odeini,
+ solution_type = solution_type,
+ outtimes = sort(unique(i_time)))
- raneff_0 <- mean_degparms(object, random = TRUE)$random$ds
- var_raneff_0 <- apply(raneff_0, 2, var)
+ out_index <- cbind(as.character(i_time), as.character(i_name))
+ out_values <- out_wide[out_index]
+ }
+ return(out_values)
+ }, mc.cores = cores)
+ res <- unlist(res_list)
+ return(res)
+ }
+ }
error.model <- switch(object[[1]]$err_mod,
const = "constant",
@@ -136,7 +193,7 @@ saemix_model <- function(object, cores = 1) {
obs = "constant")
if (object[[1]]$err_mod == "obs") {
- warning("The error model 'obs' (variance by variable) was not transferred to the saemix model")
+ warning("The error model 'obs' (variance by variable) can currently not be transferred to an saemix model")
}
error.init <- switch(object[[1]]$err_mod,
@@ -145,11 +202,15 @@ saemix_model <- function(object, cores = 1) {
b = mean(sapply(object, function(x) x$errparms[2]))),
obs = c(a = mean(sapply(object, function(x) x$errparms)), b = 1))
- res <- saemixModel(model_function, psi0,
+ psi0_matrix <- matrix(degparms_optim, nrow = 1)
+ colnames(psi0_matrix) <- names(degparms_optim)
+
+ res <- saemixModel(model_function,
+ psi0 = psi0_matrix,
"Mixed model generated from mmkin object",
- transform.par = rep(0, length(degparms_optim)),
- error.model = error.model, error.init = error.init,
- omega.init = diag(var_raneff_0)
+ transform.par = transform.par,
+ error.model = error.model,
+ error.init = error.init
)
return(res)
}
diff --git a/R/transform_odeparms.R b/R/transform_odeparms.R
index 0a25ee8c..f21d31fc 100644
--- a/R/transform_odeparms.R
+++ b/R/transform_odeparms.R
@@ -5,19 +5,19 @@
#' constants and other parameters that can only take on positive values, a
#' simple log transformation is used. For compositional parameters, such as the
#' formations fractions that should always sum up to 1 and can not be negative,
-#' the \code{\link{ilr}} transformation is used.
+#' the [ilr] transformation is used.
#'
#' The transformation of sets of formation fractions is fragile, as it supposes
#' the same ordering of the components in forward and backward transformation.
-#' This is no problem for the internal use in \code{\link{mkinfit}}.
+#' This is no problem for the internal use in [mkinfit].
#'
#' @param parms Parameters of kinetic models as used in the differential
#' equations.
#' @param transparms Transformed parameters of kinetic models as used in the
#' fitting procedure.
-#' @param mkinmod The kinetic model of class \code{\link{mkinmod}}, containing
+#' @param mkinmod The kinetic model of class [mkinmod], containing
#' the names of the model variables that are needed for grouping the
-#' formation fractions before \code{\link{ilr}} transformation, the parameter
+#' formation fractions before [ilr] transformation, the parameter
#' names and the information if the pathway to sink is included in the model.
#' @param transform_rates Boolean specifying if kinetic rate constants should
#' be transformed in the model specification used in the fitting for better
@@ -28,11 +28,15 @@
#' @param transform_fractions Boolean specifying if formation fractions
#' constants should be transformed in the model specification used in the
#' fitting for better compliance with the assumption of normal distribution
-#' of the estimator. The default (TRUE) is to do transformations. The g
-#' parameter of the DFOP and HS models are also transformed, as they can also
-#' be seen as compositional data. The transformation used for these
-#' transformations is the \code{\link{ilr}} transformation.
+#' of the estimator. The default (TRUE) is to do transformations.
+#' The g parameter of the DFOP model is also seen as a fraction.
+#' If a single fraction is transformed (g parameter of DFOP or only a single
+#' target variable e.g. a single metabolite plus a pathway to sink), a
+#' logistic transformation is used [stats::qlogis()]. In other cases, i.e. if
+#' two or more formation fractions need to be transformed whose sum cannot
+#' exceed one, the [ilr] transformation is used.
#' @return A vector of transformed or backtransformed parameters
+#' @importFrom stats plogis qlogis
#' @author Johannes Ranke
#' @examples
#'
@@ -91,8 +95,7 @@
#'
#' @export transform_odeparms
transform_odeparms <- function(parms, mkinmod,
- transform_rates = TRUE,
- transform_fractions = TRUE)
+ transform_rates = TRUE, transform_fractions = TRUE)
{
# We need the model specification for the names of the model
# variables and the information on the sink
@@ -119,8 +122,7 @@ transform_odeparms <- function(parms, mkinmod,
N <- parms[grep("^N", names(parms))]
transparms[names(N)] <- N
- # Go through state variables and apply isometric logratio transformation to
- # formation fractions if requested
+ # Go through state variables and transform formation fractions if requested
mod_vars = names(spec)
for (box in mod_vars) {
f <- parms[grep(paste("^f", box, sep = "_"), names(parms))]
@@ -128,9 +130,14 @@ transform_odeparms <- function(parms, mkinmod,
if (length(f) > 0) {
if(transform_fractions) {
if (spec[[box]]$sink) {
- trans_f <- ilr(c(f, 1 - sum(f)))
- trans_f_names <- paste("f", box, "ilr", 1:length(trans_f), sep = "_")
- transparms[trans_f_names] <- trans_f
+ if (length(f) == 1) {
+ trans_f_name <- paste("f", box, "qlogis", sep = "_")
+ transparms[trans_f_name] <- qlogis(f)
+ } else {
+ trans_f <- ilr(c(f, 1 - sum(f)))
+ trans_f_names <- paste("f", box, "ilr", 1:length(trans_f), sep = "_")
+ transparms[trans_f_names] <- trans_f
+ }
} else {
if (length(f) > 1) {
trans_f <- ilr(f)
@@ -161,7 +168,7 @@ transform_odeparms <- function(parms, mkinmod,
if (!is.na(parms["g"])) {
g <- parms["g"]
if (transform_fractions) {
- transparms["g_ilr"] <- ilr(c(g, 1 - g))
+ transparms["g_qlogis"] <- qlogis(g)
} else {
transparms["g"] <- g
}
@@ -207,20 +214,25 @@ backtransform_odeparms <- function(transparms, mkinmod,
N <- transparms[grep("^N", names(transparms))]
parms[names(N)] <- N
- # Go through state variables and apply inverse isometric logratio transformation
+ # Go through state variables and apply inverse transformations to formation fractions
mod_vars = names(spec)
for (box in mod_vars) {
# Get the names as used in the model
f_names = grep(paste("^f", box, sep = "_"), mkinmod$parms, value = TRUE)
+
# Get the formation fraction parameters
trans_f = transparms[grep(paste("^f", box, sep = "_"), names(transparms))]
if (length(trans_f) > 0) {
if(transform_fractions) {
- f <- invilr(trans_f)
- if (spec[[box]]$sink) {
- parms[f_names] <- f[1:length(f)-1]
+ if (any(grepl("qlogis", names(trans_f)))) {
+ parms[f_names] <- plogis(trans_f)
} else {
- parms[f_names] <- f
+ f <- invilr(trans_f)
+ if (spec[[box]]$sink) {
+ parms[f_names] <- f[1:length(f)-1]
+ } else {
+ parms[f_names] <- f
+ }
}
} else {
parms[names(trans_f)] <- trans_f
@@ -242,7 +254,12 @@ backtransform_odeparms <- function(transparms, mkinmod,
}
}
- # DFOP parameter g is treated as a fraction
+ # DFOP parameter g is now transformed using qlogis
+ if (!is.na(transparms["g_qlogis"])) {
+ g_qlogis <- transparms["g_qlogis"]
+ parms["g"] <- plogis(g_qlogis)
+ }
+ # In earlier times we used ilr for g, so we keep this around
if (!is.na(transparms["g_ilr"])) {
g_ilr <- transparms["g_ilr"]
parms["g"] <- invilr(g_ilr)[1]

Contact - Imprint