aboutsummaryrefslogtreecommitdiff
path: root/R/nlme.mmkin.R
blob: e58f11cb96493f1b91728b3beec306e519bcb0ea (plain) (blame)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
# Code inspired by nlme::nlme.nlsList and R/nlme_fit.R from nlmixr

# We need to assign the degradation function created in nlme.mmkin to an
# environment that is always accessible, also e.g. when evaluation is done by
# testthat or pkgdown. Therefore parent.frame() is not good enough. The
# following environment will be in the mkin namespace.
.nlme_env <- new.env(parent = emptyenv())

#' Retrieve a degradation function from the mmkin namespace
#'
#' @importFrom utils getFromNamespace
#' @return A function that was likely previously assigned from within
#'   nlme.mmkin
#' @export
get_deg_func <- function() {
  return(get("deg_func", getFromNamespace(".nlme_env", "mkin")))
}

#' Create an nlme model for an mmkin row object
#'
#' This functions sets up a nonlinear mixed effects model for an mmkin row
#' object. An mmkin row object is essentially a list of mkinfit objects that
#' have been obtained by fitting the same model to a list of datasets.
#'
#' @param model An \code{\link{mmkin}} row object.
#' @param data Ignored, data are taken from the mmkin model
#' @param fixed Ignored, all degradation parameters fitted in the
#'   mmkin model are used as fixed parameters
#' @param random If not specified, all fixed effects are complemented
#'   with uncorrelated random effects
#' @param groups See the documentation of nlme
#' @param start If not specified, mean values of the fitted degradation
#'   parameters taken from the mmkin object are used
#' @param correlation See the documentation of nlme
#' @param weights passed to nlme
#' @param subset passed to nlme
#' @param method passed to nlme
#' @param na.action passed to nlme
#' @param naPattern passed to nlme
#' @param control passed to nlme
#' @param verbose passed to nlme
#' @importFrom stats na.fail
#' @return Upon success, a fitted nlme.mmkin object, which is an nlme object
#'   with additional elements
#' @export
#' @seealso \code{\link{nlme_function}}
#' @examples
#' ds <- lapply(experimental_data_for_UBA_2019[6:10],
#'  function(x) subset(x$data[c("name", "time", "value")], name == "parent"))
#' f <- mmkin("SFO", ds, quiet = TRUE, cores = 1)
#' library(nlme)
#' endpoints(f[[1]])
#' f_nlme <- nlme(f)
#' print(f_nlme)
#' endpoints(f_nlme)
#' \dontrun{
#'   f_nlme_2 <- nlme(f, start = c(parent_0 = 100, log_k_parent_sink = 0.1))
#'   update(f_nlme_2, random = parent_0 ~ 1)
#'   # Test on some real data
#'   ds_2 <- lapply(experimental_data_for_UBA_2019[6:10],
#'    function(x) x$data[c("name", "time", "value")])
#'   m_sfo_sfo <- mkinmod(parent = mkinsub("SFO", "A1"),
#'     A1 = mkinsub("SFO"), use_of_ff = "min", quiet = TRUE)
#'   m_sfo_sfo_ff <- mkinmod(parent = mkinsub("SFO", "A1"),
#'     A1 = mkinsub("SFO"), use_of_ff = "max", quiet = TRUE)
#'   m_fomc_sfo <- mkinmod(parent = mkinsub("FOMC", "A1"),
#'     A1 = mkinsub("SFO"), quiet = TRUE)
#'   m_dfop_sfo <- mkinmod(parent = mkinsub("DFOP", "A1"),
#'     A1 = mkinsub("SFO"), quiet = TRUE)
#'
#'   f_2 <- mmkin(list("SFO-SFO" = m_sfo_sfo,
#'    "SFO-SFO-ff" = m_sfo_sfo_ff,
#'    "FOMC-SFO" = m_fomc_sfo,
#'    "DFOP-SFO" = m_dfop_sfo),
#'     ds_2, quiet = TRUE)
#'   plot(f_2["SFO-SFO", 3:4]) # Separate fits for datasets 3 and 4
#'
#'   f_nlme_sfo_sfo <- nlme(f_2["SFO-SFO", ])
#'   # plot(f_nlme_sfo_sfo) # not feasible with pkgdown figures
#'   plot(f_nlme_sfo_sfo, 3:4) # Global mixed model: Fits for datasets 3 and 4
#'
#'   # With formation fractions
#'   f_nlme_sfo_sfo_ff <- nlme(f_2["SFO-SFO-ff", ])
#'   plot(f_nlme_sfo_sfo_ff, 3:4) # chi2 different due to different df attribution
#'
#'   # For more parameters, we need to increase pnlsMaxIter and the tolerance
#'   # to get convergence
#'   f_nlme_fomc_sfo <- nlme(f_2["FOMC-SFO", ],
#'     control = list(pnlsMaxIter = 100, tolerance = 1e-4), verbose = TRUE)
#'   f_nlme_dfop_sfo <- nlme(f_2["DFOP-SFO", ],
#'     control = list(pnlsMaxIter = 120, tolerance = 5e-4), verbose = TRUE)
#'   plot(f_2["FOMC-SFO", 3:4])
#'   plot(f_nlme_fomc_sfo, 3:4)
#'
#'   plot(f_2["DFOP-SFO", 3:4])
#'   plot(f_nlme_dfop_sfo, 3:4)
#'
#'   anova(f_nlme_dfop_sfo, f_nlme_fomc_sfo, f_nlme_sfo_sfo)
#'   anova(f_nlme_dfop_sfo, f_nlme_sfo_sfo) # if we ignore FOMC
#'
#'   endpoints(f_nlme_sfo_sfo)
#'   endpoints(f_nlme_dfop_sfo)
#' }
nlme.mmkin <- function(model, data = sys.frame(sys.parent()),
  fixed, random = fixed,
  groups, start, correlation = NULL, weights = NULL,
  subset, method = c("ML", "REML"),
  na.action = na.fail, naPattern,
  control = list(), verbose= FALSE)
{
  if (nrow(model) > 1) stop("Only row objects allowed")

  thisCall <- as.list(match.call())[-1]

  # Warn in case arguments were used that are overriden
  if (any(!is.na(match(names(thisCall),
               c("fixed", "data"))))) {
    warning("'nlme.mmkin' will redefine 'fixed' and 'data'")
  }

  deg_func <- nlme_function(model)

  assign("deg_func", deg_func, getFromNamespace(".nlme_env", "mkin"))

  # For the formula, get the degradation function from the mkin namespace
  this_model_text <- paste0("value ~ mkin::get_deg_func()(",
    paste(names(formals(deg_func)), collapse = ", "), ")")
  this_model <- as.formula(this_model_text)

  thisCall[["model"]] <- this_model

  mean_dp <- mean_degparms(model)
  dp_names <- names(mean_dp)

  thisCall[["data"]] <- nlme_data(model)

  if (missing(start)) {
    thisCall[["start"]] <- mean_degparms(model, random = TRUE)
  }

  thisCall[["fixed"]] <- lapply(as.list(dp_names), function(el)
                                 eval(parse(text = paste(el, 1, sep = "~"))))

  if (missing(random)) {
    thisCall[["random"]] <- pdDiag(thisCall[["fixed"]])
  }

  val <- do.call("nlme.formula", thisCall)
  val$mmkin_orig <- model
  class(val) <- c("nlme.mmkin", "nlme", "lme")
  return(val)
}

#' @export
#' @rdname nlme.mmkin
#' @param x An nlme.mmkin object to print
#' @param data Should the data be printed?
#' @param ... Further arguments as in the generic
print.nlme.mmkin <- function(x, ...) {
  x$call$data <- "Not shown"
  NextMethod("print", x)
}

#' @export
#' @rdname nlme.mmkin
#' @param object An nlme.mmkin object to update
#' @param ... Update specifications passed to update.nlme
update.nlme.mmkin <- function(object, ...) {
  res <- NextMethod()
  res$mmkin_orig <- object$mmkin_orig
  class(res) <- c("nlme.mmkin", "nlme", "lme")
  return(res)
}

# The following is necessary as long as R bug 17761 is not fixed
# https://bugs.r-project.org/bugzilla/show_bug.cgi?id=17761
#' @export
anova.nlme.mmkin <- function(object, ...) {
  thisCall <- as.list(match.call())[-1]
  object_name <- as.character(thisCall[[1]])
  other_object_names <- sapply(thisCall[-1], as.character)

  remove_class <- function(object, classname) {
    old_class <- class(object)
    class(object) <- setdiff(old_class, classname)
    return(object)
  }
  object <- remove_class(object, "nlme.mmkin")
  other_objects <- list(...)
  other_objects <- lapply(other_objects, remove_class, "nlme.mmkin")

  env <- new.env()
  assign(object_name, object, env)
  for (i in seq_along(other_objects)) {
    assign(other_object_names[i], other_objects[[i]], env)
  }
  res <- eval(parse(text = paste0("anova.lme(", object_name, ", ",
        paste(other_object_names, collapse = ", "), ")")),
    envir = env)

  return(res)
}

Contact - Imprint