aboutsummaryrefslogtreecommitdiff
path: root/R/mkinfit.R
blob: 4ac54ce26d7a093999b7eb7ac9b2d523d7d24297 (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
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
# Copyright (C) 2010-2019 Johannes Ranke
# Portions of this code are copyright (C) 2013 Eurofins Regulatory AG
# Contact: jranke@uni-bremen.de
# The summary function is an adapted and extended version of summary.modFit
# from the FME package, v 1.1 by Soetart and Petzoldt, which was in turn
# inspired by summary.nls.lm

# This file is part of the R package mkin

# mkin is free software: you can redistribute it and/or modify it under the
# terms of the GNU General Public License as published by the Free Software
# Foundation, either version 3 of the License, or (at your option) any later
# version.

# This program is distributed in the hope that it will be useful, but WITHOUT
# ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
# FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
# details.

# You should have received a copy of the GNU General Public License along with
# this program. If not, see <http://www.gnu.org/licenses/>
if(getRversion() >= '2.15.1') utils::globalVariables(c("name", "time", "value"))

mkinfit <- function(mkinmod, observed,
  parms.ini = "auto",
  state.ini = "auto",
  fixed_parms = NULL,
  fixed_initials = names(mkinmod$diffs)[-1],
  from_max_mean = FALSE,
  solution_type = c("auto", "analytical", "eigen", "deSolve"),
  method.ode = "lsoda",
  use_compiled = "auto",
  method.modFit = c("Port", "Marq", "SANN", "Nelder-Mead", "BFGS", "CG", "L-BFGS-B"),
  maxit.modFit = "auto",
  control.modFit = list(),
  transform_rates = TRUE,
  transform_fractions = TRUE,
  plot = FALSE, quiet = FALSE,
  err = NULL,
  weight = c("none", "manual", "std", "mean", "tc"),
  tc = c(sigma_low = 0.5, rsd_high = 0.07),
  scaleVar = FALSE,
  atol = 1e-8, rtol = 1e-10, n.outtimes = 100,
  reweight.method = NULL,
  reweight.tol = 1e-8, reweight.max.iter = 10,
  trace_parms = FALSE,
  ...)
{
  # Check mkinmod and generate a model for the variable whith the highest value
  # if a suitable string is given
  parent_models_available = c("SFO", "FOMC", "DFOP", "HS", "SFORB", "IORE", "logistic")
  if (class(mkinmod) != "mkinmod") {
    presumed_parent_name = observed[which.max(observed$value), "name"]
    if (mkinmod[[1]] %in% parent_models_available) {
      speclist <- list(list(type = mkinmod, sink = TRUE))
      names(speclist) <- presumed_parent_name
      mkinmod <- mkinmod(speclist = speclist)
    } else {
      stop("Argument mkinmod must be of class mkinmod or a string containing one of\n  ",
           paste(parent_models_available, collapse = ", "))
    }
  }

  # Check optimisation method and set maximum number of iterations if specified
  method.modFit = match.arg(method.modFit)
  if (maxit.modFit != "auto") {
    if (method.modFit == "Marq") control.modFit$maxiter = maxit.modFit
    if (method.modFit == "Port") {
      control.modFit$iter.max = maxit.modFit
      control.modFit$eval.max = maxit.modFit
    }
    if (method.modFit %in% c("SANN", "Nelder-Mead", "BFGS", "CG", "L-BFGS-B")) {
        control.modFit$maxit = maxit.modFit
    }
  }

  # Get the names of the state variables in the model
  mod_vars <- names(mkinmod$diffs)

  # Get the names of observed variables
  obs_vars <- names(mkinmod$spec)

  # Subset observed data with names of observed data in the model and remove NA values
  observed <- subset(observed, name %in% obs_vars)
  observed <- subset(observed, !is.na(value))

  # Obtain data for decline from maximum mean value if requested
  if (from_max_mean) {
    # This is only used for simple decline models
    if (length(obs_vars) > 1)
      stop("Decline from maximum is only implemented for models with a single observed variable")

    means <- aggregate(value ~ time, data = observed, mean, na.rm=TRUE)
    t_of_max <- means[which.max(means$value), "time"]
    observed <- subset(observed, time >= t_of_max)
    observed$time <- observed$time - t_of_max
  }

  # Define starting values for parameters where not specified by the user
  if (parms.ini[[1]] == "auto") parms.ini = vector()

  # Warn for inital parameter specifications that are not in the model
  wrongpar.names <- setdiff(names(parms.ini), mkinmod$parms)
  if (length(wrongpar.names) > 0) {
    warning("Initial parameter(s) ", paste(wrongpar.names, collapse = ", "),
         " not used in the model")
  }

  # Warn that the sum of formation fractions may exceed one if they are not
  # fitted in the transformed way
  if (mkinmod$use_of_ff == "max" & transform_fractions == FALSE) {
    warning("The sum of formation fractions may exceed one if you do not use ",
            "transform_fractions = TRUE." )
    for (box in mod_vars) {
      # Stop if formation fractions are not transformed and we have no sink
      if (mkinmod$spec[[box]]$sink == FALSE) {
        stop("If formation fractions are not transformed during the fitting, ",
             "it is not supported to turn off pathways to sink.\n ",
             "Consider turning on the transformation of formation fractions or ",
             "setting up a model with use_of_ff = 'min'.\n")
      }
    }
  }

  # Do not allow fixing formation fractions if we are using the ilr transformation,
  # this is not supported
  if (transform_fractions == TRUE && length(fixed_parms) > 0) {
    if (any(grepl("^f_", fixed_parms))) {
      stop("Fixing formation fractions is not supported when using the ilr ",
           "transformation.")
    }
  }

  # Set initial parameter values, including a small increment (salt)
  # to avoid linear dependencies (singular matrix) in Eigenvalue based solutions
  k_salt = 0
  defaultpar.names <- setdiff(mkinmod$parms, names(parms.ini))
  for (parmname in defaultpar.names) {
    # Default values for rate constants, depending on the parameterisation
    if (grepl("^k", parmname)) {
      parms.ini[parmname] = 0.1 + k_salt
      k_salt = k_salt + 1e-4
    }
    # Default values for rate constants for reversible binding
    if (grepl("free_bound$", parmname)) parms.ini[parmname] = 0.1
    if (grepl("bound_free$", parmname)) parms.ini[parmname] = 0.02
    # Default values for IORE exponents
    if (grepl("^N", parmname)) parms.ini[parmname] = 1.1
    # Default values for the FOMC, DFOP and HS models
    if (parmname == "alpha") parms.ini[parmname] = 1
    if (parmname == "beta") parms.ini[parmname] = 10
    if (parmname == "k1") parms.ini[parmname] = 0.1
    if (parmname == "k2") parms.ini[parmname] = 0.01
    if (parmname == "tb") parms.ini[parmname] = 5
    if (parmname == "g") parms.ini[parmname] = 0.5
    if (parmname == "kmax") parms.ini[parmname] = 0.1
    if (parmname == "k0") parms.ini[parmname] = 0.0001
    if (parmname == "r") parms.ini[parmname] = 0.2
  }
  # Default values for formation fractions in case they are present
  for (box in mod_vars) {
    f_names <- mkinmod$parms[grep(paste0("^f_", box), mkinmod$parms)]
    if (length(f_names) > 0) {
      # We need to differentiate between default and specified fractions
      # and set the unspecified to 1 - sum(specified)/n_unspecified
      f_default_names <- intersect(f_names, defaultpar.names)
      f_specified_names <- setdiff(f_names, defaultpar.names)
      sum_f_specified = sum(parms.ini[f_specified_names])
      if (sum_f_specified > 1) {
        stop("Starting values for the formation fractions originating from ",
             box, " sum up to more than 1.")
      }
      if (mkinmod$spec[[box]]$sink) n_unspecified = length(f_default_names) + 1
      else {
        n_unspecified = length(f_default_names)
      }
      parms.ini[f_default_names] <- (1 - sum_f_specified) / n_unspecified
    }
  }

  # Set default for state.ini if appropriate
  parent_name = names(mkinmod$spec)[[1]]
  if (state.ini[1] == "auto") {
    parent_time_0 = subset(observed, time == 0 & name == parent_name)$value
    parent_time_0_mean = mean(parent_time_0, na.rm = TRUE)
    if (is.na(parent_time_0_mean)) {
      state.ini = c(100, rep(0, length(mkinmod$diffs) - 1))
    } else {
      state.ini = c(parent_time_0_mean, rep(0, length(mkinmod$diffs) - 1))
    }
  }

  # Name the inital state variable values if they are not named yet
  if(is.null(names(state.ini))) names(state.ini) <- mod_vars

  # Transform initial parameter values for fitting
  transparms.ini <- transform_odeparms(parms.ini, mkinmod,
                                       transform_rates = transform_rates,
                                       transform_fractions = transform_fractions)

  # Parameters to be optimised:
  # Kinetic parameters in parms.ini whose names are not in fixed_parms
  parms.fixed <- parms.ini[fixed_parms]
  parms.optim <- parms.ini[setdiff(names(parms.ini), fixed_parms)]

  transparms.fixed <- transform_odeparms(parms.fixed, mkinmod,
                                       transform_rates = transform_rates,
                                       transform_fractions = transform_fractions)
  transparms.optim <- transform_odeparms(parms.optim, mkinmod,
                                       transform_rates = transform_rates,
                                       transform_fractions = transform_fractions)

  # Inital state variables in state.ini whose names are not in fixed_initials
  state.ini.fixed <- state.ini[fixed_initials]
  state.ini.optim <- state.ini[setdiff(names(state.ini), fixed_initials)]

  # Preserve names of state variables before renaming initial state variable
  # parameters
  state.ini.optim.boxnames <- names(state.ini.optim)
  state.ini.fixed.boxnames <- names(state.ini.fixed)
  if(length(state.ini.optim) > 0) {
    names(state.ini.optim) <- paste(names(state.ini.optim), "0", sep="_")
  }
  if(length(state.ini.fixed) > 0) {
    names(state.ini.fixed) <- paste(names(state.ini.fixed), "0", sep="_")
  }

  # Decide if the solution of the model can be based on a simple analytical
  # formula, the spectral decomposition of the matrix (fundamental system)
  # or a numeric ode solver from the deSolve package
  # Prefer deSolve over eigen if a compiled model is present and use_compiled
  # is not set to FALSE
  solution_type = match.arg(solution_type)
  if (solution_type == "analytical" && length(mkinmod$spec) > 1)
     stop("Analytical solution not implemented for models with metabolites.")
  if (solution_type == "eigen" && !is.matrix(mkinmod$coefmat))
     stop("Eigenvalue based solution not possible, coefficient matrix not present.")
  if (solution_type == "auto") {
    if (length(mkinmod$spec) == 1) {
      solution_type = "analytical"
    } else {
      if (!is.null(mkinmod$cf) & use_compiled[1] != FALSE) {
        solution_type = "deSolve"
      } else {
        if (is.matrix(mkinmod$coefmat)) {
          solution_type = "eigen"
          if (max(observed$value, na.rm = TRUE) < 0.1) {
            stop("The combination of small observed values (all < 0.1) and solution_type = eigen is error-prone")
          }
        } else {
          solution_type = "deSolve"
        }
      }
    }
  }

  # Define outtimes for model solution.
  # Include time points at which observed data are available
  outtimes = sort(unique(c(observed$time, seq(min(observed$time),
                                              max(observed$time),
                                              length.out = n.outtimes))))

  weight.ini <- weight <- match.arg(weight)
  if (weight.ini == "tc") {
     observed$err = sigma_twocomp(observed$value, tc["sigma_low"], tc["rsd_high"])
     err <- "err"
  } else {
    if (!is.null(err)) weight.ini = "manual"
  }

  cost.old <- 1e100 # The first model cost should be smaller than this value
  calls <- 0 # Counter for number of model solutions
  out_predicted <- NA

  # Define the model cost function for optimisation, including (back)transformations
  cost <- function(P)
  {
    assign("calls", calls+1, inherits=TRUE) # Increase the model solution counter

    # Trace parameter values if requested
    if(trace_parms) cat(P, "\n")

    if(length(state.ini.optim) > 0) {
      odeini <- c(P[1:length(state.ini.optim)], state.ini.fixed)
      names(odeini) <- c(state.ini.optim.boxnames, state.ini.fixed.boxnames)
    } else {
      odeini <- state.ini.fixed
      names(odeini) <- state.ini.fixed.boxnames
    }

    odeparms <- c(P[(length(state.ini.optim) + 1):length(P)], transparms.fixed)

    parms <- backtransform_odeparms(odeparms, mkinmod,
                                    transform_rates = transform_rates,
                                    transform_fractions = transform_fractions)

    # Solve the system with current transformed parameter values
    out <- mkinpredict(mkinmod, parms,
                       odeini, outtimes,
                       solution_type = solution_type,
                       use_compiled = use_compiled,
                       method.ode = method.ode,
                       atol = atol, rtol = rtol, ...)

    assign("out_predicted", out, inherits=TRUE)

    mC <- modCost(out, observed, y = "value",
      err = err, weight = weight, scaleVar = scaleVar)

    # Report and/or plot if the model is improved
    if (mC$model < cost.old) {
      if(!quiet) cat("Model cost at call ", calls, ": ", mC$model, "\n")

      # Plot the data and current model output if requested
      if(plot) {
        outtimes_plot = seq(min(observed$time), max(observed$time), length.out=100)

        out_plot <- mkinpredict(mkinmod, parms,
                                odeini, outtimes_plot,
                                solution_type = solution_type,
                                use_compiled = use_compiled,
                                method.ode = method.ode,
                                atol = atol, rtol = rtol, ...)

        plot(0, type="n",
          xlim = range(observed$time), ylim = c(0, max(observed$value, na.rm=TRUE)),
          xlab = "Time", ylab = "Observed")
        col_obs <- pch_obs <- 1:length(obs_vars)
        lty_obs <- rep(1, length(obs_vars))
        names(col_obs) <- names(pch_obs) <- names(lty_obs) <- obs_vars
        for (obs_var in obs_vars) {
          points(subset(observed, name == obs_var, c(time, value)),
                 pch = pch_obs[obs_var], col = col_obs[obs_var])
        }
        matlines(out_plot$time, out_plot[-1], col = col_obs, lty = lty_obs)
        legend("topright", inset=c(0.05, 0.05), legend=obs_vars,
          col=col_obs, pch=pch_obs, lty=1:length(pch_obs))
      }

      assign("cost.old", mC$model, inherits=TRUE)
    }
    return(mC)
  }

  # Define the model cost function for the t-test, without parameter transformation
  cost_notrans <- function(P)
  {
    if(length(state.ini.optim) > 0) {
      odeini <- c(P[1:length(state.ini.optim)], state.ini.fixed)
      names(odeini) <- c(state.ini.optim.boxnames, state.ini.fixed.boxnames)
    } else {
      odeini <- state.ini.fixed
      names(odeini) <- state.ini.fixed.boxnames
    }

    odeparms <- c(P[(length(state.ini.optim) + 1):length(P)], parms.fixed)

    # Solve the system with current parameter values
    out <- mkinpredict(mkinmod, odeparms,
                       odeini, outtimes,
                       solution_type = solution_type,
                       use_compiled = use_compiled,
                       method.ode = method.ode,
                       atol = atol, rtol = rtol, ...)

    mC <- modCost(out, observed, y = "value",
      err = err, weight = weight, scaleVar = scaleVar)

    return(mC)
  }

  # Define lower and upper bounds other than -Inf and Inf for parameters
  # for which no internal transformation is requested in the call to mkinfit.
  lower <- rep(-Inf, length(c(state.ini.optim, transparms.optim)))
  upper <- rep(Inf, length(c(state.ini.optim, transparms.optim)))
  names(lower) <- names(upper) <- c(names(state.ini.optim), names(transparms.optim))

  # IORE exponentes are not transformed, but need a lower bound of zero
  index_N <- grep("^N", names(lower))
  lower[index_N] <- 0
  if (!transform_rates) {
    index_k <- grep("^k_", names(lower))
    lower[index_k] <- 0
    index_k__iore <- grep("^k__iore_", names(lower))
    lower[index_k__iore] <- 0
    other_rate_parms <- intersect(c("alpha", "beta", "k1", "k2", "tb", "r"), names(lower))
    lower[other_rate_parms] <- 0
  }

  if (!transform_fractions) {
    index_f <- grep("^f_", names(upper))
    lower[index_f] <- 0
    upper[index_f] <- 1
    other_fraction_parms <- intersect(c("g"), names(upper))
    lower[other_fraction_parms] <- 0
    upper[other_fraction_parms] <- 1
  }

  # Show parameter names if tracing is requested
  if(trace_parms) cat(names(c(state.ini.optim, transparms.optim)), "\n")

  # Do the fit and take the time
  fit_time <- system.time({
    fit <- modFit(cost, c(state.ini.optim, transparms.optim),
                  method = method.modFit, control = control.modFit,
                  lower = lower, upper = upper, ...)

    # Reiterate the fit until convergence of the variance components (IRLS)
    # if requested by the user

    if (!is.null(reweight.method)) {
      if (! reweight.method %in% c("obs", "tc")) stop("Only reweighting methods 'obs' and 'tc' are implemented")

      if (reweight.method  == "obs") {
        tc_fit <- NA
        if(!quiet) {
          cat("IRLS based on variance estimates for each observed variable\n")
          cat("Initial variance estimates are:\n")
          print(signif(fit$var_ms_unweighted, 8))
        }
      }
      if (reweight.method  == "tc") {
        tc_fit <- .fit_error_model_mad_obs(cost(fit$par)$residuals, tc, 0)

        if (is.character(tc_fit)) {
          if (!quiet) {
            cat(tc_fit, ".\n", "No reweighting will be performed.")
          }
          tc_fitted <- c(sigma_low = NA, rsd_high = NA)
        } else {
          tc_fitted <- coef(tc_fit)
          if(!quiet) {
            cat("IRLS based on variance estimates according to the two component error model\n")
            cat("Initial variance components are:\n")
            print(signif(tc_fitted))
          }
        }
      }
      reweight.diff = 1
      n.iter <- 0
      if (!is.null(err)) observed$err.ini <- observed[[err]]
      err = "err.irls"

      while (reweight.diff > reweight.tol &
             n.iter < reweight.max.iter &
             !is.character(tc_fit)) {
        n.iter <- n.iter + 1
        # Store squared residual predictors used for weighting in sr_old and define new weights
        if (reweight.method == "obs") {
          sr_old <- fit$var_ms_unweighted
          observed[err] <- sqrt(fit$var_ms_unweighted[as.character(observed$name)])
        }
        if (reweight.method == "tc") {
          sr_old <- tc_fitted

          tmp_predicted <- mkin_wide_to_long(out_predicted, time = "time")
          tmp_data <- suppressMessages(join(observed, tmp_predicted, by = c("time", "name")))

          #observed[err] <- predict(tc_fit, newdata = data.frame(mod = tmp_data[[4]]))
          observed[err] <- predict(tc_fit, newdata = data.frame(obs = observed$value))

        }
        fit <- modFit(cost, fit$par, method = method.modFit,
                      control = control.modFit, lower = lower, upper = upper, ...)

        if (reweight.method == "obs") {
          sr_new <- fit$var_ms_unweighted
        }
        if (reweight.method == "tc") {
          tc_fit <- .fit_error_model_mad_obs(cost(fit$par)$residuals, tc_fitted, n.iter)

          if (is.character(tc_fit)) {
            if (!quiet) {
              cat(tc_fit, ".\n")
            }
            break
          } else {
            tc_fitted <- coef(tc_fit)
            sr_new <- tc_fitted
          }
        }

        reweight.diff = sum((sr_new - sr_old)^2)
        if (!quiet) {
          cat("Iteration", n.iter, "yields variance estimates:\n")
          print(signif(sr_new, 8))
          cat("Sum of squared differences to last variance (component) estimates:",
              signif(reweight.diff, 2), "\n")
        }
      }
    }
  })

  # Check for convergence
  if (method.modFit == "Marq") {
    if (!fit$info %in% c(1, 2, 3)) {
      fit$warning = paste0("Optimisation by method ", method.modFit,
                           " did not converge.\n",
                           "The message returned by nls.lm is:\n",
                                    fit$message)
      warning(fit$warning)
    }
    else {
      if(!quiet) cat("Optimisation by method", method.modFit, "successfully terminated.\n")
    }
  }
  if (method.modFit == "Port") {
    if (fit$convergence != 0) {
      fit$warning = paste0("Optimisation by method ", method.modFit,
                           " did not converge:\n",
                           if(is.character(fit$counts)) fit$counts # FME bug
                           else fit$message)
      warning(fit$warning)
    } else {
      if(!quiet) cat("Optimisation by method", method.modFit, "successfully terminated.\n")
    }
  }
  if (method.modFit %in% c("SANN", "Nelder-Mead", "BFGS", "CG", "L-BFGS-B")) {
    if (fit$convergence != 0) {
      fit$warning = paste0("Optimisation by method ", method.modFit,
                           " did not converge.\n",
                           "Convergence code returned by optim() is", fit$convergence)
      warning(fit$warning)
    } else {
      if(!quiet) cat("Optimisation by method", method.modFit, "successfully terminated.\n")
    }
  }

  # Return number of iterations for SANN method and alert user to check if
  # the approximation to the optimum is sufficient
  if (method.modFit == "SANN") {
    fit$iter = maxit.modFit
    fit$warning <- paste0("Termination of the SANN algorithm does not imply convergence.\n",
      "Make sure the approximation of the optimum is adequate.")
    warning(fit$warning)
  }

  # We need to return some more data for summary and plotting
  fit$solution_type <- solution_type
  fit$transform_rates <- transform_rates
  fit$transform_fractions <- transform_fractions
  fit$method.modFit <- method.modFit
  fit$maxit.modFit <- maxit.modFit
  fit$calls <- calls
  fit$time <- fit_time

  # We also need the model for summary and plotting
  fit$mkinmod <- mkinmod

  # We need data and predictions for summary and plotting
  fit$observed <- observed
  fit$obs_vars <- obs_vars
  fit$predicted <- mkin_wide_to_long(out_predicted, time = "time")

  # Backtransform parameters
  bparms.optim = backtransform_odeparms(fit$par, fit$mkinmod,
                                        transform_rates = transform_rates,
                                        transform_fractions = transform_fractions)
  bparms.fixed = c(state.ini.fixed, parms.fixed)
  bparms.all = c(bparms.optim, parms.fixed)

  # Attach the cost functions to the fit for post-hoc parameter uncertainty analysis
  fit$cost <- cost
  fit$cost_notrans <- cost_notrans

  # Estimate the Hessian for the model cost without parameter transformations
  # to make it possible to obtain the usual t-test
  # Code ported from FME::modFit
  Jac_notrans <- try(gradient(function(p, ...) cost_notrans(p)$residuals$res,
                              bparms.optim, centered = TRUE), silent = TRUE)
  if (inherits(Jac_notrans, "try-error")) {
    warning("Calculation of the Jacobian failed for the cost function of the untransformed model.\n",
            "No t-test results will be available")
    fit$hessian_notrans <- NA
  } else {
    fit$hessian_notrans <- 2 * t(Jac_notrans) %*% Jac_notrans
  }

  # Collect initial parameter values in three dataframes
  fit$start <- data.frame(value = c(state.ini.optim,
                                    parms.optim))
  fit$start$type = c(rep("state", length(state.ini.optim)),
                     rep("deparm", length(parms.optim)))

  fit$start_transformed = data.frame(
      value = c(state.ini.optim, transparms.optim),
      lower = lower,
      upper = upper)

  fit$fixed <- data.frame(value = c(state.ini.fixed, parms.fixed))
  fit$fixed$type = c(rep("state", length(state.ini.fixed)),
                     rep("deparm", length(parms.fixed)))

  # Collect observed, predicted, residuals and weighting
  data <- merge(fit$observed, fit$predicted, by = c("time", "name"))
  data$name <- ordered(data$name, levels = obs_vars)
  data <- data[order(data$name, data$time), ]

  # Add a column named "value" holding the observed values for the case
  # that this column was selected for manual weighting, so it can be
  # shown in the summary as "err"
  data$value <- data$value.x

  fit$data <- data.frame(time = data$time,
                         variable = data$name,
                         observed = data$value.x,
                         predicted = data$value.y)

  fit$data$residual <- fit$data$observed - fit$data$predicted
  if (!is.null(data$err.ini)) fit$data$err.ini <- data$err.ini
  if (!is.null(err)) fit$data[["err"]] <- data[[err]]

  fit$atol <- atol
  fit$rtol <- rtol
  fit$weight.ini <- weight.ini
  fit$tc.ini <- tc
  fit$reweight.method <- reweight.method
  fit$reweight.tol <- reweight.tol
  fit$reweight.max.iter <- reweight.max.iter
  if (exists("tc_fitted")) fit$tc_fitted <- tc_fitted

  # Return different sets of backtransformed parameters for summary and plotting
  fit$bparms.optim <- bparms.optim
  fit$bparms.fixed <- bparms.fixed

  # Return ode and state parameters for further fitting
  fit$bparms.ode <- bparms.all[mkinmod$parms]
  fit$bparms.state <- c(bparms.all[setdiff(names(bparms.all), names(fit$bparms.ode))],
                        state.ini.fixed)
  names(fit$bparms.state) <- gsub("_0$", "", names(fit$bparms.state))

  fit$date <- date()
  fit$version <- as.character(utils::packageVersion("mkin"))
  fit$Rversion <- paste(R.version$major, R.version$minor, sep=".")

  class(fit) <- c("mkinfit", "modFit")
  return(fit)
}

summary.mkinfit <- function(object, data = TRUE, distimes = TRUE, alpha = 0.05, ...) {
  param  <- object$par
  pnames <- names(param)
  bpnames <- names(object$bparms.optim)
  p      <- length(param)
  mod_vars <- names(object$mkinmod$diffs)
  covar  <- try(solve(0.5*object$hessian), silent = TRUE)   # unscaled covariance
  covar_notrans  <- try(solve(0.5*object$hessian_notrans), silent = TRUE)   # unscaled covariance
  rdf    <- object$df.residual
  resvar <- object$ssr / rdf
  if (!is.numeric(covar) | is.na(covar[1])) {
    covar <- NULL
    se <- lci <- uci <- rep(NA, p)
  } else {
    rownames(covar) <- colnames(covar) <- pnames
    se     <- sqrt(diag(covar) * resvar)
    lci    <- param + qt(alpha/2, rdf) * se
    uci    <- param + qt(1-alpha/2, rdf) * se
  }

  if (!is.numeric(covar_notrans) | is.na(covar_notrans[1])) {
    covar_notrans <- NULL
    se_notrans <- tval <- pval <- rep(NA, p)
  } else {
    rownames(covar_notrans) <- colnames(covar_notrans) <- bpnames
    se_notrans     <- sqrt(diag(covar_notrans) * resvar)
    tval   <- object$bparms.optim/se_notrans
    pval   <- pt(abs(tval), rdf, lower.tail = FALSE)
  }

  names(se) <- pnames
  modVariance <- object$ssr / length(object$residuals)

  param <- cbind(param, se, lci, uci)
  dimnames(param) <- list(pnames, c("Estimate", "Std. Error", "Lower", "Upper"))

  bparam <- cbind(Estimate = object$bparms.optim, se_notrans, "t value" = tval, "Pr(>t)" = pval, Lower = NA, Upper = NA)

  # Transform boundaries of CI for one parameter at a time,
  # with the exception of sets of formation fractions (single fractions are OK).
  f_names_skip <- character(0)
  for (box in mod_vars) { # Figure out sets of fractions to skip
    f_names <- grep(paste("^f", box, sep = "_"), pnames, value = TRUE)
    n_paths <- length(f_names)
    if (n_paths > 1) f_names_skip <- c(f_names_skip, f_names)
  }

  for (pname in pnames) {
    if (!pname %in% f_names_skip) {
      par.lower <- param[pname, "Lower"]
      par.upper <- param[pname, "Upper"]
      names(par.lower) <- names(par.upper) <- pname
      bpl <- backtransform_odeparms(par.lower, object$mkinmod,
                                            object$transform_rates,
                                            object$transform_fractions)
      bpu <- backtransform_odeparms(par.upper, object$mkinmod,
                                            object$transform_rates,
                                            object$transform_fractions)
      bparam[names(bpl), "Lower"] <- bpl
      bparam[names(bpu), "Upper"] <- bpu
    }
  }

  ans <- list(
    version = as.character(utils::packageVersion("mkin")),
    Rversion = paste(R.version$major, R.version$minor, sep="."),
    date.fit = object$date,
    date.summary = date(),
    solution_type = object$solution_type,
    method.modFit = object$method.modFit,
    warning = object$warning,
    use_of_ff = object$mkinmod$use_of_ff,
    weight.ini = object$weight.ini,
    reweight.method = object$reweight.method,
    tc.ini = object$tc.ini,
    var_ms_unweighted = object$var_ms_unweighted,
    tc_fitted = object$tc_fitted,
    residuals = object$residuals,
    residualVariance = resvar,
    sigma = sqrt(resvar),
    modVariance = modVariance,
    df = c(p, rdf),
    cov.unscaled = covar,
    cov.scaled = covar * resvar,
    info = object$info,
    niter = object$iterations,
    calls = object$calls,
    time = object$time,
    stopmess = message,
    par = param,
    bpar = bparam)

  if (!is.null(object$version)) {
    ans$fit_version <- object$version
    ans$fit_Rversion <- object$Rversion
  }

  ans$diffs <- object$mkinmod$diffs
  if(data) ans$data <- object$data
  ans$start <- object$start
  ans$start_transformed <- object$start_transformed

  ans$fixed <- object$fixed

  ans$errmin <- mkinerrmin(object, alpha = 0.05)

  if (object$calls > 0) {
    if (!is.null(ans$cov.unscaled)){
      Corr <- cov2cor(ans$cov.unscaled)
      rownames(Corr) <- colnames(Corr) <- rownames(ans$par)
      ans$Corr <- Corr
    } else {
      warning("Could not estimate covariance matrix; singular system.")
    }
  }

  ans$bparms.ode <- object$bparms.ode
  ep <- endpoints(object)
  if (length(ep$ff) != 0)
    ans$ff <- ep$ff
  if(distimes) ans$distimes <- ep$distimes
  if(length(ep$SFORB) != 0) ans$SFORB <- ep$SFORB
  class(ans) <- c("summary.mkinfit", "summary.modFit")
  return(ans)
}

# Expanded from print.summary.modFit
print.summary.mkinfit <- function(x, digits = max(3, getOption("digits") - 3), ...) {
  if (is.null(x$fit_version)) {
    cat("mkin version:   ", x$version, "\n")
    cat("R version:      ", x$Rversion, "\n")
  } else {
    cat("mkin version used for fitting:   ", x$fit_version, "\n")
    cat("R version used for fitting:      ", x$fit_Rversion, "\n")
  }

  cat("Date of fit:    ", x$date.fit, "\n")
  cat("Date of summary:", x$date.summary, "\n")

  if (!is.null(x$warning)) cat("\n\nWarning:", x$warning, "\n\n")

  cat("\nEquations:\n")
  nice_diffs <- gsub("^(d.*) =", "\\1/dt =", x[["diffs"]])
  writeLines(strwrap(nice_diffs, exdent = 11))
  df  <- x$df
  rdf <- df[2]

  cat("\nModel predictions using solution type", x$solution_type, "\n")

  cat("\nFitted with method", x$method.modFit,
      "using", x$calls, "model solutions performed in", x$time[["elapsed"]],  "s\n")

  cat("\nWeighting:", x$weight.ini)
  if (x$weight.ini == "tc")  {
    cat(" with variance components\n")
    print(x$tc.ini)
  } else {
    cat ("\n")
  }
  if(!is.null(x$reweight.method)) {
    cat("\nIterative reweighting with method", x$reweight.method, "\n")
    if (x$reweight.method == "obs") {
      cat("Final mean squared residuals of observed variables:\n")
      print(x$var_ms_unweighted)
    }
    if (x$reweight.method == "tc") {
      cat("Final components of fitted standard deviation:\n")
      print(x$tc_fitted)
    }
  }

  cat("\nStarting values for parameters to be optimised:\n")
  print(x$start)

  cat("\nStarting values for the transformed parameters actually optimised:\n")
  print(x$start_transformed)

  cat("\nFixed parameter values:\n")
  if(length(x$fixed$value) == 0) cat("None\n")
  else print(x$fixed)

  cat("\nOptimised, transformed parameters with symmetric confidence intervals:\n")
  print(signif(x$par, digits = digits))

  if (x$calls > 0) {
    cat("\nParameter correlation:\n")
    if (!is.null(x$cov.unscaled)){
      print(x$Corr, digits = digits, ...)
    } else {
      cat("Could not estimate covariance matrix; singular system.")
    }
  }

  cat("\nResidual standard error:",
      format(signif(x$sigma, digits)), "on", rdf, "degrees of freedom\n")

  cat("\nBacktransformed parameters:\n")
  cat("Confidence intervals for internally transformed parameters are asymmetric.\n")
  if ((x$version) < "0.9-36") {
    cat("To get the usual (questionable) t-test, upgrade mkin and repeat the fit.\n")
    print(signif(x$bpar, digits = digits))
  } else {
    cat("t-test (unrealistically) based on the assumption of normal distribution\n")
    cat("for estimators of untransformed parameters.\n")
    print(signif(x$bpar[, c(1, 3, 4, 5, 6)], digits = digits))
  }

  cat("\nChi2 error levels in percent:\n")
  x$errmin$err.min <- 100 * x$errmin$err.min
  print(x$errmin, digits=digits,...)

  printSFORB <- !is.null(x$SFORB)
  if(printSFORB){
    cat("\nEstimated Eigenvalues of SFORB model(s):\n")
    print(x$SFORB, digits=digits,...)
  }

  printff <- !is.null(x$ff)
  if(printff){
    cat("\nResulting formation fractions:\n")
    print(data.frame(ff = x$ff), digits=digits,...)
  }

  printdistimes <- !is.null(x$distimes)
  if(printdistimes){
    cat("\nEstimated disappearance times:\n")
    print(x$distimes, digits=digits,...)
  }

  printdata <- !is.null(x$data)
  if (printdata){
    cat("\nData:\n")
    print(format(x$data, digits = digits, ...), row.names = FALSE)
  }

  invisible(x)
}

# Fit the median absolute deviation against the observed values,
# using the current error model for weighting
.fit_error_model_mad_obs <- function(tmp_res, tc, iteration) {
  mad_agg <- aggregate(tmp_res$res.unweighted,
                       by = list(name = tmp_res$name, time = tmp_res$x),
                       FUN = function(x) mad(x, center = 0))
  names(mad_agg) <- c("name", "time", "mad")
  error_data <- suppressMessages(
    join(data.frame(name = tmp_res$name,
                    time = tmp_res$x,
                    obs = tmp_res$obs),
         mad_agg))
  error_data_complete <- na.omit(error_data)

  tc_fit <- tryCatch(
    nls(mad ~ sigma_twocomp(obs, sigma_low, rsd_high),
      start = list(sigma_low = tc["sigma_low"], rsd_high = tc["rsd_high"]),
      weights = 1/sigma_twocomp(error_data_complete$obs,
                                tc["sigma_low"],
                                tc["rsd_high"])^2,
      data = error_data_complete,
      lower = 0,
      algorithm = "port"),
    error = function(e) paste("Fitting the error model failed in iteration", iteration))
  return(tc_fit)
}
# Alternative way to fit the error model, fitting to modelled instead of
# observed values
# .fit_error_model_mad_mod <- function(tmp_res, tc) {
#   mad_agg <- aggregate(tmp_res$res.unweighted,
#                        by = list(name = tmp_res$name, time = tmp_res$x),
#                        FUN = function(x) mad(x, center = 0))
#   names(mad_agg) <- c("name", "time", "mad")
#   mod_agg <- aggregate(tmp_res$mod,
#                        by = list(name = tmp_res$name, time = tmp_res$x),
#                        FUN = mean)
#   names(mod_agg) <- c("name", "time", "mod")
#   mod_mad <- merge(mod_agg, mad_agg)
# 
#   tc_fit <- tryCatch(
#     nls(mad ~ sigma_twocomp(mod, sigma_low, rsd_high),
#       start = list(sigma_low = tc["sigma_low"], rsd_high = tc["rsd_high"]),
#       data = mod_mad,
#       weights = 1/mod_mad$mad,
#       lower = 0,
#       algorithm = "port"),
#     error = "Fitting the error model failed in iteration")
#   return(tc_fit)
# }
# vim: set ts=2 sw=2 expandtab:

Contact - Imprint