View source: R/bootstrap.nlsfit.R
bootstrap.nlsfit | R Documentation |
Performs and bootstraps a non-linear least-squares fit to data with y and x errors.
bootstrap.nlsfit(fn, par.guess, y, x, bsamples, priors = list(param = c(), p = c(), psamples = c()), ..., lower = rep(x = -Inf, times = length(par.guess)), upper = rep(x = +Inf, times = length(par.guess)), dy, dx, CovMatrix, gr, dfn, mask, use.minpack.lm = TRUE, parallel = FALSE, error = sd, cov_fn = cov, maxiter = 500, success.infos = 1:3, relative.weights = FALSE, na.rm = FALSE)
fn |
|
par.guess |
initial guess values for the fit parameters. |
y |
the data as a one-dimensional numerical vector to be described by the fit function. |
x |
values of the explaining variable in form of a one-dimensional numerical vector. |
bsamples |
bootstrap samples of |
priors |
List possessing the elements |
... |
Additional parameters passed to |
lower |
Numeric vector of length |
upper |
Numeric vector of length |
dy, dx |
Numeric vector. Errors of the dependent and independent variable, respectively. These do not need to be specified as they can be computed from the bootstrap samples. In the case of parametric bootstrap it might would lead to a loss of information if they were computed from the pseudo-bootstrap samples. They must not be specified if a covariance matrix is given. |
CovMatrix |
complete variance-covariance matrix of dimensions
|
gr |
|
dfn |
|
mask |
logical or integer index vector. The mask is applied to select the observations from the data that are to be used in the fit. It is applied to |
use.minpack.lm |
use the |
parallel |
parallelise over bootstrap samples. The package
|
error |
Function that takes a sample vector and returns the error estimate. This is a parameter in order to support different resampling methods like jackknife. |
cov_fn |
function. Function to compute the covariance (matrix). Default is cov. |
maxiter |
integer. Maximum number of iterations that can be used in the optimization process. |
success.infos |
integer vector. When using |
relative.weights |
are the errors on y (and x) to be interpreted as relative weights instead of absolute ones? If TRUE, the covariance martix of the fit parameter results is multiplied by chi^2/dof. This is the default in many fit programs, e.g. gnuplot. |
na.rm |
logical. If set to |
returns a list of class 'bootstrapfit'. It returns all input parameters and adds in addition the following:
t0 |
the one dimensional numerical vector of length
|
t |
an array of dimensions |
bsamples |
the bootstrap samples used as an array of dimensions
|
Qval |
the p-value of the fit on the original data |
chisqr |
the residual chisqr value. |
dof |
the residual degrees of freedom of the fit. |
nx |
the number of x-values. |
tofn |
the original |
mask |
original |
Other NLS fit functions:
parametric.bootstrap.cov()
,
parametric.bootstrap()
,
parametric.nlsfit.cov()
,
parametric.nlsfit()
,
plot.bootstrapfit()
,
predict.bootstrapfit()
,
print.bootstrapfit()
,
simple.nlsfit()
,
summary.bootstrapfit()
## Declare some data. value <- c(0.1, 0.2, 0.31) dvalue <- c(0.01, 0.01, 0.015) x <- c(1, 2, 3) dx <- c(0.1, 0.1, 0.1) boot.R <- 1500 fn <- function (par, x, boot.r, ...) par[1] + par[2] * x ## Before we can use the fit with this data, we need to create bootstrap ## samples. We do not want to use the correlation matrix here. Note that you ## can simply use the parametric.nlsfit function as a convenient wrapper of ## the two steps. bsamples <- parametric.bootstrap(boot.R, c(value, x), c(dvalue, dx)) head(bsamples) fit.result <- bootstrap.nlsfit(fn, c(1, 1), value, x, bsamples) summary(fit.result) plot(fit.result, main = 'Ribbon on top') plot(fit.result, ribbon.on.top = FALSE, main = 'Ribbon below') residual_plot(fit.result, main = 'Residual Plot')
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.