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