Reference
Contents
Index
CaNNOLeS.CaNNOLeSSolver
CaNNOLeS.ParamCaNNOLeS
CaNNOLeS._check_available_method
CaNNOLeS.cannoles
CaNNOLeS.check_nan_inf
CaNNOLeS.dual_scaling
CaNNOLeS.get_nnzh
CaNNOLeS.newton_system!
CaNNOLeS.optimality_check_small_residual!
CaNNOLeS.prepare_newton_system!
CaNNOLeS.solve_ldl!
CaNNOLeS.try_to_factorize
CaNNOLeS.update_newton_hessian!
CaNNOLeS.CaNNOLeSSolver
— Typecannoles(nls)
Implementation of a solver for Nonlinear Least Squares with nonlinear constraints.
$min f(x) = ¹/₂‖F(x)‖² s.t. c(x) = 0$
For advanced usage, first define a CaNNOLeSSolver
to preallocate the memory used in the algorithm, and then call solve!
:
solver = CaNNOLeSSolver(nls; linsolve = :ma57)
solve!(solver, nls; kwargs...)
or even pre-allocate the output:
stats = GenericExecutionStats(nls)
solve!(solver, nls, stats; kwargs...)
Arguments
nls :: AbstractNLSModel{T, V}
: nonlinear least-squares model created usingNLPModels
.
Keyword arguments
x::AbstractVector = nls.meta.x0
: the initial guess;λ::AbstractVector = nls.meta.y0
: the initial Lagrange multiplier;use_initial_multiplier::Bool = false
: iftrue
useλ
for the initial stopping tests;method::Symbol = :Newton
: available methods:Newton, :LM, :Newton_noFHess
, and:Newton_vanishing
;linsolve::Symbol = :ma57
: solver to compute LDLt factorization. Available methods are::ma57
,:ldlfactorizations
;max_iter::Int = -1
: maximum number of iterations;max_eval::Real = 100000
: maximum number of evaluations computed byneval_residual(nls) + neval_cons(nls)
;max_time::Float64 = 30.0
: maximum time limit in seconds;max_inner::Int = 10000
: maximum number of inner iterations (return stalled if this limit is reached);atol::T = √eps(T)
: absolute tolerance;rtol::T = √eps(T)
: relative tolerance: the algorithm usesϵtol := atol + rtol * ‖∇F(x⁰)ᵀF(x⁰) - ∇c(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⁰)‖ and ‖c(xᵏ)‖∞ ≤ √ϵtol;verbose::Int = 0
: if > 0, display iteration details everyverbose
iteration;always_accept_extrapolation::Bool = false
: iftrue
, run even if the extrapolation step fails;δdec::Real = T(0.1)
: reducing factor on the parameterδ
.
The algorithm stops when $‖c(xᵏ)‖∞ ≤ ϵtol$ and $‖∇F(xᵏ)ᵀF(xᵏ) - ∇c(xᵏ)ᵀλᵏ‖ ≤ ϵtol * max(1, ‖λᵏ‖ / 100ncon)$.
Output
The value returned is a GenericExecutionStats
, see SolverCore.jl
.
Callback
The callback is called at each iteration. The expected signature of the callback is callback(nls, solver, stats)
, and its output is ignored. Changing any of the input arguments will affect the subsequent iterations. In particular, setting stats.status = :user
will stop the algorithm. All relevant information should be available in nlp
and solver
. Notably, you can access, and modify, the following:
solver.x
: current iterate;solver.cx
: current value of the constraints atx
;stats
: structure holding the output of the algorithm (GenericExecutionStats
), which contains, among other things:stats.solution
: current iterate;stats.multipliers
: current Lagrange multipliers wrt to the constraints;stats.primal_feas
:the primal feasibility norm atsolution
;stats.dual_feas
: the dual feasibility norm atsolution
;stats.iter
: current iteration counter;stats.objective
: current objective function value;stats.status
: current status of the algorithm. Should be:unknown
unless the algorithm has attained a stopping criterion. Changing this to anything will stop the algorithm, but you should use:user
to properly indicate the intention.stats.elapsed_time
: elapsed time in seconds.
Examples
using CaNNOLeS, ADNLPModels
nls = ADNLSModel(x -> x, ones(3), 3)
stats = cannoles(nls, linsolve = :ldlfactorizations, verbose = 0)
stats
# output
"Execution stats: first-order stationary"
using CaNNOLeS, ADNLPModels
nls = ADNLSModel(x -> x, ones(3), 3)
solver = CaNNOLeSSolver(nls, linsolve = :ldlfactorizations)
stats = solve!(solver, nls, verbose = 0)
stats
# output
"Execution stats: first-order stationary"
CaNNOLeS.ParamCaNNOLeS
— TypeParamCaNNOLeS(eig_tol,δmin,κdec,κinc,κlargeinc,ρ0,ρmax,ρmin,γA)
ParamCaNNOLeS(::Type{T})
Structure containing all the parameters used in the cannoles
call.
CaNNOLeS._check_available_method
— Method_check_available_method(method::Symbol)
Return an error if method
is not in CaNNOLeS.avail_mtds
CaNNOLeS.cannoles
— Methodcannoles(nls)
Implementation of a solver for Nonlinear Least Squares with nonlinear constraints.
$min f(x) = ¹/₂‖F(x)‖² s.t. c(x) = 0$
For advanced usage, first define a CaNNOLeSSolver
to preallocate the memory used in the algorithm, and then call solve!
:
solver = CaNNOLeSSolver(nls; linsolve = :ma57)
solve!(solver, nls; kwargs...)
or even pre-allocate the output:
stats = GenericExecutionStats(nls)
solve!(solver, nls, stats; kwargs...)
Arguments
nls :: AbstractNLSModel{T, V}
: nonlinear least-squares model created usingNLPModels
.
Keyword arguments
x::AbstractVector = nls.meta.x0
: the initial guess;λ::AbstractVector = nls.meta.y0
: the initial Lagrange multiplier;use_initial_multiplier::Bool = false
: iftrue
useλ
for the initial stopping tests;method::Symbol = :Newton
: available methods:Newton, :LM, :Newton_noFHess
, and:Newton_vanishing
;linsolve::Symbol = :ma57
: solver to compute LDLt factorization. Available methods are::ma57
,:ldlfactorizations
;max_iter::Int = -1
: maximum number of iterations;max_eval::Real = 100000
: maximum number of evaluations computed byneval_residual(nls) + neval_cons(nls)
;max_time::Float64 = 30.0
: maximum time limit in seconds;max_inner::Int = 10000
: maximum number of inner iterations (return stalled if this limit is reached);atol::T = √eps(T)
: absolute tolerance;rtol::T = √eps(T)
: relative tolerance: the algorithm usesϵtol := atol + rtol * ‖∇F(x⁰)ᵀF(x⁰) - ∇c(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⁰)‖ and ‖c(xᵏ)‖∞ ≤ √ϵtol;verbose::Int = 0
: if > 0, display iteration details everyverbose
iteration;always_accept_extrapolation::Bool = false
: iftrue
, run even if the extrapolation step fails;δdec::Real = T(0.1)
: reducing factor on the parameterδ
.
The algorithm stops when $‖c(xᵏ)‖∞ ≤ ϵtol$ and $‖∇F(xᵏ)ᵀF(xᵏ) - ∇c(xᵏ)ᵀλᵏ‖ ≤ ϵtol * max(1, ‖λᵏ‖ / 100ncon)$.
Output
The value returned is a GenericExecutionStats
, see SolverCore.jl
.
Callback
The callback is called at each iteration. The expected signature of the callback is callback(nls, solver, stats)
, and its output is ignored. Changing any of the input arguments will affect the subsequent iterations. In particular, setting stats.status = :user
will stop the algorithm. All relevant information should be available in nlp
and solver
. Notably, you can access, and modify, the following:
solver.x
: current iterate;solver.cx
: current value of the constraints atx
;stats
: structure holding the output of the algorithm (GenericExecutionStats
), which contains, among other things:stats.solution
: current iterate;stats.multipliers
: current Lagrange multipliers wrt to the constraints;stats.primal_feas
:the primal feasibility norm atsolution
;stats.dual_feas
: the dual feasibility norm atsolution
;stats.iter
: current iteration counter;stats.objective
: current objective function value;stats.status
: current status of the algorithm. Should be:unknown
unless the algorithm has attained a stopping criterion. Changing this to anything will stop the algorithm, but you should use:user
to properly indicate the intention.stats.elapsed_time
: elapsed time in seconds.
Examples
using CaNNOLeS, ADNLPModels
nls = ADNLSModel(x -> x, ones(3), 3)
stats = cannoles(nls, linsolve = :ldlfactorizations, verbose = 0)
stats
# output
"Execution stats: first-order stationary"
using CaNNOLeS, ADNLPModels
nls = ADNLSModel(x -> x, ones(3), 3)
solver = CaNNOLeSSolver(nls, linsolve = :ldlfactorizations)
stats = solve!(solver, nls, verbose = 0)
stats
# output
"Execution stats: first-order stationary"
CaNNOLeS.check_nan_inf
— MethodReturn true
if x
contains NaN
or Inf
.
CaNNOLeS.dual_scaling
— Methodsd = dual_scaling(λ::AbstractVector{T}, smax::T)
Return the dual scaling on the residual, so that the algorithm stops when max(normdual / sd, normprimal) <= ϵtol
. Return 1 if the problem has no constraints.
CaNNOLeS.get_nnzh
— Functionget_nnzh(::HessianStruct)
Return number of nonzeros in the approximatation of the Hessian.
CaNNOLeS.newton_system!
— Methodnewton_system!(d, nvar, nequ, ncon, rhs, vals, LDLT, ρold, params)
Compute an LDLt factorization of the (nvar + nequ + ncon
)-square matrix for the Newton system contained in LDLT
, i.e., sparse(LDLT.rows, LDLT.cols, LDLT.vals, N, N)
. If the factorization fails, a new factorization is attempted with an increased value for the regularization ρ as long as it is smaller than params.ρmax
. The factorization is then used to solve the linear system whose right-hand side is rhs
.
Output
d
: the solution of the linear system;solve_success
:true
if the usage of the LDLt factorization is successful;ρ
: the value of the regularization parameter used in the factorization;ρold
: the value of the regularization parameter used in the previous successful factorization, or 0 if this is the first one;nfact
: the number of factorization attempts.
CaNNOLeS.optimality_check_small_residual!
— Methodnormprimal, normdual = optimality_check_small_residual!(cgls_solver, r, λ, dual, primal, Fx, cx, Jx, Jcx, Jxtr, Jcxtλ)
Compute the norm of the primal and dual residuals. The values of r
, Jxtr
, λ
, primal
and dual
are updated.
CaNNOLeS.prepare_newton_system!
— Methodprepare_newton_system!(::HessianStruct{Ti}, vals, nls, x, λ, r, Jx_vals, Jcx_vals, δ, Fx)
Update the Newton system.
CaNNOLeS.solve_ldl!
— Functionsuccess = solve_ldl!(rhs::AbstractVector, factor::Union{Ma57, LDLFactorizations.LDLFactorization}, d::AbstractVector)
Compute the solution of LDLt d = -rhs
.
CaNNOLeS.try_to_factorize
— Functionsuccess = try_to_factorize(LDLT::LinearSolverStruct, vals::AbstractVector, nvar::Integer, nequ::Integer, ncon::Integer, eig_tol::Real)
Compute the LDLt factorization of A = sparse(LDLT.rows, LDLT.cols, LDLT.vals, N, N) where N = nvar + ncon + nequ
and return true
in case of success.
CaNNOLeS.update_newton_hessian!
— Methodupdate_newton_hessian!(::HessianStruct, nls, x, r, vals, Fx)
Update, if need for method
, the top-left block with the non-zeros values of the Hessian of the residual. For method=:Newton_vanishing
, this update is skipped if ‖F(xᵏ)‖ ≤ 1e-8
.