aboutsummaryrefslogtreecommitdiff
path: root/man/mkinpredict.Rd
blob: c6aee75f6bed947a3f6b54cf0a18acd560943f4d (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
\name{mkinpredict}
\alias{mkinpredict}
\title{
  Produce predictions from a kinetic model using specifc parameters
}
\description{
  This function produces a time series for all the observed variables in a
  kinetic model as specified by \code{\link{mkinmod}}, using a specific set of
  kinetic parameters and initial values for the state variables.
}
\usage{
  mkinpredict(mkinmod, odeparms, odeini, outtimes, solution_type = "deSolve",
	      use_compiled = "auto", method.ode = "lsoda", atol = 1e-08, rtol = 1e-10,
        map_output = TRUE, ...)
}
\arguments{
  \item{mkinmod}{
    A kinetic model as produced by \code{\link{mkinmod}}.
  }
  \item{odeparms}{
    A numeric vector specifying the parameters used in the kinetic model, which
    is generally defined as a set of ordinary differential equations.
  }
  \item{odeini}{
    A numeric vectory containing the initial values of the state variables of
    the model. Note that the state variables can differ from the observed
    variables, for example in the case of the SFORB model.
  }
  \item{outtimes}{
    A numeric vector specifying the time points for which model predictions
    should be generated.
  }
  \item{solution_type}{
    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 not applicable to some models e.g.
    using FOMC for the parent compound.
  }
  \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{use_compiled}{
    If set to \code{FALSE}, no compiled version of the \code{\link{mkinmod}} 
    model is used, even if is present. 
  }
  \item{atol}{
    Absolute error tolerance, passed to \code{\link{ode}}. Default is 1e-8,
    lower than in \code{\link{lsoda}}.
  }
  \item{rtol}{
    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.
  }
}
\value{
  A matrix in the same format as the output of \code{\link{ode}}.
}
\author{
  Johannes Ranke
}
\examples{
  SFO <- mkinmod(degradinol = list(type = "SFO"))
  # Compare solution types
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20, 
	      solution_type = "analytical")
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20, 
	      solution_type = "deSolve")
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20, 
	      solution_type = "deSolve", use_compiled = FALSE)
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20, 
	      solution_type = "eigen")


  # Compare integration methods to analytical solution
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20, 
	      solution_type = "analytical")[21,]
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20, 
	      method = "lsoda")[21,]
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20,
	      method = "ode45")[21,]
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 0:20,
	      method = "rk4")[21,]
 # rk4 is not as precise here

  # The number of output times used to make a lot of difference until the
  # default for atol was adjusted
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 
	      seq(0, 20, by = 0.1))[201,]
  mkinpredict(SFO, c(k_degradinol_sink = 0.3), c(degradinol = 100), 
	      seq(0, 20, by = 0.01))[2001,]

  # Check compiled model versions - they are faster than the eigenvalue based solutions!
  SFO_SFO = mkinmod(parent = list(type = "SFO", to = "m1"),
                    m1 = list(type = "SFO"))
  system.time(
    print(mkinpredict(SFO_SFO, c(k_parent_m1 = 0.05, k_parent_sink = 0.1, k_m1_sink = 0.01),
                c(parent = 100, m1 = 0), seq(0, 20, by = 0.1),
                solution_type = "eigen")[201,]))
  system.time(
    print(mkinpredict(SFO_SFO, c(k_parent_m1 = 0.05, k_parent_sink = 0.1, k_m1_sink = 0.01),
                c(parent = 100, m1 = 0), seq(0, 20, by = 0.1), 
                solution_type = "deSolve")[201,]))
  system.time(
    print(mkinpredict(SFO_SFO, c(k_parent_m1 = 0.05, k_parent_sink = 0.1, k_m1_sink = 0.01),
                c(parent = 100, m1 = 0), seq(0, 20, by = 0.1), 
                solution_type = "deSolve", use_compiled = FALSE)[201,]))
}
\keyword{ manip }

Contact - Imprint