aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--ChangeLog5
-rw-r--r--DESCRIPTION4
-rw-r--r--GNUmakefile27
-rw-r--r--R/mkinfit.R8
-rw-r--r--R/mkinmod.R15
-rw-r--r--R/mkinpredict.R10
-rw-r--r--man/gmkin.Rd3
-rw-r--r--man/mkinfit.Rd6
-rw-r--r--man/mkinpredict.Rd17
9 files changed, 61 insertions, 34 deletions
diff --git a/ChangeLog b/ChangeLog
index 74c161f..eba6d92 100644
--- a/ChangeLog
+++ b/ChangeLog
@@ -1,7 +1,10 @@
-2013-11-17 Johannes Ranke <jranke@uni-bremen.de> for mkin (0.9-25)
+2013-12-04 Johannes Ranke <jranke@uni-bremen.de> for mkin (0.9-25)
* Change vignette format from Sweave to knitr
* Split examples vignette to FOCUS_L and FOCUS_Z
+ * Add a starter function for the GUI: gmkin()
+ * Remove warning about constant formation fractions in mkinmod
+ as it was based on a misconception
2013-11-06 Johannes Ranke <jranke@uni-bremen.de> for mkin (0.9-24)
diff --git a/DESCRIPTION b/DESCRIPTION
index 14abd50..7b478e4 100644
--- a/DESCRIPTION
+++ b/DESCRIPTION
@@ -3,7 +3,7 @@ Type: Package
Title: Routines for fitting kinetic models with one or more state
variables to chemical degradation data
Version: 0.9-25
-Date: 2014-01-17
+Date: 2014-01-22
Authors@R: c(person("Johannes", "Ranke", role = c("aut", "cre"),
email = "jranke@uni-bremen.de"),
person("Katrin", "Lindenberger", role = "ctb"),
@@ -21,4 +21,4 @@ LazyLoad: yes
LazyData: yes
Encoding: UTF-8
VignetteBuilder: knitr
-URL: http://cran.r-project.org, http://kinfit.r-forge.r-project.org, http://github.com/jranke/mkin
+URL: http://cran.r-project.org/package=mkin, http://kinfit.r-forge.r-project.org, http://github.com/jranke/mkin
diff --git a/GNUmakefile b/GNUmakefile
index 1c124b1..cf5c845 100644
--- a/GNUmakefile
+++ b/GNUmakefile
@@ -16,11 +16,16 @@ help:
@echo ""
@echo "Development Tasks"
@echo "-----------------"
- @echo " build Create the package"
- @echo " check Invoke build and then check the package"
- @echo " install Invoke build and then install the result"
- @echo " test Install a new copy of the package and run it "
- @echo " through the testsuite"
+ @echo " build Create the package"
+ @echo " build-no-vignettes Create the package without rebuilding vignettes"
+ @echo " check Invoke build and then check the package"
+ @echo " check-no-vignettes Invoke build without rebuilding vignettes, and then check"
+ @echo " install Invoke build and then install the result"
+ @echo " install Invoke build without rebuilding vignettes and then install the result"
+ @echo " test Install a new copy of the package and run it "
+ @echo " through the testsuite"
+ @echo " test-no-vignettes Invoke build without rebuilding vignettes, and then run it"
+ @echo " through the testsuite"
@echo ""
@echo "Packaging Tasks"
@echo "---------------"
@@ -47,14 +52,26 @@ install: build
cd ..;\
"$(RBIN)/R" CMD INSTALL $(PKGNAME)_$(PKGVERS).tar.gz
+install-no-vignettes: build-no-vignettes
+ cd ..;\
+ "$(RBIN)/R" CMD INSTALL $(PKGNAME)_$(PKGVERS).tar.gz
+
check: build
cd ..;\
"$(RBIN)/R" CMD check --as-cran --no-tests $(PKGNAME)_$(PKGVERS).tar.gz
+check-no-vignettes: build-no-vignettes
+ cd ..;\
+ "$(RBIN)/R" CMD check --as-cran --no-tests $(PKGNAME)_$(PKGVERS).tar.gz
+
test: install
cd tests;\
"$(RBIN)/Rscript" doRUnit.R
+test-no-vignettes: install-no-vignettes
+ cd tests;\
+ "$(RBIN)/Rscript" doRUnit.R
+
#------------------------------------------------------------------------------
# Packaging Tasks
#------------------------------------------------------------------------------
diff --git a/R/mkinfit.R b/R/mkinfit.R
index fd53627..6c976a3 100644
--- a/R/mkinfit.R
+++ b/R/mkinfit.R
@@ -1,5 +1,3 @@
-# $Id$
-
# Copyright (C) 2010-2014 Johannes Ranke
# Contact: jranke@uni-bremen.de
# The summary function is an adapted and extended version of summary.modFit
@@ -28,6 +26,7 @@ mkinfit <- function(mkinmod, observed,
fixed_parms = NULL,
fixed_initials = names(mkinmod$diffs)[-1],
solution_type = "auto",
+ method.ode = "lsoda",
method.modFit = "Marq",
control.modFit = list(),
plot = FALSE, quiet = FALSE,
@@ -161,6 +160,7 @@ mkinfit <- function(mkinmod, observed,
# Solve the system with current transformed parameter values
out <- mkinpredict(mkinmod, parms, odeini, outtimes,
solution_type = solution_type,
+ method.ode = method.ode,
atol = atol, rtol = rtol, ...)
assign("out_predicted", out, inherits=TRUE)
@@ -177,7 +177,9 @@ mkinfit <- function(mkinmod, observed,
outtimes_plot = seq(min(observed$time), max(observed$time), length.out=100)
out_plot <- mkinpredict(mkinmod, parms, odeini, outtimes_plot,
- solution_type = solution_type, atol = atol, rtol = rtol, ...)
+ solution_type = solution_type,
+ method.ode = method.ode,
+ atol = atol, rtol = rtol, ...)
plot(0, type="n",
xlim = range(observed$time), ylim = range(observed$value, na.rm=TRUE),
diff --git a/R/mkinmod.R b/R/mkinmod.R
index c743777..2bc6ae3 100644
--- a/R/mkinmod.R
+++ b/R/mkinmod.R
@@ -1,6 +1,4 @@
-# $Id$
-
-# Copyright (C) 2010-2012 Johannes Ranke {{{
+# Copyright (C) 2010-2013 Johannes Ranke {{{
# Contact: jranke@uni-bremen.de
# This file is part of the R package mkin
@@ -46,19 +44,10 @@ mkinmod <- function(..., use_of_ff = "min", speclist = NULL)
# fractions
if(spec[[1]]$type %in% c("FOMC", "DFOP", "HS")) {
mat = FALSE
- if(!is.null(spec[[1]]$to)) {
- message <- paste(
- "Only constant formation fractions over time are implemented.",
- "Depending on the reason for the time dependence of degradation",
- "this may be unrealistic. You may want to consider using the",
- "SFORB model",
- sep="\n")
- warning(message)
- } else message <- "ok"
} else mat = TRUE
#}}}
- # Establish list of differential equations as well as map from observed {{{
+ # Establish a list of differential equations as well as a map from observed {{{
# compartments to differential equations
for (varname in obs_vars)
{
diff --git a/R/mkinpredict.R b/R/mkinpredict.R
index 7ce26ca..07c1a81 100644
--- a/R/mkinpredict.R
+++ b/R/mkinpredict.R
@@ -1,6 +1,4 @@
-# $Id$
-
-# Copyright (C) 2010-2012 Johannes Ranke
+# Copyright (C) 2010-2013 Johannes Ranke
# Contact: jranke@uni-bremen.de
# This file is part of the R package mkin
@@ -18,7 +16,10 @@
# You should have received a copy of the GNU General Public License along with
# this program. If not, see <http://www.gnu.org/licenses/>
-mkinpredict <- function(mkinmod, odeparms, odeini, outtimes, solution_type = "deSolve", map_output = TRUE, atol = 1e-8, rtol = 1e-10, ...) {
+mkinpredict <- function(mkinmod, odeparms, odeini,
+ outtimes, solution_type = "deSolve",
+ method.ode = "lsoda", atol = 1e-8, rtol = 1e-10,
+ map_output = TRUE, ...) {
# Get the names of the state variables in the model
mod_vars <- names(mkinmod$diffs)
@@ -96,6 +97,7 @@ mkinpredict <- function(mkinmod, odeparms, odeini, outtimes, solution_type = "de
times = outtimes,
func = mkindiff,
parms = odeparms,
+ method = method.ode,
atol = atol,
rtol = rtol,
...
diff --git a/man/gmkin.Rd b/man/gmkin.Rd
index 29a4f8c..8a7d132 100644
--- a/man/gmkin.Rd
+++ b/man/gmkin.Rd
@@ -8,6 +8,9 @@
\href{http://github.com/jverzani/gWidgetsWWW2}{github page of gWidgetsWWW2}
for an explanation how this toolkit works.
}
+\note{
+ This GUI is experimental and not recommended for real work.
+}
\usage{
gmkin()
}
diff --git a/man/mkinfit.Rd b/man/mkinfit.Rd
index 51e26ed..0c5d76a 100644
--- a/man/mkinfit.Rd
+++ b/man/mkinfit.Rd
@@ -20,6 +20,7 @@ mkinfit(mkinmod, observed,
state.ini = c(100, rep(0, length(mkinmod$diffs) - 1)),
fixed_parms = NULL, fixed_initials = names(mkinmod$diffs)[-1],
solution_type = "auto",
+ method.ode = "lsoda",
method.modFit = "Marq",
control.modFit = list(),
plot = FALSE, quiet = FALSE, err = NULL, weight = "none",
@@ -85,6 +86,11 @@ mkinfit(mkinmod, observed,
dependence of degradation rates and metabolites). This argument is passed
on to the helper function \code{\link{mkinpredict}}.
}
+ \item{method.ode}{
+ The solution method passed via \code{\link{mkinpredict}} to
+ \code{\link{ode}} in case the solution type is "deSolve". The default
+ "lsoda" is performant, but sometimes fails to converge.
+ }
\item{method.modFit}{
The optimisation method passed to \code{\link{modFit}}. The default "Marq"
is the Levenberg Marquardt algorithm \code{\link{nls.lm}} from the package
diff --git a/man/mkinpredict.Rd b/man/mkinpredict.Rd
index afb57e0..97db90e 100644
--- a/man/mkinpredict.Rd
+++ b/man/mkinpredict.Rd
@@ -10,7 +10,7 @@
}
\usage{
mkinpredict(mkinmod, odeparms, odeini, outtimes, solution_type = "deSolve",
- map_output = TRUE, atol = 1e-08, rtol = 1e-10, ...)
+ method.ode = "lsoda", atol = 1e-08, rtol = 1e-10, map_output = TRUE, ...)
}
\arguments{
\item{mkinmod}{
@@ -33,12 +33,13 @@
The method that should be used for producing the predictions. This should
generally be "analytical" if there is only one observed variable, and
usually "deSolve" in the case of several observed variables. The third
- possibility "eigen" is faster but produces results that the author believes
- to be less accurate.
+ possibility "eigen" is faster but not applicable to some models e.g.
+ using FOMC for the parent compound.
}
- \item{map_output}{
- Boolean to specify if the output should list values for the observed
- variables (default) or for all state variables (if set to FALSE).
+ \item{method.ode}{
+ The solution method passed via \code{\link{mkinpredict}} to
+ \code{\link{ode}} in case the solution type is "deSolve". The default
+ "lsoda" is performant, but sometimes fails to converge.
}
\item{atol}{
Absolute error tolerance, passed to \code{\link{ode}}. Default is 1e-8,
@@ -48,6 +49,10 @@
Absolute error tolerance, passed to \code{\link{ode}}. Default is 1e-10,
much lower than in \code{\link{lsoda}}.
}
+ \item{map_output}{
+ Boolean to specify if the output should list values for the observed
+ variables (default) or for all state variables (if set to FALSE).
+ }
\item{\dots}{
Further arguments passed to the ode solver in case such a solver is used.
}

Contact - Imprint