jade.cpp 39 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834
  1. ///
  2. /// @file jade.cc
  3. /// @author Ladutenko Konstantin <kostyfisik at gmail (.) com>
  4. /// @date Thu Aug 15 19:26:53 2013
  5. /// @copyright 2013 Ladutenko Konstantin
  6. /// @section LICENSE
  7. /// This file is part of JADE++.
  8. ///
  9. /// JADE++ is free software: you can redistribute it and/or modify
  10. /// it under the terms of the GNU General Public License as published by
  11. /// the Free Software Foundation, either version 3 of the License, or
  12. /// (at your option) any later version.
  13. ///
  14. /// JADE++ is distributed in the hope that it will be useful,
  15. /// but WITHOUT ANY WARRANTY; without even the implied warranty of
  16. /// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  17. /// GNU General Public License for more details.
  18. ///
  19. /// You should have received a copy of the GNU General Public License
  20. /// along with JADE++. If not, see <http://www.gnu.org/licenses/>.
  21. ///
  22. /// @brief JADE++ is a free (GPLv3+) high performance implementation of
  23. /// adaptive differential evolution optimization algorithm from
  24. /// Jingqiao Zhang and Arthur C. Sanderson book 'Adaptive Differential
  25. /// Evolution. A Robust Approach to Multimodal Problem Optimization'
  26. /// Springer, 2009. Crossover rate was patched according to PMCRADE
  27. /// approach supposed by Jie Li, Wujie Zhu, Mengjun Zhou, and Hua Wang
  28. /// in 'Power Mean Based Crossover Rate Adaptive Differential
  29. /// Evolution' in H. Deng et al. (Eds.): AICI 2011, Part II, LNAI
  30. /// 7003, pp. 34–41, 2011
  31. #include "./jade.h"
  32. #include <mpi.h>
  33. #include <random>
  34. #include <cstdio>
  35. #include <cmath>
  36. #include <algorithm>
  37. #include <map>
  38. #include <iterator>
  39. #include <string>
  40. #include <vector>
  41. #include <iostream>
  42. namespace jade {
  43. /// @todo Replace all simple kError returns with something meangfull. TODO change kError to throw exception
  44. // ********************************************************************** //
  45. // ********************************************************************** //
  46. // ********************************************************************** //
  47. int SubPopulation::Selection(std::vector<Real> crossover_u, long i) {
  48. bool is_evaluated = false;
  49. Real f_current = 0.0;
  50. for (auto f : evaluated_fitness_for_current_vectors_) {
  51. if (f.second == i) {
  52. f_current = f.first;
  53. is_evaluated = true;
  54. }
  55. } // end of searching of pre-evaluated fitness for current individual
  56. if (!is_evaluated) error_status_ = kError;
  57. Real f_best = evaluated_fitness_for_current_vectors_.front().first;
  58. Real f_crossover_u = FitnessFunction(crossover_u);
  59. bool is_success = f_crossover_u > f_current
  60. || f_crossover_u == f_best; //Selected for maxima search
  61. if (is_find_minimum_) is_success = !is_success;
  62. if (!is_success) {
  63. // Case of current x and f were new for current generation.
  64. x_vectors_next_generation_[i] = x_vectors_current_[i];
  65. for (auto &f : evaluated_fitness_for_next_generation_) {
  66. if (f.second == i) f.first = f_current;
  67. } // end of saving old fitness value in new generation.
  68. } else { // if is_success == true
  69. x_vectors_next_generation_[i] = crossover_u;
  70. for (auto &f : evaluated_fitness_for_next_generation_)
  71. if (f.second == i) f.first = f_crossover_u;
  72. to_be_archived_best_A_.push_back(x_vectors_current_[i]);
  73. successful_mutation_parameters_S_F_.push_back(mutation_F_[i]);
  74. successful_crossover_parameters_S_CR_.push_back(crossover_CR_[i]);
  75. // if (process_rank_ == kOutput)
  76. // printf("n%li f_new=%4.2f\n",i,f_crossover_u);
  77. //PrintSingleVector(crossover_u);
  78. } // end of dealing with success crossover
  79. return kDone;
  80. } // end of int SubPopulation::Selection(std::vector<Real> crossover_u);
  81. int SubPopulation::Selection(std::vector<Real> crossover_u, void* param, long i) {
  82. bool is_evaluated = false;
  83. Real f_current = 0.0;
  84. for (auto f : evaluated_fitness_for_current_vectors_) {
  85. if (f.second == i) {
  86. f_current = f.first;
  87. is_evaluated = true;
  88. }
  89. } // end of searching of pre-evaluated fitness for current individual
  90. if (!is_evaluated) error_status_ = kError;
  91. Real f_best = evaluated_fitness_for_current_vectors_.front().first;
  92. Real f_crossover_u = FitnessFunction_p(crossover_u, param);
  93. bool is_success = f_crossover_u > f_current
  94. || f_crossover_u == f_best; //Selected for maxima search
  95. if (is_find_minimum_) is_success = !is_success;
  96. if (!is_success) {
  97. // Case of current x and f were new for current generation.
  98. x_vectors_next_generation_[i] = x_vectors_current_[i];
  99. for (auto &f : evaluated_fitness_for_next_generation_) {
  100. if (f.second == i) f.first = f_current;
  101. } // end of saving old fitness value in new generation.
  102. } else { // if is_success == true
  103. x_vectors_next_generation_[i] = crossover_u;
  104. for (auto &f : evaluated_fitness_for_next_generation_)
  105. if (f.second == i) f.first = f_crossover_u;
  106. to_be_archived_best_A_.push_back(x_vectors_current_[i]);
  107. successful_mutation_parameters_S_F_.push_back(mutation_F_[i]);
  108. successful_crossover_parameters_S_CR_.push_back(crossover_CR_[i]);
  109. // if (process_rank_ == kOutput)
  110. // printf("n%li f_new=%4.2f\n",i,f_crossover_u);
  111. //PrintSingleVector(crossover_u);
  112. } // end of dealing with success crossover
  113. return kDone;
  114. } // end of int SubPopulation::Selection(std::vector<Real> crossover_u);
  115. // ********************************************************************** //
  116. // ********************************************************************** //
  117. // ********************************************************************** //
  118. int SubPopulation::ArchiveCleanUp() {
  119. auto archived_last=archived_best_A_.end();
  120. archived_best_A_.splice(archived_last, to_be_archived_best_A_);
  121. auto size_A = archived_best_A_.size();
  122. long initial_diff = size_A - subpopulation_;
  123. // if (process_rank_ == kOutput)
  124. // printf("diff = %li size_A=%li subpop=%li \n ", initial_diff, size_A, subpopulation_);
  125. // if (initial_diff < 1) return kDone;
  126. // if (process_rank_ == kOutput)
  127. // printf("diff = %li size_A=%li \n ", initial_diff, size_A);
  128. for (long i = 0; i < initial_diff; ++i) {
  129. long index_to_remove = randint(0, size_A - 1);
  130. auto element_A = archived_best_A_.begin();
  131. std::advance(element_A, index_to_remove);
  132. archived_best_A_.erase(element_A);
  133. --size_A;
  134. }
  135. const unsigned long new_size_A = archived_best_A_.size();
  136. if (new_size_A > subpopulation_) error_status_ = kError; //TODO change kError to throw exception
  137. return kDone;
  138. } // end of int SubPopulation:: ArchiveCleanUp();
  139. // ********************************************************************** //
  140. // ********************************************************************** //
  141. // ********************************************************************** //
  142. int SubPopulation::Adaption() {
  143. long elements = 0;
  144. Real sum = 0.0;
  145. for (auto CR : successful_crossover_parameters_S_CR_) {
  146. sum += CR;
  147. ++elements;
  148. } // end of collecting data for mean CR
  149. if (elements == 0) return kDone;
  150. // Should never be reached for symmetric filling of S_CR and S_F.
  151. if (successful_mutation_parameters_S_F_.size() == 0) return kDone;
  152. Real mean_a_CR = sum / static_cast<Real>(elements);
  153. // Original JADE adaption of mu_CR.
  154. if (!isPMCRADE_) {
  155. adaptor_crossover_mu_CR_ =
  156. (1 - adaptation_frequency_c_) * adaptor_crossover_mu_CR_
  157. + adaptation_frequency_c_ * mean_a_CR;
  158. } else {
  159. // PMCRADE patch for mu_CR
  160. Real std_S_CR = 0;
  161. for (auto CR : successful_crossover_parameters_S_CR_)
  162. std_S_CR += pow2(CR - mean_a_CR);
  163. std_S_CR = sqrt(std_S_CR / static_cast<Real>(elements));
  164. const Real PMCRADE_const = 0.07;
  165. if (std_S_CR < PMCRADE_const) {
  166. adaptor_crossover_mu_CR_ =
  167. (1 - adaptation_frequency_c_) * adaptor_crossover_mu_CR_
  168. + adaptation_frequency_c_ * mean_a_CR;
  169. } else {
  170. Real mean_pow2_CR = 0.0;
  171. for (auto CR : successful_crossover_parameters_S_CR_)
  172. mean_pow2_CR += pow2(CR);
  173. mean_pow2_CR = sqrt(mean_pow2_CR/static_cast<Real>(elements));
  174. adaptor_crossover_mu_CR_ =
  175. (1 - adaptation_frequency_c_) * adaptor_crossover_mu_CR_
  176. + adaptation_frequency_c_ * mean_pow2_CR;
  177. }
  178. // end of PMCRADE patch
  179. } // end of if isPMCRADE_
  180. Real sum_F = 0.0, sum_F2 = 0.0;
  181. for (auto F : successful_mutation_parameters_S_F_) {
  182. sum_F += F;
  183. sum_F2 += F*F;
  184. } // end of collection data for Lehmer mean F
  185. Real mean_l_F = sum_F2 / sum_F;
  186. adaptor_mutation_mu_F_ =
  187. (1 - adaptation_frequency_c_) * adaptor_mutation_mu_F_
  188. + adaptation_frequency_c_ * mean_l_F;
  189. return kDone;
  190. } // end of int SubPopulation::Adaption();
  191. // ********************************************************************** //
  192. // ********************************************************************** //
  193. // ********************************************************************** //
  194. int SubPopulation::RunOptimization() {
  195. //if (process_rank_ == kOutput) printf("Start optimization..\n");
  196. if (error_status_) return error_status_;
  197. adaptor_mutation_mu_F_ = 0.5;
  198. adaptor_crossover_mu_CR_ = 0.5;
  199. archived_best_A_.clear();
  200. CreateInitialPopulation();
  201. x_vectors_next_generation_ = x_vectors_current_;
  202. EvaluateCurrentVectors();
  203. evaluated_fitness_for_next_generation_ =
  204. evaluated_fitness_for_current_vectors_;
  205. for (long g = 0; g < total_generations_max_; ++g) {
  206. if (process_rank_ == kOutput && g%100 == 0) printf("%li\n",g);
  207. to_be_archived_best_A_.clear();
  208. successful_mutation_parameters_S_F_.clear();
  209. successful_crossover_parameters_S_CR_.clear();
  210. // //debug section
  211. // if (process_rank_ == kOutput)
  212. // printf("============== Generation %li =============\n", g);
  213. PrintPopulation();
  214. //PrintEvaluated();
  215. std::cin.get();
  216. //end of debug section
  217. for (unsigned long i = 0; i < subpopulation_; ++i) {
  218. SetCRiFi(i);
  219. std::vector<Real> mutated_v, crossover_u;
  220. mutated_v = Mutation(i);
  221. crossover_u = Crossover(mutated_v, i);
  222. Selection(crossover_u, i);
  223. } // end of for all individuals in subpopulation
  224. ArchiveCleanUp();
  225. Adaption();
  226. x_vectors_current_.swap(x_vectors_next_generation_);
  227. evaluated_fitness_for_current_vectors_
  228. .swap(evaluated_fitness_for_next_generation_);
  229. SortEvaluatedCurrent();
  230. if (error_status_) return error_status_;
  231. } // end of stepping generations
  232. PrintPopulation();
  233. PrintEvaluated();
  234. return kDone;
  235. } // end of int SubPopulation::RunOptimization()
  236. int SubPopulation::RunOptimization(std::ostream &out, void* param) {
  237. //if (process_rank_ == kOutput) printf("Start optimization..\n");
  238. if (error_status_) return error_status_;
  239. adaptor_mutation_mu_F_ = 0.5;
  240. adaptor_crossover_mu_CR_ = 0.5;
  241. archived_best_A_.clear();
  242. CreateInitialPopulation();
  243. x_vectors_next_generation_ = x_vectors_current_;
  244. EvaluateCurrentVectors(param);
  245. evaluated_fitness_for_next_generation_ = evaluated_fitness_for_current_vectors_;
  246. for (long g = 0; g < total_generations_max_; ++g) {
  247. if (process_rank_ == kOutput && g%100 == 0) printf("%li\n",g);
  248. to_be_archived_best_A_.clear();
  249. successful_mutation_parameters_S_F_.clear();
  250. successful_crossover_parameters_S_CR_.clear();
  251. PrintPopulation(g, out);
  252. for (unsigned long i = 0; i < subpopulation_; ++i) {
  253. SetCRiFi(i);
  254. std::vector<Real> mutated_v, crossover_u;
  255. mutated_v = Mutation(i);
  256. crossover_u = Crossover(mutated_v, i);
  257. Selection(crossover_u, param, i);
  258. } // end of for all individuals in subpopulation
  259. ArchiveCleanUp();
  260. Adaption();
  261. x_vectors_current_.swap(x_vectors_next_generation_);
  262. evaluated_fitness_for_current_vectors_
  263. .swap(evaluated_fitness_for_next_generation_);
  264. SortEvaluatedCurrent();
  265. if (error_status_) return error_status_;
  266. } // end of stepping generations
  267. PrintPopulation();
  268. PrintEvaluated();
  269. return kDone;
  270. } // end of int SubPopulation::RunOptimization()
  271. // ********************************************************************** //
  272. // ********************************************************************** //
  273. // ********************************************************************** //
  274. std::vector<Real> SubPopulation::Crossover(std::vector<Real> mutation_v, long i) {
  275. const Real CR_i = crossover_CR_[i];
  276. std::vector<Real> crossover_u, x_current;
  277. crossover_u.resize(dimension_);
  278. x_current = x_vectors_current_.at(i);
  279. unsigned long j_rand = randint(0, dimension_ - 1);
  280. for (unsigned long c = 0; c < dimension_; ++c) {
  281. if (c == j_rand || rand(0,1) < CR_i)
  282. crossover_u[c] = mutation_v[c];
  283. else
  284. crossover_u[c] = x_current[c];
  285. }
  286. // //debug section
  287. // if (process_rank_ == kOutput)
  288. // printf("x -> v -> u with CR_i=%4.2f j_rand=%li\n", CR_i, j_rand);
  289. // PrintSingleVector(mutation_v);
  290. // PrintSingleVector(x_current);
  291. // PrintSingleVector(crossover_u);
  292. // //end of debug section
  293. return crossover_u;
  294. } // end of std::vector<Real> SubPopulation::Crossover();
  295. // ********************************************************************** //
  296. // ********************************************************************** //
  297. // ********************************************************************** //
  298. std::vector<Real> SubPopulation::Mutation(long i) {
  299. std::vector<Real> mutation_v, x_best_current, x_random_current,
  300. x_random_archive_and_current, x_current;
  301. x_current = x_vectors_current_.at(i);
  302. x_best_current = GetXpBestCurrent();
  303. // //debug
  304. // if (process_rank_ == kOutput) printf("x_best: ");
  305. // PrintSingleVector(x_best_current);
  306. long index_of_random_current = -1;
  307. x_random_current = GetXRandomCurrent(&index_of_random_current, i);
  308. // //debug
  309. // if (process_rank_ == kOutput) printf("x_random: ");
  310. // PrintSingleVector(x_random_current);
  311. x_random_archive_and_current = GetXRandomArchiveAndCurrent(index_of_random_current, i);
  312. // //debug
  313. // if (process_rank_ == kOutput) printf("x_random with archive: ");
  314. // PrintSingleVector(x_random_archive_and_current);
  315. mutation_v.resize(dimension_);
  316. Real F_i = mutation_F_[i];
  317. for (unsigned long c = 0; c < dimension_; ++c) {
  318. // Mutation
  319. mutation_v[c] = x_current[c]
  320. + F_i * (x_best_current[c] - x_current[c])
  321. + F_i * (x_random_current[c]
  322. - x_random_archive_and_current[c]);
  323. // Bounds control
  324. if (mutation_v[c] > x_ubound_[c])
  325. mutation_v[c] = (x_ubound_[c] + x_current[c])/2;
  326. if (mutation_v[c] < x_lbound_[c])
  327. mutation_v[c] = (x_lbound_[c] + x_current[c])/2;
  328. }
  329. // //debug section
  330. // int isSame = 888, isSame2 = 7777, isSame3 =11111;
  331. // for (long c = 0; c < dimension_; ++c) {
  332. // Real norm = std::abs(mutation_v[c] - x_current[c]);
  333. // Real norm2 = std::abs(x_random_current[c] - x_current[c]);
  334. // Real norm3 = std::abs(x_random_current[c]
  335. // - x_random_archive_and_current[c]);
  336. // if ( norm > 0.0001) isSame = 0;
  337. // if ( norm2 > 0.0001) isSame2 = 0;
  338. // if ( norm3 > 0.0001) isSame3 = 0;
  339. // }
  340. // if (process_rank_ == kOutput) printf("mutation_v%i%i%i: ",isSame, isSame2, isSame3);
  341. // PrintSingleVector(mutation_v);
  342. // if (process_rank_ == kOutput) printf("current _v%i%i%i: ",isSame, isSame2, isSame3);
  343. // PrintSingleVector(x_current);
  344. // if (process_rank_ == kOutput)
  345. // printf(" -> f = %4.2f F_i=%4.2f\n",
  346. // FitnessFunction(mutation_v), F_i);
  347. // //end of debug section
  348. return mutation_v;
  349. } // end of std::vector<Real> SubPopulation::Mutation();
  350. // ********************************************************************** //
  351. // ********************************************************************** //
  352. // ********************************************************************** //
  353. std::vector<Real> SubPopulation::GetXpBestCurrent() {
  354. const unsigned long n_best_total = static_cast<long>
  355. (floor(subpopulation_ * best_share_p_ ));
  356. if (n_best_total == subpopulation_) error_status_ = kError; //TODO change kError to throw exception
  357. long best_n = randint(0, n_best_total);
  358. long best_n_index = -1, i = 0;
  359. for (auto x : evaluated_fitness_for_current_vectors_) {
  360. if (i == best_n) {
  361. best_n_index = x.second;
  362. break;
  363. }
  364. ++i;
  365. }
  366. return x_vectors_current_.at(best_n_index);
  367. } // end of std::vector<Real> SubPopulation::GetXpBestCurrent();
  368. // ********************************************************************** //
  369. // ********************************************************************** //
  370. // ********************************************************************** //
  371. std::vector<Real> SubPopulation::GetXRandomCurrent(long *index, long forbidden_index) {
  372. long random_n = randint(0, subpopulation_-1);
  373. while (random_n == forbidden_index) random_n = randint(0, subpopulation_-1);
  374. (*index) = random_n;
  375. return x_vectors_current_.at(random_n);
  376. } // end of std::vector<Real> SubPopulation::GetXRandomCurrent()
  377. // ********************************************************************** //
  378. // ********************************************************************** //
  379. // ********************************************************************** //
  380. std::vector<Real> SubPopulation::GetXRandomArchiveAndCurrent
  381. (unsigned long forbidden_index1, unsigned long forbidden_index2) {
  382. long archive_size = archived_best_A_.size();
  383. unsigned long random_n = randint(0, subpopulation_ + archive_size - 1);
  384. while (random_n == forbidden_index1 || random_n == forbidden_index2)
  385. random_n = randint(0, subpopulation_ + archive_size - 1);
  386. if (random_n < subpopulation_) return x_vectors_current_.at(random_n);
  387. random_n -= subpopulation_;
  388. unsigned long i = 0;
  389. for (auto x : archived_best_A_) {
  390. if (i == random_n) {
  391. // //debug
  392. // if (process_rank_ == kOutput) printf("Using Archive!!\n");
  393. return x;
  394. }
  395. ++i;
  396. } // end of selecting from archive
  397. error_status_ = kError; //TODO change kError to throw exception
  398. std::vector<Real> x;
  399. return x;
  400. } // end of std::vector<Real> SubPopulation::GetXRandomArchiveAndCurrent()
  401. // ********************************************************************** //
  402. // ********************************************************************** //
  403. // ********************************************************************** //
  404. int SubPopulation::SetCRiFi(long i) {
  405. long k = 0;
  406. while (1) {
  407. mutation_F_[i] = randc(adaptor_mutation_mu_F_, 0.1);
  408. if (mutation_F_[i] > 1) {
  409. mutation_F_[i] = 1;
  410. break;
  411. }
  412. if (mutation_F_[i] > 0) break;
  413. ++k;
  414. if (k > 10) {
  415. mutation_F_[i] = 0.001;
  416. break;
  417. }
  418. if (k > 100) printf("k");
  419. }
  420. crossover_CR_[i] = randn(adaptor_crossover_mu_CR_,0.1);
  421. if (crossover_CR_[i] > 1) crossover_CR_[i] = 1;
  422. if (crossover_CR_[i] < 0) crossover_CR_[i] = 0;
  423. return kDone;
  424. } // end of int SubPopulation::SetCRiFi(long i)
  425. // ********************************************************************** //
  426. // ********************************************************************** //
  427. // ********************************************************************** //
  428. /// %todo Change returned kError to meangfull error code.
  429. int SubPopulation::Init(long total_population, long dimension) {// NOLINT
  430. total_population_ = total_population;
  431. if (total_population_ < 1)
  432. throw std::invalid_argument("You should set population > 0!");
  433. dimension_ = dimension;
  434. if (dimension_ < 1)
  435. throw std::invalid_argument("You should set dimension > 0!");
  436. MPI_Comm_rank(MPI_COMM_WORLD, &process_rank_);
  437. MPI_Comm_size(MPI_COMM_WORLD, &number_of_processes_);
  438. if (process_rank_ < 0)
  439. throw std::invalid_argument("MPI problem: process_rank_ < 0!");
  440. if (number_of_processes_ < 1)
  441. throw std::invalid_argument("MPI problem: number_of_processes_ < 1!");
  442. std::random_device rd;
  443. generator_.seed(rd());
  444. // //debug
  445. // CheckRandom();
  446. subpopulation_ = total_population;
  447. current_generation_ = 0;
  448. x_vectors_current_.resize(subpopulation_);
  449. for (auto &x : x_vectors_current_) x.resize(dimension_);
  450. x_vectors_next_generation_.resize(subpopulation_);
  451. for (auto &x : x_vectors_next_generation_) x.resize(dimension_);
  452. evaluated_fitness_for_current_vectors_.resize(subpopulation_);
  453. // //debug
  454. // if (process_rank_ == kOutput) printf("%i, x1 size = %li \n", process_rank_, x_vectors_current_.size());
  455. x_lbound_.resize(dimension_);
  456. x_ubound_.resize(dimension_);
  457. mutation_F_.resize(subpopulation_);
  458. crossover_CR_.resize(subpopulation_);
  459. return kDone;
  460. } // end of void SubPopulation::Init()
  461. // ********************************************************************** //
  462. // ********************************************************************** //
  463. // ********************************************************************** //
  464. int SubPopulation::PrintPopulation() {
  465. if (process_rank_ == kOutput) {
  466. printf("\n");
  467. for (auto x : evaluated_fitness_for_current_vectors_) {
  468. long n = x.second;
  469. Real fitness = x.first;
  470. printf("%6.2f:% 3li||", fitness, n);
  471. for (unsigned long c = 0; c < dimension_; ++c)
  472. printf(" %+7.2f ", x_vectors_current_[n][c]);
  473. printf("\n");
  474. }
  475. // for (long i = 0; i < subpopulation_; ++i) {
  476. // printf("n%li:", i);
  477. // for (long c = 0; c < dimension_; ++c) {
  478. // printf(" %5.2f ", x_vectors_current_[i][c]);
  479. // }
  480. // printf("\n");
  481. // } // end of for each individual
  482. } // end of if output
  483. return kDone;
  484. } // end of int SubPupulation::PrintPopulation()
  485. int SubPopulation::PrintPopulation(long niter, std::ostream &out) {
  486. if (process_rank_ == kOutput) {
  487. //printf("\n");
  488. auto p = evaluated_fitness_for_current_vectors_.begin();
  489. out << niter << " " << (*p).first;
  490. for (unsigned long c = 0; c < dimension_; ++c) out << " " << x_vectors_current_[(*p).second][c];
  491. out << std::endl;
  492. } // end of if output
  493. return kDone;
  494. } // end of int SubPupulation::PrintPopulation()
  495. // ********************************************************************** //
  496. // ********************************************************************** //
  497. // ********************************************************************** //
  498. int SubPopulation::PrintEvaluated() {
  499. if (process_rank_ == kOutput) {
  500. for (auto x : evaluated_fitness_for_current_vectors_)
  501. printf("%li:%4.2f ", x.second, x.first);
  502. printf("\n");
  503. } // end of if output
  504. return kDone;
  505. } // end of int SubPopulation::PrintEvaluated()
  506. // ********************************************************************** //
  507. // ********************************************************************** //
  508. // ********************************************************************** //
  509. int SubPopulation::PrintSingleVector(std::vector<Real> x) {
  510. if (process_rank_ == kOutput) {
  511. for (auto c : x) printf("%5.2f ", c);
  512. printf("\n");
  513. } // end of output
  514. return kDone;
  515. } // end of int SubPopulation::PrintSingleVector(std::vector<Real> x)
  516. // ********************************************************************** //
  517. // ********************************************************************** //
  518. // ********************************************************************** //
  519. void SubPopulation::SetFeed(std::vector<std::vector<Real> > x_feed_vectors) {
  520. isFeed_ = true;
  521. if (dimension_ < 1)
  522. throw std::invalid_argument("You should set dimension before feed!");
  523. for (auto x : x_feed_vectors)
  524. if (x.size() != dimension_)
  525. throw std::invalid_argument
  526. ("Feed and optimization dimenstion should be the same size!");
  527. x_feed_vectors_.clear();
  528. for (auto &x : x_feed_vectors)
  529. x_feed_vectors_.push_back(x);
  530. } // end of void SubPopulation::SetFeed()
  531. // ********************************************************************** //
  532. // ********************************************************************** //
  533. // ********************************************************************** //
  534. int SubPopulation::CreateInitialPopulation() {
  535. if ((subpopulation_ - x_feed_vectors_.size()) <1)
  536. throw std::invalid_argument("Too large feed!");
  537. x_vectors_current_.resize(subpopulation_ - x_feed_vectors_.size());
  538. for (auto &x : x_vectors_current_)
  539. for (unsigned long i = 0; i < dimension_; ++i) {
  540. if (x_lbound_[i] > x_ubound_[i])
  541. throw std::invalid_argument("Wrong order of bounds!");
  542. x[i] = rand(x_lbound_[i], x_ubound_[i]); // NOLINT
  543. } // end of for each dimension
  544. // //debug
  545. // for (auto x : x_vectors_current_[0]) if (process_rank_ == kOutput) printf("%g ",x);
  546. for (auto x: x_feed_vectors_) {
  547. x_vectors_current_.push_back(x);
  548. if (process_rank_ == kOutput) {
  549. printf("--=-- Feed:\n");
  550. for (auto index:x) printf(" %+7.2f", index);
  551. printf("\n");
  552. }
  553. }
  554. if (x_vectors_current_.size() != subpopulation_)
  555. throw std::invalid_argument("Population is not full after feed!");
  556. x_feed_vectors_.clear();
  557. // if (process_rank_ == kOutput) {
  558. // printf("==--== Initial population:\n");
  559. // EvaluateCurrentVectors();
  560. // PrintPopulation();
  561. // }
  562. return kDone;
  563. } // end of int SubPopulation::CreateInitialPopulation()
  564. // ********************************************************************** //
  565. // ********************************************************************** //
  566. // ********************************************************************** //
  567. int SubPopulation::SortEvaluatedCurrent() {
  568. evaluated_fitness_for_current_vectors_
  569. .sort([=](const std::pair<Real, long>& a, // NOLINT
  570. const std::pair<Real, long>& b) { // NOLINT
  571. bool cmp = a.first < b.first;
  572. if (is_find_minimum_) return cmp;
  573. return !cmp;
  574. });
  575. return kDone;
  576. } // end of int SubPopulation::SortEvaluatedCurrent()
  577. // ********************************************************************** //
  578. // ********************************************************************** //
  579. // ********************************************************************** //
  580. int SubPopulation::EvaluateCurrentVectors() {
  581. evaluated_fitness_for_current_vectors_.clear();
  582. for (unsigned long i = 0; i < subpopulation_; ++i) { // NOLINT
  583. auto tmp = std::make_pair(FitnessFunction(x_vectors_current_[i]), i);
  584. evaluated_fitness_for_current_vectors_.push_back(tmp);
  585. }
  586. SortEvaluatedCurrent();
  587. // //debug
  588. // if (process_rank_ == kOutput) printf("\n After ");
  589. // for (auto val : evaluated_fitness_for_current_vectors_)
  590. // if (process_rank_ == kOutput) printf("%g ", val.first);
  591. return kDone;
  592. } // end of int SubPopulation::EvaluateCurrentVectors()
  593. int SubPopulation::EvaluateCurrentVectors(void* param) {
  594. evaluated_fitness_for_current_vectors_.clear();
  595. for (unsigned long i = 0; i < subpopulation_; ++i) { // NOLINT
  596. auto tmp = std::make_pair(FitnessFunction_p(x_vectors_current_[i], param), i);
  597. evaluated_fitness_for_current_vectors_.push_back(tmp);
  598. }
  599. SortEvaluatedCurrent();
  600. // //debug
  601. // if (process_rank_ == kOutput) printf("\n After ");
  602. // for (auto val : evaluated_fitness_for_current_vectors_)
  603. // if (process_rank_ == kOutput) printf("%g ", val.first);
  604. return kDone;
  605. } // end of int SubPopulation::EvaluateCurrentVectors()
  606. // ********************************************************************** //
  607. // ********************************************************************** //
  608. // ********************************************************************** //
  609. void SubPopulation::SetAllBoundsVectors
  610. (std::vector<Real> lbound, std::vector<Real> ubound) {
  611. x_lbound_.clear();
  612. x_ubound_.clear();
  613. for (auto x : lbound) x_lbound_.push_back(x);
  614. for (auto x : ubound) x_ubound_.push_back(x);
  615. }
  616. // ********************************************************************** //
  617. // ********************************************************************** //
  618. // ********************************************************************** //
  619. int SubPopulation::SetAllBounds(Real lbound, Real ubound) {
  620. if (lbound >= ubound) {
  621. error_status_ = kError;
  622. return kError; //TODO change kError to throw exception
  623. }
  624. for (auto &x : x_lbound_) x = lbound;
  625. for (auto &x : x_ubound_) x = ubound;
  626. return kDone;
  627. } // end of int SubPopulation::SetAllBounds(Real lbound, Real ubound)
  628. // ********************************************************************** //
  629. // ********************************************************************** //
  630. // ********************************************************************** //
  631. Real SubPopulation::randn(Real mean, Real stddev) {
  632. std::normal_distribution<Real> distribution(mean, stddev);
  633. return distribution(generator_);
  634. } // end of Real SubPopulation::randn(Real mean, Real stddev)
  635. // ********************************************************************** //
  636. // ********************************************************************** //
  637. // ********************************************************************** //
  638. Real SubPopulation::randc(Real location, Real scale) {
  639. std::cauchy_distribution<Real> distribution(location, scale);
  640. return distribution(generator_);
  641. } // end of Real SubPopulation::randc(Real location, Real scale)
  642. // ********************************************************************** //
  643. // ********************************************************************** //
  644. // ********************************************************************** //
  645. long SubPopulation::randint(long lbound, long ubound) { // NOLINT
  646. std::uniform_int_distribution<long> distribution(lbound, ubound); // NOLINT
  647. return distribution(generator_);
  648. }
  649. // ********************************************************************** //
  650. // ********************************************************************** //
  651. // ********************************************************************** //
  652. Real SubPopulation::rand(Real lbound, Real ubound) { // NOLINT
  653. std::uniform_real_distribution<Real> distribution(lbound, ubound);
  654. return distribution(generator_);
  655. } // end of Real rand(Real lbound, Real ubound)
  656. // ********************************************************************** //
  657. // ********************************************************************** //
  658. // ********************************************************************** //
  659. void SubPopulation::CheckRandom() {
  660. std::map<int, int> hist_n, hist_c, hist_int, hist;
  661. int factor = 1000;
  662. for (int n = 0; n < 100*factor; ++n) {
  663. ++hist_n[std::round(randn(0, 2))];
  664. ++hist_c[std::round(randc(0, 2))];
  665. ++hist_int[std::round(randint(0, 10))];
  666. ++hist[std::round(rand(0, 15))]; // NOLINT
  667. }
  668. if (process_rank_ == kOutput) {
  669. printf("Normal (0,2)\n");
  670. for (auto p : hist_n) {
  671. if (p.second > factor)
  672. printf("%i: % 4i %s\n", process_rank_, p.first,
  673. std::string(p.second/factor, '*').c_str());
  674. } // end of for p in hist_n
  675. printf("Cauchy (0,2)\n");
  676. for (auto p : hist_c) {
  677. if (p.second > factor)
  678. printf("%i: % 4i %s\n", process_rank_, p.first,
  679. std::string(p.second/factor, '*').c_str());
  680. } // end of for p in hist_c
  681. printf("Uniform int (0,10)\n");
  682. for (auto p : hist_int) {
  683. if (p.second > factor)
  684. printf("%i: % 4i %s\n", process_rank_, p.first,
  685. std::string(p.second/factor, '*').c_str());
  686. } // end of for p in hist_int
  687. printf("Uniform Real (0,15)\n");
  688. for (auto p : hist) {
  689. if (p.second > factor)
  690. printf("%i: % 4i %s\n", process_rank_, p.first,
  691. std::string(p.second/factor, '*').c_str());
  692. } // end of for p in hist
  693. } // end of if current process_rank_ == kOutput
  694. }
  695. // ********************************************************************** //
  696. // ********************************************************************** //
  697. // ********************************************************************** //
  698. int SubPopulation::SetBestShareP(Real p) {
  699. if (p < 0 || p > 1) {
  700. error_status_ = kError;
  701. return kError; //TODO change kError to throw exception
  702. }
  703. best_share_p_ = p;
  704. return kDone;
  705. } // end of int SubPopulation::SetBestShareP(Real p)
  706. // ********************************************************************** //
  707. // ********************************************************************** //
  708. // ********************************************************************** //
  709. int SubPopulation::SetAdapitonFrequencyC(Real c) {
  710. if (c < 0 || c > 1) {
  711. error_status_ = kError;
  712. return kError; //TODO change kError to throw exception
  713. }
  714. adaptation_frequency_c_ = c;
  715. return kDone;
  716. } // end of int SubPopulation::SetAdapitonFrequency(Real c)
  717. // ********************************************************************** //
  718. // ********************************************************************** //
  719. // ********************************************************************** //
  720. int SubPopulation::SetDistributionLevel(int level) {
  721. distribution_level_ = level;
  722. return kDone;
  723. }
  724. // end of int SubPopulation::SetDistributionLevel(int level)
  725. // ********************************************************************** //
  726. // ********************************************************************** //
  727. // ********************************************************************** //
  728. int SubPopulation::PrintParameters(std::string comment) {
  729. if (process_rank_ == 0) {
  730. printf("#%s dim=%li NP=%li(of %li) p=%4.2f c=%4.2f generation=%li\n",
  731. comment.c_str(),dimension_, subpopulation_, total_population_,
  732. best_share_p_, adaptation_frequency_c_, total_generations_max_
  733. );
  734. fflush(stdout);
  735. } // end of output
  736. return kDone;
  737. } // end of int SubPopulation::PrintParameters(std::string comment)
  738. // ********************************************************************** //
  739. // ********************************************************************** //
  740. // ********************************************************************** //
  741. int SubPopulation::PrintResult(std::string comment) {
  742. if (distribution_level_ == 0) {
  743. auto x = evaluated_fitness_for_current_vectors_.front();
  744. std::vector<Real> to_send {x.first};
  745. //printf("%8.5g\n ", x.first);
  746. AllGatherVectorDouble(to_send);
  747. if (process_rank_ == 0) {
  748. Real sum = 0;
  749. Real size = static_cast<Real>(recieve_Real_.size());
  750. for (auto x : recieve_Real_) sum += x;
  751. Real mean = sum/size;
  752. Real sigma = 0;
  753. for (auto x : recieve_Real_) sigma += pow2(x - mean);
  754. sigma = sqrt(sigma/size);
  755. printf("%s gen%li, mean %4.1e (stddev %4.1e = %3.2g %%) runs(%g)\n",
  756. comment.c_str(), total_generations_max_, mean,sigma,
  757. sigma*100.0/mean,size);
  758. // for (auto x : recieve_Real_)
  759. // printf("%18.15g\n", x);
  760. }
  761. }
  762. return kDone;
  763. } // end of int SubPopulation::PrintResult()
  764. // ********************************************************************** //
  765. // ********************************************************************** //
  766. // ********************************************************************** //
  767. std::vector<Real> SubPopulation::GetFinalFitness() {
  768. recieve_Real_.clear();
  769. if (distribution_level_ == 0) {
  770. auto x = evaluated_fitness_for_current_vectors_.front();
  771. std::vector<Real> to_send {x.first};
  772. AllGatherVectorDouble(to_send);
  773. }
  774. return recieve_Real_;
  775. } // end of sdt::vector<Real> SubPopulation::GetFinalFitness();
  776. // ********************************************************************** //
  777. // ********************************************************************** //
  778. // ********************************************************************** //
  779. std::vector<Real> SubPopulation::GetBest(Real *best_fitness) {
  780. auto x = evaluated_fitness_for_current_vectors_.front();
  781. (*best_fitness) = x.first;
  782. return x_vectors_current_[x.second];
  783. } // end of std::vector<Real> SubPopulation::GetBest(Real *best_fitness)
  784. // ********************************************************************** //
  785. // ********************************************************************** //
  786. // ********************************************************************** //
  787. std::vector<Real> SubPopulation::GetWorst(Real *worst_fitness) {
  788. auto x = evaluated_fitness_for_current_vectors_.back();
  789. (*worst_fitness) = x.first;
  790. return x_vectors_current_[x.second];
  791. } // end of std::vector<Real> SubPopulation::GetWorst(Real *worst_fitness)
  792. // ********************************************************************** //
  793. // ********************************************************************** //
  794. // ********************************************************************** //
  795. int SubPopulation::AllGatherVectorDouble(std::vector<Real> to_send) {
  796. long size_single = to_send.size();
  797. long size_all = size_single * number_of_processes_;
  798. recieve_Real_.clear();
  799. recieve_Real_.resize(size_all);
  800. MPI_Allgather(&to_send.front(), size_single, MPI_DOUBLE,
  801. &recieve_Real_.front(), size_single, MPI_DOUBLE,
  802. MPI_COMM_WORLD);
  803. return kDone;
  804. } // end of int SubPopulation::AllGatherVectorDouble(std::vector<Real> to_send);
  805. // ********************************************************************** //
  806. // ********************************************************************** //
  807. // ********************************************************************** //
  808. int SubPopulation::AllGatherVectorLong(std::vector<long> to_send) {
  809. long size_single = to_send.size();
  810. long size_all = size_single * number_of_processes_;
  811. recieve_long_.clear();
  812. recieve_long_.resize(size_all);
  813. MPI_Allgather(&to_send.front(), size_single, MPI_LONG,
  814. &recieve_long_.front(), size_single, MPI_LONG,
  815. MPI_COMM_WORLD);
  816. return kDone;
  817. } // end of int SubPopulation::AllGatherVectorLong(std::vector<long> to_send);
  818. // ********************************************************************** //
  819. // ********************************************************************** //
  820. // ********************************************************************** //
  821. } // end of namespace jade