Reference
Contents
Index
AdaptiveRegularization.HessDense
AdaptiveRegularization.HessGaussNewtonOp
AdaptiveRegularization.HessOp
AdaptiveRegularization.HessSparse
AdaptiveRegularization.HessSparseCOO
AdaptiveRegularization.PDataKARC
AdaptiveRegularization.PDataNLSST
AdaptiveRegularization.PDataST
AdaptiveRegularization.PDataTRK
AdaptiveRegularization.TRARCWorkspace
AdaptiveRegularization.TrustRegion
AdaptiveRegularization.TrustRegionException
AdaptiveRegularization.TRARC
AdaptiveRegularization.compute_r
AdaptiveRegularization.compute_Δq
AdaptiveRegularization.decrease
AdaptiveRegularization.hessian!
AdaptiveRegularization.increase
AdaptiveRegularization.init
AdaptiveRegularization.preprocess
AdaptiveRegularization.solve_model
AdaptiveRegularization.HessDense
— TypeHessDense(::AbstractNLPModel{T,S}, n)
Return a structure used for the evaluation of dense Hessian matrix.
AdaptiveRegularization.HessGaussNewtonOp
— TypeHessGaussNewtonOp(::AbstractNLSModel{T,S}, n)
Return a structure used for the evaluation of the Hessian matrix as an operator.
AdaptiveRegularization.HessOp
— TypeHessOp(::AbstractNLPModel{T,S}, n)
Return a structure used for the evaluation of the Hessian matrix as an operator.
AdaptiveRegularization.HessSparse
— TypeHessSparse(::AbstractNLPModel{T,S}, n)
Return a structure used for the evaluation of sparse Hessian matrix.
AdaptiveRegularization.HessSparseCOO
— TypeHessSparseCOO(::AbstractNLPModel{T,S}, n)
Return a structure used for the evaluation of sparse Hessian matrix in COO-format.
AdaptiveRegularization.PDataKARC
— TypePDataKARC(::Type{S}, ::Type{T}, n)
Return a structure used for the preprocessing of ARCqK methods.
AdaptiveRegularization.PDataNLSST
— TypePDataNLSST(::Type{S}, ::Type{T}, n)
Return a structure used for the preprocessing of Steihaug-Toint methods for Gauss-Newton approximation of nonlinear least squares.
AdaptiveRegularization.PDataST
— TypePDataST(::Type{S}, ::Type{T}, n)
Return a structure used for the preprocessing of Steihaug-Toint methods.
AdaptiveRegularization.PDataTRK
— TypePDataTRK(::Type{S}, ::Type{T}, n)
Return a structure used for the preprocessing of TRK methods.
AdaptiveRegularization.TRARCWorkspace
— TypeTRARCWorkspace(nlp, ::Type{Hess}, n)
Pre-allocate the memory used during the TRARC
call for the problem nlp
of size n
. The possible values for Hess
are: HessDense
, HessSparse
, HessSparseCOO
, HessOp
. Return a TRARCWorkspace
structure.
AdaptiveRegularization.TrustRegion
— TypeTrustRegion(α₀::T;kwargs...)
Select the main parameters used in the TRARC
algorithm with α₀
as initial TR/ARC parameter. The keyword arguments are:
max_α::T
: Maximum value forα
. DefaultT(1) / sqrt(eps(T))
.acceptance_threshold::T
: Ratio over which the step is successful. DefaultT(0.1)
.increase_threshold::T
: Ratio over which we increaseα
. DefaultT(0.75)
.reduce_threshold::T
: Ratio under which we decreaseα
. DefaultT(0.1)
.increase_factor::T
: Factor of increase ofα
. DefaultT(5.0)
.decrease_factor::T
: Factor of decrease ofα
. DefaultT(0.1)
.max_unsuccinarow::Int
: Limit on the number of successive unsucessful iterations. Default30
.
Returns a TrustRegion
structure.
This can be compared to https://github.com/JuliaSmoothOptimizers/SolverTools.jl/blob/main/src/trust-region/basic-trust-region.jl
AdaptiveRegularization.TrustRegionException
— TypeException type raised in case of error.
AdaptiveRegularization.TRARC
— FunctionTRARC(nlp; kwargs...)
Compute a local minimum of an unconstrained optimization problem using trust-region (TR)/adaptive regularization with cubics (ARC) methods.
Arguments
nlp::AbstractNLPModel
: the model solved, seeNLPModels.jl
.
The keyword arguments include
TR::TrustRegion
: structure with trust-region/ARC parameters, seeTrustRegion
. Default:TrustRegion(T(10.0))
.hess_type::Type{Hess}
: Structure used to handle the hessian. The possible values are:HessDense
,HessSparse
,HessSparseCOO
,HessOp
. Default:HessOp
.pdata_type::Type{ParamData}
Structure used for the preprocessing step. Default:PDataKARC
.solve_model::Function
Function used to solve the subproblem. Default:solve_modelKARC
.robust::Bool
:true
implements a robust evaluation of the model. Default:true
.verbose::Bool
:true
prints iteration information. Default:false
.
Additional kwargs
are used for stopping criterion, see Stopping.jl
.
Output
The returned value is a GenericExecutionStats
, see SolverCore.jl
.
This implementation uses Stopping.jl
. Therefore, it is also possible to used
TRARC(stp; kwargs...)
which returns the stp::NLPStopping
updated.
For advanced usage, the principal call to the solver uses a TRARCWorkspace
.
TRARC(stp, pdata, workspace, trust_region_parameters; kwargs...)
Some variants of TRARC are already implemented and listed in AdaptiveRegularization.ALL_solvers
.
References
This method unifies the implementation of trust-region and adaptive regularization with cubics as described in
Dussault, J.-P. (2020).
A unified efficient implementation of trust-region type algorithms for unconstrained optimization.
INFOR: Information Systems and Operational Research, 58(2), 290-309.
10.1080/03155986.2019.1624490
Examples
using AdaptiveRegularization, ADNLPModels
nlp = ADNLPModel(x -> 100 * (x[2] - x[1]^2)^2 + (x[1] - 1)^2, [-1.2; 1.0]);
stats = TRARC(nlp)
# output
"Execution stats: first-order stationary"
AdaptiveRegularization.compute_r
— Methodcompute_r(nlp, f, Δf, Δq, slope, d, xnext, gnext, robust)
Compute the actual vs predicted reduction ratio ∆f/Δq
.
Arguments:
nlp
: Current model we are trying to solvef
: current objective valueΔf
:= f - f_trial
is the actual reduction is an objective/merit/penalty function,Δq
:q - q_trial
is the reduction predicted by the model q of f.slope
: current sloped
: potential next directionxnext
: potential next iterategnext
: current gradient value, ifgood_grad
is true, then this value has been udpated.robust
: iftrue
, try to trap potential cancellation errors
Output:
r
: reduction ratio∆f/Δq
good_grad
:true
ifgnext
has been recomputedgnext
: gradient.
We assume that q
is being minimized, and therefore that
Δq > 0`.
AdaptiveRegularization.compute_Δq
— Methodcompute_Δq(Hx, d, ∇f)
Update Δq = -(∇f + 0.5 * (Hx * d)) ⋅ d
in-place.
AdaptiveRegularization.decrease
— Methoddecrease(X::TPData, α::T, TR::TrustRegion)
Return a decreased α
.
AdaptiveRegularization.hessian!
— Functionhessian!(workspace::HessDense, nlp, x)
Return the Hessian matrix of nlp
at x
in-place with memory update of workspace
.
AdaptiveRegularization.increase
— Methodincrease(X::TPData, α::T, TR::TrustRegion)
Return an increased α
.
AdaptiveRegularization.init
— Methodinit(::Type{Hess}, nlp::AbstractNLPModel{T,S}, n)
Return the hessian structure Hess
and its composite type.
AdaptiveRegularization.preprocess
— Methodpreprocess(PData::TPData, H, g, gNorm2, n1, n2, α)
Function called in the TRARC
algorithm every time a new iterate has been accepted.
Arguments
PData::TPData
: data structure used for preprocessing.H
: current Hessian matrix.g
: current gradient.gNorm2
: 2-norm of the gradient.n1
: Current count on the number of Hessian-vector products.n2
: Maximum number of Hessian-vector products accepted.α
: current value of the TR/ARC parameter.
It returns PData
.
AdaptiveRegularization.solve_model
— Methodsolve_model(PData::TPData, H, g, gNorm2, n1, n2, α)
Function called in the TRARC
algorithm to solve the subproblem.
Arguments
PData::TPData
: data structure used for preprocessing.H
: current Hessian matrix.g
: current gradient.gNorm2
: 2-norm of the gradient.n1
: Current count on the number of Hessian-vector products.n2
: Maximum number of Hessian-vector products accepted.α
: current value of the TR/ARC parameter.
It returns a couple (PData.d, PData.λ)
. Current implementations include: solve_modelKARC
, solve_modelTRK
, solve_modelST_TR
.