jade.cc 35 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734
  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. namespace jade {
  42. /// @todo Replace all simple kError returns with something meangfull. TODO change kError to throw exception
  43. // ********************************************************************** //
  44. // ********************************************************************** //
  45. // ********************************************************************** //
  46. int SubPopulation::Selection(std::vector<double> crossover_u, long i) {
  47. bool is_evaluated = false;
  48. double f_current = 0.0;
  49. for (auto f : evaluated_fitness_for_current_vectors_) {
  50. if (f.second == i) {
  51. f_current = f.first;
  52. is_evaluated = true;
  53. }
  54. } // end of searching of pre-evaluated fitness for current individual
  55. if (!is_evaluated) error_status_ = kError;
  56. double f_best = evaluated_fitness_for_current_vectors_.front().first;
  57. double f_crossover_u = FitnessFunction(crossover_u);
  58. bool is_success = f_crossover_u > f_current
  59. || f_crossover_u == f_best; //Selected for maxima search
  60. if (is_find_minimum_) is_success = !is_success;
  61. if (!is_success) {
  62. // Case of current x and f were new for current generation.
  63. x_vectors_next_generation_[i] = x_vectors_current_[i];
  64. for (auto &f : evaluated_fitness_for_next_generation_) {
  65. if (f.second == i) f.first = f_current;
  66. } // end of saving old fitness value in new generation.
  67. } else { // if is_success == true
  68. x_vectors_next_generation_[i] = crossover_u;
  69. for (auto &f : evaluated_fitness_for_next_generation_)
  70. if (f.second == i) f.first = f_crossover_u;
  71. to_be_archived_best_A_.push_back(x_vectors_current_[i]);
  72. successful_mutation_parameters_S_F_.push_back(mutation_F_[i]);
  73. successful_crossover_parameters_S_CR_.push_back(crossover_CR_[i]);
  74. // if (process_rank_ == kOutput)
  75. // printf("n%li f_new=%4.2f\n",i,f_crossover_u);
  76. //PrintSingleVector(crossover_u);
  77. } // end of dealing with success crossover
  78. return kDone;
  79. } // end of int SubPopulation::Selection(std::vector<double> crossover_u);
  80. // ********************************************************************** //
  81. // ********************************************************************** //
  82. // ********************************************************************** //
  83. int SubPopulation::ArchiveCleanUp() {
  84. auto archived_last=archived_best_A_.end();
  85. archived_best_A_.splice(archived_last, to_be_archived_best_A_);
  86. auto size_A = archived_best_A_.size();
  87. long initial_diff = size_A - subpopulation_;
  88. // if (process_rank_ == kOutput)
  89. // printf("diff = %li size_A=%li subpop=%li \n ", initial_diff, size_A, subpopulation_);
  90. // if (initial_diff < 1) return kDone;
  91. // if (process_rank_ == kOutput)
  92. // printf("diff = %li size_A=%li \n ", initial_diff, size_A);
  93. for (long i = 0; i < initial_diff; ++i) {
  94. long index_to_remove = randint(0, size_A - 1);
  95. auto element_A = archived_best_A_.begin();
  96. std::advance(element_A, index_to_remove);
  97. archived_best_A_.erase(element_A);
  98. --size_A;
  99. }
  100. const long new_size_A = archived_best_A_.size();
  101. if (new_size_A > subpopulation_) error_status_ = kError; //TODO change kError to throw exception
  102. return kDone;
  103. } // end of int SubPopulation:: ArchiveCleanUp();
  104. // ********************************************************************** //
  105. // ********************************************************************** //
  106. // ********************************************************************** //
  107. int SubPopulation::Adaption() {
  108. long elements = 0;
  109. double sum = 0.0;
  110. for (auto CR : successful_crossover_parameters_S_CR_) {
  111. sum += CR;
  112. ++elements;
  113. } // end of collecting data for mean CR
  114. if (elements == 0) return kDone;
  115. // Should never be reached for symmetric filling of S_CR and S_F.
  116. if (successful_mutation_parameters_S_F_.size() == 0) return kDone;
  117. double mean_a_CR = sum / static_cast<double>(elements);
  118. // Original JADE adaption of mu_CR.
  119. if (!isPMCRADE_) {
  120. adaptor_crossover_mu_CR_ =
  121. (1 - adaptation_frequency_c_) * adaptor_crossover_mu_CR_
  122. + adaptation_frequency_c_ * mean_a_CR;
  123. } else {
  124. // PMCRADE patch for mu_CR
  125. double std_S_CR = 0;
  126. for (auto CR : successful_crossover_parameters_S_CR_)
  127. std_S_CR += pow2(CR - mean_a_CR);
  128. std_S_CR = sqrt(std_S_CR / static_cast<double>(elements));
  129. const double PMCRADE_const = 0.07;
  130. if (std_S_CR < PMCRADE_const) {
  131. adaptor_crossover_mu_CR_ =
  132. (1 - adaptation_frequency_c_) * adaptor_crossover_mu_CR_
  133. + adaptation_frequency_c_ * mean_a_CR;
  134. } else {
  135. double mean_pow2_CR = 0.0;
  136. for (auto CR : successful_crossover_parameters_S_CR_)
  137. mean_pow2_CR += pow2(CR);
  138. mean_pow2_CR = sqrt(mean_pow2_CR/static_cast<double>(elements));
  139. adaptor_crossover_mu_CR_ =
  140. (1 - adaptation_frequency_c_) * adaptor_crossover_mu_CR_
  141. + adaptation_frequency_c_ * mean_pow2_CR;
  142. }
  143. // end of PMCRADE patch
  144. } // end of if isPMCRADE_
  145. double sum_F = 0.0, sum_F2 = 0.0;
  146. for (auto F : successful_mutation_parameters_S_F_) {
  147. sum_F += F;
  148. sum_F2 += F*F;
  149. } // end of collection data for Lehmer mean F
  150. double mean_l_F = sum_F2 / sum_F;
  151. adaptor_mutation_mu_F_ =
  152. (1 - adaptation_frequency_c_) * adaptor_mutation_mu_F_
  153. + adaptation_frequency_c_ * mean_l_F;
  154. return kDone;
  155. } // end of int SubPopulation::Adaption();
  156. // ********************************************************************** //
  157. // ********************************************************************** //
  158. // ********************************************************************** //
  159. int SubPopulation::RunOptimization() {
  160. //if (process_rank_ == kOutput) printf("Start optimization..\n");
  161. if (error_status_) return error_status_;
  162. adaptor_mutation_mu_F_ = 0.5;
  163. adaptor_crossover_mu_CR_ = 0.5;
  164. archived_best_A_.clear();
  165. CreateInitialPopulation();
  166. x_vectors_next_generation_ = x_vectors_current_;
  167. EvaluateCurrentVectors();
  168. evaluated_fitness_for_next_generation_ =
  169. evaluated_fitness_for_current_vectors_;
  170. for (long g = 0; g < total_generations_max_; ++g) {
  171. //if (process_rank_ == kOutput && g%100 == 0) {
  172. if (g%10 == 0) {
  173. double fitness = evaluated_fitness_for_current_vectors_.front().first;
  174. printf("%li -> %g\n",g,fitness);
  175. }
  176. to_be_archived_best_A_.clear();
  177. successful_mutation_parameters_S_F_.clear();
  178. successful_crossover_parameters_S_CR_.clear();
  179. // //debug section
  180. // if (process_rank_ == kOutput)
  181. // printf("============== Generation %li =============\n", g);
  182. // PrintPopulation();
  183. // PrintEvaluated();
  184. // //end of debug section
  185. for (long i = 0; i < subpopulation_; ++i) {
  186. SetCRiFi(i);
  187. std::vector<double> mutated_v, crossover_u;
  188. mutated_v = Mutation(i);
  189. crossover_u = Crossover(mutated_v, i);
  190. Selection(crossover_u, i);
  191. } // end of for all individuals in subpopulation
  192. ArchiveCleanUp();
  193. Adaption();
  194. x_vectors_current_.swap(x_vectors_next_generation_);
  195. evaluated_fitness_for_current_vectors_
  196. .swap(evaluated_fitness_for_next_generation_);
  197. SortEvaluatedCurrent();
  198. if (error_status_) return error_status_;
  199. } // end of stepping generations
  200. PrintPopulation();
  201. PrintEvaluated();
  202. return kDone;
  203. } // end of int SubPopulation::RunOptimization()
  204. // ********************************************************************** //
  205. // ********************************************************************** //
  206. // ********************************************************************** //
  207. std::vector<double> SubPopulation::Crossover(std::vector<double> mutation_v, long i) {
  208. const double CR_i = crossover_CR_[i];
  209. std::vector<double> crossover_u, x_current;
  210. crossover_u.resize(dimension_);
  211. x_current = x_vectors_current_.at(i);
  212. long j_rand = randint(0, dimension_ - 1);
  213. for (long c = 0; c < dimension_; ++c) {
  214. if (c == j_rand || rand(0,1) < CR_i)
  215. crossover_u[c] = mutation_v[c];
  216. else
  217. crossover_u[c] = x_current[c];
  218. }
  219. // //debug section
  220. // if (process_rank_ == kOutput)
  221. // printf("x -> v -> u with CR_i=%4.2f j_rand=%li\n", CR_i, j_rand);
  222. // PrintSingleVector(mutation_v);
  223. // PrintSingleVector(x_current);
  224. // PrintSingleVector(crossover_u);
  225. // //end of debug section
  226. return crossover_u;
  227. } // end of std::vector<double> SubPopulation::Crossover();
  228. // ********************************************************************** //
  229. // ********************************************************************** //
  230. // ********************************************************************** //
  231. std::vector<double> SubPopulation::Mutation(long i) {
  232. std::vector<double> mutation_v, x_best_current, x_random_current,
  233. x_random_archive_and_current, x_current;
  234. x_current = x_vectors_current_.at(i);
  235. x_best_current = GetXpBestCurrent();
  236. // //debug
  237. // if (process_rank_ == kOutput) printf("x_best: ");
  238. // PrintSingleVector(x_best_current);
  239. long index_of_random_current = -1;
  240. x_random_current = GetXRandomCurrent(&index_of_random_current, i);
  241. // //debug
  242. // if (process_rank_ == kOutput) printf("x_random: ");
  243. // PrintSingleVector(x_random_current);
  244. x_random_archive_and_current = GetXRandomArchiveAndCurrent(index_of_random_current, i);
  245. // //debug
  246. // if (process_rank_ == kOutput) printf("x_random with archive: ");
  247. // PrintSingleVector(x_random_archive_and_current);
  248. mutation_v.resize(dimension_);
  249. double F_i = mutation_F_[i];
  250. for (long c = 0; c < dimension_; ++c) {
  251. // Mutation
  252. mutation_v[c] = x_current[c]
  253. + F_i * (x_best_current[c] - x_current[c])
  254. + F_i * (x_random_current[c]
  255. - x_random_archive_and_current[c]);
  256. // Bounds control
  257. if (mutation_v[c] > x_ubound_[c])
  258. mutation_v[c] = (x_ubound_[c] + x_current[c])/2;
  259. if (mutation_v[c] < x_lbound_[c])
  260. mutation_v[c] = (x_lbound_[c] + x_current[c])/2;
  261. }
  262. // //debug section
  263. // int isSame = 888, isSame2 = 7777, isSame3 =11111;
  264. // for (long c = 0; c < dimension_; ++c) {
  265. // double norm = std::abs(mutation_v[c] - x_current[c]);
  266. // double norm2 = std::abs(x_random_current[c] - x_current[c]);
  267. // double norm3 = std::abs(x_random_current[c]
  268. // - x_random_archive_and_current[c]);
  269. // if ( norm > 0.0001) isSame = 0;
  270. // if ( norm2 > 0.0001) isSame2 = 0;
  271. // if ( norm3 > 0.0001) isSame3 = 0;
  272. // }
  273. // if (process_rank_ == kOutput) printf("mutation_v%i%i%i: ",isSame, isSame2, isSame3);
  274. // PrintSingleVector(mutation_v);
  275. // if (process_rank_ == kOutput) printf("current _v%i%i%i: ",isSame, isSame2, isSame3);
  276. // PrintSingleVector(x_current);
  277. // if (process_rank_ == kOutput)
  278. // printf(" -> f = %4.2f F_i=%4.2f\n",
  279. // FitnessFunction(mutation_v), F_i);
  280. // //end of debug section
  281. return mutation_v;
  282. } // end of std::vector<double> SubPopulation::Mutation();
  283. // ********************************************************************** //
  284. // ********************************************************************** //
  285. // ********************************************************************** //
  286. std::vector<double> SubPopulation::GetXpBestCurrent() {
  287. const long n_best_total = static_cast<long>
  288. (floor(subpopulation_ * best_share_p_ ));
  289. if (n_best_total == subpopulation_) error_status_ = kError; //TODO change kError to throw exception
  290. long best_n = randint(0, n_best_total);
  291. long best_n_index = -1, i = 0;
  292. for (auto x : evaluated_fitness_for_current_vectors_) {
  293. if (i == best_n) {
  294. best_n_index = x.second;
  295. break;
  296. }
  297. ++i;
  298. }
  299. return x_vectors_current_.at(best_n_index);
  300. } // end of std::vector<double> SubPopulation::GetXpBestCurrent();
  301. // ********************************************************************** //
  302. // ********************************************************************** //
  303. // ********************************************************************** //
  304. std::vector<double> SubPopulation::GetXRandomCurrent(long *index, long forbidden_index) {
  305. long random_n = randint(0, subpopulation_-1);
  306. while (random_n == forbidden_index) random_n = randint(0, subpopulation_-1);
  307. (*index) = random_n;
  308. return x_vectors_current_.at(random_n);
  309. } // end of std::vector<double> SubPopulation::GetXRandomCurrent()
  310. // ********************************************************************** //
  311. // ********************************************************************** //
  312. // ********************************************************************** //
  313. std::vector<double> SubPopulation::GetXRandomArchiveAndCurrent(long forbidden_index1, long forbidden_index2) {
  314. long archive_size = archived_best_A_.size();
  315. long random_n = randint(0, subpopulation_ + archive_size - 1);
  316. while (random_n == forbidden_index1 || random_n == forbidden_index2)
  317. random_n = randint(0, subpopulation_ + archive_size - 1);
  318. if (random_n < subpopulation_) return x_vectors_current_.at(random_n);
  319. random_n -= subpopulation_;
  320. long i = 0;
  321. for (auto x : archived_best_A_) {
  322. if (i == random_n) {
  323. // //debug
  324. // if (process_rank_ == kOutput) printf("Using Archive!!\n");
  325. return x;
  326. }
  327. ++i;
  328. } // end of selecting from archive
  329. error_status_ = kError; //TODO change kError to throw exception
  330. std::vector<double> x;
  331. return x;
  332. } // end of std::vector<double> SubPopulation::GetXRandomArchiveAndCurrent()
  333. // ********************************************************************** //
  334. // ********************************************************************** //
  335. // ********************************************************************** //
  336. int SubPopulation::SetCRiFi(long i) {
  337. long k = 0;
  338. while (1) {
  339. mutation_F_[i] = randc(adaptor_mutation_mu_F_, 0.1);
  340. if (mutation_F_[i] > 1) {
  341. mutation_F_[i] = 1;
  342. break;
  343. }
  344. if (mutation_F_[i] > 0) break;
  345. ++k;
  346. if (k > 10) {
  347. mutation_F_[i] = 0.001;
  348. break;
  349. }
  350. if (k > 100) printf("k");
  351. }
  352. crossover_CR_[i] = randn(adaptor_crossover_mu_CR_,0.1);
  353. if (crossover_CR_[i] > 1) crossover_CR_[i] = 1;
  354. if (crossover_CR_[i] < 0) crossover_CR_[i] = 0;
  355. return kDone;
  356. } // end of int SubPopulation::SetCRiFi(long i)
  357. // ********************************************************************** //
  358. // ********************************************************************** //
  359. // ********************************************************************** //
  360. /// %todo Change returned kError to meangfull error code.
  361. int SubPopulation::Init(long total_population, long dimension) {// NOLINT
  362. total_population_ = total_population;
  363. if (total_population_ < 1)
  364. throw std::invalid_argument("You should set population > 0!");
  365. dimension_ = dimension;
  366. if (dimension_ < 1)
  367. throw std::invalid_argument("You should set dimension > 0!");
  368. MPI_Comm_rank(MPI_COMM_WORLD, &process_rank_);
  369. MPI_Comm_size(MPI_COMM_WORLD, &number_of_processes_);
  370. if (process_rank_ < 0)
  371. throw std::invalid_argument("MPI problem: process_rank_ < 0!");
  372. if (number_of_processes_ < 1)
  373. throw std::invalid_argument("MPI problem: number_of_processes_ < 1!");
  374. std::random_device rd;
  375. generator_.seed(rd());
  376. // //debug
  377. // CheckRandom();
  378. subpopulation_ = total_population;
  379. current_generation_ = 0;
  380. x_vectors_current_.resize(subpopulation_);
  381. for (auto &x : x_vectors_current_) x.resize(dimension_);
  382. x_vectors_next_generation_.resize(subpopulation_);
  383. for (auto &x : x_vectors_next_generation_) x.resize(dimension_);
  384. evaluated_fitness_for_current_vectors_.resize(subpopulation_);
  385. // //debug
  386. // if (process_rank_ == kOutput) printf("%i, x1 size = %li \n", process_rank_, x_vectors_current_.size());
  387. x_lbound_.resize(dimension_);
  388. x_ubound_.resize(dimension_);
  389. mutation_F_.resize(subpopulation_);
  390. crossover_CR_.resize(subpopulation_);
  391. return kDone;
  392. } // end of void SubPopulation::Init()
  393. // ********************************************************************** //
  394. // ********************************************************************** //
  395. // ********************************************************************** //
  396. int SubPopulation::PrintPopulation() {
  397. if (process_rank_ == kOutput) {
  398. printf("\n");
  399. for (auto x : evaluated_fitness_for_current_vectors_) {
  400. long n = x.second;
  401. double fitness = x.first;
  402. printf("%6.2f:% 3li||", fitness, n);
  403. for (long c = 0; c < dimension_; ++c)
  404. printf(" %+7.2f ", x_vectors_current_[n][c]);
  405. printf("\n");
  406. }
  407. // for (long i = 0; i < subpopulation_; ++i) {
  408. // printf("n%li:", i);
  409. // for (long c = 0; c < dimension_; ++c) {
  410. // printf(" %5.2f ", x_vectors_current_[i][c]);
  411. // }
  412. // printf("\n");
  413. // } // end of for each individual
  414. } // end of if output
  415. return kDone;
  416. } // end of int SubPupulation::PrintPopulation()
  417. // ********************************************************************** //
  418. // ********************************************************************** //
  419. // ********************************************************************** //
  420. int SubPopulation::PrintEvaluated() {
  421. if (process_rank_ == kOutput) {
  422. for (auto x : evaluated_fitness_for_current_vectors_)
  423. printf("%li:%4.2f ", x.second, x.first);
  424. printf("\n");
  425. } // end of if output
  426. return kDone;
  427. } // end of int SubPopulation::PrintEvaluated()
  428. // ********************************************************************** //
  429. // ********************************************************************** //
  430. // ********************************************************************** //
  431. int SubPopulation::PrintSingleVector(std::vector<double> x) {
  432. if (process_rank_ == kOutput) {
  433. for (auto c : x) printf("%5.2f ", c);
  434. printf("\n");
  435. } // end of output
  436. return kDone;
  437. } // end of int SubPopulation::PrintSingleVector(std::vector<double> x)
  438. // ********************************************************************** //
  439. // ********************************************************************** //
  440. // ********************************************************************** //
  441. void SubPopulation::SetFeed(std::vector<std::vector<double> > x_feed_vectors) {
  442. isFeed_ = true;
  443. if (dimension_ == -1)
  444. throw std::invalid_argument("You should set dimension before feed!");
  445. for (auto x : x_feed_vectors)
  446. if (static_cast<long>(x.size()) != dimension_)
  447. throw std::invalid_argument
  448. ("Feed and optimization dimenstion should be the same size!");
  449. x_feed_vectors_.clear();
  450. for (auto &x : x_feed_vectors)
  451. x_feed_vectors_.push_back(x);
  452. } // end of void SubPopulation::SetFeed()
  453. // ********************************************************************** //
  454. // ********************************************************************** //
  455. // ********************************************************************** //
  456. int SubPopulation::CreateInitialPopulation() {
  457. if ((subpopulation_ - x_feed_vectors_.size()) <1)
  458. throw std::invalid_argument("Too large feed!");
  459. x_vectors_current_.resize(subpopulation_ - x_feed_vectors_.size());
  460. for (auto &x : x_vectors_current_)
  461. for (long i = 0; i < dimension_; ++i) {
  462. if (x_lbound_[i] > x_ubound_[i])
  463. throw std::invalid_argument("Wrong order of bounds!");
  464. x[i] = rand(x_lbound_[i], x_ubound_[i]); // NOLINT
  465. } // end of for each dimension
  466. // //debug
  467. // for (auto x : x_vectors_current_[0]) if (process_rank_ == kOutput) printf("%g ",x);
  468. for (auto x: x_feed_vectors_) {
  469. x_vectors_current_.push_back(x);
  470. if (process_rank_ == kOutput) {
  471. printf("--=-- Feed:\n");
  472. for (auto index:x) printf(" %+7.2f", index);
  473. printf("\n");
  474. }
  475. }
  476. if (static_cast<long>(x_vectors_current_.size()) != subpopulation_)
  477. throw std::invalid_argument("Population is not full after feed!");
  478. x_feed_vectors_.clear();
  479. // if (process_rank_ == kOutput) {
  480. // printf("==--== Initial population:\n");
  481. // EvaluateCurrentVectors();
  482. // PrintPopulation();
  483. // }
  484. return kDone;
  485. } // end of int SubPopulation::CreateInitialPopulation()
  486. // ********************************************************************** //
  487. // ********************************************************************** //
  488. // ********************************************************************** //
  489. int SubPopulation::SortEvaluatedCurrent() {
  490. evaluated_fitness_for_current_vectors_
  491. .sort([=](const std::pair<double, long>& a, // NOLINT
  492. const std::pair<double, long>& b) { // NOLINT
  493. bool cmp = a.first < b.first;
  494. if (is_find_minimum_) return cmp;
  495. return !cmp;
  496. });
  497. return kDone;
  498. } // end of int SubPopulation::SortEvaluatedCurrent()
  499. // ********************************************************************** //
  500. // ********************************************************************** //
  501. // ********************************************************************** //
  502. int SubPopulation::EvaluateCurrentVectors() {
  503. evaluated_fitness_for_current_vectors_.clear();
  504. for (long i = 0; i < subpopulation_; ++i) { // NOLINT
  505. auto tmp = std::make_pair(FitnessFunction(x_vectors_current_[i]), i);
  506. evaluated_fitness_for_current_vectors_.push_back(tmp);
  507. }
  508. SortEvaluatedCurrent();
  509. // //debug
  510. // if (process_rank_ == kOutput) printf("\n After ");
  511. // for (auto val : evaluated_fitness_for_current_vectors_)
  512. // if (process_rank_ == kOutput) printf("%g ", val.first);
  513. return kDone;
  514. } // end of int SubPopulation::EvaluateCurrentVectors()
  515. // ********************************************************************** //
  516. // ********************************************************************** //
  517. // ********************************************************************** //
  518. void SubPopulation::SetAllBoundsVectors
  519. (std::vector<double> lbound, std::vector<double> ubound) {
  520. x_lbound_.clear();
  521. x_ubound_.clear();
  522. for (auto x : lbound) x_lbound_.push_back(x);
  523. for (auto x : ubound) x_ubound_.push_back(x);
  524. }
  525. // ********************************************************************** //
  526. // ********************************************************************** //
  527. // ********************************************************************** //
  528. int SubPopulation::SetAllBounds(double lbound, double ubound) {
  529. if (lbound >= ubound) {
  530. error_status_ = kError;
  531. return kError; //TODO change kError to throw exception
  532. }
  533. for (auto &x : x_lbound_) x = lbound;
  534. for (auto &x : x_ubound_) x = ubound;
  535. return kDone;
  536. } // end of int SubPopulation::SetAllBounds(double lbound, double ubound)
  537. // ********************************************************************** //
  538. // ********************************************************************** //
  539. // ********************************************************************** //
  540. double SubPopulation::randn(double mean, double stddev) {
  541. std::normal_distribution<double> distribution(mean, stddev);
  542. return distribution(generator_);
  543. } // end of double SubPopulation::randn(double mean, double stddev)
  544. // ********************************************************************** //
  545. // ********************************************************************** //
  546. // ********************************************************************** //
  547. double SubPopulation::randc(double location, double scale) {
  548. std::cauchy_distribution<double> distribution(location, scale);
  549. return distribution(generator_);
  550. } // end of double SubPopulation::randc(double location, double scale)
  551. // ********************************************************************** //
  552. // ********************************************************************** //
  553. // ********************************************************************** //
  554. long SubPopulation::randint(long lbound, long ubound) { // NOLINT
  555. std::uniform_int_distribution<long> distribution(lbound, ubound); // NOLINT
  556. return distribution(generator_);
  557. }
  558. // ********************************************************************** //
  559. // ********************************************************************** //
  560. // ********************************************************************** //
  561. double SubPopulation::rand(double lbound, double ubound) { // NOLINT
  562. std::uniform_real_distribution<double> distribution(lbound, ubound);
  563. return distribution(generator_);
  564. } // end of double rand(double lbound, double ubound)
  565. // ********************************************************************** //
  566. // ********************************************************************** //
  567. // ********************************************************************** //
  568. void SubPopulation::CheckRandom() {
  569. std::map<int, int> hist_n, hist_c, hist_int, hist;
  570. int factor = 1000;
  571. for (int n = 0; n < 100*factor; ++n) {
  572. ++hist_n[std::round(randn(0, 2))];
  573. ++hist_c[std::round(randc(0, 2))];
  574. ++hist_int[std::round(randint(0, 10))];
  575. ++hist[std::round(rand(0, 15))]; // NOLINT
  576. }
  577. if (process_rank_ == kOutput) {
  578. printf("Normal (0,2)\n");
  579. for (auto p : hist_n) {
  580. if (p.second > factor)
  581. printf("%i: % 4i %s\n", process_rank_, p.first,
  582. std::string(p.second/factor, '*').c_str());
  583. } // end of for p in hist_n
  584. printf("Cauchy (0,2)\n");
  585. for (auto p : hist_c) {
  586. if (p.second > factor)
  587. printf("%i: % 4i %s\n", process_rank_, p.first,
  588. std::string(p.second/factor, '*').c_str());
  589. } // end of for p in hist_c
  590. printf("Uniform int (0,10)\n");
  591. for (auto p : hist_int) {
  592. if (p.second > factor)
  593. printf("%i: % 4i %s\n", process_rank_, p.first,
  594. std::string(p.second/factor, '*').c_str());
  595. } // end of for p in hist_int
  596. printf("Uniform double (0,15)\n");
  597. for (auto p : hist) {
  598. if (p.second > factor)
  599. printf("%i: % 4i %s\n", process_rank_, p.first,
  600. std::string(p.second/factor, '*').c_str());
  601. } // end of for p in hist
  602. } // end of if current process_rank_ == kOutput
  603. }
  604. // ********************************************************************** //
  605. // ********************************************************************** //
  606. // ********************************************************************** //
  607. int SubPopulation::SetBestShareP(double p) {
  608. if (p < 0 || p > 1) {
  609. error_status_ = kError;
  610. return kError; //TODO change kError to throw exception
  611. }
  612. best_share_p_ = p;
  613. return kDone;
  614. } // end of int SubPopulation::SetBestShareP(double p)
  615. // ********************************************************************** //
  616. // ********************************************************************** //
  617. // ********************************************************************** //
  618. int SubPopulation::SetAdapitonFrequencyC(double c) {
  619. if (c < 0 || c > 1) {
  620. error_status_ = kError;
  621. return kError; //TODO change kError to throw exception
  622. }
  623. adaptation_frequency_c_ = c;
  624. return kDone;
  625. } // end of int SubPopulation::SetAdapitonFrequency(double c)
  626. // ********************************************************************** //
  627. // ********************************************************************** //
  628. // ********************************************************************** //
  629. int SubPopulation::SetDistributionLevel(int level) {
  630. distribution_level_ = level;
  631. return kDone;
  632. }
  633. // end of int SubPopulation::SetDistributionLevel(int level)
  634. // ********************************************************************** //
  635. // ********************************************************************** //
  636. // ********************************************************************** //
  637. int SubPopulation::PrintParameters(std::string comment) {
  638. if (process_rank_ == 0) {
  639. printf("#%s dim=%li NP=%li(of %li) p=%4.2f c=%4.2f generation=%li\n",
  640. comment.c_str(),dimension_, subpopulation_, total_population_,
  641. best_share_p_, adaptation_frequency_c_, total_generations_max_
  642. );
  643. fflush(stdout);
  644. } // end of output
  645. return kDone;
  646. } // end of int SubPopulation::PrintParameters(std::string comment)
  647. // ********************************************************************** //
  648. // ********************************************************************** //
  649. // ********************************************************************** //
  650. int SubPopulation::PrintResult(std::string comment) {
  651. if (distribution_level_ == 0) {
  652. auto x = evaluated_fitness_for_current_vectors_.front();
  653. std::vector<double> to_send {x.first};
  654. //printf("%8.5g\n ", x.first);
  655. AllGatherVectorDouble(to_send);
  656. if (process_rank_ == 0) {
  657. double sum = 0;
  658. double size = static_cast<double>(recieve_double_.size());
  659. for (auto x : recieve_double_) sum += x;
  660. double mean = sum/size;
  661. double sigma = 0;
  662. for (auto x : recieve_double_) sigma += pow2(x - mean);
  663. sigma = sqrt(sigma/size);
  664. printf("%s gen%li, mean %4.1e (stddev %4.1e = %3.2g %%) runs(%g)\n",
  665. comment.c_str(), total_generations_max_, mean,sigma,
  666. sigma*100.0/mean,size);
  667. // for (auto x : recieve_double_)
  668. // printf("%18.15g\n", x);
  669. }
  670. }
  671. return kDone;
  672. } // end of int SubPopulation::PrintResult()
  673. // ********************************************************************** //
  674. // ********************************************************************** //
  675. // ********************************************************************** //
  676. std::vector<double> SubPopulation::GetFinalFitness() {
  677. recieve_double_.clear();
  678. if (distribution_level_ == 0) {
  679. auto x = evaluated_fitness_for_current_vectors_.front();
  680. std::vector<double> to_send {x.first};
  681. AllGatherVectorDouble(to_send);
  682. }
  683. return recieve_double_;
  684. } // end of sdt::vector<double> SubPopulation::GetFinalFitness();
  685. // ********************************************************************** //
  686. // ********************************************************************** //
  687. // ********************************************************************** //
  688. std::vector<double> SubPopulation::GetBest(double *best_fitness) {
  689. auto x = evaluated_fitness_for_current_vectors_.front();
  690. (*best_fitness) = x.first;
  691. return x_vectors_current_[x.second];
  692. } // end of std::vector<double> SubPopulation::GetBest(double *best_fitness)
  693. // ********************************************************************** //
  694. // ********************************************************************** //
  695. // ********************************************************************** //
  696. std::vector<double> SubPopulation::GetWorst(double *worst_fitness) {
  697. auto x = evaluated_fitness_for_current_vectors_.back();
  698. (*worst_fitness) = x.first;
  699. return x_vectors_current_[x.second];
  700. } // end of std::vector<double> SubPopulation::GetWorst(double *worst_fitness)
  701. // ********************************************************************** //
  702. // ********************************************************************** //
  703. // ********************************************************************** //
  704. int SubPopulation::AllGatherVectorDouble(std::vector<double> to_send) {
  705. long size_single = to_send.size();
  706. long size_all = size_single * number_of_processes_;
  707. recieve_double_.clear();
  708. recieve_double_.resize(size_all);
  709. MPI_Allgather(&to_send.front(), size_single, MPI_DOUBLE,
  710. &recieve_double_.front(), size_single, MPI_DOUBLE,
  711. MPI_COMM_WORLD);
  712. return kDone;
  713. } // end of int SubPopulation::AllGatherVectorDouble(std::vector<double> to_send);
  714. // ********************************************************************** //
  715. // ********************************************************************** //
  716. // ********************************************************************** //
  717. int SubPopulation::AllGatherVectorLong(std::vector<long> to_send) {
  718. long size_single = to_send.size();
  719. long size_all = size_single * number_of_processes_;
  720. recieve_long_.clear();
  721. recieve_long_.resize(size_all);
  722. MPI_Allgather(&to_send.front(), size_single, MPI_LONG,
  723. &recieve_long_.front(), size_single, MPI_LONG,
  724. MPI_COMM_WORLD);
  725. return kDone;
  726. } // end of int SubPopulation::AllGatherVectorLong(std::vector<long> to_send);
  727. // ********************************************************************** //
  728. // ********************************************************************** //
  729. // ********************************************************************** //
  730. } // end of namespace jade