155 int main(
int argc,
char* argv[])
163 #define EIGEN_USE_LAPACKE
164 #include "Vector/vector_dist.hpp"
165 #include "DMatrix/EMatrix.hpp"
166 #include <Eigen/Eigenvalues>
167 #include <Eigen/Jacobi>
169 #include "Vector/vector_dist.hpp"
170 #include <f15_cec_fun.hpp>
171 #include <boost/math/special_functions/sign.hpp>
178 constexpr
int dim = 10;
180 constexpr
int lambda = 7;
181 constexpr
int mu = lambda/2;
182 double psoWeight = 0.7;
185 double stopTolX = 2e-11;
186 double stopTolUpX = 2000.0;
188 size_t max_fun_eval = 30000000;
189 constexpr
int hist_size = 21;
198 double stop_fitness = 1.0;
209 constexpr
int sigma = 0;
210 constexpr
int Cov_m = 1;
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;
230 Eigen::DiagonalMatrix<double,Eigen::Dynamic>,
243 EVectorXd> > particle_type;
263 double generateGaussianNoise(
double mu,
double sigma)
265 static const double epsilon = std::numeric_limits<double>::min();
266 static const double two_pi = 2.0*3.14159265358979323846;
268 thread_local
double z1;
269 thread_local
double generate;
270 generate = !generate;
273 {
return z1 * sigma + mu;}
278 u1 = rand() * (1.0 / RAND_MAX);
279 u2 = rand() * (1.0 / RAND_MAX);
281 while ( u1 <= epsilon );
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;
289 template<
unsigned int dim>
290 EVectorXd generateGaussianVector()
295 for (
size_t i = 0 ; i < dim ; i++)
297 tmp(i) = generateGaussianNoise(0,1);
305 template<
unsigned int dim>
306 void fill_vector(
double (& f)[dim], EVectorXd & ev)
308 for (
size_t i = 0 ; i < dim ; i++)
312 void fill_vector(
const double * f, EVectorXd & ev)
314 for (
size_t i = 0 ; i < ev.size() ; i++)
323 bool operator<(
const fun_index & tmp)
const
333 for (
size_t i = 0 ; i < mu ; i++)
334 {wm[i] = log(
double(mu)+1.0) - log(
double(i)+1.0);}
338 for (
size_t i = 0 ; i < mu ; i++)
344 for (
size_t i = 0 ; i < mu ; i++)
356 double weight_sample(
int i)
375 void create_rotmat(EVectorXd & S,EVectorXd & T, EMatrixXd & R)
377 EVectorXd S_work(dim);
378 EVectorXd T_work(dim);
379 EVectorXd S_sup(dim);
380 EVectorXd T_sup(dim);
382 EMatrixXd R_tar(dim,dim);
383 EMatrixXd R_tmp(dim,dim);
384 EMatrixXd R_sup(dim,dim);
386 EMatrixXd S_tmp(2,2);
387 EMatrixXd T_tmp(2,2);
397 for (p = dim - 2; p >= 0 ; p -= 1)
400 for (q = dim - 1 ; q >= p+1 ; q-= 1)
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);
409 Eigen::JacobiRotation<double> G;
411 G.makeGivens(S_tmp(0), S_tmp(1),&z);
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();
428 S_work = R_tmp*S_work;
434 G.makeGivens(T_tmp(0), T_tmp(1),&z);
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();
448 T_work = R_tmp*T_work;
453 R = R_tar.transpose()*R;
457 EVectorXd Check(dim);
484 Eigen::DiagonalMatrix<double,Eigen::Dynamic> &
D,
487 EVectorXd best_sol_ei(dim);
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);
498 EMatrixXd R(dim,dim);
500 if (gb_vec_length > 0.0)
502 if(sigma < gb_vec_length)
504 if(sigma/gb_vec_length <= t_c*gb_vec_length)
507 {bias = sigma*gb_vec/gb_vec_length;}
513 xmean = xmean + bias;
517 EMatrixXd B_rot(dim,dim);
518 Eigen::DiagonalMatrix<double,Eigen::Dynamic> D_square(dim);
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);}
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();
529 EMatrixXd trUp = C_pso.triangularView<Eigen::Upper>();
530 EMatrixXd trDw = C_pso.triangularView<Eigen::StrictlyUpper>();
531 C_pso = trUp + trDw.transpose();
555 void broadcast_best_solution(particle_type & vd,
561 best_sol.resize(dim);
562 auto & v_cl = create_vcluster();
564 double best_old = best_sample;
565 v_cl.min(best_sample);
569 if (best < best_sample)
575 if (best_old == best_sample)
577 rank = v_cl.getProcessUnitID();
583 if (rank == v_cl.getProcessUnitID())
585 for (
size_t i = 0 ; i < dim ; i++)
586 {best_sol.get(i) = best_sample_sol.get(i);}
591 rank = std::numeric_limits<size_t>::max();
600 v_cl.Bcast(best_sol,rank);
630 double availablepercentiles[lambda];
634 for (
size_t i = 0 ; i < lambda ; i++)
636 availablepercentiles[i] = 0.0;
637 sar[i] = f_obj.get(i).f;
639 std::sort(&sar[0],&sar[lambda]);
641 for (
size_t i = 0 ; i < 2 ; i++)
643 if (perc[i] <= (100.0*0.5/lambda))
645 else if (perc[i] >= (100.0*(lambda-0.5)/lambda) )
646 {res[i] = sar[lambda-1];}
649 for (
size_t j = 0 ; j < lambda ; j++)
650 {availablepercentiles[j] = 100.0 * ((double(j)+1.0)-0.5) / lambda;}
652 for (k = 0 ; k < lambda ; k++)
653 {
if(availablepercentiles[k] >= perc[i]) {
break;}}
656 res[i] = sar[k] + (sar[k+1]-sar[k]) * (perc[i]
657 -availablepercentiles[k]) / (availablepercentiles[k+1] - availablepercentiles[k]);
664 double maxval(
double (& buf)[hist_size],
bool (& mask)[hist_size])
667 for (
size_t i = 0 ; i < hist_size ; i++)
669 if (buf[i] > max && mask[i] ==
true)
676 double minval(
double (& buf)[hist_size],
bool (& mask)[hist_size])
678 double min = std::numeric_limits<double>::max();
679 for (
size_t i = 0 ; i < hist_size ; i++)
681 if (buf[i] < min && mask[i] ==
true)
690 void cmaes_intobounds(EVectorXd & x, EVectorXd & xout,
bool (& idx)[dim],
bool & idx_any)
693 for (
size_t i = 0; i < dim ; i++)
722 EVectorXd (& arxvalid)[lambda],
723 EVectorXd (& arx)[lambda],
727 double (& weight)[dim],
728 double (& fithist)[hist_size],
730 double & validfitval,
740 bool mask[hist_size];
743 int dfitidx[hist_size];
744 double dfitsort[hist_size];
745 double prct[2] = {25.0,75.0};
748 for (
size_t i = 0 ; i < hist_size ; i++)
753 if (fithist[i] > 0.0)
759 for (
size_t i = 0 ; i < dim ; i++)
764 meandiag = C.trace()/dim;
766 cmaes_myprctile(f_obj, prct, val);
767 value = (val[1] - val[0]) / dim / meandiag / (sigma*sigma);
769 if (value >= std::numeric_limits<double>::max())
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);
775 else if(value == 0.0)
777 value = minval(fithist,mask);
779 else if (validfit == 0.0)
781 for (
size_t i = 0 ; i < hist_size ; i++)
786 for (
size_t i = 0; i < hist_size ; i++)
794 else if(i == hist_size-1)
796 for (
size_t k = 0 ; k < hist_size-1 ; k++)
797 {fithist[k] = fithist[k+1];}
803 cmaes_intobounds(xmean,tx,idx,idx_any);
810 {value = fithist[0];}
814 for (
size_t i = 0 ; i <= maxI; i++)
816 fitsort.get(i).f = fithist[i];
817 fitsort.get(i).id = i;
821 for (
size_t k = 0; k <= maxI ; k++)
822 {fitsort.get(k).f = fithist[fitsort.get(k).id];}
824 if ((maxI+1) % 2 == 0)
825 {value = (fitsort.get(maxI/2).f+fitsort.get(maxI/2+1).f)/2.0;}
827 {value = fitsort.get(maxI/2).f;}
829 for (
size_t i = 0 ; i < dim ; i++)
831 diag[i] = diag[i]/meandiag;
832 weight[i] = 2.0002 * value / diag[i];
834 if (validfitval == 1.0 && step-last_restart > 2)
844 for(
size_t i = 0 ; i < dim ; i++)
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)))) );
849 for (
size_t i = 0 ; i < dim ; i++)
853 weight[i] = pow(1.2,(std::max(1.0,mu_eff/10.0/dim)))*weight[i];
857 double arpenalty[lambda];
858 for (
size_t i = 0 ; i < lambda ; i++)
861 for (
size_t j = 0 ; j < dim ; j++)
863 arpenalty[i] += weight[j] * (arxvalid[i](j) - arx[i](j))*(arxvalid[i](j) - arx[i](j));
865 f_obj.get(i).f += arpenalty[i];
872 double adjust_sigma(
double sigma, EMatrixXd & C)
874 for (
size_t i = 0 ; i < dim ; i++)
876 if (sigma*sqrt(C(i,i)) > 5.0)
877 {sigma = 5.0/sqrt(C(i,i));}
902 void cma_step(particle_type & vd,
int step,
double & best,
907 EVectorXd xmean(dim);
908 EVectorXd mean_z(dim);
909 EVectorXd arxvalid[lambda];
910 EVectorXd arx[lambda];
912 for (
size_t i = 0 ; i < lambda ; i++)
915 arxvalid[i].resize(dim);
918 double best_sample = std::numeric_limits<double>::max();
923 int counteval = step*lambda;
925 auto it = vd.getDomainIterator();
930 if (vd.getProp<stop>(p) ==
true)
933 EVectorXd (& arz)[lambda] = vd.getProp<Zeta>(p);
937 fill_vector(vd.getPos(p),xmean);
939 for (
size_t j = 0 ; j < lambda ; j++)
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];
945 for (
size_t i = 0 ; i < dim ; i++)
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;}
952 {arxvalid[j](i) = arx[j](i);}
955 f_obj.get(j).f = hybrid_composition<dim>(arxvalid[j]);
960 if (f_obj.get(j).f < best_sample)
962 best_sample = f_obj.get(j).f;
965 for (
size_t i = 0 ; i < dim ; i++)
966 {best_sample_sol.get(i) = arxvalid[j](i);}
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));
981 for (
size_t j = 0 ; j < lambda ; j++)
982 {vd.getProp<ord>(p)[j] = f_obj.get(j).id;}
984 vd.getProp<xold>(p) = xmean;
990 for (
size_t j = 0 ; j < mu ; j++)
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]];
996 vd.getProp<xmean_st>(p) = xmean;
997 vd.getProp<meanz_st>(p) = mean_z;
1003 broadcast_best_solution(vd,best_sol,best,best_sample,best_sample_sol);
1006 bool calc_bd = counteval - eigeneval > lambda/(ccov)/dim/10;
1007 if (calc_bd ==
true)
1008 {eigeneval = counteval;}
1010 auto it2 = vd.getDomainIterator();
1011 while (it2.isNext())
1015 if (vd.getProp<stop>(p) ==
true)
1018 xmean = vd.getProp<xmean_st>(p);
1019 mean_z = vd.getProp<meanz_st>(p);
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;
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);
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);
1027 if (step % N_pso == 0)
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);
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());
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();
1041 vd.getProp<Cov_m>(p) = psoWeight*vd.getProp<Cov_m>(p) + (1.0 - psoWeight)*C_pso;
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());
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();
1057 double smaller = std::numeric_limits<double>::max();
1058 for (
size_t i = 0 ; i < dim ; i++)
1060 if (vd.getProp<sigma>(p)*sqrt(vd.getProp<
D>(p).diagonal()[i]) > 5.0)
1062 if (smaller > 5.0/sqrt(vd.getProp<
D>(p).diagonal()[i]))
1063 {smaller = 5.0/sqrt(vd.getProp<
D>(p).diagonal()[i]);}
1066 if (smaller != std::numeric_limits<double>::max())
1067 {vd.getProp<sigma>(p) = smaller;}
1070 vd.getProp<sigma>(p) = vd.getProp<sigma>(p)*exp((cs/d_amps)*(vd.getProp<path_s>(p).norm()/chiN - 1));
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();
1081 Eigen::SelfAdjointEigenSolver<EMatrixXd> eig_solver;
1083 eig_solver.compute(vd.getProp<Cov_m>(p));
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();
1090 for (
size_t i = 0 ; i < dim ; i++)
1092 if (vd.getProp<B>(p)(0,i) < 0)
1093 {vd.getProp<B>(p).col(i) = - vd.getProp<B>(p).col(i);}
1096 EMatrixXd tmp = vd.getProp<B>(p).transpose();
1100 for (
size_t i = 0 ; i < dim ; i++)
1101 {vd.getPos(p)[i] = xmean(i);}
1103 vd.getProp<sigma>(p) = adjust_sigma(vd.getProp<sigma>(p),vd.getProp<Cov_m>(p));
1106 bool stop_tol =
true;
1107 bool stop_tolX =
true;
1108 for (
size_t i = 0 ; i < dim ; i++)
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;
1114 vd.getProp<stop>(p) = stop_tol | stop_tolX;
1117 if (f_obj.get(0).f == f_obj.get(std::ceil(0.7*lambda)).f )
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";
1123 vd.getProp<stop>(p) =
true;
1126 if (vd.getProp<stop>(p) ==
true)
1127 {std::cout <<
"Stopped" << std::endl;}
1129 if (restart_cma && vd.getProp<stop>(p) ==
true)
1131 std::cout <<
"------- Restart #" << std::endl;
1133 std::cout <<
"---------------------------------" << std::endl;
1134 std::cout <<
"Best: " << best <<
" " << fun_eval << std::endl;
1135 std::cout <<
"---------------------------------" << std::endl;
1137 vd.getProp<last_restart>(p) = step;
1138 vd.getProp<xold>(p).setZero();
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;
1155 for (
size_t i = 0 ; i < dim ; i++)
1158 vd.getPos(p)[i] = 10.0*(double)rand() / RAND_MAX - 5.0;
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;
1172 auto & v_cl = create_vcluster();
1194 int main(
int argc,
char* argv[])
1197 openfpm_init(&argc,&argv);
1199 auto & v_cl = create_vcluster();
1204 for (
size_t i = 0 ; i < dim ; i++)
1206 domain.setLow(i,0.0);
1207 domain.setHigh(i,1.0);
1212 for (
size_t i = 0 ; i < dim ; i++)
1213 {bc[i] = NON_PERIODIC;};
1220 particle_type vd(16,domain,bc,g);
1224 stop_fitness = 1e-10;
1225 size_t stopeval = 1e3*dim*dim;
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;
1237 chiN = sqrt(dim)*(1.0-1.0/(4.0*dim)+1.0/(21.0*dim*dim));
1243 int seed = 24756*v_cl.rank()*v_cl.rank() + time(NULL);
1246 auto it = vd.getDomainIterator();
1252 for (
size_t i = 0 ; i < dim ; i++)
1255 vd.getPos(p)[i] = 10.0*(double)rand() / RAND_MAX - 5.0;
1258 vd.getProp<sigma>(p) = 2.0;
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;
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;
1287 if (v_cl.rank() == 0)
1288 {std::cout <<
"Starting PS-CMA-ES" << std::endl;}
1293 best = std::numeric_limits<double>::max();
1298 size_t fun_eval = 0;
1300 while (fun_eval < max_fun_eval && best > 120.000001)
1303 cma_step(vd,i+1,best,best_i,best_sol,fun_eval);
1308 if (v_cl.rank() == 0)
1310 std::cout <<
"Best solution: " << best <<
" with " << fun_eval << std::endl;
1311 std::cout <<
"at: " << std::endl;
1313 for (
size_t i = 0 ; i < best_sol.
size() ; i++)
1315 std::cout << best_sol.get(i) <<
" ";
This class represent an N-dimensional box.
Derivative second order on h (spacing)
Implementation of 1-D std::vector like structure.
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.