From a26b44d15c11ebb41083fc2efab0cc91a027b55b Mon Sep 17 00:00:00 2001 From: Johannes Ranke Date: Sat, 18 Apr 2015 00:06:38 +0200 Subject: Add a vignette showing the performance of compiled models --- .gitignore | 2 + DESCRIPTION | 2 +- GNUmakefile | 5 +- NAMESPACE | 5 +- R/mkinmod.R | 9 ++- man/mkinfit.Rd | 16 +++-- man/mkinpredict.Rd | 2 +- vignettes/compiled_models.Rmd | 87 +++++++++++++++++++++++ vignettes/compiled_models.html | 153 +++++++++++++++++++++++++++++++++++++++++ vignettes/mkin_vignettes.css | 16 +++++ 10 files changed, 282 insertions(+), 15 deletions(-) create mode 100644 vignettes/compiled_models.Rmd create mode 100644 vignettes/compiled_models.html create mode 100644 vignettes/mkin_vignettes.css diff --git a/.gitignore b/.gitignore index 39eac04e..956fa374 100644 --- a/.gitignore +++ b/.gitignore @@ -6,7 +6,9 @@ vignettes/*.blg vignettes/*.log vignettes/*.out vignettes/*.toc +vignettes/*.R vignettes/.build.timestamp vignettes/mkin.tex vignettes/FOCUS_Z.tex +vignettes/*cache/ mkin.Rcheck diff --git a/DESCRIPTION b/DESCRIPTION index 642228c0..e823ba5f 100644 --- a/DESCRIPTION +++ b/DESCRIPTION @@ -19,7 +19,7 @@ Description: Calculation routines based on the FOCUS Kinetics Report (2006). for a particular purpose. Depends: minpack.lm, rootSolve Imports: FME, deSolve -Suggests: knitr, RUnit, ccSolve +Suggests: knitr, RUnit, ccSolve, microbenchmark License: GPL LazyLoad: yes LazyData: yes diff --git a/GNUmakefile b/GNUmakefile index 126dede1..a84c9e8e 100644 --- a/GNUmakefile +++ b/GNUmakefile @@ -80,7 +80,10 @@ vignettes/FOCUS_L.html: vignettes/FOCUS_L.Rmd vignettes/FOCUS_Z.pdf: vignettes/FOCUS_Z.Rnw "$(RBIN)/Rscript" -e "tools::buildVignette(file = 'vignettes/FOCUS_Z.Rnw', dir = 'vignettes')" -vignettes: vignettes/mkin.pdf vignettes/FOCUS_L.html vignettes/FOCUS_Z.pdf +vignettes/compiled_models.html: vignettes/compiled_models.Rmd + "$(RBIN)/Rscript" -e "tools::buildVignette(file = 'vignettes/compiled_models.Rmd', dir = 'vignettes')" + +vignettes: vignettes/mkin.pdf vignettes/FOCUS_L.html vignettes/FOCUS_Z.pdf vignettes/compiled_models.html sd: "$(RBIN)/Rscript" -e "library(staticdocs); build_site()" diff --git a/NAMESPACE b/NAMESPACE index 92f0f271..49be428e 100644 --- a/NAMESPACE +++ b/NAMESPACE @@ -4,11 +4,10 @@ S3method(plot, mkinfit) S3method(summary, mkinfit) S3method(print, summary.mkinfit) -# Import all packages listed as Imports or Depends +# Import packages listed as Imports or Depends import( FME, deSolve, minpack.lm, - rootSolve, - ccSolve + rootSolve ) diff --git a/R/mkinmod.R b/R/mkinmod.R index 3fff0878..149d2beb 100644 --- a/R/mkinmod.R +++ b/R/mkinmod.R @@ -273,7 +273,7 @@ mkinmod <- function(..., use_of_ff = "min", speclist = NULL) }#}}} # Create a function compiled from C code if possible #{{{ - if (require(ccSolve)) { + if (requireNamespace("ccSolve", quietly = TRUE)) { diffs.C <- paste(diffs, collapse = ";\n") diffs.C <- paste0(diffs.C, ";") for (i in seq_along(diffs)) { @@ -295,9 +295,12 @@ mkinmod <- function(..., use_of_ff = "min", speclist = NULL) replacement <- paste0("* y[", i - 1, "]") diffs.C <- gsub(pattern, replacement, diffs.C) } - if (sum(sapply(spec, function(x) x$type %in% c("SFO", "FOMC", "DFOP", "SFORB"))) == length(spec)) { + if (sum(sapply(spec, function(x) x$type %in% + c("SFO", "FOMC", "DFOP", "SFORB"))) == length(spec)) { message("Compiling differential equation model from auto-generated C code...") - model$compiled <- compile.ode(diffs.C, language = "C", parms = parms, declaration = "double time = *t;") + model$compiled <- ccSolve::compile.ode(diffs.C, language = "C", + parms = parms, + declaration = "double time = *t;") } } # }}} diff --git a/man/mkinfit.Rd b/man/mkinfit.Rd index 5f43e145..801f540c 100644 --- a/man/mkinfit.Rd +++ b/man/mkinfit.Rd @@ -21,7 +21,7 @@ mkinfit(mkinmod, observed, parms.ini = "auto", state.ini = "auto", fixed_parms = NULL, fixed_initials = names(mkinmod$diffs)[-1], - solution_type = c("auto", "analytical", "eigen", "deSolve") + 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"), @@ -241,11 +241,15 @@ SFO_SFO <- mkinmod( parent = list(type = "SFO", to = "m1", sink = TRUE), m1 = list(type = "SFO")) # Fit the model to the FOCUS example dataset D using defaults -system.time(fit <- mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "eigen")) -system.time(fit.deSolve <- mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "deSolve", use_compiled = FALSE)) -system.time(fit.deSolve.compiled <- mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "deSolve", use_compiled = TRUE)) -# The performance benefit of using the compiled version is immense compared with deSolve without compilation, and -# about a factor of two compared with the eigenvalue based solution +system.time(fit <- mkinfit(SFO_SFO, FOCUS_2006_D, + solution_type = "eigen")) +system.time(fit.deSolve <- mkinfit(SFO_SFO, FOCUS_2006_D, + solution_type = "deSolve", use_compiled = FALSE)) +system.time(fit.deSolve.compiled <- mkinfit(SFO_SFO, FOCUS_2006_D, + solution_type = "deSolve", use_compiled = TRUE)) +# The performance benefit of using the compiled version is immense compared +# with deSolve without compilation, and about a factor of two compared with the +# eigenvalue based solution coef(fit) coef(fit.deSolve) coef(fit.deSolve.compiled) diff --git a/man/mkinpredict.Rd b/man/mkinpredict.Rd index 7450a7ce..c6aee75f 100644 --- a/man/mkinpredict.Rd +++ b/man/mkinpredict.Rd @@ -10,7 +10,7 @@ } \usage{ mkinpredict(mkinmod, odeparms, odeini, outtimes, solution_type = "deSolve", - method.ode = "lsoda", use_compiled = "auto", atol = 1e-08, rtol = 1e-10, + use_compiled = "auto", method.ode = "lsoda", atol = 1e-08, rtol = 1e-10, map_output = TRUE, ...) } \arguments{ diff --git a/vignettes/compiled_models.Rmd b/vignettes/compiled_models.Rmd new file mode 100644 index 00000000..bac284c5 --- /dev/null +++ b/vignettes/compiled_models.Rmd @@ -0,0 +1,87 @@ +--- +title: "Performance benefit by using compiled model definitions in mkin" +output: + html_document: + css: mkin_vignettes.css + toc: true + mathjax: null + theme: united +--- + + +```{r, include = FALSE} +library(knitr) +opts_chunk$set(tidy = FALSE, cache = TRUE) +if (!require("ccSolve")) + message("Please install the ccSolve package for this vignette to produce sensible output") + +``` + +# Benchmark for a model that can also be solved with Eigenvalues + +This evaluation is taken from the example section of mkinfit. When using an mkin version +greater than 0.9-36 and the ccSolve package is installed and functional, you will get a +message that the model is being compiled when defining a model using mkinmod. + +```{r create_SFO_SFO} +library("mkin") +SFO_SFO <- mkinmod( + parent = list(type = "SFO", to = "m1", sink = TRUE), + m1 = list(type = "SFO")) +``` + +We can compare the performance of the Eigenvalue based solution against the +compiled version and the R implementation of the differential equations using +the microbenchmark package. + + +```{r benchmark_SFO_SFO, echo=-(1:2)} +# Redefining the model, in order not to confuse the knitr cache which leads to segfaults +suppressMessages(SFO_SFO <- mkinmod( + parent = list(type = "SFO", to = "m1", sink = TRUE), + m1 = list(type = "SFO"))) +library("microbenchmark") +mb.1 <- microbenchmark( + mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "deSolve", use_compiled = FALSE, + quiet = TRUE), + mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "eigen", quiet = TRUE), + mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "deSolve", quiet = TRUE), + times = 3, control = list(warmup = 1)) +smb.1 <- summary(mb.1)[-1] +rownames(smb.1) <- c("deSolve, not compiled", "Eigenvalue based", "deSolve, compiled") +print(smb.1) +``` + +We see that using the compiled model is almost a factor of 8 faster than using the R version +with the default ode solver, and it is even faster than the Eigenvalue based solution implemented +in R which does not need iterative solution of the ODEs: + +```{r} +smb.1["median"]/smb.1["deSolve, compiled", "median"] +``` + +# Benchmark for a model that can not be solved with Eigenvalues + +This evaluation is also taken from the example section of mkinfit. + +```{r benchmark_FOMC_SFO} +FOMC_SFO <- mkinmod( + parent = list(type = "FOMC", to = "m1", sink = TRUE), + m1 = list(type = "SFO")) + +mb.2 <- microbenchmark( + mkinfit(FOMC_SFO, FOCUS_2006_D, use_compiled = FALSE, quiet = TRUE), + mkinfit(FOMC_SFO, FOCUS_2006_D, quiet = TRUE), + times = 3, control = list(warmup = 1)) +smb.2 <- summary(mb.2)[-1] +rownames(smb.2) <- c("deSolve, not compiled", "deSolve, compiled") +print(smb.2) +smb.2["median"]/smb.2["deSolve, compiled", "median"] + +``` + +Here we get a performance benefit of more than a factor of 8 using the version +of the differential equation model compiled from C code using the ccSolve package! diff --git a/vignettes/compiled_models.html b/vignettes/compiled_models.html new file mode 100644 index 00000000..2f2a6edb --- /dev/null +++ b/vignettes/compiled_models.html @@ -0,0 +1,153 @@ + + + + + + + + + + + + +Performance benefit by using compiled model definitions in mkin + + + + + + + + + + + + + + + + + + + + + + +
+ + + + + + + +
+

Benchmark for a model that can also be solved with Eigenvalues

+

This evaluation is taken from the example section of mkinfit. When using an mkin version greater than 0.9-36 and the ccSolve package is installed and functional, you will get a message that the model is being compiled when defining a model using mkinmod.

+
library("mkin")
+SFO_SFO <- mkinmod(
+  parent = list(type = "SFO", to = "m1", sink = TRUE),
+  m1 = list(type = "SFO"))
+
## Compiling differential equation model from auto-generated C code...
+

We can compare the performance of the Eigenvalue based solution against the compiled version and the R implementation of the differential equations using the microbenchmark package.

+
library("microbenchmark")
+mb.1 <- microbenchmark(
+  mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "deSolve", use_compiled = FALSE, 
+          quiet = TRUE),
+  mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "eigen", quiet = TRUE),
+  mkinfit(SFO_SFO, FOCUS_2006_D, solution_type = "deSolve", quiet = TRUE),
+  times = 3, control = list(warmup = 1))
+smb.1 <- summary(mb.1)[-1]
+rownames(smb.1) <- c("deSolve, not compiled", "Eigenvalue based", "deSolve, compiled")
+print(smb.1)
+
##                             min        lq      mean    median        uq
+## deSolve, not compiled 6192.0125 6195.3470 6211.0309 6198.6816 6220.5401
+## Eigenvalue based       956.7604 1008.7224 1026.2572 1060.6844 1061.0055
+## deSolve, compiled      869.6880  871.9315  883.4929  874.1751  890.3953
+##                             max neval
+## deSolve, not compiled 6242.3986     3
+## Eigenvalue based      1061.3266     3
+## deSolve, compiled      906.6155     3
+

We see that using the compiled model is almost a factor of 8 faster than using the R version with the default ode solver, and it is even faster than the Eigenvalue based solution implemented in R which does not need iterative solution of the ODEs:

+
smb.1["median"]/smb.1["deSolve, compiled", "median"]
+
##                         median
+## deSolve, not compiled 7.120877
+## Eigenvalue based      1.205328
+## deSolve, compiled     1.000000
+
+
+

Benchmark for a model that can not be solved with Eigenvalues

+

This evaluation is also taken from the example section of mkinfit.

+
FOMC_SFO <- mkinmod(
+  parent = list(type = "FOMC", to = "m1", sink = TRUE),
+  m1 = list(type = "SFO"))
+
## Compiling differential equation model from auto-generated C code...
+
mb.2 <- microbenchmark(
+  mkinfit(FOMC_SFO, FOCUS_2006_D, use_compiled = FALSE, quiet = TRUE),
+  mkinfit(FOMC_SFO, FOCUS_2006_D, quiet = TRUE),
+  times = 3, control = list(warmup = 1))
+smb.2 <- summary(mb.2)[-1]
+rownames(smb.2) <- c("deSolve, not compiled", "deSolve, compiled")
+print(smb.2)
+
##                             min        lq      mean    median        uq
+## deSolve, not compiled 13.297283 13.427702 13.481155 13.558121 13.573092
+## deSolve, compiled      1.486926  1.526887  1.546851  1.566848  1.576813
+##                             max neval
+## deSolve, not compiled 13.588063     3
+## deSolve, compiled      1.586778     3
+
smb.2["median"]/smb.2["deSolve, compiled", "median"]
+
##                         median
+## deSolve, not compiled 8.653119
+## deSolve, compiled     1.000000
+

Here we get a performance benefit of more than a factor of 8 using the version of the differential equation model compiled from C code using the ccSolve package!

+
+ + +
+ + + + + + diff --git a/vignettes/mkin_vignettes.css b/vignettes/mkin_vignettes.css new file mode 100644 index 00000000..3bd91ab1 --- /dev/null +++ b/vignettes/mkin_vignettes.css @@ -0,0 +1,16 @@ +/* Thanks to Steve Powell for http://rpubs.com/stevepowell99/floating-css */ +#TOC { + position: fixed; + left: 0; + top: 0; + width: 200px; + height: 100%; + overflow:auto; +} + +body { + max-width: 800px; + margin: auto; + margin-left:210px; + line-height: 20px; +} -- cgit v1.2.1