Processing math: 0%

Nelder-Mead is a derivative-free simplex method.


Definition and Syntax

1
2
bool nm(arma::vec& init_out_vals, std::function<double (const arma::vec& vals_inp, arma::vec* grad_out, void* opt_data)> opt_objfn, void* opt_data);
bool nm(arma::vec& init_out_vals, std::function<double (const arma::vec& vals_inp, arma::vec* grad_out, void* opt_data)> opt_objfn, void* opt_data, algo_settings_t& settings);

Function arguments:

  • init_out_vals a column vector of initial values; will be replaced by the solution values.
  • opt_objfn the function to be minimized, taking three arguments:
    • vals_inp a vector of inputs;
    • grad_out an empty vector, as Nelder-Mead does not require the gradient to be known/exist; and
    • opt_data additional parameters passed to the function.
  • opt_data additional parameters passed to the function.
  • settings parameters controlling the optimization routine; see below.

Optimization control parameters:

  • double err_tol the value controlling how small should be before 'convergence' is declared.
  • int iter_max the maximum number of iterations/updates before the algorithm exits.
  • bool vals_bound whether the search space is bounded. If true, then
    • arma::vec lower_bounds this defines the lower bounds.
    • arma::vec upper_bounds this defines the upper bounds.

  • double nm_par_alpha reflection parameter.
  • double nm_par_gamma expansion parameter.
  • double nm_par_beta contraction parameter.
  • double nm_par_delta shrikage parameter.

Details

Let x^{(i)} denote the simplex values at stage i of the algorithm. The Nelder-Mead updating rule is as follows.

  • Step 1. Sort the vertices in order of function values, from smallest to largest:
  • f(x^{(i)}(1,:)) \leq f(x^{(i)}(2,:)) \leq \cdots \leq f(x^{(i)}(n+1,:))
  • Step 2. Calculate the centroid value up to the nth vertex: \bar{x} = \frac{1}{n} \sum_{j=1}^n x^{(i)}(j,:) and compute the reflection point: x^r = \bar{x} + \alpha (\bar{x} - x^{(i)}(n+1,:)) where \alpha is set by nm_par_alpha.

    If f(x^r) \geq f(x^{(i)}(1,:)) and f(x^r) < f(x^{(i)}(n,:)), then x^{(i+1)}(n+1,:) = x^r, \ \ \textbf{ and go to Step 1.} Otherwise continue to Step 3.
  • Step 3. If f(x^r) \geq f(x^{(i)}(1,:)) then go to Step 4, otherwise compute the expansion point: x^e = \bar{x} + \gamma (x^r - \bar{x}) where where \gamma is set by nm_par_gamma.

    Set x^{(i+1)}(n+1,:) = \begin{cases} x^e & \text{ if } f(x^e) < f(x^r) \\ x^r & \text{ else } \end{cases} and go to Step 1.
  • Steps 4 & 5. If f(x^r) < f(x^{(i)}(n,:)), then compute the outside or inside contraction: x^{c} = \begin{cases} \bar{x} + \beta(x^r - \bar{x}) & \text{ if } f(x^r) < f(x^{(i)}(n+1,:)) \\ \bar{x} - \beta(x^r - \bar{x}) & \text{ else} \end{cases} where \beta is set by nm_par_beta.

    If f(x^c) < f(x^{(i)}(n+1,:)), then x^{(i+1)}(n+1,:) = x^c, \ \ \textbf{ and go to Step 1.} Otherwise go to Step 6.
  • Step 6. Shrink the simplex toward x^{(i)}(1,:): x^{(i+1)}(j,:) = x^{(i)}(1,:) + \delta (x^{(i)}(j,:) - x^{(i)}(1,:)), \ \ j = 2, \ldots, n+1 where \delta is set by nm_par_delta. Go to Step 1.