OptimizationLBFGSB.jl

OptimizationLBFGSB.jl is a package that wraps the L-BFGS-B fortran routine via the LBFGSB.jl package.

Installation

To use this package, install the OptimizationLBFGSB package:

using Pkg
Pkg.add("OptimizationLBFGSB")

Methods

  • LBFGSB: The popular quasi-Newton method that leverages limited memory BFGS approximation of the inverse of the Hessian. It directly supports box-constraints.

    This can also handle arbitrary non-linear constraints through an Augmented Lagrangian method with bounds constraints described in 17.4 of Numerical Optimization by Nocedal and Wright. Thus serving as a general-purpose nonlinear optimization solver.

OptimizationLBFGSB.LBFGSBType
struct LBFGSB

L-BFGS-B Nonlinear Optimization Code from LBFGSB. It is a quasi-Newton optimization algorithm that supports bounds.

References

  • R. H. Byrd, P. Lu and J. Nocedal. A Limited Memory Algorithm for Bound Constrained Optimization, (1995), SIAM Journal on Scientific and Statistical Computing , 16, 5, pp. 1190-1208.
  • C. Zhu, R. H. Byrd and J. Nocedal. L-BFGS-B: Algorithm 778: L-BFGS-B, FORTRAN routines for large scale bound constrained optimization (1997), ACM Transactions on Mathematical Software, Vol 23, Num. 4, pp. 550 - 560.
  • J.L. Morales and J. Nocedal. L-BFGS-B: Remark on Algorithm 778: L-BFGS-B, FORTRAN routines for large scale bound constrained optimization (2011), to appear in ACM Transactions on Mathematical Software.
source

Examples

Unconstrained rosenbrock problem

using OptimizationBase, OptimizationLBFGSB, ADTypes, Zygote

rosenbrock(x, p) = (p[1] - x[1])^2 + p[2] * (x[2] - x[1]^2)^2
x0 = zeros(2)
p = [1.0, 100.0]

optf = OptimizationFunction(rosenbrock, ADTypes.AutoZygote())
prob = OptimizationProblem(optf, x0, p)
sol = solve(prob, LBFGSB())
retcode: Success
u: 2-element Vector{Float64}:
 0.9999997057368228
 0.999999398151528

With nonlinear and bounds constraints

function con2_c(res, x, p)
    res .= [x[1]^2 + x[2]^2, (x[2] * sin(x[1]) + x[1]) - 5]
end

optf = OptimizationFunction(rosenbrock, ADTypes.AutoZygote(), cons = con2_c)
prob = OptimizationProblem(optf, x0, p, lcons = [1.0, -Inf],
    ucons = [1.0, 0.0], lb = [-1.0, -1.0],
    ub = [1.0, 1.0])
res = solve(prob, LBFGSB(), maxiters = 100)
retcode: Success
u: 2-element Vector{Float64}:
 0.783397417853095
 0.6215211044097776