Title: | Unconstrained Optimization using the Subplex Algorithm |
---|---|
Description: | The subplex algorithm for unconstrained optimization, developed by Tom Rowan. |
Authors: | Aaron A. King [aut, trl, cre] , Tom Rowan [aut] |
Maintainer: | Aaron A. King <[email protected]> |
License: | GPL-3 |
Version: | 1.9 |
Built: | 2024-11-06 15:19:17 UTC |
Source: | https://github.com/kingaa/subplex |
The subplex package implements Tom Rowan's subspace-searching simplex algorithm for unconstrained minimization of a function.
Subplex is a subspace-searching simplex method for the unconstrained optimization of general multivariate functions. Like the Nelder-Mead simplex method it generalizes, the subplex method is well suited for optimizing noisy objective functions. The number of function evaluations required for convergence typically increases only linearly with the problem size, so for most applications the subplex method is much more efficient than the simplex method.
Subplex was written in FORTRAN by Tom Rowan (Oak Ridge National Laboratory). The FORTRAN source code is maintained on the netlib repository at https://www.netlib.org/opt/subplex.tgz.
Aaron A. King (kingaa at umich dot edu)
T. Rowan, "Functional Stability Analysis of Numerical Algorithms", Ph.D. thesis, Department of Computer Sciences, University of Texas at Austin, 1990.
subplex
minimizes a function.
subplex(par, fn, control = list(), hessian = FALSE, ...)
subplex(par, fn, control = list(), hessian = FALSE, ...)
par |
Initial guess of the parameters to be optimized over. |
fn |
The function to be minimized. Its first argument must be the vector of parameters to be optimized over. It should return a scalar result. |
control |
A list of control parameters, consisting of some or all of the following:
|
hessian |
If |
... |
Additional arguments to be passed to the function |
The convergence codes are as follows:
invalid input
number of function evaluations needed exceeds maxnfe
success: tolerance tol
satisfied
limit of machine precision reached
For more details, see the source code.
subplex
returns a list containing the following:
par |
Estimated parameters that minimize the function. |
value |
Minimized value of the function. |
count |
Number of function evaluations required. |
convergence |
Convergence code (see Details). |
message |
A character string giving a diagnostic message from the optimizer, or 'NULL'. |
hessian |
Hessian matrix. |
Aaron A. King [email protected]
T. Rowan, “Functional Stability Analysis of Numerical Algorithms”, Ph.D. thesis, Department of Computer Sciences, University of Texas at Austin, 1990.
ripple <- function (x) { r <- sqrt(sum(x^2)) 1-exp(-r^2)*cos(10*r)^2 } subplex(par=c(1),fn=ripple,hessian=TRUE) subplex(par=c(0.1,3),fn=ripple,hessian=TRUE) subplex(par=c(0.1,3,2),fn=ripple,hessian=TRUE) ## Rosenbrock Banana function rosen <- function (x) { x1 <- x[1] x2 <- x[2] 100*(x2-x1*x1)^2+(1-x1)^2 } subplex(par=c(11,-33),fn=rosen) ## Rosenbrock Banana function (using names) rosen <- function (x, g = 0, h = 0) { x1 <- x['a'] x2 <- x['b']-h 100*(x2-x1*x1)^2+(1-x1)^2+g } subplex(par=c(b=11,a=-33),fn=rosen,h=22,control=list(abstol=1e-9,parscale=5),hessian=TRUE)
ripple <- function (x) { r <- sqrt(sum(x^2)) 1-exp(-r^2)*cos(10*r)^2 } subplex(par=c(1),fn=ripple,hessian=TRUE) subplex(par=c(0.1,3),fn=ripple,hessian=TRUE) subplex(par=c(0.1,3,2),fn=ripple,hessian=TRUE) ## Rosenbrock Banana function rosen <- function (x) { x1 <- x[1] x2 <- x[2] 100*(x2-x1*x1)^2+(1-x1)^2 } subplex(par=c(11,-33),fn=rosen) ## Rosenbrock Banana function (using names) rosen <- function (x, g = 0, h = 0) { x1 <- x['a'] x2 <- x['b']-h 100*(x2-x1*x1)^2+(1-x1)^2+g } subplex(par=c(b=11,a=-33),fn=rosen,h=22,control=list(abstol=1e-9,parscale=5),hessian=TRUE)