API
FidesProblem
To solve an optimization problem with Fides, a FidesProblem must first be created:
Fides.FidesProblem — Type
FidesProblem(f, grad!, x0; hess! = nothing, lb = nothing, ub = nothing)Optimization problem to be minimized with the Fides Newton Trust Region optimizer.
Arguments
f: The objective function to minimize. Accepts aVectororComponentVectoras input and return a scalar.grad!: In-place function to compute the gradient offon the formgrad!(g, x).x0: Initial starting point for the optimization. Can be aVectororComponentVectorfrom ComponentArrays.jl.hess!: (Optional) In-place function to compute the Hessian offon the formhess!(H, x). If not provided, a Hessian approximation method must be selected when callingsolve.lb: Lower parameter bounds. Defaults to-Infif not specified.ub: Upper parameter bounds. Defaults toInfif not specified.
See also solve and FidesOptions.
FidesProblem(fides_obj, x0; lb = nothing, ub = nothing)Optimization problem created from a function that computes:
- Objective and gradient;
fides_obj(x) -> (obj, g). - Objective, gradient and Hessian;
fides_obj(x) -> (obj, g, H).
Internally, Fides computes the objective function and derivatives simultaneously. Therefore, this constructor is the most runtime-efficient option when intermediate quantities can be reused between the objective and derivative computations.
Description of Fides method
Fides implements an Interior Trust Region Reflective method for boundary-constrained optimization, as described in [1, 2]. Optimization problems on the following form are targeted:
\[\min_x f(x) \quad \text{s.t.} \quad lb \leq x \leq ub\]
At each iteration, the Fides approximates the objective function by a second-order model:
\[\min_x m_k(x) = f(x_k) + \nabla f(x_k)^T (x - x_k) + 0.5 (x - x_k)^T B_k (x - x_k) \quad \text{s.t.} \quad |x - x_k| \leq \Delta_k\]
Where, Δₖ is the trust region radius reflecting the confidence in the second-order approximation, ∇f(xₖ) is the gradient of f at the current iteration xₖ, and Bₖ is a symmetric positive-semidefinite matrix, that is either the exact Hessian (if hess! is provided) or an approximation.
References
- Coleman, T. F., & Li, Y. (1994). On the convergence of interior-reflective Newton methods for nonlinear minimization subject to bounds. Mathematical programming, 67(1), 189-224.
- Coleman, T. F., & Li, Y. (1996). An interior trust region approach for nonlinear minimization subject to bounds. SIAM Journal on optimization, 6(2), 418-445.
Thereafter, the FidesProblem is solved using the solve function, which accepts numerous tuning options:
Fides.solve — Function
solve(prob::FidesProblem, hess_update; options = FidesOptions())Solve the given FidesProblem using the Fides Trust Region method, with the specified hess_update method for computing the Hessian matrix.
In case a custom Hessian is provided to prob, use hess_update = Fides.CustomHessian. Otherwise, a Hessian approximation must be provided, and a complete list of available approximations can be found in the API documentation.
See also FidesOptions.
sourceFides.FidesOptions — Type
FidesOptions(; kwargs...)Options for the Fides Optimizer.
Keyword arguments
maxiter = 1000: Maximum number of allowed iterationsfatol = 1e-8: Absolute tolerance for convergence based on objective (f) valuefrtol = 1e-8: Relative tolerance for convergence based on objective (f) valuegatol = 1e-6: Absolute tolerance for convergence based on the gradientgrtol = 0.0: Relative tolerance for convergence based on the gradientxtol = 0.0: Tolerance for convergence based onx(parameter vector)maxtime = Inf: Maximum amount of wall-time in secondsverbose: The logging (verbosity) level of the optimizer. Allowed values are:warning(default): Only warnings are printed.info: Information is printed for each iterations.error: Only errors are printed.debug: Detailed information is printed, typically only of interest for developers.
stepback_strategy: Refinement method if proposed step reaches optimization boundary. Allowed options are:reflect(default): Recursive reflections at boundaryrefine: Perform optimization to refine stepreflect_single: Single reflection at boundarymixed: Mix reflections and truncationstrunace: Truncate step at boundary and re-solve
subspace_solver: Subspace dimension in which the Trust region subproblem is solved. Allowed options are:2D(default): Two dimensional Newton/Gradient subspacescg: CG subspace via Steihaug’s methodfull: Full on R^n
delta_init = 1.0: Initial trust region radiusmu = 0.25: Acceptance threshold for trust region ratioeta = 0.75: Trust region increase threshold for trust region ratiotheta_max = 0.95: Maximal fraction of step that would hit boundsgamma1 = 0.25: Factor by which trust region radius will be decreasedgamma2 = 2.0: Factor by which trust region radius will be increasedhistory_file = nothing: Records statistics when set
Results are stored in a FidesSolution struct:
Fides.FidesSolution — Type
FidesSolutionSolution information from a Fides optmization run.
Fields
xmin: Minimizing parameter vector found by the optimization rungfmin: Minimum objective value found by the optimization runniterations: Number of iterations for the optimization runruntime: Runtime in seconds for the optimization runretcode: Return code from the optimization run. Possible values are:DELTA_TOO_SMALL: Trust Region Radius too small to proceedDID_NOT_RUN: Optimizer did not runEXCEEDED_BOUNDARY: Exceeded specified boundariesFTOL: Converged according to fval differenceGTOL: Converged according to gradient normXTOL: Converged according to x differenceMAXITER: Reached maximum number of allowed iterationsMAXTIME: Reached maximum runtimeNOT_FINITE: Encountered non-finite fval/grad/hess
Hessian Options
Multiple Hessian options and approximation methods are available. When the Hessian is too costly or difficult to compute, the BFGS method is often performant:
Fides.CustomHessian — Type
CustomHessian()User‑provided Hessian function.
The Hessian function should be provided when creating a FidesProblem.
See also: FidesProblem
sourceFides.BFGS — Type
BFGS(; init_hess = nothing, enforce_curv_cond::Bool = true)The Broyden-Fletcher-Goldfarb-Shanno (BFGS) update strategy is a rank-2 update method that preserves both symmetry and positive-semidefiniteness [1].
Keyword arguments
init_hess = nothing: Initial Hessian for the update scheme. If provided as aMatrix, the given matrix is used; if set tonothing(default), the identity matrix is used.enforce_curv_cond = true: Whether the update should attempt to preserve positive definiteness. Iftrue, updates from steps that violate the curvature condition are discarded.
References
- Nocedal, Jorge, and Stephen J. Wright, eds. Numerical optimization. New York, NY: Springer New York, 1999.
Fides.SR1 — Type
SR1(; init_hess = nothing)The Symmetric Rank 1 update strategy as described in [1]. This is a rank 1 update strategy that preserves symmetry but does not preserve positive-semidefiniteness.
Keyword arguments
init_hess = nothing: Initial Hessian for the update scheme. If provided as aMatrix, the given matrix is used; if set tonothing(default), the identity matrix is used.
References
- Nocedal, Jorge, and Stephen J. Wright, eds. Numerical optimization (Chapter 6.2). New York, NY: Springer New York, 1999.
Fides.DFP — Type
DFP(; init_hess = nothing, enforce_curv_cond::Bool = true)The Davidon-Fletcher-Powell update strategy [1]. This is a rank 2 update strategy that preserves symmetry and positive-semidefiniteness.
Keyword arguments
init_hess = nothing: Initial Hessian for the update scheme. If provided as aMatrix, the given matrix is used; if set tonothing(default), the identity matrix is used.enforce_curv_cond = true: Whether the update should attempt to preserve positive definiteness. Iftrue, updates from steps that violate the curvature condition are discarded.
References
- Avriel, M. (2003). Nonlinear programming: analysis and methods. Courier Corporation.
Fides.Broyden — Type
Broyden(phi; init_hess = nothing, enforce_curv_cond::Bool = true)The update scheme, as described in [1], which is a generalization of the BFGS/DFP methods where phi controls the convex combination between the two. This rank-2 update strategy preserves both symmetry and positive-semidefiniteness when 0 ≤ phi ≤ 1.
Arguments
phi::AbstractFloat: The convex combination parameter interpolating between BFGS (phi=0) and DFP (phi=1).
Keyword arguments
init_hess = nothing: Initial Hessian for the update scheme. If provided as aMatrix, the given matrix is used; if set tonothing(default), the identity matrix is used.enforce_curv_cond = true: Whether the update should attempt to preserve positive definiteness. Iftrue, updates from steps that violate the curvature condition are discarded.
References
- Nocedal, Jorge, and Stephen J. Wright, eds. Numerical optimization. New York, NY: Springer New York, 1999.
Fides.BG — Type
BG(; init_hess = nothing}Broydens “good” method as introduced in [1]. This is a rank 1 update strategy that does not preserve symmetry or positive definiteness.
Keyword arguments
init_hess = nothing: Initial Hessian for the update scheme. If provided as aMatrix, the given matrix is used; if set tonothing(default), the identity matrix is used.
References
- Broyden, C. G. (1965). A class of methods for solving nonlinear simultaneous equations. Mathematics of computation, 19(92), 577-593.
Fides.BB — Type
BB(; init_hess = nothing}The Broydens “bad” method as introduced in [1]. This is a rank 1 update strategy that does not preserve symmetry or positive definiteness.
Keyword arguments
init_hess = nothing: Initial Hessian for the update scheme. If provided as aMatrix, the given matrix is used; if set tonothing(default), the identity matrix is used.
References
- Broyden, C. G. (1965). A class of methods for solving nonlinear simultaneous equations. Mathematics of computation, 19(92), 577-593.