Table of contents
Differential Evolution (DE) is a stochastic genetic search algorithm for global optimization of problems of the form
\min_{x \in X} f(x)
where f is potentially ill-behaved in one or more ways, such as non-convexity and/or non-differentiability. The updating rule for DE is described below.
Let X^{(i)} denote the N \times d dimensional array of input values at stage i of the algorithm, where each row corresponds to a different vector of candidate solutions.
The Mutation Step. For unique random indices a,b,c \in \{1, \ldots, d\}, set the mutation proposal X^{(*)} as follows.
If
de_mutation_method = 1
, use the 'rand' method:X^{(*)} = X^{(i)}(c,:) + F \times \left( X^{(i)}(a,:) - X^{(i)}(b,:) \right)
where F is the mutation parameter, set via
de_par_F
.If
de_mutation_method = 2
, use the 'best' method:X^{(*)} = X^{(i)}(\text{best},:) + F \times ( X^{(i)}(a,:) - X^{(i)}(b,:) )
where
X^{(i)} (\text{best},:) := \arg \min \left\{ f(X^{(i)}(1,:)), \ldots, f(X^{(i)}(N,:)) \right\}
The Crossover Step.
Choose a random integer r_k \in \{1, \ldots, d \}.
Draw a vector u of independent uniform random variables of length d
For each j \in \{ 1, \ldots, N \} and k \in \{ 1, \ldots, d \}, set
X_c^{(*)} (j,k) = \begin{cases} X^*(j,k) & \text{ if } u_k \leq CR \text{ or } k = r_k \\ X^{(i)} (j,k) & \text{ else } \end{cases}
where CR \in [0,1] is the crossover parameter, set via
de_par_CR
.
The Update Step.
X^{(i+1)} (j,:) = \begin{cases} X_c^*(j,:) & \text{ if } f(X_c^*(j,:)) < f(X^{(i)}(j,:)) \\ X^{(i)} (j,:) & \text{ else } \end{cases}
The algorithm stops when at least one of the following conditions are met:
- the relative improvement in the objective function from the best candidate solution is less than
rel_objfn_change_tol
betweende_settings.check_freq
number of generations;- the total number of generations exceeds
de_settings.n_gen
.
TBW.
.. doxygenfunction:: de(ColVec_t& init_out_vals, std::function<fp_t (const ColVec_t& vals_inp, ColVec_t* grad_out, void* opt_data)> opt_objfn, void* opt_data) :project: optimlib
.. doxygenfunction:: de(ColVec_t& init_out_vals, std::function<fp_t (const ColVec_t& vals_inp, ColVec_t* grad_out, void* opt_data)> opt_objfn, void* opt_data, algo_settings_t& settings) :project: optimlib
.. doxygenfunction:: de_prmm(ColVec_t& init_out_vals, std::function<fp_t (const ColVec_t& vals_inp, ColVec_t* grad_out, void* opt_data)> opt_objfn, void* opt_data) :project: optimlib
.. doxygenfunction:: de_prmm(ColVec_t& init_out_vals, std::function<fp_t (const ColVec_t& vals_inp, ColVec_t* grad_out, void* opt_data)> opt_objfn, void* opt_data, algo_settings_t& settings) :project: optimlib
The basic control parameters are:
fp_t rel_objfn_change_tol
: the error tolerance value controlling how small the relative change in best candidate solution should be before 'convergence' is declared.size_t iter_max
: the maximum number of iterations/updates before the algorithm exits.bool vals_bound
: whether the search space of the algorithm is bounded. Iftrue
, thenColVec_t lower_bounds
: defines the lower bounds of the search space.ColVec_t upper_bounds
: defines the upper bounds of the search space.
In addition to these:
int print_level
: Set print level.- Level 1: Print iteration count and error value.
- Level 2: Level 1 and print best input values and corresponding objective function value.
- Level 3: Level 2 and print full population matrix, X.
Code to run this example is given below.
.. toggle-header:: :header: **Armadillo (Click to show/hide)** .. code:: cpp #define OPTIM_ENABLE_ARMA_WRAPPERS #include "optim.hpp" #define OPTIM_PI 3.14159265358979 double ackley_fn(const arma::vec& vals_inp, arma::vec* grad_out, void* opt_data) { const double x = vals_inp(0); const double y = vals_inp(1); const double obj_val = 20 + std::exp(1) - 20*std::exp( -0.2*std::sqrt(0.5*(x*x + y*y)) ) - std::exp( 0.5*(std::cos(2 * OPTIM_PI * x) + std::cos(2 * OPTIM_PI * y)) ); return obj_val; } int main() { arma::vec x = arma::ones(2,1) + 1.0; // initial values: (2,2) bool success = optim::de(x, ackley_fn, nullptr); if (success) { std::cout << "de: Ackley test completed successfully." << std::endl; } else { std::cout << "de: Ackley test completed unsuccessfully." << std::endl; } arma::cout << "de: solution to Ackley test:\n" << x << arma::endl; return 0; }
.. toggle-header:: :header: **Eigen (Click to show/hide)** .. code:: cpp #define OPTIM_ENABLE_EIGEN_WRAPPERS #include "optim.hpp" #define OPTIM_PI 3.14159265358979 double ackley_fn(const Eigen::VectorXd& vals_inp, Eigen::VectorXd* grad_out, void* opt_data) { const double x = vals_inp(0); const double y = vals_inp(1); const double obj_val = 20 + std::exp(1) - 20*std::exp( -0.2*std::sqrt(0.5*(x*x + y*y)) ) - std::exp( 0.5*(std::cos(2 * OPTIM_PI * x) + std::cos(2 * OPTIM_PI * y)) ); return obj_val; } int main() { Eigen::VectorXd x = 2.0 * Eigen::VectorXd::Ones(2); // initial values: (2,2) bool success = optim::de(x, ackley_fn, nullptr); if (success) { std::cout << "de: Ackley test completed successfully." << std::endl; } else { std::cout << "de: Ackley test completed unsuccessfully." << std::endl; } arma::cout << "de: solution to Ackley test:\n" << x << arma::endl; return 0; }
Code to run this example is given below.
.. toggle-header:: :header: **Armadillo Code (Click to show/hide)** .. code:: cpp #define OPTIM_ENABLE_ARMA_WRAPPERS #include "optim.hpp" #define OPTIM_PI 3.14159265358979 struct rastrigin_fn_data { double A; }; double rastrigin_fn(const arma::vec& vals_inp, arma::vec* grad_out, void* opt_data) { const int n = vals_inp.n_elem; rastrigin_fn_data* objfn_data = reinterpret_cast<rastrigin_fn_data*>(opt_data); const double A = objfn_data->A; double obj_val = A*n + arma::accu( arma::pow(vals_inp,2) - A*arma::cos(2 * OPTIM_PI * vals_inp) ); return obj_val; } int main() { rastrigin_fn_data test_data; test_data.A = 10; arma::vec x = arma::ones(2,1) + 1.0; // initial values: (2,2) bool success = optim::de(x, rastrigin_fn, &test_data); if (success) { std::cout << "de: Rastrigin test completed successfully." << std::endl; } else { std::cout << "de: Rastrigin test completed unsuccessfully." << std::endl; } arma::cout << "de: solution to Rastrigin test:\n" << x << arma::endl; return 0; }
.. toggle-header:: :header: **Eigen Code (Click to show/hide)** .. code:: cpp #define OPTIM_ENABLE_EIGEN_WRAPPERS #include "optim.hpp" #define OPTIM_PI 3.14159265358979 struct rastrigin_fn_data { double A; }; double rastrigin_fn(const Eigen::VectorXd& vals_inp, Eigen::VectorXd* grad_out, void* opt_data) { const int n = vals_inp.n_elem; rastrigin_fn_data* objfn_data = reinterpret_cast<rastrigin_fn_data*>(opt_data); const double A = objfn_data->A; double obj_val = A*n + vals_inp.array().pow(2).sum() - A * (2 * OPTIM_PI * vals_inp).array().cos().sum(); return obj_val; } int main() { rastrigin_fn_data test_data; test_data.A = 10; Eigen::VectorXd x = 2.0 * Eigen::VectorXd::Ones(2); // initial values: (2,2) bool success = optim::de(x, rastrigin_fn, &test_data); if (success) { std::cout << "de: Rastrigin test completed successfully." << std::endl; } else { std::cout << "de: Rastrigin test completed unsuccessfully." << std::endl; } arma::cout << "de: solution to Rastrigin test:\n" << x << arma::endl; return 0; }