differential_evolution.hpp 9.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236
  1. /*
  2. * Copyright Nick Thompson, 2024
  3. * Use, modification and distribution are subject to the
  4. * Boost Software License, Version 1.0. (See accompanying file
  5. * LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  6. */
  7. #ifndef BOOST_MATH_OPTIMIZATION_DIFFERENTIAL_EVOLUTION_HPP
  8. #define BOOST_MATH_OPTIMIZATION_DIFFERENTIAL_EVOLUTION_HPP
  9. #include <atomic>
  10. #include <boost/math/optimization/detail/common.hpp>
  11. #include <cmath>
  12. #include <limits>
  13. #include <mutex>
  14. #include <random>
  15. #include <sstream>
  16. #include <stdexcept>
  17. #include <thread>
  18. #include <utility>
  19. #include <vector>
  20. namespace boost::math::optimization {
  21. // Storn, R., Price, K. (1997). Differential evolution-a simple and efficient heuristic for global optimization over
  22. // continuous spaces.
  23. // Journal of global optimization, 11, 341-359.
  24. // See:
  25. // https://www.cp.eng.chula.ac.th/~prabhas//teaching/ec/ec2012/storn_price_de.pdf
  26. // We provide the parameters in a struct-there are too many of them and they are too unwieldy to pass individually:
  27. template <typename ArgumentContainer> struct differential_evolution_parameters {
  28. using Real = typename ArgumentContainer::value_type;
  29. using DimensionlessReal = decltype(Real()/Real());
  30. ArgumentContainer lower_bounds;
  31. ArgumentContainer upper_bounds;
  32. // mutation factor is also called scale factor or just F in the literature:
  33. DimensionlessReal mutation_factor = static_cast<DimensionlessReal>(0.65);
  34. DimensionlessReal crossover_probability = static_cast<DimensionlessReal>(0.5);
  35. // Population in each generation:
  36. size_t NP = 500;
  37. size_t max_generations = 1000;
  38. ArgumentContainer const *initial_guess = nullptr;
  39. unsigned threads = std::thread::hardware_concurrency();
  40. };
  41. template <typename ArgumentContainer>
  42. void validate_differential_evolution_parameters(differential_evolution_parameters<ArgumentContainer> const &de_params) {
  43. using std::isfinite;
  44. using std::isnan;
  45. std::ostringstream oss;
  46. detail::validate_bounds(de_params.lower_bounds, de_params.upper_bounds);
  47. if (de_params.NP < 4) {
  48. oss << __FILE__ << ":" << __LINE__ << ":" << __func__;
  49. oss << ": The population size must be at least 4, but requested population size of " << de_params.NP << ".";
  50. throw std::invalid_argument(oss.str());
  51. }
  52. // From: "Differential Evolution: A Practical Approach to Global Optimization (Natural Computing Series)"
  53. // > The scale factor, F in (0,1+), is a positive real number that controls the rate at which the population evolves.
  54. // > While there is no upper limit on F, effective values are seldom greater than 1.0.
  55. // ...
  56. // Also see "Limits on F", Section 2.5.1:
  57. // > This discontinuity at F = 1 reduces the number of mutants by half and can result in erratic convergence...
  58. auto F = de_params.mutation_factor;
  59. if (isnan(F) || F >= 1 || F <= 0) {
  60. oss << __FILE__ << ":" << __LINE__ << ":" << __func__;
  61. oss << ": F in (0, 1) is required, but got F=" << F << ".";
  62. throw std::domain_error(oss.str());
  63. }
  64. if (de_params.max_generations < 1) {
  65. oss << __FILE__ << ":" << __LINE__ << ":" << __func__;
  66. oss << ": There must be at least one generation.";
  67. throw std::invalid_argument(oss.str());
  68. }
  69. if (de_params.initial_guess) {
  70. detail::validate_initial_guess(*de_params.initial_guess, de_params.lower_bounds, de_params.upper_bounds);
  71. }
  72. if (de_params.threads == 0) {
  73. oss << __FILE__ << ":" << __LINE__ << ":" << __func__;
  74. oss << ": There must be at least one thread.";
  75. throw std::invalid_argument(oss.str());
  76. }
  77. }
  78. template <typename ArgumentContainer, class Func, class URBG>
  79. ArgumentContainer differential_evolution(
  80. const Func cost_function, differential_evolution_parameters<ArgumentContainer> const &de_params, URBG &gen,
  81. std::invoke_result_t<Func, ArgumentContainer> target_value =
  82. std::numeric_limits<std::invoke_result_t<Func, ArgumentContainer>>::quiet_NaN(),
  83. std::atomic<bool> *cancellation = nullptr,
  84. std::vector<std::pair<ArgumentContainer, std::invoke_result_t<Func, ArgumentContainer>>> *queries = nullptr,
  85. std::atomic<std::invoke_result_t<Func, ArgumentContainer>> *current_minimum_cost = nullptr) {
  86. using Real = typename ArgumentContainer::value_type;
  87. using DimensionlessReal = decltype(Real()/Real());
  88. using ResultType = std::invoke_result_t<Func, ArgumentContainer>;
  89. using std::clamp;
  90. using std::isnan;
  91. using std::round;
  92. using std::uniform_real_distribution;
  93. validate_differential_evolution_parameters(de_params);
  94. const size_t dimension = de_params.lower_bounds.size();
  95. auto NP = de_params.NP;
  96. auto population = detail::random_initial_population(de_params.lower_bounds, de_params.upper_bounds, NP, gen);
  97. if (de_params.initial_guess) {
  98. population[0] = *de_params.initial_guess;
  99. }
  100. std::vector<ResultType> cost(NP, std::numeric_limits<ResultType>::quiet_NaN());
  101. std::atomic<bool> target_attained = false;
  102. // This mutex is only used if the queries are stored:
  103. std::mutex mt;
  104. std::vector<std::thread> thread_pool;
  105. auto const threads = de_params.threads;
  106. for (size_t j = 0; j < threads; ++j) {
  107. // Note that if some members of the population take way longer to compute,
  108. // then this parallelization strategy is very suboptimal.
  109. // However, we tried using std::async (which should be robust to this particular problem),
  110. // but the overhead was just totally unacceptable on ARM Macs (the only platform tested).
  111. // As the economists say "there are no solutions, only tradeoffs".
  112. thread_pool.emplace_back([&, j]() {
  113. for (size_t i = j; i < cost.size(); i += threads) {
  114. cost[i] = cost_function(population[i]);
  115. if (current_minimum_cost && cost[i] < *current_minimum_cost) {
  116. *current_minimum_cost = cost[i];
  117. }
  118. if (queries) {
  119. std::scoped_lock lock(mt);
  120. queries->push_back(std::make_pair(population[i], cost[i]));
  121. }
  122. if (!isnan(target_value) && cost[i] <= target_value) {
  123. target_attained = true;
  124. }
  125. }
  126. });
  127. }
  128. for (auto &thread : thread_pool) {
  129. thread.join();
  130. }
  131. std::vector<ArgumentContainer> trial_vectors(NP);
  132. for (size_t i = 0; i < NP; ++i) {
  133. if constexpr (detail::has_resize_v<ArgumentContainer>) {
  134. trial_vectors[i].resize(dimension);
  135. }
  136. }
  137. std::vector<URBG> thread_generators(threads);
  138. for (size_t j = 0; j < threads; ++j) {
  139. thread_generators[j].seed(gen());
  140. }
  141. // std::vector<bool> isn't threadsafe!
  142. std::vector<int> updated_indices(NP, 0);
  143. for (size_t generation = 0; generation < de_params.max_generations; ++generation) {
  144. if (cancellation && *cancellation) {
  145. break;
  146. }
  147. if (target_attained) {
  148. break;
  149. }
  150. thread_pool.resize(0);
  151. for (size_t j = 0; j < threads; ++j) {
  152. thread_pool.emplace_back([&, j]() {
  153. auto& tlg = thread_generators[j];
  154. uniform_real_distribution<DimensionlessReal> unif01(DimensionlessReal(0), DimensionlessReal(1));
  155. for (size_t i = j; i < cost.size(); i += threads) {
  156. if (target_attained) {
  157. return;
  158. }
  159. if (cancellation && *cancellation) {
  160. return;
  161. }
  162. size_t r1, r2, r3;
  163. do {
  164. r1 = tlg() % NP;
  165. } while (r1 == i);
  166. do {
  167. r2 = tlg() % NP;
  168. } while (r2 == i || r2 == r1);
  169. do {
  170. r3 = tlg() % NP;
  171. } while (r3 == i || r3 == r2 || r3 == r1);
  172. for (size_t k = 0; k < dimension; ++k) {
  173. // See equation (4) of the reference:
  174. auto guaranteed_changed_idx = tlg() % dimension;
  175. if (unif01(tlg) < de_params.crossover_probability || k == guaranteed_changed_idx) {
  176. auto tmp = population[r1][k] + de_params.mutation_factor * (population[r2][k] - population[r3][k]);
  177. auto const &lb = de_params.lower_bounds[k];
  178. auto const &ub = de_params.upper_bounds[k];
  179. // Some others recommend regenerating the indices rather than clamping;
  180. // I dunno seems like it could get stuck regenerating . . .
  181. trial_vectors[i][k] = clamp(tmp, lb, ub);
  182. } else {
  183. trial_vectors[i][k] = population[i][k];
  184. }
  185. }
  186. auto const trial_cost = cost_function(trial_vectors[i]);
  187. if (isnan(trial_cost)) {
  188. continue;
  189. }
  190. if (queries) {
  191. std::scoped_lock lock(mt);
  192. queries->push_back(std::make_pair(trial_vectors[i], trial_cost));
  193. }
  194. if (trial_cost < cost[i] || isnan(cost[i])) {
  195. cost[i] = trial_cost;
  196. if (!isnan(target_value) && cost[i] <= target_value) {
  197. target_attained = true;
  198. }
  199. if (current_minimum_cost && cost[i] < *current_minimum_cost) {
  200. *current_minimum_cost = cost[i];
  201. }
  202. // Can't do this! It's a race condition!
  203. //population[i] = trial_vectors[i];
  204. // Instead mark all the indices that need to be updated:
  205. updated_indices[i] = 1;
  206. }
  207. }
  208. });
  209. }
  210. for (auto &thread : thread_pool) {
  211. thread.join();
  212. }
  213. for (size_t i = 0; i < NP; ++i) {
  214. if (updated_indices[i]) {
  215. population[i] = trial_vectors[i];
  216. updated_indices[i] = 0;
  217. }
  218. }
  219. }
  220. auto it = std::min_element(cost.begin(), cost.end());
  221. return population[std::distance(cost.begin(), it)];
  222. }
  223. } // namespace boost::math::optimization
  224. #endif