OpenFPM  5.2.0
Project that contain the implementation of distributed structures
main.cpp
1 
148 
149 // This example need EIGEN, if we do not have NOP the full example
150 
151 #include "config.h"
152 
153 #ifndef HAVE_EIGEN
154 
155 int main(int argc, char* argv[])
156 {
157  return 0;
158 }
159 
160 #else
161 
162 
163 #define EIGEN_USE_LAPACKE
164 #include "Vector/vector_dist.hpp"
165 #include "DMatrix/EMatrix.hpp"
166 #include <Eigen/Eigenvalues>
167 #include <Eigen/Jacobi>
168 #include <limits>
169 #include "Vector/vector_dist.hpp"
170 #include <f15_cec_fun.hpp>
171 #include <boost/math/special_functions/sign.hpp>
172 
174 
176 
177 // PARAMETERS
178 constexpr int dim = 10;
179 // when you set dim set also lambda to to 4+std::floor(3*log(dim))
180 constexpr int lambda = 7;
181 constexpr int mu = lambda/2;
182 double psoWeight = 0.7;
183 // number of cma-step before pso step
184 int N_pso = 200;
185 double stopTolX = 2e-11;
186 double stopTolUpX = 2000.0;
187 int restart_cma = 1;
188 size_t max_fun_eval = 30000000;
189 constexpr int hist_size = 21;
190 
191 // Convenient global variables (Their value is set after)
192 double mu_eff = 1.0;
193 double cs = 1.0;
194 double cc = 1.0;
195 double ccov = 1.0;
196 double chiN;
197 double d_amps = 1.0;
198 double stop_fitness = 1.0;
199 int eigeneval = 0;
200 double t_c = 0.1;
201 double b = 0.1;
202 
204 
206 
208 
209 constexpr int sigma = 0;
210 constexpr int Cov_m = 1;
211 constexpr int B = 2;
212 constexpr int D = 3;
213 constexpr int Zeta = 4;
214 constexpr int path_s = 5;
215 constexpr int path_c = 6;
216 constexpr int ord = 7;
217 constexpr int stop = 8;
218 constexpr int fithist = 9;
219 constexpr int weight = 10;
220 constexpr int validfit = 11;
221 constexpr int xold = 12;
222 constexpr int last_restart = 13;
223 constexpr int iniphase = 14;
224 constexpr int xmean_st = 15;
225 constexpr int meanz_st = 16;
226 
227 typedef vector_dist<dim,double, aggregate<double,
228  EMatrixXd,
229  EMatrixXd,
230  Eigen::DiagonalMatrix<double,Eigen::Dynamic>,
231  EVectorXd[lambda],
232  EVectorXd,
233  EVectorXd,
234  int[lambda],
235  int,
236  double [hist_size],
237  double [dim],
238  double,
239  EVectorXd,
240  int,
241  bool,
242  EVectorXd,
243  EVectorXd> > particle_type;
244 
246 
262 
263 double generateGaussianNoise(double mu, double sigma)
264 {
265  static const double epsilon = std::numeric_limits<double>::min();
266  static const double two_pi = 2.0*3.14159265358979323846;
267 
268  thread_local double z1;
269  thread_local double generate;
270  generate = !generate;
271 
272  if (!generate)
273  {return z1 * sigma + mu;}
274 
275  double u1, u2;
276  do
277  {
278  u1 = rand() * (1.0 / RAND_MAX);
279  u2 = rand() * (1.0 / RAND_MAX);
280  }
281  while ( u1 <= epsilon );
282 
283  double z0;
284  z0 = sqrt(-2.0 * log(u2)) * cos(two_pi * u1);
285  z1 = sqrt(-2.0 * log(u2)) * sin(two_pi * u1);
286  return z0 * sigma + mu;
287 }
288 
289 template<unsigned int dim>
290 EVectorXd generateGaussianVector()
291 {
292  EVectorXd tmp;
293  tmp.resize(dim);
294 
295  for (size_t i = 0 ; i < dim ; i++)
296  {
297  tmp(i) = generateGaussianNoise(0,1);
298  }
299 
300  return tmp;
301 }
302 
304 
305 template<unsigned int dim>
306 void fill_vector(double (& f)[dim], EVectorXd & ev)
307 {
308  for (size_t i = 0 ; i < dim ; i++)
309  {ev(i) = f[i];}
310 }
311 
312 void fill_vector(const double * f, EVectorXd & ev)
313 {
314  for (size_t i = 0 ; i < ev.size() ; i++)
315  {ev(i) = f[i];}
316 }
317 
318 struct fun_index
319 {
320  double f;
321  int id;
322 
323  bool operator<(const fun_index & tmp) const
324  {
325  return f < tmp.f;
326  }
327 };
328 
329 double wm[mu];
330 
331 void init_weight()
332 {
333  for (size_t i = 0 ; i < mu ; i++)
334  {wm[i] = log(double(mu)+1.0) - log(double(i)+1.0);}
335 
336  double tot = 0.0;
337 
338  for (size_t i = 0 ; i < mu ; i++)
339  {tot += wm[i];}
340 
341  double sum = 0.0;
342  double sum2 = 0.0;
343 
344  for (size_t i = 0 ; i < mu ; i++)
345  {
346  wm[i] /= tot;
347  sum += wm[i];
348  sum2 += wm[i]*wm[i];
349  }
350 
351  // also set mu_eff
352  mu_eff=sum*sum/sum2;
353 
354 }
355 
356 double weight_sample(int i)
357 {
358  return wm[i];
359 }
360 
374 
375 void create_rotmat(EVectorXd & S,EVectorXd & T, EMatrixXd & R)
376 {
377  EVectorXd S_work(dim);
378  EVectorXd T_work(dim);
379  EVectorXd S_sup(dim);
380  EVectorXd T_sup(dim);
381 
382  EMatrixXd R_tar(dim,dim);
383  EMatrixXd R_tmp(dim,dim);
384  EMatrixXd R_sup(dim,dim);
385  double G_S,G_C;
386  EMatrixXd S_tmp(2,2);
387  EMatrixXd T_tmp(2,2);
388  int p,q,i;
389 
390  S_work = S;
391  T_work = T;
392 
393  R.setIdentity();
394  R_tar = R;
395  R_tmp = R;
396 
397  for (p = dim - 2; p >= 0 ; p -= 1)
398  {
399 
400  for (q = dim - 1 ; q >= p+1 ; q-= 1)
401  {
402  T_tmp(0) = T_work(p);
403  T_tmp(1) = T_work(q);
404  S_tmp(0) = S_work(p);
405  S_tmp(1) = S_work(q);
406 
407  // Perform Givens Rotation on start vector
408 
409  Eigen::JacobiRotation<double> G;
410  double z;
411  G.makeGivens(S_tmp(0), S_tmp(1),&z);
412 
413  // Check direction of rotation
414  double sign = 1.0;
415  if (z < 0.0)
416  {sign = -1.0;}
417 
418  // Build a Rotation Matrix out of G_C and G_S
419  R_tmp.setIdentity();
420  R_tmp(p,p) = sign*G.c();
421  R_tmp(q,q) = sign*G.c();
422  R_tmp(p,q) = sign*-G.s();
423  R_tmp(q,p) = sign*G.s();
424 
425  // Rotate start vector and update R
426  // S_work = R_tmp*S_work
427 
428  S_work = R_tmp*S_work;
429  // R = R_tmp*R
430  R = R_tmp*R;
431 
432  // Perform Givens Rotation on target vector
433 
434  G.makeGivens(T_tmp(0), T_tmp(1),&z);
435 
436  sign = 1.0;
437  if (z < 0.0)
438  {sign = -1.0;}
439 
440  R_tmp.setIdentity();
441  R_tmp(p,p) = sign*G.c();
442  R_tmp(q,q) = sign*G.c();
443  R_tmp(p,q) = sign*-G.s();
444  R_tmp(q,p) = sign*G.s();
445 
446  // Rotate target vector and update R_tar
447 
448  T_work = R_tmp*T_work;
449  R_tar = R_tmp*R_tar;
450  }
451  }
452 
453  R = R_tar.transpose()*R;
454 
455  // Check the rotation
456 
457  EVectorXd Check(dim);
458  Check = R*S;
459 }
460 
462 
478 
479 void updatePso(openfpm::vector<double> & best_sol,
480  double sigma,
481  EVectorXd & xmean,
482  EVectorXd & xold,
483  EMatrixXd & B,
484  Eigen::DiagonalMatrix<double,Eigen::Dynamic> & D,
485  EMatrixXd & C_pso)
486 {
487  EVectorXd best_sol_ei(dim);
488 
489  double bias_weight = psoWeight;
490  fill_vector(&best_sol.get(0),best_sol_ei);
491  EVectorXd gb_vec = best_sol_ei-xmean;
492  double gb_vec_length = sqrt(gb_vec.transpose() * gb_vec);
493  EVectorXd b_main = B.col(dim-1);
494  EVectorXd bias(dim);
495  bias.setZero();
496 
497  // Rotation Matrix
498  EMatrixXd R(dim,dim);
499 
500  if (gb_vec_length > 0.0)
501  {
502  if(sigma < gb_vec_length)
503  {
504  if(sigma/gb_vec_length <= t_c*gb_vec_length)
505  {bias = 0.5*gb_vec;}
506  else
507  {bias = sigma*gb_vec/gb_vec_length;}
508  }
509  else
510  {bias.setZero();}
511  }
512 
513  xmean = xmean + bias;
514 
515  if (psoWeight < 1.0)
516  {
517  EMatrixXd B_rot(dim,dim);
518  Eigen::DiagonalMatrix<double,Eigen::Dynamic> D_square(dim);
519 
520  EVectorXd gb_vec_old = best_sol_ei - xold;
521  create_rotmat(b_main,gb_vec_old,R);
522  for (size_t i = 0 ; i < dim ; i++)
523  {B_rot.col(i) = R*B.col(i);}
524 
525  for (size_t i = 0 ; i < dim ; i++)
526  {D_square.diagonal()[i] = D.diagonal()[i] * D.diagonal()[i];}
527  C_pso = B_rot * D_square * B_rot.transpose();
528 
529  EMatrixXd trUp = C_pso.triangularView<Eigen::Upper>();
530  EMatrixXd trDw = C_pso.triangularView<Eigen::StrictlyUpper>();
531  C_pso = trUp + trDw.transpose();
532  }
533 }
534 
536 
554 
555 void broadcast_best_solution(particle_type & vd,
556  openfpm::vector<double> & best_sol,
557  double & best,
558  double best_sample,
559  openfpm::vector<double> & best_sample_sol)
560 {
561  best_sol.resize(dim);
562  auto & v_cl = create_vcluster();
563 
564  double best_old = best_sample;
565  v_cl.min(best_sample);
566  v_cl.execute();
567 
568  // The old solution remain the best
569  if (best < best_sample)
570  {return;}
571 
572  best = best_sample;
573 
574  size_t rank;
575  if (best_old == best_sample)
576  {
577  rank = v_cl.getProcessUnitID();
578 
579  // we own the minimum and we decide who broad cast
580  v_cl.min(rank);
581  v_cl.execute();
582 
583  if (rank == v_cl.getProcessUnitID())
584  {
585  for (size_t i = 0 ; i < dim ; i++)
586  {best_sol.get(i) = best_sample_sol.get(i);}
587  }
588  }
589  else
590  {
591  rank = std::numeric_limits<size_t>::max();
592 
593  // we do not own decide who broad cast
594  v_cl.min(rank);
595  v_cl.execute();
596  }
597 
598  // now we broad cast the best solution across processors
599 
600  v_cl.Bcast(best_sol,rank);
601  v_cl.execute();
602 }
603 
605 
626 
627 void cmaes_myprctile(openfpm::vector<fun_index> & f_obj, double (& perc)[2], double (& res)[2])
628 {
629  double sar[lambda];
630  double availablepercentiles[lambda];
631  int idx[hist_size];
632  int i,k;
633 
634  for (size_t i = 0 ; i < lambda ; i++)
635  {
636  availablepercentiles[i] = 0.0;
637  sar[i] = f_obj.get(i).f;
638  }
639  std::sort(&sar[0],&sar[lambda]);
640 
641  for (size_t i = 0 ; i < 2 ; i++)
642  {
643  if (perc[i] <= (100.0*0.5/lambda))
644  {res[i] = sar[0];}
645  else if (perc[i] >= (100.0*(lambda-0.5)/lambda) )
646  {res[i] = sar[lambda-1];}
647  else
648  {
649  for (size_t j = 0 ; j < lambda ; j++)
650  {availablepercentiles[j] = 100.0 * ((double(j)+1.0)-0.5) / lambda;}
651 
652  for (k = 0 ; k < lambda ; k++)
653  {if(availablepercentiles[k] >= perc[i]) {break;}}
654  k-=1;
655 
656  res[i] = sar[k] + (sar[k+1]-sar[k]) * (perc[i]
657  -availablepercentiles[k]) / (availablepercentiles[k+1] - availablepercentiles[k]);
658  }
659  }
660 }
661 
663 
664 double maxval(double (& buf)[hist_size], bool (& mask)[hist_size])
665 {
666  double max = 0.0;
667  for (size_t i = 0 ; i < hist_size ; i++)
668  {
669  if (buf[i] > max && mask[i] == true)
670  {max = buf[i];}
671  }
672 
673  return max;
674 }
675 
676 double minval(double (& buf)[hist_size], bool (& mask)[hist_size])
677 {
678  double min = std::numeric_limits<double>::max();
679  for (size_t i = 0 ; i < hist_size ; i++)
680  {
681  if (buf[i] < min && mask[i] == true)
682  {min = buf[i];}
683  }
684 
685  return min;
686 }
687 
689 
690 void cmaes_intobounds(EVectorXd & x, EVectorXd & xout,bool (& idx)[dim], bool & idx_any)
691 {
692  idx_any = false;
693  for (size_t i = 0; i < dim ; i++)
694  {
695  if(x(i) < -5.0)
696  {
697  xout(i) = -5.0;
698  idx[i] = true;
699  idx_any = true;
700  }
701  else if (x(i) > 5.0)
702  {
703  xout(i) = 5.0;
704  idx[i] = true;
705  idx_any = true;
706  }
707  else
708  {
709  xout(i) = x(i);
710  idx[i] = false;
711  }
712  }
713 }
714 
716 
718 
719 void cmaes_handlebounds(openfpm::vector<fun_index> & f_obj,
720  double sigma,
721  double & validfit,
722  EVectorXd (& arxvalid)[lambda],
723  EVectorXd (& arx)[lambda],
724  EMatrixXd & C,
725  EVectorXd & xmean,
726  EVectorXd & xold,
727  double (& weight)[dim],
728  double (& fithist)[hist_size],
729  bool & iniphase,
730  double & validfitval,
731  double mu_eff,
732  int step,
733  int last_restart)
734 {
735  double val[2];
736  double value;
737  double diag[dim];
738  double meandiag;
739  int i,k,maxI;
740  bool mask[hist_size];
741  bool idx[dim];
742  EVectorXd tx(dim);
743  int dfitidx[hist_size];
744  double dfitsort[hist_size];
745  double prct[2] = {25.0,75.0};
746  bool idx_any;
747 
748  for (size_t i = 0 ; i < hist_size ; i++)
749  {
750  dfitsort[i] = 0.0;
751  dfitidx[i] = 0;
752 
753  if (fithist[i] > 0.0)
754  {mask[i] = true;}
755  else
756  {mask[i] = false;}
757  }
758 
759  for (size_t i = 0 ; i < dim ; i++)
760  {diag[i] = C(i,i);}
761 
762  maxI = 0;
763 
764  meandiag = C.trace()/dim;
765 
766  cmaes_myprctile(f_obj, prct, val);
767  value = (val[1] - val[0]) / dim / meandiag / (sigma*sigma);
768 
769  if (value >= std::numeric_limits<double>::max())
770  {
771  auto & v_cl = create_vcluster();
772  std::cout << "Process " << v_cl.rank() << " warning: Non-finite fitness range" << std::endl;
773  value = maxval(fithist,mask);
774  }
775  else if(value == 0.0)
776  {
777  value = minval(fithist,mask);
778  }
779  else if (validfit == 0.0)
780  {
781  for (size_t i = 0 ; i < hist_size ; i++)
782  {fithist[i] = -1.0;}
783  validfit = 1;
784  }
785 
786  for (size_t i = 0; i < hist_size ; i++)
787  {
788  if(fithist[i] < 0.0)
789  {
790  fithist[i] = value;
791  maxI = i;
792  break;
793  }
794  else if(i == hist_size-1)
795  {
796  for (size_t k = 0 ; k < hist_size-1 ; k++)
797  {fithist[k] = fithist[k+1];}
798  fithist[i] = value;
799  maxI = i;
800  }
801  }
802 
803  cmaes_intobounds(xmean,tx,idx,idx_any);
804 
805  if (iniphase)
806  {
807  if (idx_any)
808  {
809  if(maxI == 0)
810  {value = fithist[0];}
811  else
812  {
813  openfpm::vector<fun_index> fitsort(maxI+1);
814  for (size_t i = 0 ; i <= maxI; i++)
815  {
816  fitsort.get(i).f = fithist[i];
817  fitsort.get(i).id = i;
818  }
819 
820  fitsort.sort();
821  for (size_t k = 0; k <= maxI ; k++)
822  {fitsort.get(k).f = fithist[fitsort.get(k).id];}
823 
824  if ((maxI+1) % 2 == 0)
825  {value = (fitsort.get(maxI/2).f+fitsort.get(maxI/2+1).f)/2.0;}
826  else
827  {value = fitsort.get(maxI/2).f;}
828  }
829  for (size_t i = 0 ; i < dim ; i++)
830  {
831  diag[i] = diag[i]/meandiag;
832  weight[i] = 2.0002 * value / diag[i];
833  }
834  if (validfitval == 1.0 && step-last_restart > 2)
835  {
836  iniphase = false;
837  }
838  }
839  }
840 
841  if(idx_any)
842  {
843  tx = xmean - tx;
844  for(size_t i = 0 ; i < dim ; i++)
845  {
846  idx[i] = (idx[i] && (fabs(tx(i)) > 3.0*std::max(1.0,sqrt(dim)/mu_eff) * sigma * sqrt(diag[i])));
847  idx[i] = (idx[i] && (std::copysign(1.0,tx(i)) == std::copysign(1.0,(xmean(i)-xold(i)))) );
848  }
849  for (size_t i = 0 ; i < dim ; i++)
850  {
851  if (idx[i] == true)
852  {
853  weight[i] = pow(1.2,(std::max(1.0,mu_eff/10.0/dim)))*weight[i];
854  }
855  }
856  }
857  double arpenalty[lambda];
858  for (size_t i = 0 ; i < lambda ; i++)
859  {
860  arpenalty[i] = 0.0;
861  for (size_t j = 0 ; j < dim ; j++)
862  {
863  arpenalty[i] += weight[j] * (arxvalid[i](j) - arx[i](j))*(arxvalid[i](j) - arx[i](j));
864  }
865  f_obj.get(i).f += arpenalty[i];
866  }
867 // fitness%sel = fitness%raw + bnd%arpenalty;
868 }
869 
871 
872 double adjust_sigma(double sigma, EMatrixXd & C)
873 {
874  for (size_t i = 0 ; i < dim ; i++)
875  {
876  if (sigma*sqrt(C(i,i)) > 5.0)
877  {sigma = 5.0/sqrt(C(i,i));}
878  }
879 
880  return sigma;
881 }
882 
901 
902 void cma_step(particle_type & vd, int step, double & best,
903  int & best_i, openfpm::vector<double> & best_sol,
904  size_t & fun_eval)
905 {
906  size_t fe = 0;
907  EVectorXd xmean(dim);
908  EVectorXd mean_z(dim);
909  EVectorXd arxvalid[lambda];
910  EVectorXd arx[lambda];
911 
912  for (size_t i = 0 ; i < lambda ; i++)
913  {
914  arx[i].resize(dim);
915  arxvalid[i].resize(dim);
916  }
917 
918  double best_sample = std::numeric_limits<double>::max();
919  openfpm::vector<double> best_sample_sol(dim);
920 
921  openfpm::vector<fun_index> f_obj(lambda);
922 
923  int counteval = step*lambda;
924 
925  auto it = vd.getDomainIterator();
926  while (it.isNext())
927  {
928  auto p = it.get();
929 
930  if (vd.getProp<stop>(p) == true)
931  {++it;continue;}
932 
933  EVectorXd (& arz)[lambda] = vd.getProp<Zeta>(p);
934 
935  // fill the mean vector;
936 
937  fill_vector(vd.getPos(p),xmean);
938 
939  for (size_t j = 0 ; j < lambda ; j++)
940  {
941  vd.getProp<Zeta>(p)[j] = generateGaussianVector<dim>();
942  arx[j] = xmean + vd.getProp<sigma>(p)*vd.getProp<B>(p)*vd.getProp<D>(p)*vd.getProp<Zeta>(p)[j];
943 
944  // sample point has to be inside -5.0 and 5.0
945  for (size_t i = 0 ; i < dim ; i++)
946  {
947  if (arx[j](i) < -5.0)
948  {arxvalid[j](i) = -5.0;}
949  else if (arx[j](i) > 5.0)
950  {arxvalid[j](i) = 5.0;}
951  else
952  {arxvalid[j](i) = arx[j](i);}
953  }
954 
955  f_obj.get(j).f = hybrid_composition<dim>(arxvalid[j]);
956  f_obj.get(j).id = j;
957  fe++;
958 
959  // Get the best ever
960  if (f_obj.get(j).f < best_sample)
961  {
962  best_sample = f_obj.get(j).f;
963 
964  // Copy the new mean as position of the particle
965  for (size_t i = 0 ; i < dim ; i++)
966  {best_sample_sol.get(i) = arxvalid[j](i);}
967  }
968  }
969 
970  // Add penalities for out of bound points
971  cmaes_handlebounds(f_obj,vd.getProp<sigma>(p),
972  vd.getProp<validfit>(p),arxvalid,
973  arx,vd.getProp<Cov_m>(p),
974  xmean,vd.getProp<xold>(p),vd.getProp<weight>(p),
975  vd.getProp<fithist>(p),vd.getProp<iniphase>(p),
976  vd.getProp<validfit>(p),mu_eff,
977  step,vd.getProp<last_restart>(p));
978 
979  f_obj.sort();
980 
981  for (size_t j = 0 ; j < lambda ; j++)
982  {vd.getProp<ord>(p)[j] = f_obj.get(j).id;}
983 
984  vd.getProp<xold>(p) = xmean;
985 
986  // Calculate weighted mean
987 
988  xmean.setZero();
989  mean_z.setZero();
990  for (size_t j = 0 ; j < mu ; j++)
991  {
992  xmean += weight_sample(j)*arx[vd.getProp<ord>(p)[j]];
993  mean_z += weight_sample(j)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[j]];
994  }
995 
996  vd.getProp<xmean_st>(p) = xmean;
997  vd.getProp<meanz_st>(p) = mean_z;
998 
999  ++it;
1000  }
1001 
1002  // Find the best point across processors
1003  broadcast_best_solution(vd,best_sol,best,best_sample,best_sample_sol);
1004 
1005  // bool calculate B and D
1006  bool calc_bd = counteval - eigeneval > lambda/(ccov)/dim/10;
1007  if (calc_bd == true)
1008  {eigeneval = counteval;}
1009 
1010  auto it2 = vd.getDomainIterator();
1011  while (it2.isNext())
1012  {
1013  auto p = it2.get();
1014 
1015  if (vd.getProp<stop>(p) == true)
1016  {++it2;continue;}
1017 
1018  xmean = vd.getProp<xmean_st>(p);
1019  mean_z = vd.getProp<meanz_st>(p);
1020 
1021  vd.getProp<path_s>(p) = vd.getProp<path_s>(p)*(1.0 - cs) + sqrt(cs*(2.0-cs)*mu_eff)*vd.getProp<B>(p)*mean_z;
1022 
1023  double hsig = vd.getProp<path_s>(p).norm()/sqrt(1.0-pow((1.0-cs),(2.0*double((step-vd.getProp<last_restart>(p))))))/chiN < 1.4 + 2.0/(dim+1);
1024 
1025  vd.getProp<path_c>(p) = (1-cc)*vd.getProp<path_c>(p) + hsig * sqrt(cc*(2-cc)*mu_eff)*(vd.getProp<B>(p)*vd.getProp<D>(p)*mean_z);
1026 
1027  if (step % N_pso == 0)
1028  {
1029  EMatrixXd C_pso(dim,dim);
1030  updatePso(best_sol,vd.getProp<sigma>(p),xmean,vd.getProp<xold>(p),vd.getProp<B>(p),vd.getProp<D>(p),C_pso);
1031 
1032  // Adapt covariance matrix C
1033  vd.getProp<Cov_m>(p) = (1.0-ccov+(1.0-hsig)*ccov*cc*(2.0-cc)/mu_eff)*vd.getProp<Cov_m>(p) +
1034  ccov*(1.0/mu_eff)*(vd.getProp<path_c>(p)*vd.getProp<path_c>(p).transpose());
1035 
1036  for (size_t i = 0 ; i < mu ; i++)
1037  {vd.getProp<Cov_m>(p) += ccov*(1.0-1.0/mu_eff)*(vd.getProp<B>(p)*vd.getProp<D>(p)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[i]])*weight_sample(i)*
1038  (vd.getProp<B>(p)*vd.getProp<D>(p)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[i]]).transpose();
1039  }
1040 
1041  vd.getProp<Cov_m>(p) = psoWeight*vd.getProp<Cov_m>(p) + (1.0 - psoWeight)*C_pso;
1042  }
1043  else
1044  {
1045  // Adapt covariance matrix C
1046  vd.getProp<Cov_m>(p) = (1.0-ccov+(1.0-hsig)*ccov*cc*(2.0-cc)/mu_eff)*vd.getProp<Cov_m>(p) +
1047  ccov*(1.0/mu_eff)*(vd.getProp<path_c>(p)*vd.getProp<path_c>(p).transpose());
1048 
1049  for (size_t i = 0 ; i < mu ; i++)
1050  {vd.getProp<Cov_m>(p) += ccov*(1.0-1.0/mu_eff)*(vd.getProp<B>(p)*vd.getProp<D>(p)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[i]])*weight_sample(i)*
1051  (vd.getProp<B>(p)*vd.getProp<D>(p)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[i]]).transpose();
1052  }
1053  }
1054 
1055  // Numeric error
1056 
1057  double smaller = std::numeric_limits<double>::max();
1058  for (size_t i = 0 ; i < dim ; i++)
1059  {
1060  if (vd.getProp<sigma>(p)*sqrt(vd.getProp<D>(p).diagonal()[i]) > 5.0)
1061  {
1062  if (smaller > 5.0/sqrt(vd.getProp<D>(p).diagonal()[i]))
1063  {smaller = 5.0/sqrt(vd.getProp<D>(p).diagonal()[i]);}
1064  }
1065  }
1066  if (smaller != std::numeric_limits<double>::max())
1067  {vd.getProp<sigma>(p) = smaller;}
1068 
1069  //Adapt step-size sigma
1070  vd.getProp<sigma>(p) = vd.getProp<sigma>(p)*exp((cs/d_amps)*(vd.getProp<path_s>(p).norm()/chiN - 1));
1071 
1072  // Update B and D from C
1073 
1074  if (calc_bd)
1075  {
1076  EMatrixXd trUp = vd.getProp<Cov_m>(p).triangularView<Eigen::Upper>();
1077  EMatrixXd trDw = vd.getProp<Cov_m>(p).triangularView<Eigen::StrictlyUpper>();
1078  vd.getProp<Cov_m>(p) = trUp + trDw.transpose();
1079 
1080  // Eigen decomposition
1081  Eigen::SelfAdjointEigenSolver<EMatrixXd> eig_solver;
1082 
1083  eig_solver.compute(vd.getProp<Cov_m>(p));
1084 
1085  for (size_t i = 0 ; i < eig_solver.eigenvalues().size() ; i++)
1086  {vd.getProp<D>(p).diagonal()[i] = sqrt(eig_solver.eigenvalues()[i]);}
1087  vd.getProp<B>(p) = eig_solver.eigenvectors();
1088 
1089  // Make first component always positive
1090  for (size_t i = 0 ; i < dim ; i++)
1091  {
1092  if (vd.getProp<B>(p)(0,i) < 0)
1093  {vd.getProp<B>(p).col(i) = - vd.getProp<B>(p).col(i);}
1094  }
1095 
1096  EMatrixXd tmp = vd.getProp<B>(p).transpose();
1097  }
1098 
1099  // Copy the new mean as position of the particle
1100  for (size_t i = 0 ; i < dim ; i++)
1101  {vd.getPos(p)[i] = xmean(i);}
1102 
1103  vd.getProp<sigma>(p) = adjust_sigma(vd.getProp<sigma>(p),vd.getProp<Cov_m>(p));
1104 
1105  // Stop conditions
1106  bool stop_tol = true;
1107  bool stop_tolX = true;
1108  for (size_t i = 0 ; i < dim ; i++)
1109  {
1110  stop_tol &= (vd.getProp<sigma>(p)*std::max(fabs(vd.getProp<path_c>(p)(i)),sqrt(vd.getProp<Cov_m>(p)(i,i)))) < stopTolX;
1111  stop_tolX &= vd.getProp<sigma>(p)*sqrt(vd.getProp<D>(p).diagonal()[i]) > stopTolUpX;
1112  }
1113 
1114  vd.getProp<stop>(p) = stop_tol | stop_tolX;
1115 
1116  // Escape flat fitness, or better terminate?
1117  if (f_obj.get(0).f == f_obj.get(std::ceil(0.7*lambda)).f )
1118  {
1119  vd.getProp<sigma>(p) = vd.getProp<sigma>(p)*exp(0.2+cs/d_amps);
1120  std::cout << "warning: flat fitness, consider reformulating the objective";
1121 
1122  // Stop it
1123  vd.getProp<stop>(p) = true;
1124  }
1125 
1126  if (vd.getProp<stop>(p) == true)
1127  {std::cout << "Stopped" << std::endl;}
1128 
1129  if (restart_cma && vd.getProp<stop>(p) == true)
1130  {
1131  std::cout << "------- Restart #" << std::endl;
1132 
1133  std::cout << "---------------------------------" << std::endl;
1134  std::cout << "Best: " << best << " " << fun_eval << std::endl;
1135  std::cout << "---------------------------------" << std::endl;
1136 
1137  vd.getProp<last_restart>(p) = step;
1138  vd.getProp<xold>(p).setZero();
1139 
1140  for (size_t i = 0 ; i < vd.getProp<D>(p).diagonal().size() ; i++)
1141  {vd.getProp<D>(p).diagonal()[i] = 1.0;}
1142  vd.getProp<B>(p).resize(dim,dim);
1143  vd.getProp<B>(p).setIdentity();
1144  vd.getProp<Cov_m>(p) = vd.getProp<B>(p)*vd.getProp<D>(p)*vd.getProp<D>(p)*vd.getProp<B>(p);
1145  vd.getProp<path_s>(p).resize(dim);
1146  vd.getProp<path_s>(p).setZero(dim);
1147  vd.getProp<path_c>(p).resize(dim);
1148  vd.getProp<path_c>(p).setZero(dim);
1149  vd.getProp<stop>(p) = false;
1150  vd.getProp<iniphase>(p) = true;
1151  vd.getProp<last_restart>(p) = 0;
1152  vd.getProp<sigma>(p) = 2.0;
1153 
1154  // a different point in space
1155  for (size_t i = 0 ; i < dim ; i++)
1156  {
1157  // we define x, assign a random position between 0.0 and 1.0
1158  vd.getPos(p)[i] = 10.0*(double)rand() / RAND_MAX - 5.0;
1159  }
1160 
1161  // Initialize the bound history
1162 
1163  for (size_t i = 0 ; i < hist_size ; i++)
1164  {vd.getProp<fithist>(p)[i] = -1.0;}
1165  vd.getProp<fithist>(p)[0] = 1.0;
1166  vd.getProp<validfit>(p) = 0.0;
1167  }
1168 
1169  ++it2;
1170  }
1171 
1172  auto & v_cl = create_vcluster();
1173  v_cl.sum(fe);
1174  v_cl.execute();
1175 
1176  fun_eval += fe;
1177 }
1178 
1180 
1194 int main(int argc, char* argv[])
1195 {
1196  // initialize the library
1197  openfpm_init(&argc,&argv);
1198 
1199  auto & v_cl = create_vcluster();
1200 
1201  // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
1202  Box<dim,double> domain;
1203 
1204  for (size_t i = 0 ; i < dim ; i++)
1205  {
1206  domain.setLow(i,0.0);
1207  domain.setHigh(i,1.0);
1208  }
1209 
1210  // Here we define the boundary conditions of our problem
1211  size_t bc[dim];
1212  for (size_t i = 0 ; i < dim ; i++)
1213  {bc[i] = NON_PERIODIC;};
1214 
1215  prepare_f15<dim>();
1216 
1217  // extended boundary around the domain, and the processor domain
1218  Ghost<dim,double> g(0.0);
1219 
1220  particle_type vd(16,domain,bc,g);
1221 
1222  // Initialize constants
1223 
1224  stop_fitness = 1e-10;
1225  size_t stopeval = 1e3*dim*dim;
1226 
1227  // Strategy parameter setting: Selection
1228  init_weight();
1229 
1230  // Strategy parameter setting: Adaptation
1231  cc = 4.0 / (dim+4.0);
1232  cs = (mu_eff+2.0) / (double(dim)+mu_eff+3.0);
1233  ccov = (1.0/mu_eff) * 2.0/((dim+1.41)*(dim+1.41)) +
1234  (1.0 - 1.0/mu_eff)* std::min(1.0,(2.0*mu_eff-1.0)/((dim+2.0)*(dim+2.0) + mu_eff));
1235  d_amps = 1 + 2*std::max(0.0, sqrt((mu_eff-1.0)/(dim+1))-1) + cs;
1236 
1237  chiN = sqrt(dim)*(1.0-1.0/(4.0*dim)+1.0/(21.0*dim*dim));
1238 
1240 
1241 
1242  // initialize the srand
1243  int seed = 24756*v_cl.rank()*v_cl.rank() + time(NULL);
1244  srand(seed);
1245 
1246  auto it = vd.getDomainIterator();
1247 
1248  while (it.isNext())
1249  {
1250  auto p = it.get();
1251 
1252  for (size_t i = 0 ; i < dim ; i++)
1253  {
1254  // we define x, assign a random position between 0.0 and 1.0
1255  vd.getPos(p)[i] = 10.0*(double)rand() / RAND_MAX - 5.0;
1256  }
1257 
1258  vd.getProp<sigma>(p) = 2.0;
1259 
1260  // Initialize the covariant Matrix,B and D to identity
1261 
1262  vd.getProp<D>(p).resize(dim);
1263  for (size_t i = 0 ; i < vd.getProp<D>(p).diagonal().size() ; i++)
1264  {vd.getProp<D>(p).diagonal()[i] = 1.0;}
1265  vd.getProp<B>(p).resize(dim,dim);
1266  vd.getProp<B>(p).setIdentity();
1267  vd.getProp<Cov_m>(p) = vd.getProp<B>(p)*vd.getProp<D>(p)*vd.getProp<D>(p)*vd.getProp<B>(p);
1268  vd.getProp<path_s>(p).resize(dim);
1269  vd.getProp<path_s>(p).setZero(dim);
1270  vd.getProp<path_c>(p).resize(dim);
1271  vd.getProp<path_c>(p).setZero(dim);
1272  vd.getProp<stop>(p) = false;
1273  vd.getProp<iniphase>(p) = true;
1274  vd.getProp<last_restart>(p) = 0;
1275 
1276  // Initialize the bound history
1277 
1278  for (size_t i = 0 ; i < hist_size ; i++)
1279  {vd.getProp<fithist>(p)[i] = -1.0;}
1280  vd.getProp<fithist>(p)[0] = 1.0;
1281  vd.getProp<validfit>(p) = 0.0;
1282 
1283  // next particle
1284  ++it;
1285  }
1286 
1287  if (v_cl.rank() == 0)
1288  {std::cout << "Starting PS-CMA-ES" << std::endl;}
1289 
1290  double best = 0.0;
1291  int best_i = 0;
1292 
1293  best = std::numeric_limits<double>::max();
1294  openfpm::vector<double> best_sol(dim);
1295  // now do several iteration
1296 
1297  int stop_cond = 0;
1298  size_t fun_eval = 0;
1299  int i = 0;
1300  while (fun_eval < max_fun_eval && best > 120.000001)
1301  {
1302  // sample offspring
1303  cma_step(vd,i+1,best,best_i,best_sol,fun_eval);
1304 
1305  i++;
1306  }
1307 
1308  if (v_cl.rank() == 0)
1309  {
1310  std::cout << "Best solution: " << best << " with " << fun_eval << std::endl;
1311  std::cout << "at: " << std::endl;
1312 
1313  for (size_t i = 0 ; i < best_sol.size() ; i++)
1314  {
1315  std::cout << best_sol.get(i) << " ";
1316  }
1317  }
1318 
1319  openfpm_finalize();
1320 
1322 
1332 }
1333 
1334 #endif
1335 
This class represent an N-dimensional box.
Definition: Box.hpp:60
Derivative second order on h (spacing)
Definition: Derivative.hpp:29
Definition: Ghost.hpp:40
Implementation of 1-D std::vector like structure.
Definition: map_vector.hpp:204
size_t size()
Stub size.
Definition: map_vector.hpp:212
Distributed vector.
KeyT const ValueT ValueT OffsetIteratorT OffsetIteratorT int
[in] The number of segments that comprise the sorting data
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...
Definition: aggregate.hpp:221
It model an expression expr1 + ... exprn.
Definition: sum.hpp:93