nloptr-package

R interface to NLopt

R interface to NLopt

nloptr is an R interface to NLopt, a free/open-source library for nonlinear optimization started by Steven G. Johnson, providing a common interface for a number of different free optimization routines available online as well as original implementations of various other algorithms. The NLopt library is available under the GNU Lesser General Public License (LGPL), and the copyrights are owned by a variety of authors. Most of the information here has been taken from the NLopt website, where more details are available. package

Details

NLopt addresses general nonlinear optimization problems of the form:

minf(x)xRnminf(x)xinRn \min f(x)\quad x\in R^nmin f(x) x in R^n s.t. g(x)0h(x)=0lbxubs.t.g(x)<=0h(x)=0lb<=x<=ub \textrm{s.t. }\\ g(x) \leq 0\\ h(x) = 0\\ lb \leq x \leq ubs.t. g(x) <= 0 h(x) = 0 lb <= x <= ub

where f(x)f(x) is the objective function to be minimized and xx

represents the nn optimization parameters. This problem may optionally be subject to the bound constraints (also called box constraints), lblb

and ubub. For partially or totally unconstrained problems the bounds can take -Inf or Inf. One may also optionally have mm

nonlinear inequality constraints (sometimes called a nonlinear programming problem), which can be specified in g(x)g(x), and equality constraints that can be specified in h(x)h(x). Note that not all of the algorithms in NLopt can handle constraints.

An optimization problem can be solved with the general nloptr

interface, or using one of the wrapper functions for the separate algorithms; auglag, bobyqa, ccsaq, cobyla, crs2lm, direct, directL, isres, lbfgs, mlsl, mma, neldermead, newuoa, sbplx, slsqp, stogo, tnewton, varmetric.

Package:nloptr
Type:Package
Version:2.0.3
Date:2022-05-26
License:L-GPL >= 3

Note

See ?nloptr for more examples.

Examples

# Example problem, number 71 from the Hock-Schittkowsky test suite. # # \min_{x} x1 * x4 * (x1 + x2 + x3) + x3 # s.t. # x1 * x2 * x3 * x4 >= 25 # x1 ^ 2 + x2 ^ 2 + x3 ^ 2 + x4 ^ 2 = 40 # 1 <= x1, x2, x3, x4 <= 5 # # we re-write the inequality as # 25 - x1 * x2 * x3 * x4 <= 0 # # and the equality as # x1 ^ 2 + x2 ^ 2 + x3 ^ 2 + x4 ^ 2 - 40 = 0 # # x0 = (1, 5, 5, 1) # # optimal solution = (1.000000, 4.742999, 3.821151, 1.379408) library('nloptr') # # f(x) = x1 * x4 * (x1 + x2 + x3) + x3 # eval_f <- function(x) { list("objective" = x[1] * x[4] * (x[1] + x[2] + x[3]) + x[3], "gradient" = c(x[1] * x[4] + x[4] * (x[1] + x[2] + x[3]), x[1] * x[4], x[1] * x[4] + 1.0, x[1] * (x[1] + x[2] + x[3]))) } # constraint functions # inequalities eval_g_ineq <- function(x) { constr <- c(25 - x[1] * x[2] * x[3] * x[4]) grad <- c(-x[2] * x[3] * x[4], -x[1] * x[3] * x[4], -x[1] * x[2] * x[4], -x[1] * x[2] * x[3] ) list("constraints" = constr, "jacobian" = grad) } # equalities eval_g_eq <- function(x) { constr <- c(x[1] ^ 2 + x[2] ^ 2 + x[3] ^ 2 + x[4] ^ 2 - 40) grad <- c(2.0 * x[1], 2.0 * x[2], 2.0 * x[3], 2.0 * x[4]) list("constraints" = constr, "jacobian" = grad) } # initial values x0 <- c(1, 5, 5, 1) # lower and upper bounds of control lb <- c(1, 1, 1, 1) ub <- c(5, 5, 5, 5) local_opts <- list("algorithm" = "NLOPT_LD_MMA", "xtol_rel" = 1.0e-7) opts <- list("algorithm" = "NLOPT_LD_AUGLAG", "xtol_rel" = 1.0e-7, "maxeval" = 1000, "local_opts" = local_opts) res <- nloptr(x0 = x0, eval_f = eval_f, lb = lb, ub = ub, eval_g_ineq = eval_g_ineq, eval_g_eq = eval_g_eq, opts = opts) print(res)

References

Steven G. Johnson, The NLopt nonlinear-optimization package, https://nlopt.readthedocs.io/en/latest/

See Also

optim nlm nlminb

Rsolnp::Rsolnp Rsolnp::solnp nloptr

auglag bobyqa ccsaq

cobyla crs2lm direct

directL isres lbfgs

mlsl mma neldermead

newuoa sbplx slsqp

stogo tnewton varmetric

Author(s)

Steven G. Johnson and others (C code)

Jelmer Ypma (R interface)

Hans W. Borchers (wrappers)