<!-- Generated by pkgdown: do not edit by hand -->
<!DOCTYPE html>
<html lang="en">
<head>
<meta charset="utf-8">
<meta http-equiv="X-UA-Compatible" content="IE=edge">
<meta name="viewport" content="width=device-width, initial-scale=1.0">
<title>Fit a kinetic model to data with one or more state variables — mkinfit • mkin</title>
<!-- jquery -->
<script src="https://cdnjs.cloudflare.com/ajax/libs/jquery/3.3.1/jquery.min.js" integrity="sha256-FgpCb/KJQlLNfOu91ta32o/NMZxltwRo8QtmkMRdAu8=" crossorigin="anonymous"></script>
<!-- Bootstrap -->
<link rel="stylesheet" href="https://cdnjs.cloudflare.com/ajax/libs/twitter-bootstrap/3.3.7/css/bootstrap.min.css" integrity="sha256-916EbMg70RQy9LHiGkXzG8hSg9EdNy97GazNG/aiY1w=" crossorigin="anonymous" />
<script src="https://cdnjs.cloudflare.com/ajax/libs/twitter-bootstrap/3.3.7/js/bootstrap.min.js" integrity="sha256-U5ZEeKfGNOja007MMD3YBI0A3OSZOQbeG6z2f2Y0hu8=" crossorigin="anonymous"></script>
<!-- Font Awesome icons -->
<link rel="stylesheet" href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/4.7.0/css/font-awesome.min.css" integrity="sha256-eZrrJcwDc/3uDhsdt61sL2oOBY362qM3lon1gyExkL0=" crossorigin="anonymous" />
<!-- clipboard.js -->
<script src="https://cdnjs.cloudflare.com/ajax/libs/clipboard.js/2.0.4/clipboard.min.js" integrity="sha256-FiZwavyI2V6+EXO1U+xzLG3IKldpiTFf3153ea9zikQ=" crossorigin="anonymous"></script>
<!-- sticky kit -->
<script src="https://cdnjs.cloudflare.com/ajax/libs/sticky-kit/1.1.3/sticky-kit.min.js" integrity="sha256-c4Rlo1ZozqTPE2RLuvbusY3+SU1pQaJC0TjuhygMipw=" crossorigin="anonymous"></script>
<!-- pkgdown -->
<link href="../pkgdown.css" rel="stylesheet">
<script src="../pkgdown.js"></script>
<meta property="og:title" content="Fit a kinetic model to data with one or more state variables — mkinfit" />
<meta property="og:description" content="This function uses the Flexible Modelling Environment package
FME to create a function calculating the model cost, i.e. the
deviation between the kinetic model and the observed data. This model cost is
then minimised using the Port algorithm nlminb,
using the specified initial or fixed parameters and starting values.
Per default, parameters in the kinetic models are internally transformed in order
to better satisfy the assumption of a normal distribution of their estimators.
In each step of the optimsation, the kinetic model is solved using the
function mkinpredict. The variance of the residuals for each
observed variable can optionally be iteratively reweighted until convergence
using the argument reweight.method = "obs"." />
<meta name="twitter:card" content="summary" />
<!-- mathjax -->
<script src="https://cdnjs.cloudflare.com/ajax/libs/mathjax/2.7.5/MathJax.js" integrity="sha256-nvJJv9wWKEm88qvoQl9ekL2J+k/RWIsaSScxxlsrv8k=" crossorigin="anonymous"></script>
<script src="https://cdnjs.cloudflare.com/ajax/libs/mathjax/2.7.5/config/TeX-AMS-MML_HTMLorMML.js" integrity="sha256-84DKXVJXs0/F8OTMzX4UR909+jtl4G7SPypPavF+GfA=" crossorigin="anonymous"></script>
<!--[if lt IE 9]>
<script src="https://oss.maxcdn.com/html5shiv/3.7.3/html5shiv.min.js"></script>
<script src="https://oss.maxcdn.com/respond/1.4.2/respond.min.js"></script>
<![endif]-->
</head>
<body>
<div class="container template-reference-topic">
<header>
<div class="navbar navbar-default navbar-fixed-top" role="navigation">
<div class="container">
<div class="navbar-header">
<button type="button" class="navbar-toggle collapsed" data-toggle="collapse" data-target="#navbar" aria-expanded="false">
<span class="sr-only">Toggle navigation</span>
<span class="icon-bar"></span>
<span class="icon-bar"></span>
<span class="icon-bar"></span>
</button>
<span class="navbar-brand">
<a class="navbar-link" href="../index.html">mkin</a>
<span class="version label label-default" data-toggle="tooltip" data-placement="bottom" title="Released version">0.9.48.1</span>
</span>
</div>
<div id="navbar" class="navbar-collapse collapse">
<ul class="nav navbar-nav">
<li>
<a href="../reference/index.html">Functions and data</a>
</li>
<li class="dropdown">
<a href="#" class="dropdown-toggle" data-toggle="dropdown" role="button" aria-expanded="false">
Articles
<span class="caret"></span>
</a>
<ul class="dropdown-menu" role="menu">
<li>
<a href="../articles/mkin.html">Introduction to mkin</a>
</li>
<li>
<a href="../articles/FOCUS_D.html">Example evaluation of FOCUS Example Dataset D</a>
</li>
<li>
<a href="../articles/FOCUS_L.html">Example evaluation of FOCUS Laboratory Data L1 to L3</a>
</li>
<li>
<a href="../articles/web_only/FOCUS_Z.html">Example evaluation of FOCUS Example Dataset Z</a>
</li>
<li>
<a href="../articles/web_only/compiled_models.html">Performance benefit by using compiled model definitions in mkin</a>
</li>
<li>
<a href="../articles/twa.html">Calculation of time weighted average concentrations with mkin</a>
</li>
</ul>
</li>
<li>
<a href="../news/index.html">News</a>
</li>
</ul>
<ul class="nav navbar-nav navbar-right">
</ul>
</div><!--/.nav-collapse -->
</div><!--/.container -->
</div><!--/.navbar -->
</header>
<div class="row">
<div class="col-md-9 contents">
<div class="page-header">
<h1>Fit a kinetic model to data with one or more state variables</h1>
<div class="hidden name"><code>mkinfit.Rd</code></div>
</div>
<div class="ref-description">
<p>This function uses the Flexible Modelling Environment package
<code>FME</code> to create a function calculating the model cost, i.e. the
deviation between the kinetic model and the observed data. This model cost is
then minimised using the Port algorithm <code><a href='https://www.rdocumentation.org/packages/stats/topics/nlminb'>nlminb</a></code>,
using the specified initial or fixed parameters and starting values.
Per default, parameters in the kinetic models are internally transformed in order
to better satisfy the assumption of a normal distribution of their estimators.
In each step of the optimsation, the kinetic model is solved using the
function <code><a href='mkinpredict.html'>mkinpredict</a></code>. The variance of the residuals for each
observed variable can optionally be iteratively reweighted until convergence
using the argument <code>reweight.method = "obs"</code>.</p>
</div>
<pre class="usage"><span class='fu'>mkinfit</span>(<span class='no'>mkinmod</span>, <span class='no'>observed</span>,
<span class='kw'>parms.ini</span> <span class='kw'>=</span> <span class='st'>"auto"</span>,
<span class='kw'>state.ini</span> <span class='kw'>=</span> <span class='st'>"auto"</span>,
<span class='kw'>fixed_parms</span> <span class='kw'>=</span> <span class='kw'>NULL</span>, <span class='kw'>fixed_initials</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/names'>names</a></span>(<span class='no'>mkinmod</span>$<span class='no'>diffs</span>)[-<span class='fl'>1</span>],
<span class='kw'>from_max_mean</span> <span class='kw'>=</span> <span class='fl'>FALSE</span>,
<span class='kw'>solution_type</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/c'>c</a></span>(<span class='st'>"auto"</span>, <span class='st'>"analytical"</span>, <span class='st'>"eigen"</span>, <span class='st'>"deSolve"</span>),
<span class='kw'>method.ode</span> <span class='kw'>=</span> <span class='st'>"lsoda"</span>,
<span class='kw'>use_compiled</span> <span class='kw'>=</span> <span class='st'>"auto"</span>,
<span class='kw'>method.modFit</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/c'>c</a></span>(<span class='st'>"Port"</span>, <span class='st'>"Marq"</span>, <span class='st'>"SANN"</span>, <span class='st'>"Nelder-Mead"</span>, <span class='st'>"BFGS"</span>, <span class='st'>"CG"</span>, <span class='st'>"L-BFGS-B"</span>),
<span class='kw'>maxit.modFit</span> <span class='kw'>=</span> <span class='st'>"auto"</span>,
<span class='kw'>control.modFit</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/list'>list</a></span>(),
<span class='kw'>transform_rates</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>,
<span class='kw'>transform_fractions</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>,
<span class='kw'>plot</span> <span class='kw'>=</span> <span class='fl'>FALSE</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>FALSE</span>, <span class='kw'>err</span> <span class='kw'>=</span> <span class='kw'>NULL</span>,
<span class='kw'>weight</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/c'>c</a></span>(<span class='st'>"none"</span>, <span class='st'>"manual"</span>, <span class='st'>"std"</span>, <span class='st'>"mean"</span>, <span class='st'>"tc"</span>),
<span class='kw'>tc</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/c'>c</a></span>(<span class='kw'>sigma_low</span> <span class='kw'>=</span> <span class='fl'>0.5</span>, <span class='kw'>rsd_high</span> <span class='kw'>=</span> <span class='fl'>0.07</span>),
<span class='kw'>scaleVar</span> <span class='kw'>=</span> <span class='fl'>FALSE</span>,
<span class='kw'>atol</span> <span class='kw'>=</span> <span class='fl'>1e-8</span>, <span class='kw'>rtol</span> <span class='kw'>=</span> <span class='fl'>1e-10</span>, <span class='kw'>n.outtimes</span> <span class='kw'>=</span> <span class='fl'>100</span>,
<span class='kw'>reweight.method</span> <span class='kw'>=</span> <span class='kw'>NULL</span>,
<span class='kw'>reweight.tol</span> <span class='kw'>=</span> <span class='fl'>1e-8</span>, <span class='kw'>reweight.max.iter</span> <span class='kw'>=</span> <span class='fl'>10</span>,
<span class='kw'>trace_parms</span> <span class='kw'>=</span> <span class='fl'>FALSE</span>, <span class='no'>...</span>)</pre>
<h2 class="hasAnchor" id="arguments"><a class="anchor" href="#arguments"></a>Arguments</h2>
<table class="ref-arguments">
<colgroup><col class="name" /><col class="desc" /></colgroup>
<tr>
<th>mkinmod</th>
<td><p>A list of class <code><a href='mkinmod.html'>mkinmod</a></code>, containing the kinetic model to be
fitted to the data, or one of the shorthand names ("SFO", "FOMC", "DFOP",
"HS", "SFORB"). If a shorthand name is given, a parent only degradation
model is generated for the variable with the highest value in
<code>observed</code>.</p></td>
</tr>
<tr>
<th>observed</th>
<td><p>The observed data. It has to be in the long format as described in
<code>modFit</code>, i.e. the first column called "name" must contain the
name of the observed variable for each data point. The second column must
contain the times of observation, named "time". The third column must be
named "value" and contain the observed values. Optionally, a further column
can contain weights for each data point. Its name must be passed as a
further argument named <code>err</code> which is then passed on to
<code>modFit</code>.</p></td>
</tr>
<tr>
<th>parms.ini</th>
<td><p>A named vector of initial values for the parameters, including parameters
to be optimised and potentially also fixed parameters as indicated by
<code>fixed_parms</code>. If set to "auto", initial values for rate constants
are set to default values. Using parameter names that are not in the model
gives an error.</p>
<p>It is possible to only specify a subset of the parameters that the model
needs. You can use the parameter lists "bparms.ode" from a previously
fitted model, which contains the differential equation parameters from this
model. This works nicely if the models are nested. An example is given
below.</p></td>
</tr>
<tr>
<th>state.ini</th>
<td><p>A named vector of initial values for the state variables of the model. In
case the observed variables are represented by more than one model
variable, the names will differ from the names of the observed variables
(see <code>map</code> component of <code><a href='mkinmod.html'>mkinmod</a></code>). The default is to set
the initial value of the first model variable to the mean of the time zero
values for the variable with the maximum observed value, and all others to 0.
If this variable has no time zero observations, its initial value is set to 100.</p></td>
</tr>
<tr>
<th>fixed_parms</th>
<td><p>The names of parameters that should not be optimised but rather kept at the
values specified in <code>parms.ini</code>.</p></td>
</tr>
<tr>
<th>fixed_initials</th>
<td><p>The names of model variables for which the initial state at time 0 should
be excluded from the optimisation. Defaults to all state variables except
for the first one.</p></td>
</tr>
<tr>
<th>from_max_mean</th>
<td><p>If this is set to TRUE, and the model has only one observed variable, then
data before the time of the maximum observed value (after averaging for each
sampling time) are discarded, and this time is subtracted from all
remaining time values, so the time of the maximum observed mean value is
the new time zero.</p></td>
</tr>
<tr>
<th>solution_type</th>
<td><p>If set to "eigen", the solution of the system of differential equations is
based on the spectral decomposition of the coefficient matrix in cases that
this is possible. If set to "deSolve", a numerical ode solver from package
<code>deSolve</code> is used. If set to "analytical", an analytical
solution of the model is used. This is only implemented for simple
degradation experiments with only one state variable, i.e. with no
metabolites. The default is "auto", which uses "analytical" if possible,
otherwise "eigen" if the model can be expressed using eigenvalues and
eigenvectors, and finally "deSolve" for the remaining models (time
dependence of degradation rates and metabolites). This argument is passed
on to the helper function <code><a href='mkinpredict.html'>mkinpredict</a></code>.</p></td>
</tr>
<tr>
<th>method.ode</th>
<td><p>The solution method passed via <code><a href='mkinpredict.html'>mkinpredict</a></code> to
<code>ode</code> in case the solution type is "deSolve". The default
"lsoda" is performant, but sometimes fails to converge.</p></td>
</tr>
<tr>
<th>use_compiled</th>
<td><p>If set to <code>FALSE</code>, no compiled version of the <code><a href='mkinmod.html'>mkinmod</a></code>
model is used, in the calls to <code><a href='mkinpredict.html'>mkinpredict</a></code> even if
a compiled verion is present.</p></td>
</tr>
<tr>
<th>method.modFit</th>
<td><p>The optimisation method passed to <code>modFit</code>.</p>
<p>In order to optimally deal with problems where local minima occur, the
"Port" algorithm is now used per default as it is less prone to get trapped
in local minima and depends less on starting values for parameters than
the Levenberg Marquardt variant selected by "Marq". However, "Port" needs
more iterations.</p>
<p>The former default "Marq" is the Levenberg Marquardt algorithm
<code>nls.lm</code> from the package <code>minpack.lm</code> and usually needs
the least number of iterations.</p>
<p>The "Pseudo" algorithm is not included because it needs finite parameter bounds
which are currently not supported.</p>
<p>The "Newton" algorithm is not included because its number of iterations
can not be controlled by <code>control.modFit</code> and it does not appear
to provide advantages over the other algorithms.</p></td>
</tr>
<tr>
<th>maxit.modFit</th>
<td><p>Maximum number of iterations in the optimisation. If not "auto", this will
be passed to the method called by <code>modFit</code>, overriding
what may be specified in the next argument <code>control.modFit</code>.</p></td>
</tr>
<tr>
<th>control.modFit</th>
<td><p>Additional arguments passed to the optimisation method used by
<code>modFit</code>.</p></td>
</tr>
<tr>
<th>transform_rates</th>
<td><p>Boolean specifying if kinetic rate constants should be transformed in the
model specification used in the fitting for better compliance with the
assumption of normal distribution of the estimator. If TRUE, also
alpha and beta parameters of the FOMC model are log-transformed, as well
as k1 and k2 rate constants for the DFOP and HS models and the break point
tb of the HS model.
If FALSE, zero is used as a lower bound for the rates in the optimisation.</p></td>
</tr>
<tr>
<th>transform_fractions</th>
<td><p>Boolean specifying if formation fractions constants should be transformed in the
model specification used in the fitting for better compliance with the
assumption of normal distribution of the estimator. The default (TRUE) is
to do transformations. If TRUE, the g parameter of the DFOP and HS
models are also transformed, as they can also be seen as compositional
data. The transformation used for these transformations is the
<code><a href='ilr.html'>ilr</a></code> transformation.</p></td>
</tr>
<tr>
<th>plot</th>
<td><p>Should the observed values and the numerical solutions be plotted at each
stage of the optimisation?</p></td>
</tr>
<tr>
<th>quiet</th>
<td><p>Suppress printing out the current model cost after each improvement?</p></td>
</tr>
<tr>
<th>err </th>
<td><p>either <code>NULL</code>, or the name of the column with the
<em>error</em> estimates, used to weigh the residuals (see details of
<code>modCost</code>); if <code>NULL</code>, then the residuals are not weighed.</p></td>
</tr>
<tr>
<th>weight</th>
<td><p>only if <code>err</code>=<code>NULL</code>: how to weight the residuals, one of "none",
"std", "mean", see details of <code>modCost</code>, or "tc" for the
two component error model. The option "manual" is available for
the case that <code>err</code>!=<code>NULL</code>, but it is not necessary to specify it.</p></td>
</tr>
<tr>
<th>tc</th>
<td><p>The two components of the error model as used for (initial)
weighting</p></td>
</tr>
<tr>
<th>scaleVar</th>
<td><p>Will be passed to <code>modCost</code>. Default is not to scale Variables
according to the number of observations.</p></td>
</tr>
<tr>
<th>atol</th>
<td><p>Absolute error tolerance, passed to <code>ode</code>. Default is 1e-8,
lower than in <code>lsoda</code>.</p></td>
</tr>
<tr>
<th>rtol</th>
<td><p>Absolute error tolerance, passed to <code>ode</code>. Default is 1e-10,
much lower than in <code>lsoda</code>.</p></td>
</tr>
<tr>
<th>n.outtimes</th>
<td><p>The length of the dataseries that is produced by the model prediction
function <code><a href='mkinpredict.html'>mkinpredict</a></code>. This impacts the accuracy of
the numerical solver if that is used (see <code>solution_type</code> argument.
The default value is 100.</p></td>
</tr>
<tr>
<th>reweight.method</th>
<td><p>The method used for iteratively reweighting residuals, also known
as iteratively reweighted least squares (IRLS). Default is NULL,
i.e. no iterative weighting.
The first reweighting method is called "obs", meaning that each
observed variable is assumed to have its own variance. This variance
is estimated from the fit (mean squared residuals) and used for weighting
the residuals in each iteration until convergence of this estimate up to
<code>reweight.tol</code> or up to the maximum number of iterations
specified by <code>reweight.max.iter</code>.
The second reweighting method is called "tc" (two-component error model).
When using this method, the two components of an error model similar to
the one described by
Rocke and Lorenzato (1995) are estimated from the fit and the resulting
variances are used for weighting the residuals in each iteration until
convergence of these components or up to the maximum number of iterations
specified. Note that this method deviates from the model by Rocke and
Lorenzato, as their model implies that the errors follow a lognormal
distribution for large values, not a normal distribution as assumed by this
method.</p></td>
</tr>
<tr>
<th>reweight.tol</th>
<td><p>Tolerance for convergence criterion for the variance components
in IRLS fits.</p></td>
</tr>
<tr>
<th>reweight.max.iter</th>
<td><p>Maximum iterations in IRLS fits.</p></td>
</tr>
<tr>
<th>trace_parms</th>
<td><p>Should a trace of the parameter values be listed?</p></td>
</tr>
<tr>
<th>…</th>
<td><p>Further arguments that will be passed to <code>modFit</code>.</p></td>
</tr>
</table>
<h2 class="hasAnchor" id="value"><a class="anchor" href="#value"></a>Value</h2>
<p>A list with "mkinfit" and "modFit" in the class attribute.
A summary can be obtained by <code><a href='summary.mkinfit.html'>summary.mkinfit</a></code>.</p>
<h2 class="hasAnchor" id="see-also"><a class="anchor" href="#see-also"></a>See also</h2>
<div class='dont-index'><p>Plotting methods <code><a href='plot.mkinfit.html'>plot.mkinfit</a></code> and
<code><a href='mkinparplot.html'>mkinparplot</a></code>.</p>
<p>Comparisons of models fitted to the same data can be made using <code><a href='https://www.rdocumentation.org/packages/stats/topics/AIC'>AIC</a></code>
by virtue of the method <code><a href='logLik.mkinfit.html'>logLik.mkinfit</a></code>.</p>
<p>Fitting of several models to several datasets in a single call to
<code><a href='mmkin.html'>mmkin</a></code>.</p></div>
<h2 class="hasAnchor" id="note"><a class="anchor" href="#note"></a>Note</h2>
<p>The implementation of iteratively reweighted least squares is inspired by the
work of the KinGUII team at Bayer Crop Science (Walter Schmitt and Zhenglei
Gao). A similar implemention can also be found in CAKE 2.0, which is the
other GUI derivative of mkin, sponsored by Syngenta.</p>
<h2 class="hasAnchor" id="note"><a class="anchor" href="#note"></a>Note</h2>
<p>When using the "IORE" submodel for metabolites, fitting with
"transform_rates = TRUE" (the default) often leads to failures of the
numerical ODE solver. In this situation it may help to switch off the
internal rate transformation.</p>
<h2 class="hasAnchor" id="source"><a class="anchor" href="#source"></a>Source</h2>
<p>Rocke, David M. und Lorenzato, Stefan (1995) A two-component model for
measurement error in analytical chemistry. Technometrics 37(2), 176-184.</p>
<h2 class="hasAnchor" id="examples"><a class="anchor" href="#examples"></a>Examples</h2>
<pre class="examples"><div class='input'><span class='co'># Use shorthand notation for parent only degradation</span>
<span class='no'>fit</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='st'>"FOMC"</span>, <span class='no'>FOCUS_2006_C</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/summary'>summary</a></span>(<span class='no'>fit</span>)</div><div class='output co'>#> mkin version used for fitting: 0.9.48.1
#> R version used for fitting: 3.5.2
#> Date of fit: Mon Feb 25 14:53:09 2019
#> Date of summary: Mon Feb 25 14:53:09 2019
#>
#> Equations:
#> d_parent/dt = - (alpha/beta) * 1/((time/beta) + 1) * parent
#>
#> Model predictions using solution type analytical
#>
#> Fitted with method Port using 64 model solutions performed in 0.168 s
#>
#> Weighting: none
#>
#> Starting values for parameters to be optimised:
#> value type
#> parent_0 85.1 state
#> alpha 1.0 deparm
#> beta 10.0 deparm
#>
#> Starting values for the transformed parameters actually optimised:
#> value lower upper
#> parent_0 85.100000 -Inf Inf
#> log_alpha 0.000000 -Inf Inf
#> log_beta 2.302585 -Inf Inf
#>
#> Fixed parameter values:
#> None
#>
#> Optimised, transformed parameters with symmetric confidence intervals:
#> Estimate Std. Error Lower Upper
#> parent_0 85.87000 2.2460 80.38000 91.3700
#> log_alpha 0.05192 0.1605 -0.34080 0.4446
#> log_beta 0.65100 0.2801 -0.03452 1.3360
#>
#> Parameter correlation:
#> parent_0 log_alpha log_beta
#> parent_0 1.0000 -0.2033 -0.3624
#> log_alpha -0.2033 1.0000 0.9547
#> log_beta -0.3624 0.9547 1.0000
#>
#> Residual standard error: 2.275 on 6 degrees of freedom
#>
#> Backtransformed parameters:
#> Confidence intervals for internally transformed parameters are asymmetric.
#> t-test (unrealistically) based on the assumption of normal distribution
#> for estimators of untransformed parameters.
#> Estimate t value Pr(>t) Lower Upper
#> parent_0 85.870 38.230 1.069e-08 80.3800 91.370
#> alpha 1.053 6.231 3.953e-04 0.7112 1.560
#> beta 1.917 3.570 5.895e-03 0.9661 3.806
#>
#> Chi2 error levels in percent:
#> err.min n.optim df
#> All data 6.657 3 6
#> parent 6.657 3 6
#>
#> Estimated disappearance times:
#> DT50 DT90 DT50back
#> parent 1.785 15.15 4.56
#>
#> Data:
#> time variable observed predicted residual
#> 0 parent 85.1 85.875 -0.7749
#> 1 parent 57.9 55.191 2.7091
#> 3 parent 29.9 31.845 -1.9452
#> 7 parent 14.6 17.012 -2.4124
#> 14 parent 9.7 9.241 0.4590
#> 28 parent 6.6 4.754 1.8460
#> 63 parent 4.0 2.102 1.8977
#> 91 parent 3.9 1.441 2.4590
#> 119 parent 0.6 1.092 -0.4919</div><div class='input'>
<span class='co'># One parent compound, one metabolite, both single first order.</span>
<span class='co'># Use mkinsub for convenience in model formulation. Pathway to sink included per default.</span>
<span class='no'>SFO_SFO</span> <span class='kw'><-</span> <span class='fu'><a href='mkinmod.html'>mkinmod</a></span>(
<span class='kw'>parent</span> <span class='kw'>=</span> <span class='fu'><a href='mkinsub.html'>mkinsub</a></span>(<span class='st'>"SFO"</span>, <span class='st'>"m1"</span>),
<span class='kw'>m1</span> <span class='kw'>=</span> <span class='fu'><a href='mkinsub.html'>mkinsub</a></span>(<span class='st'>"SFO"</span>))</div><div class='output co'>#> <span class='message'>Successfully compiled differential equation model from auto-generated C code.</span></div><div class='input'><span class='co'># Fit the model to the FOCUS example dataset D using defaults</span>
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/print'>print</a></span>(<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/system.time'>system.time</a></span>(<span class='no'>fit</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO</span>, <span class='no'>FOCUS_2006_D</span>,
<span class='kw'>solution_type</span> <span class='kw'>=</span> <span class='st'>"eigen"</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)))</div><div class='output co'>#> User System verstrichen
#> 1.038 0.000 1.039 </div><div class='input'><span class='fu'><a href='https://www.rdocumentation.org/packages/stats/topics/coef'>coef</a></span>(<span class='no'>fit</span>)</div><div class='output co'>#> parent_0 log_k_parent_sink log_k_parent_m1 log_k_m1_sink
#> 99.59848 -3.03822 -2.98030 -5.24750 </div><div class='input'><span class='fu'><a href='endpoints.html'>endpoints</a></span>(<span class='no'>fit</span>)</div><div class='output co'>#> $ff
#> parent_sink parent_m1 m1_sink
#> 0.485524 0.514476 1.000000
#>
#> $SFORB
#> logical(0)
#>
#> $distimes
#> DT50 DT90
#> parent 7.022929 23.32967
#> m1 131.760712 437.69961
#> </div><div class='input'><span class='co'># deSolve is slower when no C compiler (gcc) was available during model generation</span>
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/print'>print</a></span>(<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/system.time'>system.time</a></span>(<span class='no'>fit.deSolve</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO</span>, <span class='no'>FOCUS_2006_D</span>,
<span class='kw'>solution_type</span> <span class='kw'>=</span> <span class='st'>"deSolve"</span>)))</div><div class='output co'>#> Model cost at call 1 : 18915.53
#> Model cost at call 2 : 18915.53
#> Model cost at call 6 : 11424.02
#> Model cost at call 10 : 11424
#> Model cost at call 12 : 4094.396
#> Model cost at call 16 : 4094.396
#> Model cost at call 19 : 1340.595
#> Model cost at call 20 : 1340.593
#> Model cost at call 25 : 1072.239
#> Model cost at call 28 : 1072.236
#> Model cost at call 30 : 874.2615
#> Model cost at call 33 : 874.2611
#> Model cost at call 35 : 616.2377
#> Model cost at call 37 : 616.2372
#> Model cost at call 40 : 467.4386
#> Model cost at call 42 : 467.4381
#> Model cost at call 46 : 398.2914
#> Model cost at call 48 : 398.2914
#> Model cost at call 49 : 398.2913
#> Model cost at call 51 : 395.0712
#> Model cost at call 54 : 395.0711
#> Model cost at call 56 : 378.3298
#> Model cost at call 59 : 378.3298
#> Model cost at call 62 : 376.9812
#> Model cost at call 64 : 376.9811
#> Model cost at call 67 : 375.2085
#> Model cost at call 69 : 375.2085
#> Model cost at call 70 : 375.2085
#> Model cost at call 71 : 375.2085
#> Model cost at call 72 : 374.5723
#> Model cost at call 74 : 374.5723
#> Model cost at call 77 : 374.0075
#> Model cost at call 79 : 374.0075
#> Model cost at call 80 : 374.0075
#> Model cost at call 82 : 373.1711
#> Model cost at call 84 : 373.1711
#> Model cost at call 87 : 372.6445
#> Model cost at call 88 : 372.1614
#> Model cost at call 90 : 372.1614
#> Model cost at call 91 : 372.1614
#> Model cost at call 94 : 371.6464
#> Model cost at call 99 : 371.4299
#> Model cost at call 101 : 371.4299
#> Model cost at call 104 : 371.4071
#> Model cost at call 106 : 371.4071
#> Model cost at call 107 : 371.4071
#> Model cost at call 109 : 371.2524
#> Model cost at call 113 : 371.2524
#> Model cost at call 114 : 371.2136
#> Model cost at call 115 : 371.2136
#> Model cost at call 116 : 371.2136
#> Model cost at call 119 : 371.2134
#> Model cost at call 120 : 371.2134
#> Model cost at call 122 : 371.2134
#> Model cost at call 123 : 371.2134
#> Model cost at call 125 : 371.2134
#> Model cost at call 126 : 371.2134
#> Model cost at call 135 : 371.2134
#> Model cost at call 146 : 371.2134
#> Optimisation by method Port successfully terminated.
#> User System verstrichen
#> 0.827 0.000 0.827 </div><div class='input'><span class='fu'><a href='https://www.rdocumentation.org/packages/stats/topics/coef'>coef</a></span>(<span class='no'>fit.deSolve</span>)</div><div class='output co'>#> parent_0 log_k_parent_sink log_k_parent_m1 log_k_m1_sink
#> 99.59848 -3.03822 -2.98030 -5.24750 </div><div class='input'><span class='fu'><a href='endpoints.html'>endpoints</a></span>(<span class='no'>fit.deSolve</span>)</div><div class='output co'>#> $ff
#> parent_sink parent_m1 m1_sink
#> 0.485524 0.514476 1.000000
#>
#> $SFORB
#> logical(0)
#>
#> $distimes
#> DT50 DT90
#> parent 7.022929 23.32967
#> m1 131.760711 437.69961
#> </div><div class='input'>
# Use stepwise fitting, using optimised parameters from parent only fit, FOMC
</div><div class='input'><span class='no'>FOMC_SFO</span> <span class='kw'><-</span> <span class='fu'><a href='mkinmod.html'>mkinmod</a></span>(
<span class='kw'>parent</span> <span class='kw'>=</span> <span class='fu'><a href='mkinsub.html'>mkinsub</a></span>(<span class='st'>"FOMC"</span>, <span class='st'>"m1"</span>),
<span class='kw'>m1</span> <span class='kw'>=</span> <span class='fu'><a href='mkinsub.html'>mkinsub</a></span>(<span class='st'>"SFO"</span>))</div><div class='output co'>#> <span class='message'>Successfully compiled differential equation model from auto-generated C code.</span></div><div class='input'><span class='co'># Fit the model to the FOCUS example dataset D using defaults</span>
<span class='no'>fit.FOMC_SFO</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>FOMC_SFO</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='co'># Use starting parameters from parent only FOMC fit</span>
<span class='no'>fit.FOMC</span> <span class='kw'>=</span> <span class='fu'>mkinfit</span>(<span class='st'>"FOMC"</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='no'>fit.FOMC_SFO</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>FOMC_SFO</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>,
<span class='kw'>parms.ini</span> <span class='kw'>=</span> <span class='no'>fit.FOMC</span>$<span class='no'>bparms.ode</span>)
<span class='co'># Use stepwise fitting, using optimised parameters from parent only fit, SFORB</span>
<span class='no'>SFORB_SFO</span> <span class='kw'><-</span> <span class='fu'><a href='mkinmod.html'>mkinmod</a></span>(
<span class='kw'>parent</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/list'>list</a></span>(<span class='kw'>type</span> <span class='kw'>=</span> <span class='st'>"SFORB"</span>, <span class='kw'>to</span> <span class='kw'>=</span> <span class='st'>"m1"</span>, <span class='kw'>sink</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>),
<span class='kw'>m1</span> <span class='kw'>=</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/list'>list</a></span>(<span class='kw'>type</span> <span class='kw'>=</span> <span class='st'>"SFO"</span>))</div><div class='output co'>#> <span class='message'>Successfully compiled differential equation model from auto-generated C code.</span></div><div class='input'><span class='co'># Fit the model to the FOCUS example dataset D using defaults</span>
<span class='no'>fit.SFORB_SFO</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFORB_SFO</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='no'>fit.SFORB_SFO.deSolve</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFORB_SFO</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>solution_type</span> <span class='kw'>=</span> <span class='st'>"deSolve"</span>,
<span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='co'># Use starting parameters from parent only SFORB fit (not really needed in this case)</span>
<span class='no'>fit.SFORB</span> <span class='kw'>=</span> <span class='fu'>mkinfit</span>(<span class='st'>"SFORB"</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='no'>fit.SFORB_SFO</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFORB_SFO</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>parms.ini</span> <span class='kw'>=</span> <span class='no'>fit.SFORB</span>$<span class='no'>bparms.ode</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)</div><div class='input'>
</div><div class='input'><span class='co'># Weighted fits, including IRLS</span>
<span class='no'>SFO_SFO.ff</span> <span class='kw'><-</span> <span class='fu'><a href='mkinmod.html'>mkinmod</a></span>(<span class='kw'>parent</span> <span class='kw'>=</span> <span class='fu'><a href='mkinsub.html'>mkinsub</a></span>(<span class='st'>"SFO"</span>, <span class='st'>"m1"</span>),
<span class='kw'>m1</span> <span class='kw'>=</span> <span class='fu'><a href='mkinsub.html'>mkinsub</a></span>(<span class='st'>"SFO"</span>), <span class='kw'>use_of_ff</span> <span class='kw'>=</span> <span class='st'>"max"</span>)</div><div class='output co'>#> <span class='message'>Successfully compiled differential equation model from auto-generated C code.</span></div><div class='input'><span class='no'>f.noweight</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO.ff</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/summary'>summary</a></span>(<span class='no'>f.noweight</span>)</div><div class='output co'>#> mkin version used for fitting: 0.9.48.1
#> R version used for fitting: 3.5.2
#> Date of fit: Mon Feb 25 14:53:22 2019
#> Date of summary: Mon Feb 25 14:53:22 2019
#>
#> Equations:
#> d_parent/dt = - k_parent * parent
#> d_m1/dt = + f_parent_to_m1 * k_parent * parent - k_m1 * m1
#>
#> Model predictions using solution type deSolve
#>
#> Fitted with method Port using 186 model solutions performed in 0.898 s
#>
#> Weighting: none
#>
#> Starting values for parameters to be optimised:
#> value type
#> parent_0 100.7500 state
#> k_parent 0.1000 deparm
#> k_m1 0.1001 deparm
#> f_parent_to_m1 0.5000 deparm
#>
#> Starting values for the transformed parameters actually optimised:
#> value lower upper
#> parent_0 100.750000 -Inf Inf
#> log_k_parent -2.302585 -Inf Inf
#> log_k_m1 -2.301586 -Inf Inf
#> f_parent_ilr_1 0.000000 -Inf Inf
#>
#> Fixed parameter values:
#> value type
#> m1_0 0 state
#>
#> Optimised, transformed parameters with symmetric confidence intervals:
#> Estimate Std. Error Lower Upper
#> parent_0 99.60000 1.61400 96.3300 102.9000
#> log_k_parent -2.31600 0.04187 -2.4010 -2.2310
#> log_k_m1 -5.24800 0.13610 -5.5230 -4.9720
#> f_parent_ilr_1 0.04096 0.06477 -0.0904 0.1723
#>
#> Parameter correlation:
#> parent_0 log_k_parent log_k_m1 f_parent_ilr_1
#> parent_0 1.0000 0.5178 -0.1701 -0.5489
#> log_k_parent 0.5178 1.0000 -0.3285 -0.5451
#> log_k_m1 -0.1701 -0.3285 1.0000 0.7466
#> f_parent_ilr_1 -0.5489 -0.5451 0.7466 1.0000
#>
#> Residual standard error: 3.211 on 36 degrees of freedom
#>
#> Backtransformed parameters:
#> Confidence intervals for internally transformed parameters are asymmetric.
#> t-test (unrealistically) based on the assumption of normal distribution
#> for estimators of untransformed parameters.
#> Estimate t value Pr(>t) Lower Upper
#> parent_0 99.600000 61.720 2.024e-38 96.330000 1.029e+02
#> k_parent 0.098700 23.880 5.700e-24 0.090660 1.074e-01
#> k_m1 0.005261 7.349 5.758e-09 0.003992 6.933e-03
#> f_parent_to_m1 0.514500 22.490 4.375e-23 0.468100 5.606e-01
#>
#> Chi2 error levels in percent:
#> err.min n.optim df
#> All data 6.398 4 15
#> parent 6.459 2 7
#> m1 4.690 2 8
#>
#> Resulting formation fractions:
#> ff
#> parent_m1 0.5145
#> parent_sink 0.4855
#>
#> Estimated disappearance times:
#> DT50 DT90
#> parent 7.023 23.33
#> m1 131.761 437.70
#>
#> Data:
#> time variable observed predicted residual
#> 0 parent 99.46 99.59848 -1.385e-01
#> 0 parent 102.04 99.59848 2.442e+00
#> 1 parent 93.50 90.23787 3.262e+00
#> 1 parent 92.50 90.23787 2.262e+00
#> 3 parent 63.23 74.07319 -1.084e+01
#> 3 parent 68.99 74.07319 -5.083e+00
#> 7 parent 52.32 49.91206 2.408e+00
#> 7 parent 55.13 49.91206 5.218e+00
#> 14 parent 27.27 25.01257 2.257e+00
#> 14 parent 26.64 25.01257 1.627e+00
#> 21 parent 11.50 12.53462 -1.035e+00
#> 21 parent 11.64 12.53462 -8.946e-01
#> 35 parent 2.85 3.14787 -2.979e-01
#> 35 parent 2.91 3.14787 -2.379e-01
#> 50 parent 0.69 0.71624 -2.624e-02
#> 50 parent 0.63 0.71624 -8.624e-02
#> 75 parent 0.05 0.06074 -1.074e-02
#> 75 parent 0.06 0.06074 -7.381e-04
#> 0 m1 0.00 0.00000 0.000e+00
#> 0 m1 0.00 0.00000 0.000e+00
#> 1 m1 4.84 4.80296 3.704e-02
#> 1 m1 5.64 4.80296 8.370e-01
#> 3 m1 12.91 13.02400 -1.140e-01
#> 3 m1 12.96 13.02400 -6.400e-02
#> 7 m1 22.97 25.04476 -2.075e+00
#> 7 m1 24.47 25.04476 -5.748e-01
#> 14 m1 41.69 36.69002 5.000e+00
#> 14 m1 33.21 36.69002 -3.480e+00
#> 21 m1 44.37 41.65310 2.717e+00
#> 21 m1 46.44 41.65310 4.787e+00
#> 35 m1 41.22 43.31312 -2.093e+00
#> 35 m1 37.95 43.31312 -5.363e+00
#> 50 m1 41.19 41.21831 -2.831e-02
#> 50 m1 40.01 41.21831 -1.208e+00
#> 75 m1 40.09 36.44703 3.643e+00
#> 75 m1 33.85 36.44703 -2.597e+00
#> 100 m1 31.04 31.98163 -9.416e-01
#> 100 m1 33.13 31.98163 1.148e+00
#> 120 m1 25.15 28.78984 -3.640e+00
#> 120 m1 33.31 28.78984 4.520e+00</div><div class='input'><span class='no'>f.irls</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO.ff</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>reweight.method</span> <span class='kw'>=</span> <span class='st'>"obs"</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/summary'>summary</a></span>(<span class='no'>f.irls</span>)</div><div class='output co'>#> mkin version used for fitting: 0.9.48.1
#> R version used for fitting: 3.5.2
#> Date of fit: Mon Feb 25 14:53:24 2019
#> Date of summary: Mon Feb 25 14:53:24 2019
#>
#> Equations:
#> d_parent/dt = - k_parent * parent
#> d_m1/dt = + f_parent_to_m1 * k_parent * parent - k_m1 * m1
#>
#> Model predictions using solution type deSolve
#>
#> Fitted with method Port using 551 model solutions performed in 2.572 s
#>
#> Weighting: none
#>
#> Iterative reweighting with method obs
#> Final mean squared residuals of observed variables:
#> parent m1
#> 11.573407 7.407845
#>
#> Starting values for parameters to be optimised:
#> value type
#> parent_0 100.7500 state
#> k_parent 0.1000 deparm
#> k_m1 0.1001 deparm
#> f_parent_to_m1 0.5000 deparm
#>
#> Starting values for the transformed parameters actually optimised:
#> value lower upper
#> parent_0 100.750000 -Inf Inf
#> log_k_parent -2.302585 -Inf Inf
#> log_k_m1 -2.301586 -Inf Inf
#> f_parent_ilr_1 0.000000 -Inf Inf
#>
#> Fixed parameter values:
#> value type
#> m1_0 0 state
#>
#> Optimised, transformed parameters with symmetric confidence intervals:
#> Estimate Std. Error Lower Upper
#> parent_0 99.67000 1.79200 96.04000 103.300
#> log_k_parent -2.31200 0.04560 -2.40400 -2.219
#> log_k_m1 -5.25100 0.12510 -5.50500 -4.998
#> f_parent_ilr_1 0.03785 0.06318 -0.09027 0.166
#>
#> Parameter correlation:
#> parent_0 log_k_parent log_k_m1 f_parent_ilr_1
#> parent_0 1.0000 0.5083 -0.1979 -0.6148
#> log_k_parent 0.5083 1.0000 -0.3894 -0.6062
#> log_k_m1 -0.1979 -0.3894 1.0000 0.7417
#> f_parent_ilr_1 -0.6148 -0.6062 0.7417 1.0000
#>
#> Residual standard error: 1.054 on 36 degrees of freedom
#>
#> Backtransformed parameters:
#> Confidence intervals for internally transformed parameters are asymmetric.
#> t-test (unrealistically) based on the assumption of normal distribution
#> for estimators of untransformed parameters.
#> Estimate t value Pr(>t) Lower Upper
#> parent_0 99.67000 55.630 8.185e-37 96.040000 1.033e+02
#> k_parent 0.09906 21.930 1.016e-22 0.090310 1.087e-01
#> k_m1 0.00524 7.996 8.486e-10 0.004066 6.753e-03
#> f_parent_to_m1 0.51340 23.000 2.038e-23 0.468100 5.584e-01
#>
#> Chi2 error levels in percent:
#> err.min n.optim df
#> All data 6.399 4 15
#> parent 6.466 2 7
#> m1 4.679 2 8
#>
#> Resulting formation fractions:
#> ff
#> parent_m1 0.5134
#> parent_sink 0.4866
#>
#> Estimated disappearance times:
#> DT50 DT90
#> parent 6.997 23.24
#> m1 132.282 439.43
#>
#> Data:
#> time variable observed predicted residual err
#> 0 parent 99.46 99.67218 -2.122e-01 3.402
#> 0 parent 102.04 99.67218 2.368e+00 3.402
#> 1 parent 93.50 90.27153 3.228e+00 3.402
#> 1 parent 92.50 90.27153 2.228e+00 3.402
#> 3 parent 63.23 74.04648 -1.082e+01 3.402
#> 3 parent 68.99 74.04648 -5.056e+00 3.402
#> 7 parent 52.32 49.82092 2.499e+00 3.402
#> 7 parent 55.13 49.82092 5.309e+00 3.402
#> 14 parent 27.27 24.90288 2.367e+00 3.402
#> 14 parent 26.64 24.90288 1.737e+00 3.402
#> 21 parent 11.50 12.44765 -9.476e-01 3.402
#> 21 parent 11.64 12.44765 -8.076e-01 3.402
#> 35 parent 2.85 3.11002 -2.600e-01 3.402
#> 35 parent 2.91 3.11002 -2.000e-01 3.402
#> 50 parent 0.69 0.70374 -1.374e-02 3.402
#> 50 parent 0.63 0.70374 -7.374e-02 3.402
#> 75 parent 0.05 0.05913 -9.134e-03 3.402
#> 75 parent 0.06 0.05913 8.662e-04 3.402
#> 0 m1 0.00 0.00000 0.000e+00 2.722
#> 0 m1 0.00 0.00000 0.000e+00 2.722
#> 1 m1 4.84 4.81328 2.672e-02 2.722
#> 1 m1 5.64 4.81328 8.267e-01 2.722
#> 3 m1 12.91 13.04779 -1.378e-01 2.722
#> 3 m1 12.96 13.04779 -8.779e-02 2.722
#> 7 m1 22.97 25.07615 -2.106e+00 2.722
#> 7 m1 24.47 25.07615 -6.062e-01 2.722
#> 14 m1 41.69 36.70729 4.983e+00 2.722
#> 14 m1 33.21 36.70729 -3.497e+00 2.722
#> 21 m1 44.37 41.65050 2.720e+00 2.722
#> 21 m1 46.44 41.65050 4.790e+00 2.722
#> 35 m1 41.22 43.28866 -2.069e+00 2.722
#> 35 m1 37.95 43.28866 -5.339e+00 2.722
#> 50 m1 41.19 41.19339 -3.386e-03 2.722
#> 50 m1 40.01 41.19339 -1.183e+00 2.722
#> 75 m1 40.09 36.43820 3.652e+00 2.722
#> 75 m1 33.85 36.43820 -2.588e+00 2.722
#> 100 m1 31.04 31.98971 -9.497e-01 2.722
#> 100 m1 33.13 31.98971 1.140e+00 2.722
#> 120 m1 25.15 28.80898 -3.659e+00 2.722
#> 120 m1 33.31 28.80898 4.501e+00 2.722</div><div class='input'><span class='no'>f.w.mean</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO.ff</span>, <span class='no'>FOCUS_2006_D</span>, <span class='kw'>weight</span> <span class='kw'>=</span> <span class='st'>"mean"</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/summary'>summary</a></span>(<span class='no'>f.w.mean</span>)</div><div class='output co'>#> mkin version used for fitting: 0.9.48.1
#> R version used for fitting: 3.5.2
#> Date of fit: Mon Feb 25 14:53:25 2019
#> Date of summary: Mon Feb 25 14:53:25 2019
#>
#> Equations:
#> d_parent/dt = - k_parent * parent
#> d_m1/dt = + f_parent_to_m1 * k_parent * parent - k_m1 * m1
#>
#> Model predictions using solution type deSolve
#>
#> Fitted with method Port using 155 model solutions performed in 0.709 s
#>
#> Weighting: mean
#>
#> Starting values for parameters to be optimised:
#> value type
#> parent_0 100.7500 state
#> k_parent 0.1000 deparm
#> k_m1 0.1001 deparm
#> f_parent_to_m1 0.5000 deparm
#>
#> Starting values for the transformed parameters actually optimised:
#> value lower upper
#> parent_0 100.750000 -Inf Inf
#> log_k_parent -2.302585 -Inf Inf
#> log_k_m1 -2.301586 -Inf Inf
#> f_parent_ilr_1 0.000000 -Inf Inf
#>
#> Fixed parameter values:
#> value type
#> m1_0 0 state
#>
#> Optimised, transformed parameters with symmetric confidence intervals:
#> Estimate Std. Error Lower Upper
#> parent_0 99.7300 1.93200 95.81000 103.6000
#> log_k_parent -2.3090 0.04837 -2.40700 -2.2110
#> log_k_m1 -5.2550 0.12070 -5.49900 -5.0100
#> f_parent_ilr_1 0.0354 0.06344 -0.09327 0.1641
#>
#> Parameter correlation:
#> parent_0 log_k_parent log_k_m1 f_parent_ilr_1
#> parent_0 1.0000 0.5004 -0.2143 -0.6514
#> log_k_parent 0.5004 1.0000 -0.4282 -0.6383
#> log_k_m1 -0.2143 -0.4282 1.0000 0.7390
#> f_parent_ilr_1 -0.6514 -0.6383 0.7390 1.0000
#>
#> Residual standard error: 0.09829 on 36 degrees of freedom
#>
#> Backtransformed parameters:
#> Confidence intervals for internally transformed parameters are asymmetric.
#> t-test (unrealistically) based on the assumption of normal distribution
#> for estimators of untransformed parameters.
#> Estimate t value Pr(>t) Lower Upper
#> parent_0 99.730000 51.630 1.166e-35 95.81000 1.036e+02
#> k_parent 0.099360 20.670 7.304e-22 0.09007 1.096e-01
#> k_m1 0.005224 8.287 3.649e-10 0.00409 6.672e-03
#> f_parent_to_m1 0.512500 22.860 2.497e-23 0.46710 5.578e-01
#>
#> Chi2 error levels in percent:
#> err.min n.optim df
#> All data 6.401 4 15
#> parent 6.473 2 7
#> m1 4.671 2 8
#>
#> Resulting formation fractions:
#> ff
#> parent_m1 0.5125
#> parent_sink 0.4875
#>
#> Estimated disappearance times:
#> DT50 DT90
#> parent 6.976 23.18
#> m1 132.696 440.81
#>
#> Data:
#> time variable observed predicted residual
#> 0 parent 99.46 99.73057 -0.270570
#> 0 parent 102.04 99.73057 2.309430
#> 1 parent 93.50 90.29805 3.201945
#> 1 parent 92.50 90.29805 2.201945
#> 3 parent 63.23 74.02503 -10.795028
#> 3 parent 68.99 74.02503 -5.035028
#> 7 parent 52.32 49.74838 2.571618
#> 7 parent 55.13 49.74838 5.381618
#> 14 parent 27.27 24.81588 2.454124
#> 14 parent 26.64 24.81588 1.824124
#> 21 parent 11.50 12.37885 -0.878849
#> 21 parent 11.64 12.37885 -0.738849
#> 35 parent 2.85 3.08022 -0.230219
#> 35 parent 2.91 3.08022 -0.170219
#> 50 parent 0.69 0.69396 -0.003958
#> 50 parent 0.63 0.69396 -0.063958
#> 75 parent 0.05 0.05789 -0.007888
#> 75 parent 0.06 0.05789 0.002112
#> 0 m1 0.00 0.00000 0.000000
#> 0 m1 0.00 0.00000 0.000000
#> 1 m1 4.84 4.82149 0.018512
#> 1 m1 5.64 4.82149 0.818512
#> 3 m1 12.91 13.06669 -0.156692
#> 3 m1 12.96 13.06669 -0.106692
#> 7 m1 22.97 25.10106 -2.131058
#> 7 m1 24.47 25.10106 -0.631058
#> 14 m1 41.69 36.72092 4.969077
#> 14 m1 33.21 36.72092 -3.510923
#> 21 m1 44.37 41.64835 2.721647
#> 21 m1 46.44 41.64835 4.791647
#> 35 m1 41.22 43.26923 -2.049225
#> 35 m1 37.95 43.26923 -5.319225
#> 50 m1 41.19 41.17364 0.016361
#> 50 m1 40.01 41.17364 -1.163639
#> 75 m1 40.09 36.43122 3.658776
#> 75 m1 33.85 36.43122 -2.581224
#> 100 m1 31.04 31.99612 -0.956124
#> 100 m1 33.13 31.99612 1.133876
#> 120 m1 25.15 28.82413 -3.674128
#> 120 m1 33.31 28.82413 4.485872</div><div class='input'><span class='no'>f.w.value</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO.ff</span>, <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/subset'>subset</a></span>(<span class='no'>FOCUS_2006_D</span>, <span class='no'>value</span> <span class='kw'>!=</span> <span class='fl'>0</span>), <span class='kw'>err</span> <span class='kw'>=</span> <span class='st'>"value"</span>,
<span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/summary'>summary</a></span>(<span class='no'>f.w.value</span>)</div><div class='output co'>#> mkin version used for fitting: 0.9.48.1
#> R version used for fitting: 3.5.2
#> Date of fit: Mon Feb 25 14:53:26 2019
#> Date of summary: Mon Feb 25 14:53:26 2019
#>
#> Equations:
#> d_parent/dt = - k_parent * parent
#> d_m1/dt = + f_parent_to_m1 * k_parent * parent - k_m1 * m1
#>
#> Model predictions using solution type deSolve
#>
#> Fitted with method Port using 174 model solutions performed in 0.81 s
#>
#> Weighting: manual
#>
#> Starting values for parameters to be optimised:
#> value type
#> parent_0 100.7500 state
#> k_parent 0.1000 deparm
#> k_m1 0.1001 deparm
#> f_parent_to_m1 0.5000 deparm
#>
#> Starting values for the transformed parameters actually optimised:
#> value lower upper
#> parent_0 100.750000 -Inf Inf
#> log_k_parent -2.302585 -Inf Inf
#> log_k_m1 -2.301586 -Inf Inf
#> f_parent_ilr_1 0.000000 -Inf Inf
#>
#> Fixed parameter values:
#> value type
#> m1_0 0 state
#>
#> Optimised, transformed parameters with symmetric confidence intervals:
#> Estimate Std. Error Lower Upper
#> parent_0 99.6600 2.712000 94.14000 105.2000
#> log_k_parent -2.2980 0.008118 -2.31500 -2.2820
#> log_k_m1 -5.2410 0.096690 -5.43800 -5.0450
#> f_parent_ilr_1 0.0231 0.057990 -0.09474 0.1409
#>
#> Parameter correlation:
#> parent_0 log_k_parent log_k_m1 f_parent_ilr_1
#> parent_0 1.00000 0.6843 -0.08687 -0.7564
#> log_k_parent 0.68435 1.0000 -0.12695 -0.5812
#> log_k_m1 -0.08687 -0.1269 1.00000 0.5195
#> f_parent_ilr_1 -0.75644 -0.5812 0.51952 1.0000
#>
#> Residual standard error: 0.08396 on 34 degrees of freedom
#>
#> Backtransformed parameters:
#> Confidence intervals for internally transformed parameters are asymmetric.
#> t-test (unrealistically) based on the assumption of normal distribution
#> for estimators of untransformed parameters.
#> Estimate t value Pr(>t) Lower Upper
#> parent_0 99.660000 36.75 2.957e-29 94.14000 1.052e+02
#> k_parent 0.100400 123.20 5.927e-47 0.09878 1.021e-01
#> k_m1 0.005295 10.34 2.447e-12 0.00435 6.444e-03
#> f_parent_to_m1 0.508200 24.79 1.184e-23 0.46660 5.497e-01
#>
#> Chi2 error levels in percent:
#> err.min n.optim df
#> All data 6.461 4 15
#> parent 6.520 2 7
#> m1 4.744 2 8
#>
#> Resulting formation fractions:
#> ff
#> parent_m1 0.5082
#> parent_sink 0.4918
#>
#> Estimated disappearance times:
#> DT50 DT90
#> parent 6.902 22.93
#> m1 130.916 434.89
#>
#> Data:
#> time variable observed predicted residual err
#> 0 parent 99.46 99.65571 -0.195715 99.46
#> 0 parent 102.04 99.65571 2.384285 102.04
#> 1 parent 93.50 90.13383 3.366170 93.50
#> 1 parent 92.50 90.13383 2.366170 92.50
#> 3 parent 63.23 73.73252 -10.502518 63.23
#> 3 parent 68.99 73.73252 -4.742518 68.99
#> 7 parent 52.32 49.34027 2.979728 52.32
#> 7 parent 55.13 49.34027 5.789728 55.13
#> 14 parent 27.27 24.42873 2.841271 27.27
#> 14 parent 26.64 24.42873 2.211271 26.64
#> 21 parent 11.50 12.09484 -0.594842 11.50
#> 21 parent 11.64 12.09484 -0.454842 11.64
#> 35 parent 2.85 2.96482 -0.114824 2.85
#> 35 parent 2.91 2.96482 -0.054824 2.91
#> 50 parent 0.69 0.65733 0.032670 0.69
#> 50 parent 0.63 0.65733 -0.027330 0.63
#> 75 parent 0.05 0.05339 -0.003386 0.05
#> 75 parent 0.06 0.05339 0.006614 0.06
#> 1 m1 4.84 4.82570 0.014301 4.84
#> 1 m1 5.64 4.82570 0.814301 5.64
#> 3 m1 12.91 13.06402 -0.154020 12.91
#> 3 m1 12.96 13.06402 -0.104020 12.96
#> 7 m1 22.97 25.04656 -2.076564 22.97
#> 7 m1 24.47 25.04656 -0.576564 24.47
#> 14 m1 41.69 36.53601 5.153988 41.69
#> 14 m1 33.21 36.53601 -3.326012 33.21
#> 21 m1 44.37 41.34639 3.023609 44.37
#> 21 m1 46.44 41.34639 5.093609 46.44
#> 35 m1 41.22 42.82669 -1.606690 41.22
#> 35 m1 37.95 42.82669 -4.876690 37.95
#> 50 m1 41.19 40.67342 0.516578 41.19
#> 50 m1 40.01 40.67342 -0.663422 40.01
#> 75 m1 40.09 35.91105 4.178947 40.09
#> 75 m1 33.85 35.91105 -2.061053 33.85
#> 100 m1 31.04 31.48161 -0.441612 31.04
#> 100 m1 33.13 31.48161 1.648388 33.13
#> 120 m1 25.15 28.32018 -3.170181 25.15
#> 120 m1 33.31 28.32018 4.989819 33.31</div><div class='input'>
</div><div class='input'><span class='co'># Manual weighting</span>
<span class='no'>dw</span> <span class='kw'><-</span> <span class='no'>FOCUS_2006_D</span>
<span class='no'>errors</span> <span class='kw'><-</span> <span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/c'>c</a></span>(<span class='kw'>parent</span> <span class='kw'>=</span> <span class='fl'>2</span>, <span class='kw'>m1</span> <span class='kw'>=</span> <span class='fl'>1</span>)
<span class='no'>dw</span>$<span class='no'>err.man</span> <span class='kw'><-</span> <span class='no'>errors</span>[<span class='no'>FOCUS_2006_D</span>$<span class='no'>name</span>]
<span class='no'>f.w.man</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO.ff</span>, <span class='no'>dw</span>, <span class='kw'>err</span> <span class='kw'>=</span> <span class='st'>"err.man"</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>)
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/summary'>summary</a></span>(<span class='no'>f.w.man</span>)</div><div class='output co'>#> mkin version used for fitting: 0.9.48.1
#> R version used for fitting: 3.5.2
#> Date of fit: Mon Feb 25 14:53:28 2019
#> Date of summary: Mon Feb 25 14:53:28 2019
#>
#> Equations:
#> d_parent/dt = - k_parent * parent
#> d_m1/dt = + f_parent_to_m1 * k_parent * parent - k_m1 * m1
#>
#> Model predictions using solution type deSolve
#>
#> Fitted with method Port using 270 model solutions performed in 1.25 s
#>
#> Weighting: manual
#>
#> Starting values for parameters to be optimised:
#> value type
#> parent_0 100.7500 state
#> k_parent 0.1000 deparm
#> k_m1 0.1001 deparm
#> f_parent_to_m1 0.5000 deparm
#>
#> Starting values for the transformed parameters actually optimised:
#> value lower upper
#> parent_0 100.750000 -Inf Inf
#> log_k_parent -2.302585 -Inf Inf
#> log_k_m1 -2.301586 -Inf Inf
#> f_parent_ilr_1 0.000000 -Inf Inf
#>
#> Fixed parameter values:
#> value type
#> m1_0 0 state
#>
#> Optimised, transformed parameters with symmetric confidence intervals:
#> Estimate Std. Error Lower Upper
#> parent_0 99.49000 1.33200 96.7800 102.2000
#> log_k_parent -2.32100 0.03550 -2.3930 -2.2490
#> log_k_m1 -5.24100 0.21280 -5.6730 -4.8100
#> f_parent_ilr_1 0.04571 0.08966 -0.1361 0.2275
#>
#> Parameter correlation:
#> parent_0 log_k_parent log_k_m1 f_parent_ilr_1
#> parent_0 1.00000 0.5312 -0.09456 -0.3351
#> log_k_parent 0.53123 1.0000 -0.17800 -0.3360
#> log_k_m1 -0.09456 -0.1780 1.00000 0.7616
#> f_parent_ilr_1 -0.33514 -0.3360 0.76156 1.0000
#>
#> Residual standard error: 2.628 on 36 degrees of freedom
#>
#> Backtransformed parameters:
#> Confidence intervals for internally transformed parameters are asymmetric.
#> t-test (unrealistically) based on the assumption of normal distribution
#> for estimators of untransformed parameters.
#> Estimate t value Pr(>t) Lower Upper
#> parent_0 99.490000 74.69 2.221e-41 96.780000 1.022e+02
#> k_parent 0.098140 28.17 2.012e-26 0.091320 1.055e-01
#> k_m1 0.005292 4.70 1.873e-05 0.003437 8.148e-03
#> f_parent_to_m1 0.516200 16.30 1.686e-18 0.452000 5.798e-01
#>
#> Chi2 error levels in percent:
#> err.min n.optim df
#> All data 6.400 4 15
#> parent 6.454 2 7
#> m1 4.708 2 8
#>
#> Resulting formation fractions:
#> ff
#> parent_m1 0.5162
#> parent_sink 0.4838
#>
#> Estimated disappearance times:
#> DT50 DT90
#> parent 7.063 23.46
#> m1 130.971 435.08
#>
#> Data:
#> time variable observed predicted residual err
#> 0 parent 99.46 99.48598 -0.025979 1
#> 0 parent 102.04 99.48598 2.554021 1
#> 1 parent 93.50 90.18612 3.313880 1
#> 1 parent 92.50 90.18612 2.313880 1
#> 3 parent 63.23 74.11316 -10.883163 1
#> 3 parent 68.99 74.11316 -5.123163 1
#> 7 parent 52.32 50.05030 2.269705 1
#> 7 parent 55.13 50.05030 5.079705 1
#> 14 parent 27.27 25.17975 2.090250 1
#> 14 parent 26.64 25.17975 1.460250 1
#> 21 parent 11.50 12.66765 -1.167654 1
#> 21 parent 11.64 12.66765 -1.027654 1
#> 35 parent 2.85 3.20616 -0.356164 1
#> 35 parent 2.91 3.20616 -0.296164 1
#> 50 parent 0.69 0.73562 -0.045619 1
#> 50 parent 0.63 0.73562 -0.105619 1
#> 75 parent 0.05 0.06326 -0.013256 1
#> 75 parent 0.06 0.06326 -0.003256 1
#> 0 m1 0.00 0.00000 0.000000 2
#> 0 m1 0.00 0.00000 0.000000 2
#> 1 m1 4.84 4.78729 0.052713 2
#> 1 m1 5.64 4.78729 0.852713 2
#> 3 m1 12.91 12.98785 -0.077848 2
#> 3 m1 12.96 12.98785 -0.027848 2
#> 7 m1 22.97 24.99695 -2.026946 2
#> 7 m1 24.47 24.99695 -0.526946 2
#> 14 m1 41.69 36.66353 5.026472 2
#> 14 m1 33.21 36.66353 -3.453528 2
#> 21 m1 44.37 41.65681 2.713186 2
#> 21 m1 46.44 41.65681 4.783186 2
#> 35 m1 41.22 43.35031 -2.130314 2
#> 35 m1 37.95 43.35031 -5.400314 2
#> 50 m1 41.19 41.25637 -0.066368 2
#> 50 m1 40.01 41.25637 -1.246368 2
#> 75 m1 40.09 36.46057 3.629429 2
#> 75 m1 33.85 36.46057 -2.610571 2
#> 100 m1 31.04 31.96929 -0.929293 2
#> 100 m1 33.13 31.96929 1.160707 2
#> 120 m1 25.15 28.76062 -3.610621 2
#> 120 m1 33.31 28.76062 4.549379 2</div><div class='input'><span class='no'>f.w.man.irls</span> <span class='kw'><-</span> <span class='fu'>mkinfit</span>(<span class='no'>SFO_SFO.ff</span>, <span class='no'>dw</span>, <span class='kw'>err</span> <span class='kw'>=</span> <span class='st'>"err.man"</span>, <span class='kw'>quiet</span> <span class='kw'>=</span> <span class='fl'>TRUE</span>,
<span class='kw'>reweight.method</span> <span class='kw'>=</span> <span class='st'>"obs"</span>)
<span class='fu'><a href='https://www.rdocumentation.org/packages/base/topics/summary'>summary</a></span>(<span class='no'>f.w.man.irls</span>)</div><div class='output co'>#> mkin version used for fitting: 0.9.48.1
#> R version used for fitting: 3.5.2
#> Date of fit: Mon Feb 25 14:53:31 2019
#> Date of summary: Mon Feb 25 14:53:31 2019
#>
#> Equations:
#> d_parent/dt = - k_parent * parent
#> d_m1/dt = + f_parent_to_m1 * k_parent * parent - k_m1 * m1
#>
#> Model predictions using solution type deSolve
#>
#> Fitted with method Port using 692 model solutions performed in 3.302 s
#>
#> Weighting: manual
#>
#> Iterative reweighting with method obs
#> Final mean squared residuals of observed variables:
#> parent m1
#> 11.573406 7.407846
#>
#> Starting values for parameters to be optimised:
#> value type
#> parent_0 100.7500 state
#> k_parent 0.1000 deparm
#> k_m1 0.1001 deparm
#> f_parent_to_m1 0.5000 deparm
#>
#> Starting values for the transformed parameters actually optimised:
#> value lower upper
#> parent_0 100.750000 -Inf Inf
#> log_k_parent -2.302585 -Inf Inf
#> log_k_m1 -2.301586 -Inf Inf
#> f_parent_ilr_1 0.000000 -Inf Inf
#>
#> Fixed parameter values:
#> value type
#> m1_0 0 state
#>
#> Optimised, transformed parameters with symmetric confidence intervals:
#> Estimate Std. Error Lower Upper
#> parent_0 99.67000 1.79200 96.04000 103.300
#> log_k_parent -2.31200 0.04560 -2.40400 -2.220
#> log_k_m1 -5.25100 0.12510 -5.50500 -4.998
#> f_parent_ilr_1 0.03785 0.06318 -0.09027 0.166
#>
#> Parameter correlation:
#> parent_0 log_k_parent log_k_m1 f_parent_ilr_1
#> parent_0 1.0000 0.5083 -0.1979 -0.6148
#> log_k_parent 0.5083 1.0000 -0.3894 -0.6062
#> log_k_m1 -0.1979 -0.3894 1.0000 0.7417
#> f_parent_ilr_1 -0.6148 -0.6062 0.7417 1.0000
#>
#> Residual standard error: 1.054 on 36 degrees of freedom
#>
#> Backtransformed parameters:
#> Confidence intervals for internally transformed parameters are asymmetric.
#> t-test (unrealistically) based on the assumption of normal distribution
#> for estimators of untransformed parameters.
#> Estimate t value Pr(>t) Lower Upper
#> parent_0 99.67000 55.630 8.185e-37 96.040000 1.033e+02
#> k_parent 0.09906 21.930 1.016e-22 0.090310 1.087e-01
#> k_m1 0.00524 7.996 8.486e-10 0.004066 6.753e-03
#> f_parent_to_m1 0.51340 23.000 2.039e-23 0.468100 5.584e-01
#>
#> Chi2 error levels in percent:
#> err.min n.optim df
#> All data 6.399 4 15
#> parent 6.466 2 7
#> m1 4.679 2 8
#>
#> Resulting formation fractions:
#> ff
#> parent_m1 0.5134
#> parent_sink 0.4866
#>
#> Estimated disappearance times:
#> DT50 DT90
#> parent 6.997 23.24
#> m1 132.282 439.43
#>
#> Data:
#> time variable observed predicted residual err.ini err
#> 0 parent 99.46 99.67217 -2.122e-01 1 3.402
#> 0 parent 102.04 99.67217 2.368e+00 1 3.402
#> 1 parent 93.50 90.27152 3.228e+00 1 3.402
#> 1 parent 92.50 90.27152 2.228e+00 1 3.402
#> 3 parent 63.23 74.04648 -1.082e+01 1 3.402
#> 3 parent 68.99 74.04648 -5.056e+00 1 3.402
#> 7 parent 52.32 49.82092 2.499e+00 1 3.402
#> 7 parent 55.13 49.82092 5.309e+00 1 3.402
#> 14 parent 27.27 24.90288 2.367e+00 1 3.402
#> 14 parent 26.64 24.90288 1.737e+00 1 3.402
#> 21 parent 11.50 12.44765 -9.477e-01 1 3.402
#> 21 parent 11.64 12.44765 -8.077e-01 1 3.402
#> 35 parent 2.85 3.11002 -2.600e-01 1 3.402
#> 35 parent 2.91 3.11002 -2.000e-01 1 3.402
#> 50 parent 0.69 0.70375 -1.375e-02 1 3.402
#> 50 parent 0.63 0.70375 -7.375e-02 1 3.402
#> 75 parent 0.05 0.05913 -9.134e-03 1 3.402
#> 75 parent 0.06 0.05913 8.661e-04 1 3.402
#> 0 m1 0.00 0.00000 0.000e+00 2 2.722
#> 0 m1 0.00 0.00000 0.000e+00 2 2.722
#> 1 m1 4.84 4.81328 2.672e-02 2 2.722
#> 1 m1 5.64 4.81328 8.267e-01 2 2.722
#> 3 m1 12.91 13.04779 -1.378e-01 2 2.722
#> 3 m1 12.96 13.04779 -8.779e-02 2 2.722
#> 7 m1 22.97 25.07615 -2.106e+00 2 2.722
#> 7 m1 24.47 25.07615 -6.062e-01 2 2.722
#> 14 m1 41.69 36.70729 4.983e+00 2 2.722
#> 14 m1 33.21 36.70729 -3.497e+00 2 2.722
#> 21 m1 44.37 41.65050 2.719e+00 2 2.722
#> 21 m1 46.44 41.65050 4.789e+00 2 2.722
#> 35 m1 41.22 43.28866 -2.069e+00 2 2.722
#> 35 m1 37.95 43.28866 -5.339e+00 2 2.722
#> 50 m1 41.19 41.19339 -3.387e-03 2 2.722
#> 50 m1 40.01 41.19339 -1.183e+00 2 2.722
#> 75 m1 40.09 36.43820 3.652e+00 2 2.722
#> 75 m1 33.85 36.43820 -2.588e+00 2 2.722
#> 100 m1 31.04 31.98971 -9.497e-01 2 2.722
#> 100 m1 33.13 31.98971 1.140e+00 2 2.722
#> 120 m1 25.15 28.80897 -3.659e+00 2 2.722
#> 120 m1 33.31 28.80897 4.501e+00 2 2.722</div></pre>
</div>
<div class="col-md-3 hidden-xs hidden-sm" id="sidebar">
<h2>Contents</h2>
<ul class="nav nav-pills nav-stacked">
<li><a href="#arguments">Arguments</a></li>
<li><a href="#value">Value</a></li>
<li><a href="#see-also">See also</a></li>
<li><a href="#note">Note</a></li>
<li><a href="#note">Note</a></li>
<li><a href="#source">Source</a></li>
<li><a href="#examples">Examples</a></li>
</ul>
<h2>Author</h2>
<p>Johannes Ranke</p>
</div>
</div>
<footer>
<div class="copyright">
<p>Developed by Johannes Ranke.</p>
</div>
<div class="pkgdown">
<p>Site built with <a href="https://pkgdown.r-lib.org/">pkgdown</a> 1.3.0.</p>
</div>
</footer>
</div>
</body>
</html>