optim {stats} | R Documentation |
General-purpose Optimization
Description
General-purpose optimization based on Nelder–Mead, quasi-Newton and conjugate-gradient algorithms. It includes an option for box-constrained optimization and simulated annealing.
Usage
optim(par, fn, gr = NULL, ...,
method = c("Nelder-Mead", "BFGS", "CG", "L-BFGS-B", "SANN",
"Brent"),
lower = -Inf, upper = Inf,
control = list(), hessian = FALSE)
optimHess(par, fn, gr = NULL, ..., control = list())
Arguments
par |
Initial values for the parameters to be optimized over. |
fn |
A function to be minimized (or maximized), with first argument the vector of parameters over which minimization is to take place. It should return a scalar result. |
gr |
A function to return the gradient for the For the |
... |
Further arguments to be passed to |
method |
The method to be used. See ‘Details’. Can be abbreviated. |
lower , upper |
Bounds on the variables for the |
control |
a |
hessian |
Logical. Should a numerically differentiated Hessian matrix be returned? |
Details
Note that arguments after ...
must be matched exactly.
By default optim
performs minimization, but it will maximize
if control$fnscale
is negative. optimHess
is an
auxiliary function to compute the Hessian at a later stage if
hessian = TRUE
was forgotten.
The default method is an implementation of that of Nelder and Mead (1965), that uses only function values and is robust but relatively slow. It will work reasonably well for non-differentiable functions.
Method "BFGS"
is a quasi-Newton method (also known as a variable
metric algorithm), specifically that published simultaneously in 1970
by Broyden, Fletcher, Goldfarb and Shanno.
This uses function values
and gradients to build up a picture of the surface to be optimized.
Method "CG"
is a conjugate gradients method based on that by
Fletcher and Reeves (1964) (but with the option of
Polak–Ribiere or Beale–Sorenson updates).
Conjugate gradient methods will generally
be more fragile than the BFGS method, but as they do not store a
matrix they may be successful in much larger optimization problems.
Method "L-BFGS-B"
is that of Byrd et al. (1995) which
allows box constraints, that is each variable can be given a lower
and/or upper bound. The initial value must satisfy the constraints.
This uses a limited-memory modification of the BFGS quasi-Newton
method. If non-trivial bounds are supplied, this method will be
selected, with a warning.
Nocedal and Wright (1999) is a comprehensive reference for the previous three methods.
Method "SANN"
is by default a variant of simulated annealing
given in Belisle (1992). Simulated-annealing belongs to the class of
stochastic global optimization methods. It uses only function values
but is relatively slow. It will also work for non-differentiable
functions. This implementation uses the Metropolis function for the
acceptance probability. By default the next candidate point is
generated from a Gaussian Markov kernel with scale proportional to the
actual temperature. If a function to generate a new candidate point is
given, method "SANN"
can also be used to solve combinatorial
optimization problems. Temperatures are decreased according to the
logarithmic cooling schedule as given in
Belisle (1992, p. 890);
specifically, the temperature is set to
temp / log(((t-1) %/% tmax)*tmax + exp(1))
, where t
is
the current iteration step and temp
and tmax
are
specifiable via control
, see below. Note that the
"SANN"
method depends critically on the settings of the control
parameters. It is not a general-purpose method but can be very useful
in getting to a good value on a very rough surface.
Method "Brent"
is for one-dimensional problems only, using
optimize()
. It can be useful in cases where
optim()
is used inside other functions where only method
can be specified, such as in mle
from package stats4.
Function fn
can return NA
or Inf
if the function
cannot be evaluated at the supplied value, but the initial value must
have a computable finite value of fn
.
(Except for method "L-BFGS-B"
where the values should always be
finite.)
optim
can be used recursively, and for a single parameter
as well as many. It also accepts a zero-length par
, and just
evaluates the function with that argument.
The control
argument is a list that can supply any of the
following components:
trace
Non-negative integer. If positive, tracing information on the progress of the optimization is produced. Higher values may produce more tracing information: for method
"L-BFGS-B"
there are six levels of tracing. (To understand exactly what these do see the source code: higher levels give more detail.)fnscale
An overall scaling to be applied to the value of
fn
andgr
during optimization. If negative, turns the problem into a maximization problem. Optimization is performed onfn(par)/fnscale
.parscale
A vector of scaling values for the parameters. Optimization is performed on
par/parscale
and these should be comparable in the sense that a unit change in any element produces about a unit change in the scaled value. Not used (nor needed) formethod = "Brent"
.ndeps
A vector of step sizes for the finite-difference approximation to the gradient, on
par/parscale
scale. Defaults to1e-3
.maxit
The maximum number of iterations. Defaults to
100
for the derivative-based methods, and500
for"Nelder-Mead"
.For
"SANN"
maxit
gives the total number of function evaluations: there is no other stopping criterion. Defaults to10000
.abstol
The absolute convergence tolerance. Only useful for non-negative functions, as a tolerance for reaching zero.
reltol
Relative convergence tolerance. The algorithm stops if it is unable to reduce the value by a factor of
reltol * (abs(val) + reltol)
at a step. Defaults tosqrt(.Machine$double.eps)
, typically about1e-8
.alpha
,beta
,gamma
Scaling parameters for the
"Nelder-Mead"
method.alpha
is the reflection factor (default 1.0),beta
the contraction factor (0.5) andgamma
the expansion factor (2.0).REPORT
The frequency of reports for the
"BFGS"
,"L-BFGS-B"
and"SANN"
methods ifcontrol$trace
is positive. Defaults to every 10 iterations for"BFGS"
and"L-BFGS-B"
, or every 100 temperatures for"SANN"
.warn.1d.NelderMead
a
logical
indicating if the (default)"Nelder-Mead"
method should signal a warning when used for one-dimensional minimization. As the warning is sometimes inappropriate, you can suppress it by setting this option to false.type
for the conjugate-gradients method. Takes value
1
for the Fletcher–Reeves update,2
for Polak–Ribiere and3
for Beale–Sorenson.lmm
is an integer giving the number of BFGS updates retained in the
"L-BFGS-B"
method, It defaults to5
.factr
controls the convergence of the
"L-BFGS-B"
method. Convergence occurs when the reduction in the objective is within this factor of the machine tolerance. Default is1e7
, that is a tolerance of about1e-8
.pgtol
helps control the convergence of the
"L-BFGS-B"
method. It is a tolerance on the projected gradient in the current search direction. This defaults to zero, when the check is suppressed.temp
controls the
"SANN"
method. It is the starting temperature for the cooling schedule. Defaults to10
.tmax
is the number of function evaluations at each temperature for the
"SANN"
method. Defaults to10
.
Any names given to par
will be copied to the vectors passed to
fn
and gr
. Note that no other attributes of par
are copied over.
The parameter vector passed to fn
has special semantics and may
be shared between calls: the function should not change or copy it.
Value
For optim
, a list with components:
par |
The best set of parameters found. |
value |
The value of |
counts |
A two-element integer vector giving the number of calls
to |
convergence |
An integer code.
|
message |
A character string giving any additional information
returned by the optimizer, or |
hessian |
Only if argument |
For optimHess
, the description of the hessian
component
applies.
Note
optim
will work with one-dimensional par
s, but the
default method does not work well (and will warn). Method
"Brent"
uses optimize
and needs bounds to be available;
"BFGS"
often works well enough if not.
Source
The code for methods "Nelder-Mead"
, "BFGS"
and
"CG"
was based originally on Pascal code in Nash (1990) that was
translated by p2c
and then hand-optimized. Dr Nash has agreed
that the code can be made freely available.
The code for method "L-BFGS-B"
is based on Fortran code by Zhu,
Byrd, Lu-Chen and Nocedal obtained from Netlib (file
‘opt/lbfgs_bcm.shar’: another version is in ‘toms/778’).
The code for method "SANN"
was contributed by A. Trapletti.
References
Belisle, C. J. P. (1992).
Convergence theorems for a class of simulated annealing algorithms on
R^d
.
Journal of Applied Probability, 29, 885–895.
doi:10.2307/3214721.
Byrd, R. H., Lu, P., Nocedal, J. and Zhu, C. (1995). A limited memory algorithm for bound constrained optimization. SIAM Journal on Scientific Computing, 16, 1190–1208. doi:10.1137/0916069.
Fletcher, R. and Reeves, C. M. (1964). Function minimization by conjugate gradients. Computer Journal 7, 148–154. doi:10.1093/comjnl/7.2.149.
Nash, J. C. (1990). Compact Numerical Methods for Computers. Linear Algebra and Function Minimisation. Adam Hilger.
Nelder, J. A. and Mead, R. (1965). A simplex algorithm for function minimization. Computer Journal, 7, 308–313. doi:10.1093/comjnl/7.4.308.
Nocedal, J. and Wright, S. J. (1999). Numerical Optimization. Springer.
See Also
optimize
for one-dimensional minimization and
constrOptim
for constrained optimization.
Examples
require(graphics)
fr <- function(x) { ## Rosenbrock Banana function
x1 <- x[1]
x2 <- x[2]
100 * (x2 - x1 * x1)^2 + (1 - x1)^2
}
grr <- function(x) { ## Gradient of 'fr'
x1 <- x[1]
x2 <- x[2]
c(-400 * x1 * (x2 - x1 * x1) - 2 * (1 - x1),
200 * (x2 - x1 * x1))
}
optim(c(-1.2,1), fr)
(res <- optim(c(-1.2,1), fr, grr, method = "BFGS"))
optimHess(res$par, fr, grr)
optim(c(-1.2,1), fr, NULL, method = "BFGS", hessian = TRUE)
## These do not converge in the default number of steps
optim(c(-1.2,1), fr, grr, method = "CG")
optim(c(-1.2,1), fr, grr, method = "CG", control = list(type = 2))
optim(c(-1.2,1), fr, grr, method = "L-BFGS-B")
flb <- function(x)
{ p <- length(x); sum(c(1, rep(4, p-1)) * (x - c(1, x[-p])^2)^2) }
## 25-dimensional box constrained
optim(rep(3, 25), flb, NULL, method = "L-BFGS-B",
lower = rep(2, 25), upper = rep(4, 25)) # par[24] is *not* at boundary
## "wild" function , global minimum at about -15.81515
fw <- function (x)
10*sin(0.3*x)*sin(1.3*x^2) + 0.00001*x^4 + 0.2*x+80
plot(fw, -50, 50, n = 1000, main = "optim() minimising 'wild function'")
res <- optim(50, fw, method = "SANN",
control = list(maxit = 20000, temp = 20, parscale = 20))
res
## Now improve locally {typically only by a small bit}:
(r2 <- optim(res$par, fw, method = "BFGS"))
points(r2$par, r2$value, pch = 8, col = "red", cex = 2)
## Combinatorial optimization: Traveling salesman problem
library(stats) # normally loaded
eurodistmat <- as.matrix(eurodist)
distance <- function(sq) { # Target function
sq2 <- embed(sq, 2)
sum(eurodistmat[cbind(sq2[,2], sq2[,1])])
}
genseq <- function(sq) { # Generate new candidate sequence
idx <- seq(2, NROW(eurodistmat)-1)
changepoints <- sample(idx, size = 2, replace = FALSE)
tmp <- sq[changepoints[1]]
sq[changepoints[1]] <- sq[changepoints[2]]
sq[changepoints[2]] <- tmp
sq
}
sq <- c(1:nrow(eurodistmat), 1) # Initial sequence: alphabetic
distance(sq)
# rotate for conventional orientation
loc <- -cmdscale(eurodist, add = TRUE)$points
x <- loc[,1]; y <- loc[,2]
s <- seq_len(nrow(eurodistmat))
tspinit <- loc[sq,]
plot(x, y, type = "n", asp = 1, xlab = "", ylab = "",
main = "initial solution of traveling salesman problem", axes = FALSE)
arrows(tspinit[s,1], tspinit[s,2], tspinit[s+1,1], tspinit[s+1,2],
angle = 10, col = "green")
text(x, y, labels(eurodist), cex = 0.8)
set.seed(123) # chosen to get a good soln relatively quickly
res <- optim(sq, distance, genseq, method = "SANN",
control = list(maxit = 30000, temp = 2000, trace = TRUE,
REPORT = 500))
res # Near optimum distance around 12842
tspres <- loc[res$par,]
plot(x, y, type = "n", asp = 1, xlab = "", ylab = "",
main = "optim() 'solving' traveling salesman problem", axes = FALSE)
arrows(tspres[s,1], tspres[s,2], tspres[s+1,1], tspres[s+1,2],
angle = 10, col = "red")
text(x, y, labels(eurodist), cex = 0.8)
## 1-D minimization: "Brent" or optimize() being preferred.. but NM may be ok and "unavoidable",
## ---------------- so we can suppress the check+warning :
system.time(rO <- optimize(function(x) (x-pi)^2, c(0, 10)))
system.time(ro <- optim(1, function(x) (x-pi)^2, control=list(warn.1d.NelderMead = FALSE)))
rO$minimum - pi # 0 (perfect), on one platform
ro$par - pi # ~= 1.9e-4 on one platform
utils::str(ro)