Optimization

This page documents library components that attempt to find the minimum or maximum of a user supplied function. An introduction to the general purpose non-linear optimizers in this section can be found here.

Main Routines find_min find_min_single_variable find_min_using_approximate_derivatives find_min_bobyqa solve_qp_using_smo oca find_max find_max_single_variable find_max_using_approximate_derivatives find_max_bobyqa
Strategies cg_search_strategy bfgs_search_strategy lbfgs_search_strategy objective_delta_stop_strategy gradient_norm_stop_strategy
Helper Routines derivative negate_function make_line_search_function poly_min_extrap lagrange_poly_min_extrap line_search
derivative dlib/optimization.h dlib/optimization/optimization_abstract.h This is a function that takes another function as input and returns a function object that numerically computes the derivative of the input function. negate_function dlib/optimization.h dlib/optimization/optimization_abstract.h This is a function that takes another function as input and returns a function object that computes the negation of the input function. make_line_search_function dlib/optimization.h dlib/optimization/optimization_line_search_abstract.h This is a function that takes another function f(x) as input and returns a function object l(z) = f(start + z*direction). It is useful for turning multi-variable functions into single-variable functions for use with the line_search routine. poly_min_extrap dlib/optimization.h dlib/optimization/optimization_line_search_abstract.h This function finds the 3rd degree polynomial that interpolates a set of points and returns you the minimum of that polynomial. lagrange_poly_min_extrap dlib/optimization.h dlib/optimization/optimization_line_search_abstract.h This function finds the second order polynomial that interpolates a set of points and returns you the minimum of that polynomial. line_search dlib/optimization.h dlib/optimization/optimization_line_search_abstract.h Performs a gradient based line search on a given function and returns the input that makes the function significantly smaller. cg_search_strategy dlib/optimization.h dlib/optimization/optimization_search_strategies_abstract.h This object represents a strategy for determining which direction a line search should be carried out along. This particular object is an implementation of the Polak-Ribiere conjugate gradient method for determining this direction.

This method uses an amount of memory that is linear in the number of variables to be optimized. So it is capable of handling problems with a very large number of variables. However, it is generally not as good as the L-BFGS algorithm (see the lbfgs_search_strategy class).

optimization_ex.cpp.html
bfgs_search_strategy dlib/optimization.h dlib/optimization/optimization_search_strategies_abstract.h This object represents a strategy for determining which direction a line search should be carried out along. This particular object is an implementation of the BFGS quasi-newton method for determining this direction.

This method uses an amount of memory that is quadratic in the number of variables to be optimized. It is generally very effective but if your problem has a very large number of variables then it isn't appropriate. Instead You should try the lbfgs_search_strategy.

optimization_ex.cpp.html
lbfgs_search_strategy dlib/optimization.h dlib/optimization/optimization_search_strategies_abstract.h This object represents a strategy for determining which direction a line search should be carried out along. This particular object is an implementation of the L-BFGS quasi-newton method for determining this direction.

This method uses an amount of memory that is linear in the number of variables to be optimized. This makes it an excellent method to use when an optimization problem has a large number of variables.

optimization_ex.cpp.html
objective_delta_stop_strategy dlib/optimization.h dlib/optimization/optimization_stop_strategies_abstract.h This object represents a strategy for deciding if an optimization algorithm should terminate. This particular object looks at the change in the objective function from one iteration to the next and bases its decision on how large this change is. If the change is below a user given threshold then the search stops. optimization_ex.cpp.html gradient_norm_stop_strategy dlib/optimization.h dlib/optimization/optimization_stop_strategies_abstract.h This object represents a strategy for deciding if an optimization algorithm should terminate. This particular object looks at the norm (i.e. the length) of the current gradient vector and stops if it is smaller than a user given threshold. find_min dlib/optimization.h dlib/optimization/optimization_abstract.h Performs an unconstrained minimization of a nonlinear function using some search strategy (e.g. bfgs_search_strategy). optimization_ex.cpp.html find_min_single_variable dlib/optimization.h dlib/optimization/optimization_line_search_abstract.h Performs a bound constrained minimization of a nonlinear function. The function must be of a single variable. Derivatives are not required. solve_qp_using_smo dlib/optimization.h dlib/optimization/optimization_solve_qp_using_smo_abstract.h This function solves the following quadratic program:
   Minimize: f(alpha) == 0.5*trans(alpha)*Q*alpha - trans(alpha)*b
   subject to the following constraints:
      sum(alpha) == C 
      min(alpha) >= 0 
oca dlib/optimization.h dlib/optimization/optimization_oca_abstract.h This object is a tool for solving the following optimization problem:
   Minimize: f(w) == 0.5*dot(w,w) + C*R(w)

   Where R(w) is a convex function and C > 0


For a detailed discussion you should consult the following papers from the Journal of Machine Learning Research:
Optimized Cutting Plane Algorithm for Large-Scale Risk Minimization by Vojtech Franc, Soren Sonnenburg; 10(Oct):2157--2192, 2009.
Bundle Methods for Regularized Risk Minimization by Choon Hui Teo, S.V.N. Vishwanthan, Alex J. Smola, Quoc V. Le; 11(Jan):311-365, 2010.
find_min_bobyqa dlib/optimization.h dlib/optimization/optimization_bobyqa_abstract.h This function defines the dlib interface to the BOBYQA software developed by M.J.D Powell. BOBYQA is a method for optimizing a function in the absence of derivative information. Powell described it as a method that seeks the least value of a function of many variables, by applying a trust region method that forms quadratic models by interpolation. There is usually some freedom in the interpolation conditions, which is taken up by minimizing the Frobenius norm of the change to the second derivative of the model, beginning with the zero matrix. The values of the variables are constrained by upper and lower bounds.

The following paper, published in 2009 by Powell, describes the detailed working of the BOBYQA algorithm.

The BOBYQA algorithm for bound constrained optimization without derivatives by M.J.D. Powell

Note that BOBYQA only works on functions of two or more variables. So if you need to perform derivative-free optimization on a function of a single variable then you should use the find_min_single_variable function.

optimization_ex.cpp.html model_selection_ex.cpp.html
find_max_bobyqa dlib/optimization.h dlib/optimization/optimization_bobyqa_abstract.h This function is identical to the find_min_bobyqa routine except that it negates the objective function before performing optimization. Thus this function will attempt to find the maximizer of the objective rather than the minimizer.

Note that BOBYQA only works on functions of two or more variables. So if you need to perform derivative-free optimization on a function of a single variable then you should use the find_max_single_variable function.

optimization_ex.cpp.html model_selection_ex.cpp.html
find_min_using_approximate_derivatives dlib/optimization.h dlib/optimization/optimization_abstract.h Performs an unconstrained minimization of a nonlinear function using some search strategy (e.g. bfgs_search_strategy). This version doesn't take a gradient function but instead numerically approximates the gradient. optimization_ex.cpp.html find_max dlib/optimization.h dlib/optimization/optimization_abstract.h Performs an unconstrained maximization of a nonlinear function using some search strategy (e.g. bfgs_search_strategy). find_max_single_variable dlib/optimization.h dlib/optimization/optimization_line_search_abstract.h Performs a bound constrained maximization of a nonlinear function. The function must be of a single variable. Derivatives are not required. find_max_using_approximate_derivatives dlib/optimization.h dlib/optimization/optimization_abstract.h Performs an unconstrained maximization of a nonlinear function using some search strategy (e.g. bfgs_search_strategy). This version doesn't take a gradient function but instead numerically approximates the gradient.