Title: | GSL Multi-Start Nonlinear Least-Squares Fitting |
---|---|
Description: | An R interface to nonlinear least-squares optimization with the GNU Scientific Library (GSL), see M. Galassi et al. (2009, ISBN:0954612078). The available trust region methods include the Levenberg-Marquardt algorithm with and without geodesic acceleration, the Steihaug-Toint conjugate gradient algorithm for large systems and several variants of Powell's dogleg algorithm. Multi-start optimization based on quasi-random samples is implemented using a modified version of the algorithm in Hickernell and Yuan (1997, OR Transactions). Bindings are provided to tune a number of parameters affecting the low-level aspects of the trust region algorithms. The interface mimics R's nls() function and returns model objects inheriting from the same class. |
Authors: | Joris Chau [aut, cre] |
Maintainer: | Joris Chau <[email protected]> |
License: | LGPL-3 |
Version: | 1.3.2 |
Built: | 2024-11-28 06:27:43 UTC |
Source: | CRAN |
Returns the analysis of variance (or deviance) tables for two or
more fitted "gsl_nls"
objects.
## S3 method for class 'gsl_nls' anova(object, ...)
## S3 method for class 'gsl_nls' anova(object, ...)
object |
An object inheriting from class |
... |
Additional objects inheriting from class |
A data.frame object of class "anova"
similar to anova
representing the
analysis-of-variance table of the fitted model objects when printed.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + 1 + rnorm(n, sd = 0.1) ) ## model obj1 <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) obj2 <- gsl_nls(fn = y ~ A * exp(-lam * x) + b, data = xy, start = c(A = 1, lam = 1, b = 0)) anova(obj1, obj2)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + 1 + rnorm(n, sd = 0.1) ) ## model obj1 <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) obj2 <- gsl_nls(fn = y ~ A * exp(-lam * x) + b, data = xy, start = c(A = 1, lam = 1, b = 0)) anova(obj1, obj2)
Returns the fitted model coefficients from a "gsl_nls"
object.
coefficients
can also be used as an alias.
## S3 method for class 'gsl_nls' coef(object, ...)
## S3 method for class 'gsl_nls' coef(object, ...)
object |
An object inheriting from class |
... |
At present no optional arguments are used. |
Named numeric vector of fitted coefficients similar to coef
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) coef(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) coef(obj)
Returns asymptotic or profile likelihood confidence intervals for the parameters in a
fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' confint(object, parm, level = 0.95, method = c("asymptotic", "profile"), ...)
## S3 method for class 'gsl_nls' confint(object, parm, level = 0.95, method = c("asymptotic", "profile"), ...)
object |
An object inheriting from class |
parm |
A character vector of parameter names for which to evaluate confidence intervals, defaults to all parameters. |
level |
A numeric scalar between 0 and 1 giving the level of the parameter confidence intervals. |
method |
Method to be used, either |
... |
At present no optional arguments are used. |
Method "asymptotic"
assumes (approximate) normality of the errors in the model and calculates
standard asymptotic confidence intervals based on the quantiles of a t-distribution. Method "profile"
calculates profile likelihood confidence intervals using the confint.nls
method
in the MASS package and for this reason is only available for "gsl_nls"
objects that
are also of class "nls"
.
A matrix with columns giving the lower and upper confidence limits for each parameter.
confint
, confint.nls
in package MASS.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) ## asymptotic ci's confint(obj) ## Not run: ## profile ci's (requires MASS) confint(obj, method = "profile") ## End(Not run)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) ## asymptotic ci's confint(obj) ## Not run: ## profile ci's (requires MASS) confint(obj, method = "profile") ## End(Not run)
confintd
is a generic function to compute confidence intervals for continuous functions
of the parameters in a fitted model. The function invokes particular methods which depend on the
class
of the first argument.
confintd(object, expr, level = 0.95, ...)
confintd(object, expr, level = 0.95, ...)
object |
A fitted model object. |
expr |
An expression or character vector that can be transformed to an |
level |
A numeric scalar between 0 and 1 giving the level of the derived parameter confidence intervals. |
... |
Additional argument(s) for methods |
A matrix with columns giving the fitted values and lower and upper confidence limits for each derived parameter. The row names list the individual derived parameter expressions.
Returns fitted values and confidence intervals for continuous functions of parameters
from a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' confintd(object, expr, level = 0.95, dtype = "symbolic", ...)
## S3 method for class 'gsl_nls' confintd(object, expr, level = 0.95, dtype = "symbolic", ...)
object |
A fitted model object. |
expr |
An expression or character vector that can be transformed to an |
level |
A numeric scalar between 0 and 1 giving the level of the derived parameter confidence intervals. |
dtype |
A character string equal to |
... |
Additional argument(s) for methods |
This method assumes (approximate) normality of the errors in the model and confidence intervals are
calculated using the delta method, i.e. a first-order Taylor approximation of the (continuous)
function of the parameters. If dtype = "symbolic"
(the default), expr
is differentiated
with respect to the parameters using symbolic differentiation with deriv
. As such,
each expression in expr
must contain only operators that are known to deriv
.
If dtype = "numeric"
, expr
is differentiated using numeric differentiation with
numericDeriv
, which should be used if expr
cannot be derived symbolically
with deriv
.
A matrix with columns giving the fitted values and lower and upper confidence limits for each derived parameter. The row names list the individual derived parameter expressions.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) ## delta method ci's confintd(obj, expr = c("log(lam)", "A / lam"))
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) ## delta method ci's confintd(obj, expr = c("log(lam)", "A / lam"))
Returns the deviance of a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' deviance(object, ...)
## S3 method for class 'gsl_nls' deviance(object, ...)
object |
An object inheriting from class |
... |
At present no optional arguments are used. |
Numeric deviance value similar to deviance
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) deviance(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) deviance(obj)
Returns the residual degrees-of-freedom from a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' df.residual(object, ...)
## S3 method for class 'gsl_nls' df.residual(object, ...)
object |
An object inheriting from class |
... |
At present no optional arguments are used. |
Integer residual degrees-of-freedom similar to df.residual
.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) df.residual(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) df.residual(obj)
Returns the fitted responses from a "gsl_nls"
object. fitted.values
can also be used as an alias.
## S3 method for class 'gsl_nls' fitted(object, ...)
## S3 method for class 'gsl_nls' fitted(object, ...)
object |
An object inheriting from class |
... |
At present no optional arguments are used. |
Numeric vector of fitted responses similar to fitted
.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) fitted(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) fitted(obj)
Returns the model formula from a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' formula(x, ...)
## S3 method for class 'gsl_nls' formula(x, ...)
x |
An object inheriting from class |
... |
At present no optional arguments are used. |
If the object inherits from class "nls"
returns the fitted model as a formula similar
to formula
. Otherwise returns the fitted model as a function.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) formula(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) formula(obj)
Determine the nonlinear least-squares estimates of the parameters of a
nonlinear model using the gsl_multifit_nlinear
module present in
the GNU Scientific Library (GSL).
gsl_nls(fn, ...) ## S3 method for class 'formula' gsl_nls( fn, data = parent.frame(), start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D"), control = gsl_nls_control(), lower, upper, jac = NULL, fvv = NULL, trace = FALSE, subset, weights, na.action, model = FALSE, ... ) ## S3 method for class 'function' gsl_nls( fn, y, start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D"), control = gsl_nls_control(), lower, upper, jac = NULL, fvv = NULL, trace = FALSE, weights, ... )
gsl_nls(fn, ...) ## S3 method for class 'formula' gsl_nls( fn, data = parent.frame(), start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D"), control = gsl_nls_control(), lower, upper, jac = NULL, fvv = NULL, trace = FALSE, subset, weights, na.action, model = FALSE, ... ) ## S3 method for class 'function' gsl_nls( fn, y, start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D"), control = gsl_nls_control(), lower, upper, jac = NULL, fvv = NULL, trace = FALSE, weights, ... )
fn |
a nonlinear model defined either as a two-sided formula including variables and parameters, or as a function returning a numeric vector, with first argument the vector of parameters to be estimated. See the individual method descriptions below. |
data |
an optional data frame in which to evaluate the variables in |
y |
numeric response vector if |
start |
a vector, list or matrix of initial parameter values or parameter ranges.
|
algorithm |
character string specifying the algorithm to use. The following choices are supported:
|
control |
an optional list of control parameters to tune the least squares iterations and multistart algorithm.
See |
lower |
a named list or named numeric vector of parameter lower bounds. If missing (default), the parameters are unconstrained from below. |
upper |
a named list or named numeric vector of parameter upper bounds. If missing (default), the parameters are unconstrained from above. |
jac |
either |
fvv |
either |
trace |
logical value indicating if a trace of the iteration progress should be printed.
Default is |
subset |
an optional vector specifying a subset of observations to be used in the fitting process.
This argument is only used if |
weights |
an optional numeric vector of (fixed) weights. When present, the objective function is weighted least squares. |
na.action |
a function which indicates what should happen when the data contain |
model |
a logical value. If |
... |
additional arguments passed to the calls of |
If fn
is a formula
returns a list object of class nls
.
If fn
is a function
returns a list object of class gsl_nls
.
See the individual method descriptions for the structures of the returned lists and the generic functions
applicable to objects of both classes.
gsl_nls(formula)
: If fn
is a formula
, the returned list object is of classes gsl_nls
and nls
.
Therefore, all generic functions applicable to objects of class nls
, such as anova
, coef
, confint
,
deviance
, df.residual
, fitted
, formula
, logLik
, nobs
, predict
, print
, profile
,
residuals
, summary
, vcov
and weights
are also applicable to the returned list object.
In addition, a method confintd
is available for inference of derived parameters.
gsl_nls(function)
: If fn
is a function
, the first argument must be the vector of parameters and
the function should return a numeric vector containing the nonlinear model evaluations at
the provided parameter and predictor or covariate vectors. In addition, the argument y
needs to contain the numeric vector of observed responses, equal in length to the numeric
vector returned by fn
. The returned list object is (only) of class gsl_nls
.
Although the returned object is not of class nls
, the following generic functions remain
applicable for an object of class gsl_nls
: anova
, coef
, confint
, deviance
,
df.residual
, fitted
, formula
, logLik
, nobs
, predict
, print
,
residuals
, summary
, vcov
and weights
. In addition, a method confintd
is available for inference of derived parameters.
If start
is a list or matrix of parameter ranges, or contains any missing values, a modified version of the multi-start algorithm described in
Hickernell and Yuan (1997) is applied. Note that the start
parameter ranges are only used to bound the domain for the
starting values, i.e. the resulting parameter estimates are not constrained to lie within these bounds, use lower
and/or upper
for
this purpose instead. Quasi-random starting values are sampled in the unit hypercube from a Sobol sequence if p < 41
and from a Halton sequence (up to p = 1229
) otherwise.
The initial starting values are scaled to the specified parameter ranges using an inverse (scaled) logistic function favoring starting values near the center of the
(scaled) domain. The trust region algorithm as specified by algorithm
used for the inexpensive and expensive local search (see Algorithm 2.1 of Hickernell
and Yuan (1997)) are the same, only differing in the number of search iterations mstart_p
versus mstart_maxiter
, where
mstart_p
is typically much smaller than mstart_maxiter
. When a new stationary point is detected, the scaling step from the unit hypercube to
the starting value domain is updated using the diagonal of the estimated trust method's scaling matrix D
, which improves optimization performance
especially when the parameters live on very different scales. The multi-start algorithm terminates when NSP (number of stationary points)
is larger than or equal to mstart_minsp
and NWSP (number of worse stationary points) is larger than or equal to mstart_r
times NSP,
or when the maximum number of major iterations mstart_maxstart
is reached. After termination of the multi-start algorithm, a full
single-start optimization is executed starting from the best multi-start solution.
If start
contains missing (or infinite) values, the multi-start algorithm is executed without fixed parameter ranges for the missing parameters.
The ranges for the missing parameters are initialized to the unit interval and dynamically increased or decreased in each major iteration
of the multi-start algorithm. The decision to increase or decrease a parameter range is driven by the minimum and maximum parameter values
attained by the first mstart_q
inexpensive local searches ordered by their squared loss, which typically provide a decent indication of the
order of magnitude of the parameter range in which to search for the optimal solution. Note that this procedure is not expected to always
return a global minimum of the nonlinear least-squares objective. Especially when the objective function contains many local optima,
the algorithm may be unable to select parameter ranges that include the global minimizing solution. In this case, it may help to increase
the values of mstart_n
, mstart_r
or mstart_minsp
to avoid early termination of the algorithm at the cost of
increased computational effort.
M. Galassi et al., GNU Scientific Library Reference Manual (3rd Ed.), ISBN 0954612078.
Hickernell, F.J. and Yuan, Y. (1997) “A simple multistart algorithm for global optimization”, OR Transactions, Vol. 1 (2).
https://www.gnu.org/software/gsl/doc/html/nls.html
# Example 1: exponential model # (https://www.gnu.org/software/gsl/doc/html/nls.html#exponential-fitting-example) ## data set.seed(1) n <- 25 x <- (seq_len(n) - 1) * 3 / (n - 1) f <- function(A, lam, b, x) A * exp(-lam * x) + b y <- f(A = 5, lam = 1.5, b = 1, x) + rnorm(n, sd = 0.25) ## model fit ex1_fit <- gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = 0, lam = 0, b = 0) ## starting values ) summary(ex1_fit) ## model summary predict(ex1_fit, interval = "prediction") ## prediction intervals ## multi-start gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = list(A = c(0, 100), lam = c(0, 10), b = c(-10, 10)) ## starting ranges ) ## missing starting values gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = NA, lam = NA, b = NA) ## unknown start ) ## analytic Jacobian 1 gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = 0, lam = 0, b = 0), ## starting values jac = function(par) with(as.list(par), ## jacobian cbind(A = exp(-lam * x), lam = -A * x * exp(-lam * x), b = 1) ) ) ## analytic Jacobian 2 gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = 0, lam = 0, b = 0), ## starting values jac = TRUE ## automatic derivation ) ## self-starting model gsl_nls( fn = y ~ SSasymp(x, Asym, R0, lrc), ## model formula data = data.frame(x = x, y = y) ## model fit data ) # Example 2: Gaussian function # (https://www.gnu.org/software/gsl/doc/html/nls.html#geodesic-acceleration-example-2) ## data set.seed(1) n <- 100 x <- seq_len(n) / n f <- function(a, b, c, x) a * exp(-(x - b)^2 / (2 * c^2)) y <- f(a = 5, b = 0.4, c = 0.15, x) * rnorm(n, mean = 1, sd = 0.1) ## Levenberg-Marquardt (default) gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values trace = TRUE ## verbose output ) ## Levenberg-Marquardt w/ geodesic acceleration 1 gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values algorithm = "lmaccel", ## algorithm trace = TRUE ## verbose output ) ## Levenberg-Marquardt w/ geodesic acceleration 2 ## second directional derivative fvv <- function(par, v, x) { with(as.list(par), { zi <- (x - b) / c ei <- exp(-zi^2 / 2) 2 * v[["a"]] * v[["b"]] * zi / c * ei + 2 * v[["a"]] * v[["c"]] * zi^2 / c * ei - v[["b"]]^2 * a / c^2 * (1 - zi^2) * ei - 2 * v[["b"]] * v[["c"]] * a / c^2 * zi * (2 - zi^2) * ei - v[["c"]]^2 * a / c^2 * zi^2 * (3 - zi^2) * ei }) } ## analytic fvv 1 gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values algorithm = "lmaccel", ## algorithm trace = TRUE, ## verbose output fvv = fvv, ## analytic fvv x = x ## argument passed to fvv ) ## analytic fvv 2 gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values algorithm = "lmaccel", ## algorithm trace = TRUE, ## verbose output fvv = TRUE ## automatic derivation ) # Example 3: Branin function # (https://www.gnu.org/software/gsl/doc/html/nls.html#comparing-trs-methods-example) ## Branin model function branin <- function(x) { a <- c(-5.1 / (4 * pi^2), 5 / pi, -6, 10, 1 / (8 * pi)) f1 <- x[2] + a[1] * x[1]^2 + a[2] * x[1] + a[3] f2 <- sqrt(a[4] * (1 + (1 - a[5]) * cos(x[1]))) c(f1, f2) } ## Dogleg minimization w/ model as function gsl_nls( fn = branin, ## model function y = c(0, 0), ## response vector start = c(x1 = 6, x2 = 14.5), ## starting values algorithm = "dogleg" ## algorithm ) # Available example problems nls_test_list() ## BOD regression ## (https://www.itl.nist.gov/div898/strd/nls/nls_main.shtml) (boxbod <- nls_test_problem(name = "BoxBOD")) with(boxbod, gsl_nls( fn = fn, data = data, start = list(b1 = NA, b2 = NA) ) ) ## Rosenbrock function (rosenbrock <- nls_test_problem(name = "Rosenbrock")) with(rosenbrock, gsl_nls( fn = fn, y = y, start = c(x1 = NA, x2 = NA), jac = jac ) )
# Example 1: exponential model # (https://www.gnu.org/software/gsl/doc/html/nls.html#exponential-fitting-example) ## data set.seed(1) n <- 25 x <- (seq_len(n) - 1) * 3 / (n - 1) f <- function(A, lam, b, x) A * exp(-lam * x) + b y <- f(A = 5, lam = 1.5, b = 1, x) + rnorm(n, sd = 0.25) ## model fit ex1_fit <- gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = 0, lam = 0, b = 0) ## starting values ) summary(ex1_fit) ## model summary predict(ex1_fit, interval = "prediction") ## prediction intervals ## multi-start gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = list(A = c(0, 100), lam = c(0, 10), b = c(-10, 10)) ## starting ranges ) ## missing starting values gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = NA, lam = NA, b = NA) ## unknown start ) ## analytic Jacobian 1 gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = 0, lam = 0, b = 0), ## starting values jac = function(par) with(as.list(par), ## jacobian cbind(A = exp(-lam * x), lam = -A * x * exp(-lam * x), b = 1) ) ) ## analytic Jacobian 2 gsl_nls( fn = y ~ A * exp(-lam * x) + b, ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(A = 0, lam = 0, b = 0), ## starting values jac = TRUE ## automatic derivation ) ## self-starting model gsl_nls( fn = y ~ SSasymp(x, Asym, R0, lrc), ## model formula data = data.frame(x = x, y = y) ## model fit data ) # Example 2: Gaussian function # (https://www.gnu.org/software/gsl/doc/html/nls.html#geodesic-acceleration-example-2) ## data set.seed(1) n <- 100 x <- seq_len(n) / n f <- function(a, b, c, x) a * exp(-(x - b)^2 / (2 * c^2)) y <- f(a = 5, b = 0.4, c = 0.15, x) * rnorm(n, mean = 1, sd = 0.1) ## Levenberg-Marquardt (default) gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values trace = TRUE ## verbose output ) ## Levenberg-Marquardt w/ geodesic acceleration 1 gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values algorithm = "lmaccel", ## algorithm trace = TRUE ## verbose output ) ## Levenberg-Marquardt w/ geodesic acceleration 2 ## second directional derivative fvv <- function(par, v, x) { with(as.list(par), { zi <- (x - b) / c ei <- exp(-zi^2 / 2) 2 * v[["a"]] * v[["b"]] * zi / c * ei + 2 * v[["a"]] * v[["c"]] * zi^2 / c * ei - v[["b"]]^2 * a / c^2 * (1 - zi^2) * ei - 2 * v[["b"]] * v[["c"]] * a / c^2 * zi * (2 - zi^2) * ei - v[["c"]]^2 * a / c^2 * zi^2 * (3 - zi^2) * ei }) } ## analytic fvv 1 gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values algorithm = "lmaccel", ## algorithm trace = TRUE, ## verbose output fvv = fvv, ## analytic fvv x = x ## argument passed to fvv ) ## analytic fvv 2 gsl_nls( fn = y ~ a * exp(-(x - b)^2 / (2 * c^2)), ## model formula data = data.frame(x = x, y = y), ## model fit data start = c(a = 1, b = 0, c = 1), ## starting values algorithm = "lmaccel", ## algorithm trace = TRUE, ## verbose output fvv = TRUE ## automatic derivation ) # Example 3: Branin function # (https://www.gnu.org/software/gsl/doc/html/nls.html#comparing-trs-methods-example) ## Branin model function branin <- function(x) { a <- c(-5.1 / (4 * pi^2), 5 / pi, -6, 10, 1 / (8 * pi)) f1 <- x[2] + a[1] * x[1]^2 + a[2] * x[1] + a[3] f2 <- sqrt(a[4] * (1 + (1 - a[5]) * cos(x[1]))) c(f1, f2) } ## Dogleg minimization w/ model as function gsl_nls( fn = branin, ## model function y = c(0, 0), ## response vector start = c(x1 = 6, x2 = 14.5), ## starting values algorithm = "dogleg" ## algorithm ) # Available example problems nls_test_list() ## BOD regression ## (https://www.itl.nist.gov/div898/strd/nls/nls_main.shtml) (boxbod <- nls_test_problem(name = "BoxBOD")) with(boxbod, gsl_nls( fn = fn, data = data, start = list(b1 = NA, b2 = NA) ) ) ## Rosenbrock function (rosenbrock <- nls_test_problem(name = "Rosenbrock")) with(rosenbrock, gsl_nls( fn = fn, y = y, start = c(x1 = NA, x2 = NA), jac = jac ) )
Allow the user to tune the characteristics of the gsl_nls
and gsl_nls_large
nonlinear least squares algorithms.
gsl_nls_control( maxiter = 100, scale = "more", solver = "qr", fdtype = "forward", factor_up = 2, factor_down = 3, avmax = 0.75, h_df = sqrt(.Machine$double.eps), h_fvv = 0.02, xtol = sqrt(.Machine$double.eps), ftol = sqrt(.Machine$double.eps), gtol = sqrt(.Machine$double.eps), mstart_n = 30, mstart_p = 5, mstart_q = mstart_n%/%10, mstart_r = 4, mstart_s = 2, mstart_tol = 0.25, mstart_maxiter = 10, mstart_maxstart = 250, mstart_minsp = 1, ... )
gsl_nls_control( maxiter = 100, scale = "more", solver = "qr", fdtype = "forward", factor_up = 2, factor_down = 3, avmax = 0.75, h_df = sqrt(.Machine$double.eps), h_fvv = 0.02, xtol = sqrt(.Machine$double.eps), ftol = sqrt(.Machine$double.eps), gtol = sqrt(.Machine$double.eps), mstart_n = 30, mstart_p = 5, mstart_q = mstart_n%/%10, mstart_r = 4, mstart_s = 2, mstart_tol = 0.25, mstart_maxiter = 10, mstart_maxstart = 250, mstart_minsp = 1, ... )
maxiter |
positive integer, termination occurs when the number of iterations reaches |
scale |
character, scaling method or damping strategy determining the diagonal scaling matrix D. The following options are supported:
|
solver |
character, method used to solve the linear least squares system resulting as a subproblem in each iteration.
For large-scale problems fitted with
|
fdtype |
character, method used to numerically approximate the Jacobian and/or second-order derivatives
when geodesic acceleration is used. Either |
factor_up |
numeric factor by which to increase the trust region radius when a search step is accepted. Too large values may destabilize the search, too small values slow down the search, defaults to 2. |
factor_down |
numeric factor by which to decrease the trust region radius when a search step is rejected. Too large values may destabilize the search, too small values slow down the search, defaults to 3. |
avmax |
numeric value, the ratio of the acceleration term to the velocity term when using geodesic acceleration to
solve the nonlinear least squares problem. Any steps with a ratio larger than |
h_df |
numeric value, the step size for approximating the Jacobian matrix with finite differences, defaults to |
h_fvv |
numeric value, the step size for approximating the second directional derivative when geodesic acceleration
is used to solve the nonlinear least squares problem, defaults to 0.02. This is only used if no analytic second
directional derivative ( |
xtol |
numeric value, termination occurs when the relative change in parameters between iterations is |
ftol |
numeric value, termination occurs when the relative change in sum of squared residuals between iterations is |
gtol |
numeric value, termination occurs when the relative size of the gradient of the sum of squared residuals is |
mstart_n |
positive integer, number of quasi-random points drawn in each major iteration, parameter |
mstart_p |
positive integer, number of iterations of inexpensive local search to concentrate the sample, parameter |
mstart_q |
positive integer, number of points retained in the concentrated sample, parameter |
mstart_r |
positive integer, scaling factor of number of stationary points determining when the multi-start algorithm terminates, parameter |
mstart_s |
positive integer, minimum number of iterations a point needs to be retained before starting an efficient local search, parameter |
mstart_tol |
numeric value, multiplicative tolerance |
mstart_maxiter |
positive integer, maximum number of iterations in the efficient local search algorithm (Algorithm B, Hickernell and Yuan (1997)), defaults to 10. |
mstart_maxstart |
positive integer, minimum number of major iterations (Algorithm 2.1, Hickernell and Yuan (1997)) before the multi-start algorithm terminates, defaults to 250. |
mstart_minsp |
positive integer, minimum number of detected stationary points before the multi-start algorithm terminates, defaults to 1. |
... |
any additional arguments (currently not used). |
A list
with exactly twenty-one components:
maxiter
scale
solver
fdtype
factor_up
factor_down
avmax
h_df
h_fvv
xtol
ftol
gtol
mstart_n
mstart_p
mstart_q
mstart_r
mstart_s
mstart_tol
mstart_maxiter
mstart_maxstart
mstart_minsp
with meanings as explained under 'Arguments'.
ftol
is disabled in some versions of the GSL library.
M. Galassi et al., GNU Scientific Library Reference Manual (3rd Ed.), ISBN 0954612078.
Hickernell, F.J. and Yuan, Y. (1997) “A simple multistart algorithm for global optimization”, OR Transactions, Vol. 1 (2).
https://www.gnu.org/software/gsl/doc/html/nls.html#tunable-parameters
## default tuning parameters gsl_nls_control()
## default tuning parameters gsl_nls_control()
Determine the nonlinear least-squares estimates of the parameters of a large
nonlinear model system using the gsl_multilarge_nlinear
module present in
the GNU Scientific Library (GSL).
gsl_nls_large(fn, ...) ## S3 method for class 'formula' gsl_nls_large( fn, data = parent.frame(), start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D", "cgst"), control = gsl_nls_control(), jac, fvv, trace = FALSE, subset, weights, na.action, model = FALSE, ... ) ## S3 method for class 'function' gsl_nls_large( fn, y, start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D", "cgst"), control = gsl_nls_control(), jac, fvv, trace = FALSE, weights, ... )
gsl_nls_large(fn, ...) ## S3 method for class 'formula' gsl_nls_large( fn, data = parent.frame(), start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D", "cgst"), control = gsl_nls_control(), jac, fvv, trace = FALSE, subset, weights, na.action, model = FALSE, ... ) ## S3 method for class 'function' gsl_nls_large( fn, y, start, algorithm = c("lm", "lmaccel", "dogleg", "ddogleg", "subspace2D", "cgst"), control = gsl_nls_control(), jac, fvv, trace = FALSE, weights, ... )
fn |
a nonlinear model defined either as a two-sided formula including variables and parameters, or as a function returning a numeric vector, with first argument the vector of parameters to be estimated. See the individual method descriptions below. |
data |
an optional data frame in which to evaluate the variables in |
y |
numeric response vector if |
start |
a named list or named numeric vector of starting estimates. |
algorithm |
character string specifying the algorithm to use. The following choices are supported:
|
control |
an optional list of control parameters to tune the least squares iterations and multistart algorithm.
See |
jac |
a function returning the |
fvv |
a function returning an |
trace |
logical value indicating if a trace of the iteration progress should be printed.
Default is |
subset |
an optional vector specifying a subset of observations to be used in the fitting process.
This argument is only used if |
weights |
an optional numeric vector of (fixed) weights. When present, the objective function is weighted least squares. |
na.action |
a function which indicates what should happen when the data contain |
model |
a logical value. If |
... |
additional arguments passed to the calls of |
If fn
is a formula
returns a list object of class nls
.
If fn
is a function
returns a list object of class gsl_nls
.
See the individual method descriptions for the structures of the returned lists and the generic functions
applicable to objects of both classes.
gsl_nls_large(formula)
: If fn
is a formula
, the returned list object is of classes gsl_nls
and nls
.
Therefore, all generic functions applicable to objects of class nls
, such as anova
, coef
, confint
,
deviance
, df.residual
, fitted
, formula
, logLik
, nobs
, predict
, print
, profile
,
residuals
, summary
, vcov
and weights
are also applicable to the returned list object.
In addition, a method confintd
is available for inference of derived parameters.
gsl_nls_large(function)
: If fn
is a function
, the first argument must be the vector of parameters and
the function should return a numeric vector containing the nonlinear model evaluations at
the provided parameter and predictor or covariate vectors. In addition, the argument y
needs to contain the numeric vector of observed responses, equal in length to the numeric
vector returned by fn
. The returned list object is (only) of class gsl_nls
.
Although the returned object is not of class nls
, the following generic functions remain
applicable for an object of class gsl_nls
: anova
, coef
, confint
, deviance
,
df.residual
, fitted
, formula
, logLik
, nobs
, predict
, print
,
residuals
, summary
, vcov
and weights
. In addition, a method confintd
is available for inference of derived parameters.
M. Galassi et al., GNU Scientific Library Reference Manual (3rd Ed.), ISBN 0954612078.
https://www.gnu.org/software/gsl/doc/html/nls.html
# Large NLS example # (https://www.gnu.org/software/gsl/doc/html/nls.html#large-nonlinear-least-squares-example) ## number of parameters p <- 250 ## model function f <- function(theta) { c(sqrt(1e-5) * (theta - 1), sum(theta^2) - 0.25) } ## jacobian function jac <- function(theta) { rbind(diag(sqrt(1e-5), nrow = length(theta)), 2 * t(theta)) } ## dense Levenberg-Marquardt gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values algorithm = "lm", ## algorithm jac = jac, ## jacobian control = list(maxiter = 250) ) ## dense Steihaug-Toint conjugate gradient gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values jac = jac, ## jacobian algorithm = "cgst" ## algorithm ) ## sparse Jacobian function jacsp <- function(theta) { rbind(Matrix::Diagonal(x = sqrt(1e-5), n = length(theta)), 2 * t(theta)) } ## sparse Levenberg-Marquardt gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values algorithm = "lm", ## algorithm jac = jacsp, ## sparse jacobian control = list(maxiter = 250) ) ## sparse Steihaug-Toint conjugate gradient gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values jac = jacsp, ## sparse jacobian algorithm = "cgst" ## algorithm )
# Large NLS example # (https://www.gnu.org/software/gsl/doc/html/nls.html#large-nonlinear-least-squares-example) ## number of parameters p <- 250 ## model function f <- function(theta) { c(sqrt(1e-5) * (theta - 1), sum(theta^2) - 0.25) } ## jacobian function jac <- function(theta) { rbind(diag(sqrt(1e-5), nrow = length(theta)), 2 * t(theta)) } ## dense Levenberg-Marquardt gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values algorithm = "lm", ## algorithm jac = jac, ## jacobian control = list(maxiter = 250) ) ## dense Steihaug-Toint conjugate gradient gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values jac = jac, ## jacobian algorithm = "cgst" ## algorithm ) ## sparse Jacobian function jacsp <- function(theta) { rbind(Matrix::Diagonal(x = sqrt(1e-5), n = length(theta)), 2 * t(theta)) } ## sparse Levenberg-Marquardt gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values algorithm = "lm", ## algorithm jac = jacsp, ## sparse jacobian control = list(maxiter = 250) ) ## sparse Steihaug-Toint conjugate gradient gsl_nls_large( fn = f, ## model y = rep(0, p + 1), ## (dummy) responses start = 1:p, ## start values jac = jacsp, ## sparse jacobian algorithm = "cgst" ## algorithm )
Returns leverage values (hat values) from a fitted "gsl_nls"
object based on the estimated
variance-covariance matrix of the model parameters.
## S3 method for class 'gsl_nls' hatvalues(model, ...)
## S3 method for class 'gsl_nls' hatvalues(model, ...)
model |
An object inheriting from class |
... |
At present no optional arguments are used. |
Numeric vector of leverage values similar to hatvalues
.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) hatvalues(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) hatvalues(obj)
Returns the model log-likelihood of a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' logLik(object, REML = FALSE, ...)
## S3 method for class 'gsl_nls' logLik(object, REML = FALSE, ...)
object |
An object inheriting from class |
REML |
logical value; included for compatibility reasons only, should not be used. |
... |
At present no optional arguments are used. |
Numeric object of class "logLik"
identical to logLik
.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) logLik(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) logLik(obj)
Returns an overview of 59 NLS test problems originating primarily from the NIST Statistical Reference Datasets (StRD) archive; Bates and Watts (1988); and More, Garbow and Hillstrom (1981).
nls_test_list(fields = c("name", "class", "p", "n", "check"))
nls_test_list(fields = c("name", "class", "p", "n", "check"))
fields |
optional character vector to return a subset of columns in the data.frame. |
a data.frame with high-level information about the available test problems. The following columns are returned by default:
"name"
Name of the test problem for use in nls_test_problem
.
"class"
Either "formula"
if the model is defined as a formula or "function"
if defined as a function.
"p"
Default number of parameters in the test problem.
"n"
Default number of residuals in the test problem.
"check"
One of the following three options: (1) "p, n fixed"
if the listed values for p
and n
are the only ones possible;
(2) "p <= n free"
if the values for p
and n
are not fixed, but p
must be smaller or equal to n
;
(3) "p == n free"
if the values for p
and n
are not fixed, but p
must be equal to n
.
D.M. Bates and Watts, D.G. (1988). Nonlinear Regression Analysis and Its Applications, Wiley, ISBN: 0471816434.
J.J. Moré, Garbow, B.S. and Hillstrom, K.E. (1981). Testing unconstrained optimization software, ACM Transactions on Mathematical Software, 7(1), 17-41.
https://www.itl.nist.gov/div898/strd/nls/nls_main.shtml
https://people.math.sc.edu/Burkardt/f_src/test_nls/test_nls.html
## available test problems nls_test_list()
## available test problems nls_test_list()
Fetches the model definition and model data required to solve a single NLS test problem with gsl_nls
(or nls
if the model is defined as a formula). Use nls_test_list
to
list the names of the available NLS test problems.
nls_test_problem(name, p = NA, n = NA)
nls_test_problem(name, p = NA, n = NA)
name |
Name of the NLS test problem, as returned in the |
p |
The number of parameters in the NLS test problem constrained by the |
n |
The number of residuals in the NLS test problem constrained by the |
If the model is defined as a formula, a list of class "nls_test_formula"
with elements:
data
a data.frame with n
rows containing the data (predictor and response values) used in the regression problem.
fn
a formula defining the test problem model.
start
a named vector of length p
with suggested starting values for the parameters.
target
a named vector of length p
with the certified target values for the parameters corresponding to the
best-available solutions.
If the model is defined as a function, a list of class "nls_test_function"
with elements:
fn
a function defining the test problem model. fn
takes a vector of parameters
of length p
as its first argument and returns a numeric vector of length n
.
fn
y
a numeric vector of length n
containing the response values.
start
a numeric named vector of length p
with suggested starting values for the parameters.
jac
a function defining the analytic Jacobian matrix of the model fn
. jac
takes a vector of parameters of length p
as its first argument and returns an n
by p
dimensional matrix.
target
a numeric named vector of length p
with the certified target values for the parameters, or a vector of
NA
's if no target solution is available.
For several problems the optimal least-squares objective of the target solution can be obtained at multiple different parameter locations.
D.M. Bates and Watts, D.G. (1988). Nonlinear Regression Analysis and Its Applications, Wiley, ISBN: 0471816434.
J.J. Moré, Garbow, B.S. and Hillstrom, K.E. (1981). Testing unconstrained optimization software, ACM Transactions on Mathematical Software, 7(1), 17-41.
https://www.itl.nist.gov/div898/strd/nls/nls_main.shtml
https://people.math.sc.edu/Burkardt/f_src/test_nls/test_nls.html
## example regression problem ratkowsky2 <- nls_test_problem(name = "Ratkowsky2") with(ratkowsky2, gsl_nls( fn = fn, data = data, start = start ) ) ## example optimization problem rosenbrock <- nls_test_problem(name = "Rosenbrock") with(rosenbrock, gsl_nls( fn = fn, y = y, start = start, jac = jac ) )
## example regression problem ratkowsky2 <- nls_test_problem(name = "Ratkowsky2") with(ratkowsky2, gsl_nls( fn = fn, data = data, start = start ) ) ## example optimization problem rosenbrock <- nls_test_problem(name = "Rosenbrock") with(rosenbrock, gsl_nls( fn = fn, y = y, start = start, jac = jac ) )
Returns the number of observations from a "gsl_nls"
object.
## S3 method for class 'gsl_nls' nobs(object, ...)
## S3 method for class 'gsl_nls' nobs(object, ...)
object |
An object inheriting from class |
... |
At present no optional arguments are used. |
Integer number of observations similar to nobs
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) nobs(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) nobs(obj)
Returns predicted values for the expected response from a fitted "gsl_nls"
object.
Asymptotic confidence or prediction (tolerance) intervals at a given level
can be evaluated
by specifying the appropriate interval
argument.
## S3 method for class 'gsl_nls' predict( object, newdata, scale = NULL, interval = c("none", "confidence", "prediction"), level = 0.95, ... )
## S3 method for class 'gsl_nls' predict( object, newdata, scale = NULL, interval = c("none", "confidence", "prediction"), level = 0.95, ... )
object |
An object inheriting from class |
newdata |
A named list or data.frame in which to look for variables with which to predict. If
|
scale |
A numeric scalar or vector. If it is set, it is used as the residual standard deviation (or vector of residual standard deviations) in the computation of the standard errors, otherwise this information is extracted from the model fit. |
interval |
A character string indicating if confidence or prediction (tolerance) intervals at the specified level should be returned. |
level |
A numeric scalar between 0 and 1 giving the confidence level for the intervals (if any) to be calculated. |
... |
At present no optional arguments are used. |
If interval = "none"
(default), a vector of predictions for the mean response. Otherwise,
a matrix with columns fit
, lwr
and upr
. The first column (fit
) contains
predictions for the mean response. The other two columns contain lower (lwr
) and upper (upr
)
confidence or prediction bounds at the specified level
.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) predict(obj) predict(obj, newdata = data.frame(x = 1:(2 * n) / n)) predict(obj, interval = "confidence") predict(obj, interval = "prediction", level = 0.99)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) predict(obj) predict(obj, newdata = data.frame(x = 1:(2 * n) / n)) predict(obj, interval = "confidence") predict(obj, interval = "prediction", level = 0.99)
Returns the model residuals from a fitted "gsl_nls"
object.
resid
can also be used as an alias.
## S3 method for class 'gsl_nls' residuals(object, type = c("response", "pearson"), ...)
## S3 method for class 'gsl_nls' residuals(object, type = c("response", "pearson"), ...)
object |
An object inheriting from class |
type |
character; if |
... |
At present no optional arguments are used. |
Numeric vector of model residuals similar to residuals
.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) residuals(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) residuals(obj)
Returns the estimated (unweighted) residual standard deviation of a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' sigma(object, ...)
## S3 method for class 'gsl_nls' sigma(object, ...)
object |
An object inheriting from class |
... |
At present no optional arguments are used. |
Numeric residual standard deviation value similar to sigma
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) sigma(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) sigma(obj)
Returns the model summary for a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' summary(object, correlation = FALSE, symbolic.cor = FALSE, ...)
## S3 method for class 'gsl_nls' summary(object, correlation = FALSE, symbolic.cor = FALSE, ...)
object |
An object inheriting from class |
correlation |
logical; if |
symbolic.cor |
logical; if |
... |
At present no optional arguments are used. |
List object of class "summary.nls"
identical to summary.nls
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) summary(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) summary(obj)
Returns the estimated variance-covariance matrix of the model parameters
from a fitted "gsl_nls"
object.
## S3 method for class 'gsl_nls' vcov(object, ...)
## S3 method for class 'gsl_nls' vcov(object, ...)
object |
An object inheriting from class |
... |
At present no optional arguments are used. |
A matrix containing the estimated covariances between the parameter estimates similar
to vcov
with row and column names corresponding to the parameter names given by coef.gsl_nls
.
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) vcov(obj)
## data set.seed(1) n <- 25 xy <- data.frame( x = (1:n) / n, y = 2.5 * exp(-1.5 * (1:n) / n) + rnorm(n, sd = 0.1) ) ## model obj <- gsl_nls(fn = y ~ A * exp(-lam * x), data = xy, start = c(A = 1, lam = 1)) vcov(obj)