Solvers

Solver list

Problem typeSolvers
Unconstrained Nonlinear Optimization Problemlbfgs, tron, trunk, R2, fomo
Unconstrained Nonlinear Least Squarestrunk, tron
Bound-constrained Nonlinear Optimization Problemtron
Bound-constrained Nonlinear Least Squarestron

Solver list

JSOSolvers.lbfgsFunction
lbfgs(nlp; kwargs...)

An implementation of a limited memory BFGS line-search method for unconstrained minimization.

For advanced usage, first define a LBFGSSolver to preallocate the memory used in the algorithm, and then call solve!.

solver = LBFGSSolver(nlp; mem::Int = 5)
solve!(solver, nlp; kwargs...)

Arguments

  • nlp::AbstractNLPModel{T, V} represents the model to solve, see NLPModels.jl.

The keyword arguments may include

  • x::V = nlp.meta.x0: the initial guess.
  • mem::Int = 5: algorithm parameter, see LBFGSParameterSet.
  • atol::T = √eps(T): absolute tolerance.
  • rtol::T = √eps(T): relative tolerance, the algorithm stops when ‖∇f(xᵏ)‖ ≤ atol + rtol * ‖∇f(x⁰)‖.
  • callback: function called at each iteration, see Callback section.
  • max_eval::Int = -1: maximum number of objective function evaluations.
  • max_time::Float64 = 30.0: maximum time limit in seconds.
  • max_iter::Int = typemax(Int): maximum number of iterations.
  • τ₁::T = T(0.9999): algorithm parameter, see LBFGSParameterSet.
  • bk_max:: Int = 25: algorithm parameter, see LBFGSParameterSet.
  • verbose::Int = 0: if > 0, display iteration details every verbose iteration.
  • verbose_subsolver::Int = 0: if > 0, display iteration information every verbose_subsolver iteration of the subsolver.

Output

The returned value is a GenericExecutionStats, see SolverCore.jl.

Examples

using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3));
stats = lbfgs(nlp)

# output

"Execution stats: first-order stationary"
using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3));
solver = LBFGSSolver(nlp; mem = 5);
stats = solve!(solver, nlp)

# output

"Execution stats: first-order stationary"
source
JSOSolvers.tronFunction
tron(nlp; kwargs...)

A pure Julia implementation of a trust-region solver for bound-constrained optimization:

    min f(x)    s.t.    ℓ ≦ x ≦ u

For advanced usage, first define a TronSolver to preallocate the memory used in the algorithm, and then call solve!:

solver = TronSolver(nlp; kwargs...)
solve!(solver, nlp; kwargs...)

Arguments

  • nlp::AbstractNLPModel{T, V} represents the model to solve, see NLPModels.jl.

The keyword arguments may include

  • x::V = nlp.meta.x0: the initial guess.
  • μ₀::T = T(1 / 100): algorithm parameter, see TRONParameterSet.
  • μ₁::T = T(1): algorithm parameter, see TRONParameterSet.
  • σ::T = T(10): algorithm parameter, see TRONParameterSet.
  • max_eval::Int = -1: maximum number of objective function evaluations.
  • max_time::Float64 = 30.0: maximum time limit in seconds.
  • max_iter::Int = typemax(Int): maximum number of iterations.
  • max_cgiter::Int = 50: subproblem's iteration limit.
  • use_only_objgrad::Bool = false: If true, the algorithm uses only the function objgrad instead of obj and grad.
  • cgtol::T = T(0.1): subproblem tolerance.
  • atol::T = √eps(T): absolute tolerance.
  • rtol::T = √eps(T): relative tolerance, the algorithm stops when ‖x - Proj(x - ∇f(xᵏ))‖ ≤ atol + rtol * ‖∇f(x⁰)‖. Proj denotes here the projection over the bounds.
  • callback: function called at each iteration, see Callback section.
  • verbose::Int = 0: if > 0, display iteration details every verbose iteration.
  • subsolver_verbose::Int = 0: if > 0, display iteration information every subsolver_verbose iteration of the subsolver.

The keyword arguments of TronSolver are passed to the TRONTrustRegion constructor.

Output

The value returned is a GenericExecutionStats, see SolverCore.jl.

References

TRON is described in

Chih-Jen Lin and Jorge J. Moré, *Newton's Method for Large Bound-Constrained
Optimization Problems*, SIAM J. Optim., 9(4), 1100–1127, 1999.
DOI: 10.1137/S1052623498345075

Examples

using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x), ones(3), zeros(3), 2 * ones(3));
stats = tron(nlp)
using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x), ones(3), zeros(3), 2 * ones(3));
solver = TronSolver(nlp);
stats = solve!(solver, nlp)
source
tron(nls; kwargs...)

A pure Julia implementation of a trust-region solver for bound-constrained nonlinear least-squares problems:

min ½‖F(x)‖²    s.t.    ℓ ≦ x ≦ u

For advanced usage, first define a TronSolverNLS to preallocate the memory used in the algorithm, and then call solve!: solver = TronSolverNLS(nls, subsolver::Symbol=:lsmr; kwargs...) solve!(solver, nls; kwargs...)

Arguments

  • nls::AbstractNLSModel{T, V} represents the model to solve, see NLPModels.jl.

The keyword arguments may include

  • x::V = nlp.meta.x0: the initial guess.
  • subsolver::Symbol = :lsmr: Krylov.jl method used as subproblem solver, see JSOSolvers.tronls_allowed_subsolvers for a list.
  • μ₀::T = T(1 / 100): algorithm parameter, see TRONLSParameterSet.
  • μ₁::T = T(1): algorithm parameter, see TRONLSParameterSet.
  • σ::T = T(10): algorithm parameter, see TRONLSParameterSet.
  • max_eval::Int = -1: maximum number of objective function evaluations.
  • max_time::Float64 = 30.0: maximum time limit in seconds.
  • max_iter::Int = typemax(Int): maximum number of iterations.
  • max_cgiter::Int = 50: subproblem iteration limit.
  • cgtol::T = T(0.1): subproblem tolerance.
  • atol::T = √eps(T): absolute tolerance.
  • rtol::T = √eps(T): relative tolerance, the algorithm stops when ‖x - Proj(x - ∇f(xᵏ))‖ ≤ atol + rtol * ‖∇f(x⁰)‖. Proj denotes here the projection over the bounds.
  • Fatol::T = √eps(T): absolute tolerance on the residual.
  • Frtol::T = eps(T): relative tolerance on the residual, the algorithm stops when ‖F(xᵏ)‖ ≤ Fatol + Frtol * ‖F(x⁰)‖.
  • callback: function called at each iteration, see Callback section.
  • verbose::Int = 0: if > 0, display iteration details every verbose iteration.
  • subsolver_verbose::Int = 0: if > 0, display iteration information every subsolver_verbose iteration of the subsolver.

The keyword arguments of TronSolverNLS are passed to the TRONTrustRegion constructor.

Output

The value returned is a GenericExecutionStats, see SolverCore.jl.

References

This is an adaptation for bound-constrained nonlinear least-squares problems of the TRON method described in

Chih-Jen Lin and Jorge J. Moré, *Newton's Method for Large Bound-Constrained
Optimization Problems*, SIAM J. Optim., 9(4), 1100–1127, 1999.
DOI: 10.1137/S1052623498345075

Examples

using JSOSolvers, ADNLPModels
F(x) = [x[1] - 1.0; 10 * (x[2] - x[1]^2)]
x0 = [-1.2; 1.0]
nls = ADNLSModel(F, x0, 2, zeros(2), 0.5 * ones(2))
stats = tron(nls)
using JSOSolvers, ADNLPModels
F(x) = [x[1] - 1.0; 10 * (x[2] - x[1]^2)]
x0 = [-1.2; 1.0]
nls = ADNLSModel(F, x0, 2, zeros(2), 0.5 * ones(2))
solver = TronSolverNLS(nls)
stats = solve!(solver, nls)
source
JSOSolvers.trunkFunction
trunk(nlp; kwargs...)

A trust-region solver for unconstrained optimization using exact second derivatives.

For advanced usage, first define a TrunkSolver to preallocate the memory used in the algorithm, and then call solve!:

solver = TrunkSolver(nlp, subsolver::Symbol = :cg)
solve!(solver, nlp; kwargs...)

Arguments

  • nlp::AbstractNLPModel{T, V} represents the model to solve, see NLPModels.jl.

The keyword arguments may include

  • subsolver_logger::AbstractLogger = NullLogger(): subproblem's logger.
  • x::V = nlp.meta.x0: the initial guess.
  • atol::T = √eps(T): absolute tolerance.
  • rtol::T = √eps(T): relative tolerance, the algorithm stops when ‖∇f(xᵏ)‖ ≤ atol + rtol * ‖∇f(x⁰)‖.
  • callback: function called at each iteration, see Callback section.
  • max_eval::Int = -1: maximum number of objective function evaluations.
  • max_time::Float64 = 30.0: maximum time limit in seconds.
  • max_iter::Int = typemax(Int): maximum number of iterations.
  • bk_max::Int = 10: algorithm parameter, see TRUNKParameterSet.
  • monotone::Bool = true: algorithm parameter, see TRUNKParameterSet.
  • nm_itmax::Int = 25: algorithm parameter, see TRUNKParameterSet.
  • verbose::Int = 0: if > 0, display iteration information every verbose iteration.
  • subsolver_verbose::Int = 0: if > 0, display iteration information every subsolver_verbose iteration of the subsolver.
  • M: linear operator that models a Hermitian positive-definite matrix of size n; passed to Krylov subsolvers.

Output

The returned value is a GenericExecutionStats, see SolverCore.jl.

References

This implementation follows the description given in

A. R. Conn, N. I. M. Gould, and Ph. L. Toint,
Trust-Region Methods, volume 1 of MPS/SIAM Series on Optimization.
SIAM, Philadelphia, USA, 2000.
DOI: 10.1137/1.9780898719857

The main algorithm follows the basic trust-region method described in Section 6. The backtracking linesearch follows Section 10.3.2. The nonmonotone strategy follows Section 10.1.3, Algorithm 10.1.2.

Examples

using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3))
stats = trunk(nlp)
using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3))
solver = TrunkSolver(nlp)
stats = solve!(solver, nlp)
source
trunk(nls; kwargs...)

A pure Julia implementation of a trust-region solver for nonlinear least-squares problems:

min ½‖F(x)‖²

For advanced usage, first define a TrunkSolverNLS to preallocate the memory used in the algorithm, and then call solve!:

solver = TrunkSolverNLS(nls, subsolver::Symbol = :lsmr)
solve!(solver, nls; kwargs...)

Arguments

  • nls::AbstractNLSModel{T, V} represents the model to solve, see NLPModels.jl.

The keyword arguments may include

  • x::V = nlp.meta.x0: the initial guess.
  • atol::T = √eps(T): absolute tolerance.
  • rtol::T = √eps(T): relative tolerance, the algorithm stops when ‖∇f(xᵏ)‖ ≤ atol + rtol * ‖∇f(x⁰)‖.
  • Fatol::T = √eps(T): absolute tolerance on the residual.
  • Frtol::T = eps(T): relative tolerance on the residual, the algorithm stops when ‖F(xᵏ)‖ ≤ Fatol + Frtol * ‖F(x⁰)‖.
  • callback: function called at each iteration, see Callback section.
  • max_eval::Int = -1: maximum number of objective function evaluations.
  • max_time::Float64 = 30.0: maximum time limit in seconds.
  • max_iter::Int = typemax(Int): maximum number of iterations.
  • bk_max::Int = 10: algorithm parameter, see TRUNKLSParameterSet.
  • monotone::Bool = true: algorithm parameter, see TRUNKLSParameterSet.
  • nm_itmax::Int = 25: algorithm parameter, see TRUNKLSParameterSet.
  • verbose::Int = 0: if > 0, display iteration details every verbose iteration.
  • subsolver_verbose::Int = 0: if > 0, display iteration information every subsolver_verbose iteration of the subsolver.

See JSOSolvers.trunkls_allowed_subsolvers for a list of available Krylov solvers.

Output

The value returned is a GenericExecutionStats, see SolverCore.jl.

References

This implementation follows the description given in

A. R. Conn, N. I. M. Gould, and Ph. L. Toint,
Trust-Region Methods, volume 1 of MPS/SIAM Series on Optimization.
SIAM, Philadelphia, USA, 2000.
DOI: 10.1137/1.9780898719857

The main algorithm follows the basic trust-region method described in Section 6. The backtracking linesearch follows Section 10.3.2. The nonmonotone strategy follows Section 10.1.3, Algorithm 10.1.2.

Examples

using JSOSolvers, ADNLPModels
F(x) = [x[1] - 1.0; 10 * (x[2] - x[1]^2)]
x0 = [-1.2; 1.0]
nls = ADNLSModel(F, x0, 2)
stats = trunk(nls)
using JSOSolvers, ADNLPModels
F(x) = [x[1] - 1.0; 10 * (x[2] - x[1]^2)]
x0 = [-1.2; 1.0]
nls = ADNLSModel(F, x0, 2)
solver = TrunkSolverNLS(nls)
stats = solve!(solver, nls)
source
JSOSolvers.R2Function
fo(nlp; kwargs...)
R2(nlp; kwargs...)
TR(nlp; kwargs...)

A First-Order (FO) model-based method for unconstrained optimization. Supports quadratic regularization and trust region method with linear model.

For advanced usage, first define a FomoSolver to preallocate the memory used in the algorithm, and then call solve!:

solver = FoSolver(nlp)
solve!(solver, nlp; kwargs...)

R2 and TR runs fo with the dedicated step_backend keyword argument.

Arguments

  • nlp::AbstractNLPModel{T, V} is the model to solve, see NLPModels.jl.

Keyword arguments

  • x::V = nlp.meta.x0: the initial guess.
  • atol::T = √eps(T): absolute tolerance.
  • rtol::T = √eps(T): relative tolerance: algorithm stops when ‖∇f(xᵏ)‖ ≤ atol + rtol * ‖∇f(x⁰)‖.
  • η1 = eps(T)^(1 // 4): algorithm parameter, see FOMOParameterSet.
  • η2 = T(95/100): algorithm parameter, see FOMOParameterSet.
  • γ1 = T(1/2): algorithm parameter, see FOMOParameterSet.
  • γ2 = T(2): algorithm parameter, see FOMOParameterSet.
  • αmax = 1/eps(T): algorithm parameter, see FOMOParameterSet.
  • max_eval::Int = -1: maximum number of evaluation of the objective function.
  • max_time::Float64 = 30.0: maximum time limit in seconds.
  • max_iter::Int = typemax(Int): maximum number of iterations.
  • M = 1 : algorithm parameter, see FOMOParameterSet.
  • verbose::Int = 0: if > 0, display iteration details every verbose iteration.
  • step_backend = r2_step(): algorithm parameter, see FOMOParameterSet.

Output

The value returned is a GenericExecutionStats, see SolverCore.jl.

Examples

using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3))
stats = fo(nlp) # run with step_backend = r2_step(), equivalent to R2(nlp)

# output

"Execution stats: first-order stationary"
using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3))
solver = FoSolver(nlp);
stats = solve!(solver, nlp)

# output

"Execution stats: first-order stationary"
source
JSOSolvers.fomoFunction
fomo(nlp; kwargs...)

A First-Order with MOmentum (FOMO) model-based method for unconstrained optimization. Supports quadratic regularization and trust region method with linear model.

Algorithm description

The step is computed along d = - (1-βmax) .* ∇f(xk) - βmax .* mk with mk the memory of past gradients (initialized at 0), and updated at each successful iteration as mk .= ∇f(xk) .* (1 - βmax) .+ mk .* βmax and βmax ∈ [0,β] chosen as to ensure d is gradient-related, i.e., the following 2 conditions are satisfied: (1-βmax) .* ∇f(xk) + βmax .* ∇f(xk)ᵀmk ≥ θ1 * ‖∇f(xk)‖² (1) ‖∇f(xk)‖ ≥ θ2 * ‖(1-βmax) . ∇f(xk) + βmax . mk‖ (2) In the nonmonotone case, (1) rewrites (1-βmax) .* ∇f(xk) + βmax .* ∇f(xk)ᵀmk + (fm - fk)/μk ≥ θ1 * ‖∇f(xk)‖², with fm the largest objective value over the last M successful iterations, and fk = f(xk).

Advanced usage

For advanced usage, first define a FomoSolver to preallocate the memory used in the algorithm, and then call solve!:

solver = FomoSolver(nlp)
solve!(solver, nlp; kwargs...)

No momentum: if the user does not whish to use momentum (β = 0), it is recommended to use the memory-optimized fo method.

Arguments

  • nlp::AbstractNLPModel{T, V} is the model to solve, see NLPModels.jl.

Keyword arguments

  • x::V = nlp.meta.x0: the initial guess.
  • atol::T = √eps(T): absolute tolerance.
  • rtol::T = √eps(T): relative tolerance: algorithm stops when ‖∇f(xᵏ)‖ ≤ atol + rtol * ‖∇f(x⁰)‖.
  • callback: function called at each iteration, see Callback section.
  • η1 = eps(T)^(1 // 4), η2 = T(95/100): step acceptance parameters.
  • γ1 = T(1/2), γ2 = T(2): regularization update parameters.
  • γ3 = T(1/2) : momentum factor βmax update parameter in case of unsuccessful iteration.
  • αmax = 1/eps(T): maximum step parameter for fomo algorithm.
  • max_eval::Int = -1: maximum number of objective evaluations.
  • max_time::Float64 = 30.0: maximum time limit in seconds.
  • max_iter::Int = typemax(Int): maximum number of iterations.
  • β = T(9/10) ∈ [0,1): target decay rate for the momentum.
  • θ1 = T(1/10): momentum contribution parameter for convergence condition (1).
  • θ2 = eps(T)^(1/3): momentum contribution parameter for convergence condition (2).
  • M = 1 : requires objective decrease over the M last iterates (nonmonotone context). M=1 implies monotone behaviour.
  • verbose::Int = 0: if > 0, display iteration details every verbose iteration.
  • step_backend = r2_step(): step computation mode. Options are r2_step() for quadratic regulation step and tr_step() for first-order trust-region.

Output

The value returned is a GenericExecutionStats, see SolverCore.jl.

Examples

fomo

using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3))
stats = fomo(nlp)

# output

"Execution stats: first-order stationary"
using JSOSolvers, ADNLPModels
nlp = ADNLPModel(x -> sum(x.^2), ones(3))
solver = FomoSolver(nlp);
stats = solve!(solver, nlp)

# output

"Execution stats: first-order stationary"
source