aboutsummaryrefslogtreecommitdiff
path: root/R/mhmkin.R
blob: 6265a59e5ffdee63e39f7804560cf2ab8e045246 (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
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
#' Fit nonlinear mixed-effects models built from one or more kinetic
#' degradation models and one or more error models
#'
#' The name of the methods expresses that (**m**ultiple) **h**ierarchichal
#' (also known as multilevel) **m**ulticompartment **kin**etic models are
#' fitted. Our kinetic models are nonlinear, so we can use various nonlinear
#' mixed-effects model fitting functions.
#'
#' @param objects A list of [mmkin] objects containing fits of the same
#' degradation models to the same data, but using different error models.
#' Alternatively, a single [mmkin] object containing fits of several
#' degradation models to the same data
#' @param backend The backend to be used for fitting. Currently, only saemix is
#' supported
#' @param no_random_effect Default is NULL and will be passed to [saem]. If a
#' character vector is supplied, it will be passed to all calls to [saem],
#' which will exclude random effects for all matching parameters. Alternatively,
#' a list of character vectors or an object of class [illparms.mhmkin] can be
#' specified. They have to have the same dimensions that the return object of
#' the current call will have, i.e. the number of rows must match the number
#' of degradation models in the mmkin object(s), and the number of columns must
#' match the number of error models used in the mmkin object(s).
#' @param algorithm The algorithm to be used for fitting (currently not used)
#' @param \dots Further arguments that will be passed to the nonlinear mixed-effects
#' model fitting function.
#' @param cores The number of cores to be used for multicore processing. This
#' is only used when the \code{cluster} argument is \code{NULL}. On Windows
#' machines, cores > 1 is not supported, you need to use the \code{cluster}
#' argument to use multiple logical processors. Per default, all cores detected
#' by [parallel::detectCores()] are used, except on Windows where the default
#' is 1.
#' @param cluster A cluster as returned by [makeCluster] to be used for
#' parallel execution.
#' @importFrom parallel mclapply parLapply detectCores
#' @return A two-dimensional [array] of fit objects and/or try-errors that can
#' be indexed using the degradation model names for the first index (row index)
#' and the error model names for the second index (column index), with class
#' attribute 'mhmkin'.
#' @author Johannes Ranke
#' @seealso \code{\link{[.mhmkin}} for subsetting [mhmkin] objects
#' @export
mhmkin <- function(objects, ...) {
  UseMethod("mhmkin")
}

#' @export
#' @rdname mhmkin
mhmkin.mmkin <- function(objects, ...) {
  mhmkin(list(objects), ...)
}

#' @export
#' @rdname mhmkin
#' @examples
#' \dontrun{
#' # We start with separate evaluations of all the first six datasets with two
#' # degradation models and two error models
#' f_sep_const <- mmkin(c("SFO", "FOMC"), ds_fomc[1:6], cores = 2, quiet = TRUE)
#' f_sep_tc <- update(f_sep_const, error_model = "tc")
#' # The mhmkin function sets up hierarchical degradation models aka
#' # nonlinear mixed-effects models for all four combinations, specifying
#' # uncorrelated random effects for all degradation parameters
#' f_saem_1 <- mhmkin(list(f_sep_const, f_sep_tc), cores = 2)
#' status(f_saem_1)
#' # The 'illparms' function shows that in all hierarchical fits, at least
#' # one random effect is ill-defined (the confidence interval for the
#' # random effect expressed as standard deviation includes zero)
#' illparms(f_saem_1)
#' # Therefore we repeat the fits, excluding the ill-defined random effects
#' f_saem_2 <- update(f_saem_1, no_random_effect = illparms(f_saem_1))
#' status(f_saem_2)
#' illparms(f_saem_2)
#' # Model comparisons show that FOMC with two-component error is preferable,
#' # and confirms our reduction of the default parameter model
#' anova(f_saem_1)
#' anova(f_saem_2)
#' # The convergence plot for the selected model looks fine
#' saemix::plot(f_saem_2[["FOMC", "tc"]]$so, plot.type = "convergence")
#' # The plot of predictions versus data shows that we have a pretty data-rich
#' # situation with homogeneous distribution of residuals, because we used the
#' # same degradation model, error model and parameter distribution model that
#' # was used in the data generation.
#' plot(f_saem_2[["FOMC", "tc"]])
#' # We can specify the same parameter model reductions manually
#' no_ranef <- list("parent_0", "log_beta", "parent_0", c("parent_0", "log_beta"))
#' dim(no_ranef) <- c(2, 2)
#' f_saem_2m <- update(f_saem_1, no_random_effect = no_ranef)
#' anova(f_saem_2m)
#' }
mhmkin.list <- function(objects, backend = "saemix", algorithm = "saem",
  no_random_effect = NULL,
  ...,
  cores = if (Sys.info()["sysname"] == "Windows") 1 else parallel::detectCores(), cluster = NULL)
{
  call <- match.call()
  dot_args <- list(...)
  backend_function <- switch(backend,
    saemix = "saem"
  )

  deg_models <- lapply(objects[[1]][, 1], function(x) x$mkinmod)
  names(deg_models) <- dimnames(objects[[1]])$model
  n.deg <- length(deg_models)

  ds <- lapply(objects[[1]][1, ], function(x) x$data)

  for (other in objects[-1]) {
    # Check if the degradation models in all objects are the same
    for (deg_model_name in names(deg_models)) {
      if (!all.equal(other[[deg_model_name, 1]]$mkinmod$spec,
            deg_models[[deg_model_name]]$spec))
      {
        stop("The mmkin objects have to be based on the same degradation models")
      }
    }
    # Check if they have been fitted to the same dataset
    other_object_ds <- lapply(other[1, ], function(x) x$data)
    for (i in 1:length(ds)) {
      if (!all.equal(ds[[i]][c("time", "variable", "observed")],
          other_object_ds[[i]][c("time", "variable", "observed")]))
      {
        stop("The mmkin objects have to be fitted to the same datasets")
      }
    }
  }

  n.o <- length(objects)

  error_models = sapply(objects, function(x) x[[1]]$err_mod)
  n.e <- length(error_models)

  n.fits <- n.deg * n.e
  fit_indices <- matrix(1:n.fits, ncol = n.e)
  dimnames(fit_indices) <- list(degradation = names(deg_models),
                                error = error_models)

  if (is.null(no_random_effect) || is.null(dim(no_random_effect))) {
    no_ranef <- rep(list(no_random_effect), n.fits)
    dim(no_ranef) <- dim(fit_indices)
  } else {
    if (!identical(dim(no_random_effect), dim(fit_indices))) {
      stop("Dimensions of argument 'no_random_effect' are not suitable")
    }
    if (is(no_random_effect, "illparms.mhmkin")) {
      no_ranef_dim <- dim(no_random_effect)
      no_ranef <- lapply(no_random_effect, function(x) {
        no_ranef_split <- strsplit(x, ", ")
        ret <- sapply(no_ranef_split, function(y) {
          gsub("sd\\((.*)\\)", "\\1", y)
        })
        return(ret)
      })
      dim(no_ranef) <- no_ranef_dim
    } else {
      no_ranef <- no_random_effect
    }
  }

  fit_function <- function(fit_index) {
    w <- which(fit_indices == fit_index, arr.ind = TRUE)
    deg_index <- w[1]
    error_index <- w[2]
    mmkin_row <- objects[[error_index]][deg_index, ]
    res <- try(do.call(backend_function,
        args = c(
          list(mmkin_row),
          dot_args,
          list(no_random_effect = no_ranef[[deg_index, error_index]]))))
    return(res)
  }


  fit_time <- system.time({
    if (is.null(cluster)) {
      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
  attr(results, "time") <- fit_time
  class(results) <- switch(backend,
    saemix = c("mhmkin.saem.mmkin", "mhmkin")
  )
  return(results)
}

#' Subsetting method for mhmkin objects
#'
#' @param x An [mhmkin] object.
#' @param i Row index selecting the fits for specific models
#' @param j Column index selecting the fits to specific datasets
#' @param drop If FALSE, the method always returns an mhmkin object, otherwise
#'   either a list of fit objects or a single fit object.
#' @return An object inheriting from \code{\link{mhmkin}}.
#' @rdname mhmkin
#' @export
`[.mhmkin` <- function(x, i, j, ..., drop = FALSE) {
  original_class <- class(x)
  class(x) <- NULL
  x_sub <- x[i, j, drop = drop]

  if (!drop) {
    class(x_sub) <- original_class
  }
  return(x_sub)
}

#' Print method for mhmkin objects
#'
#' @rdname mhmkin
#' @export
print.mhmkin <- function(x, ...) {
  cat("<mhmkin> object\n")
  cat("Status of individual fits:\n\n")
  print(status(x))
}

#' @export
AIC.mhmkin <- function(object, ..., k = 2) {
  if (inherits(object[[1]], "saem.mmkin")) {
    check_failed <- function(x) if (inherits(x$so, "try-error")) TRUE else FALSE
  }
  res <- sapply(object, function(x) {
    if (check_failed(x)) return(NA)
    else return(AIC(x$so, k = k))
  })
  dim(res) <- dim(object)
  dimnames(res) <- dimnames(object)
  return(res)
}

#' @export
BIC.mhmkin <- function(object, ...) {
  if (inherits(object[[1]], "saem.mmkin")) {
    check_failed <- function(x) if (inherits(x$so, "try-error")) TRUE else FALSE
  }
  res <- sapply(object, function(x) {
    if (check_failed(x)) return(NA)
    else return(BIC(x$so))
  })
  dim(res) <- dim(object)
  dimnames(res) <- dimnames(object)
  return(res)
}

#' @export
update.mhmkin <- function(object, ..., evaluate = TRUE) {
  call <- attr(object, "call")
  # For some reason we get mhkin.list in call[[1]] when using mhmkin from the
  # loaded package so we need to fix this so we do not have to export
  # mhmkin.list in addition to the S3 method mhmkin
  call[[1]] <- mhmkin

  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
}

#' @export
anova.mhmkin <- function(object, ...,
  method = c("is", "lin", "gq"), test = FALSE, model.names = "auto") {
  if (identical(model.names, "auto")) {
    model.names <- outer(rownames(object), colnames(object), paste)
  }
  rlang::inject(anova(!!!(object), method = method, test = test,
      model.names = model.names))
}

Contact - Imprint