152 int main(
int argc,
char* argv[])
160 #define EIGEN_USE_LAPACKE 161 #include "Vector/vector_dist.hpp" 162 #include "DMatrix/EMatrix.hpp" 163 #include <Eigen/Eigenvalues> 164 #include <Eigen/Jacobi> 166 #include "Vector/vector_dist.hpp" 167 #include <f15_cec_fun.hpp> 168 #include <boost/math/special_functions/sign.hpp> 175 constexpr
int dim = 10;
177 constexpr
int lambda = 7;
178 constexpr
int mu = lambda/2;
179 double psoWeight = 0.7;
182 double stopTolX = 2e-11;
183 double stopTolUpX = 2000.0;
185 size_t max_fun_eval = 30000000;
186 constexpr
int hist_size = 21;
195 double stop_fitness = 1.0;
206 constexpr
int sigma = 0;
207 constexpr
int Cov_m = 1;
210 constexpr
int Zeta = 4;
211 constexpr
int path_s = 5;
212 constexpr
int path_c = 6;
213 constexpr
int ord = 7;
214 constexpr
int stop = 8;
215 constexpr
int fithist = 9;
216 constexpr
int weight = 10;
217 constexpr
int validfit = 11;
218 constexpr
int xold = 12;
219 constexpr
int last_restart = 13;
220 constexpr
int iniphase = 14;
221 constexpr
int xmean_st = 15;
222 constexpr
int meanz_st = 16;
227 Eigen::DiagonalMatrix<double,Eigen::Dynamic>,
240 EVectorXd> > particle_type;
260 double generateGaussianNoise(
double mu,
double sigma)
262 static const double epsilon = std::numeric_limits<double>::min();
263 static const double two_pi = 2.0*3.14159265358979323846;
265 thread_local
double z1;
266 thread_local
double generate;
267 generate = !generate;
270 {
return z1 * sigma + mu;}
275 u1 = rand() * (1.0 / RAND_MAX);
276 u2 = rand() * (1.0 / RAND_MAX);
278 while ( u1 <= epsilon );
281 z0 = sqrt(-2.0 * log(u2)) * cos(two_pi * u1);
282 z1 = sqrt(-2.0 * log(u2)) * sin(two_pi * u1);
283 return z0 * sigma + mu;
286 template<
unsigned int dim>
287 EVectorXd generateGaussianVector()
292 for (
size_t i = 0 ; i < dim ; i++)
294 tmp(i) = generateGaussianNoise(0,1);
302 template<
unsigned int dim>
303 void fill_vector(
double (& f)[dim], EVectorXd & ev)
305 for (
size_t i = 0 ; i < dim ; i++)
309 void fill_vector(
const double * f, EVectorXd & ev)
311 for (
size_t i = 0 ; i < ev.size() ; i++)
320 bool operator<(
const fun_index & tmp)
const 330 for (
size_t i = 0 ; i < mu ; i++)
331 {wm[i] = log(
double(mu)+1.0) - log(
double(i)+1.0);}
335 for (
size_t i = 0 ; i < mu ; i++)
341 for (
size_t i = 0 ; i < mu ; i++)
353 double weight_sample(
int i)
372 void create_rotmat(EVectorXd & S,EVectorXd & T, EMatrixXd & R)
374 EVectorXd S_work(dim);
375 EVectorXd T_work(dim);
376 EVectorXd S_sup(dim);
377 EVectorXd T_sup(dim);
379 EMatrixXd R_tar(dim,dim);
380 EMatrixXd R_tmp(dim,dim);
381 EMatrixXd R_sup(dim,dim);
383 EMatrixXd S_tmp(2,2);
384 EMatrixXd T_tmp(2,2);
394 for (p = dim - 2; p >= 0 ; p -= 1)
397 for (q = dim - 1 ; q >= p+1 ; q-= 1)
399 T_tmp(0) = T_work(p);
400 T_tmp(1) = T_work(q);
401 S_tmp(0) = S_work(p);
402 S_tmp(1) = S_work(q);
406 Eigen::JacobiRotation<double> G;
408 G.makeGivens(S_tmp(0), S_tmp(1),&z);
417 R_tmp(p,p) = sign*G.c();
418 R_tmp(q,q) = sign*G.c();
419 R_tmp(p,q) = sign*-G.s();
420 R_tmp(q,p) = sign*G.s();
425 S_work = R_tmp*S_work;
431 G.makeGivens(T_tmp(0), T_tmp(1),&z);
438 R_tmp(p,p) = sign*G.c();
439 R_tmp(q,q) = sign*G.c();
440 R_tmp(p,q) = sign*-G.s();
441 R_tmp(q,p) = sign*G.s();
445 T_work = R_tmp*T_work;
450 R = R_tar.transpose()*R;
454 EVectorXd Check(dim);
481 Eigen::DiagonalMatrix<double,Eigen::Dynamic> &
D,
484 EVectorXd best_sol_ei(dim);
486 double bias_weight = psoWeight;
487 fill_vector(&best_sol.get(0),best_sol_ei);
488 EVectorXd gb_vec = best_sol_ei-xmean;
489 double gb_vec_length = sqrt(gb_vec.transpose() * gb_vec);
490 EVectorXd b_main = B.col(dim-1);
495 EMatrixXd R(dim,dim);
497 if (gb_vec_length > 0.0)
499 if(sigma < gb_vec_length)
501 if(sigma/gb_vec_length <= t_c*gb_vec_length)
504 {bias = sigma*gb_vec/gb_vec_length;}
510 xmean = xmean + bias;
514 EMatrixXd B_rot(dim,dim);
515 Eigen::DiagonalMatrix<double,Eigen::Dynamic> D_square(dim);
517 EVectorXd gb_vec_old = best_sol_ei - xold;
518 create_rotmat(b_main,gb_vec_old,R);
519 for (
size_t i = 0 ; i < dim ; i++)
520 {B_rot.col(i) = R*B.col(i);}
522 for (
size_t i = 0 ; i < dim ; i++)
523 {D_square.diagonal()[i] =
D.diagonal()[i] *
D.diagonal()[i];}
524 C_pso = B_rot * D_square * B_rot.transpose();
526 EMatrixXd trUp = C_pso.triangularView<Eigen::Upper>();
527 EMatrixXd trDw = C_pso.triangularView<Eigen::StrictlyUpper>();
528 C_pso = trUp + trDw.transpose();
552 void broadcast_best_solution(particle_type & vd,
558 best_sol.resize(dim);
559 auto & v_cl = create_vcluster();
561 double best_old = best_sample;
562 v_cl.min(best_sample);
566 if (best < best_sample)
572 if (best_old == best_sample)
574 rank = v_cl.getProcessUnitID();
580 if (rank == v_cl.getProcessUnitID())
582 for (
size_t i = 0 ; i < dim ; i++)
583 {best_sol.get(i) = best_sample_sol.get(i);}
588 rank = std::numeric_limits<size_t>::max();
597 v_cl.Bcast(best_sol,rank);
627 double availablepercentiles[lambda];
631 for (
size_t i = 0 ; i < lambda ; i++)
633 availablepercentiles[i] = 0.0;
634 sar[i] = f_obj.get(i).f;
636 std::sort(&sar[0],&sar[lambda]);
638 for (
size_t i = 0 ; i < 2 ; i++)
640 if (perc[i] <= (100.0*0.5/lambda))
642 else if (perc[i] >= (100.0*(lambda-0.5)/lambda) )
643 {res[i] = sar[lambda-1];}
646 for (
size_t j = 0 ; j < lambda ; j++)
647 {availablepercentiles[j] = 100.0 * ((double(j)+1.0)-0.5) / lambda;}
649 for (k = 0 ; k < lambda ; k++)
650 {
if(availablepercentiles[k] >= perc[i]) {
break;}}
653 res[i] = sar[k] + (sar[k+1]-sar[k]) * (perc[i]
654 -availablepercentiles[k]) / (availablepercentiles[k+1] - availablepercentiles[k]);
661 double maxval(
double (& buf)[hist_size],
bool (& mask)[hist_size])
664 for (
size_t i = 0 ; i < hist_size ; i++)
666 if (buf[i] > max && mask[i] ==
true)
673 double minval(
double (& buf)[hist_size],
bool (& mask)[hist_size])
675 double min = std::numeric_limits<double>::max();
676 for (
size_t i = 0 ; i < hist_size ; i++)
678 if (buf[i] < min && mask[i] ==
true)
687 void cmaes_intobounds(EVectorXd & x, EVectorXd & xout,
bool (& idx)[dim],
bool & idx_any)
690 for (
size_t i = 0; i < dim ; i++)
719 EVectorXd (& arxvalid)[lambda],
720 EVectorXd (& arx)[lambda],
724 double (& weight)[dim],
725 double (& fithist)[hist_size],
727 double & validfitval,
737 bool mask[hist_size];
740 int dfitidx[hist_size];
741 double dfitsort[hist_size];
742 double prct[2] = {25.0,75.0};
745 for (
size_t i = 0 ; i < hist_size ; i++)
750 if (fithist[i] > 0.0)
756 for (
size_t i = 0 ; i < dim ; i++)
761 meandiag = C.trace()/dim;
763 cmaes_myprctile(f_obj, prct, val);
764 value = (val[1] - val[0]) / dim / meandiag / (sigma*sigma);
766 if (value >= std::numeric_limits<double>::max())
768 auto & v_cl = create_vcluster();
769 std::cout <<
"Process " << v_cl.rank() <<
" warning: Non-finite fitness range" << std::endl;
770 value = maxval(fithist,mask);
772 else if(value == 0.0)
774 value = minval(fithist,mask);
776 else if (validfit == 0.0)
778 for (
size_t i = 0 ; i < hist_size ; i++)
783 for (
size_t i = 0; i < hist_size ; i++)
791 else if(i == hist_size-1)
793 for (
size_t k = 0 ; k < hist_size-1 ; k++)
794 {fithist[k] = fithist[k+1];}
800 cmaes_intobounds(xmean,tx,idx,idx_any);
807 {value = fithist[0];}
811 for (
size_t i = 0 ; i <= maxI; i++)
813 fitsort.get(i).f = fithist[i];
814 fitsort.get(i).id = i;
818 for (
size_t k = 0; k <= maxI ; k++)
819 {fitsort.get(k).f = fithist[fitsort.get(k).id];}
821 if ((maxI+1) % 2 == 0)
822 {value = (fitsort.get(maxI/2).f+fitsort.get(maxI/2+1).f)/2.0;}
824 {value = fitsort.get(maxI/2).f;}
826 for (
size_t i = 0 ; i < dim ; i++)
828 diag[i] = diag[i]/meandiag;
829 weight[i] = 2.0002 * value / diag[i];
831 if (validfitval == 1.0 && step-last_restart > 2)
841 for(
size_t i = 0 ; i < dim ; i++)
843 idx[i] = (idx[i] && (fabs(tx(i)) > 3.0*std::max(1.0,sqrt(dim)/mu_eff) * sigma * sqrt(diag[i])));
844 idx[i] = (idx[i] && (std::copysign(1.0,tx(i)) == std::copysign(1.0,(xmean(i)-xold(i)))) );
846 for (
size_t i = 0 ; i < dim ; i++)
850 weight[i] = pow(1.2,(std::max(1.0,mu_eff/10.0/dim)))*weight[i];
854 double arpenalty[lambda];
855 for (
size_t i = 0 ; i < lambda ; i++)
858 for (
size_t j = 0 ; j < dim ; j++)
860 arpenalty[i] += weight[j] * (arxvalid[i](j) - arx[i](j))*(arxvalid[i](j) - arx[i](j));
862 f_obj.get(i).f += arpenalty[i];
869 double adjust_sigma(
double sigma, EMatrixXd & C)
871 for (
size_t i = 0 ; i < dim ; i++)
873 if (sigma*sqrt(C(i,i)) > 5.0)
874 {sigma = 5.0/sqrt(C(i,i));}
899 void cma_step(particle_type & vd,
int step,
double & best,
904 EVectorXd xmean(dim);
905 EVectorXd mean_z(dim);
906 EVectorXd arxvalid[lambda];
907 EVectorXd arx[lambda];
909 for (
size_t i = 0 ; i < lambda ; i++)
912 arxvalid[i].resize(dim);
915 double best_sample = std::numeric_limits<double>::max();
920 int counteval = step*lambda;
922 auto it = vd.getDomainIterator();
927 if (vd.getProp<stop>(p) ==
true)
930 EVectorXd (& arz)[lambda] = vd.getProp<Zeta>(p);
934 fill_vector(vd.getPos(p),xmean);
936 for (
size_t j = 0 ; j < lambda ; j++)
938 vd.getProp<Zeta>(p)[j] = generateGaussianVector<dim>();
939 arx[j] = xmean + vd.getProp<sigma>(p)*vd.getProp<B>(p)*vd.getProp<
D>(p)*vd.getProp<Zeta>(p)[j];
942 for (
size_t i = 0 ; i < dim ; i++)
944 if (arx[j](i) < -5.0)
945 {arxvalid[j](i) = -5.0;}
946 else if (arx[j](i) > 5.0)
947 {arxvalid[j](i) = 5.0;}
949 {arxvalid[j](i) = arx[j](i);}
952 f_obj.get(j).f = hybrid_composition<dim>(arxvalid[j]);
957 if (f_obj.get(j).f < best_sample)
959 best_sample = f_obj.get(j).f;
962 for (
size_t i = 0 ; i < dim ; i++)
963 {best_sample_sol.get(i) = arxvalid[j](i);}
968 cmaes_handlebounds(f_obj,vd.getProp<sigma>(p),
969 vd.getProp<validfit>(p),arxvalid,
970 arx,vd.getProp<Cov_m>(p),
971 xmean,vd.getProp<xold>(p),vd.getProp<weight>(p),
972 vd.getProp<fithist>(p),vd.getProp<iniphase>(p),
973 vd.getProp<validfit>(p),mu_eff,
974 step,vd.getProp<last_restart>(p));
978 for (
size_t j = 0 ; j < lambda ; j++)
979 {vd.getProp<ord>(p)[j] = f_obj.get(j).id;}
981 vd.getProp<xold>(p) = xmean;
987 for (
size_t j = 0 ; j < mu ; j++)
989 xmean += weight_sample(j)*arx[vd.getProp<ord>(p)[j]];
990 mean_z += weight_sample(j)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[j]];
993 vd.getProp<xmean_st>(p) = xmean;
994 vd.getProp<meanz_st>(p) = mean_z;
1000 broadcast_best_solution(vd,best_sol,best,best_sample,best_sample_sol);
1003 bool calc_bd = counteval - eigeneval > lambda/(ccov)/dim/10;
1004 if (calc_bd ==
true)
1005 {eigeneval = counteval;}
1007 auto it2 = vd.getDomainIterator();
1008 while (it2.isNext())
1012 if (vd.getProp<stop>(p) ==
true)
1015 xmean = vd.getProp<xmean_st>(p);
1016 mean_z = vd.getProp<meanz_st>(p);
1018 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 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 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 if (step % N_pso == 0)
1026 EMatrixXd C_pso(dim,dim);
1027 updatePso(best_sol,vd.getProp<sigma>(p),xmean,vd.getProp<xold>(p),vd.getProp<B>(p),vd.getProp<
D>(p),C_pso);
1030 vd.getProp<Cov_m>(p) = (1.0-ccov+(1.0-hsig)*ccov*cc*(2.0-cc)/mu_eff)*vd.getProp<Cov_m>(p) +
1031 ccov*(1.0/mu_eff)*(vd.getProp<path_c>(p)*vd.getProp<path_c>(p).transpose());
1033 for (
size_t i = 0 ; i < mu ; i++)
1034 {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)*
1035 (vd.getProp<B>(p)*vd.getProp<
D>(p)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[i]]).transpose();
1038 vd.getProp<Cov_m>(p) = psoWeight*vd.getProp<Cov_m>(p) + (1.0 - psoWeight)*C_pso;
1043 vd.getProp<Cov_m>(p) = (1.0-ccov+(1.0-hsig)*ccov*cc*(2.0-cc)/mu_eff)*vd.getProp<Cov_m>(p) +
1044 ccov*(1.0/mu_eff)*(vd.getProp<path_c>(p)*vd.getProp<path_c>(p).transpose());
1046 for (
size_t i = 0 ; i < mu ; i++)
1047 {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)*
1048 (vd.getProp<B>(p)*vd.getProp<
D>(p)*vd.getProp<Zeta>(p)[vd.getProp<ord>(p)[i]]).transpose();
1054 double smaller = std::numeric_limits<double>::max();
1055 for (
size_t i = 0 ; i < dim ; i++)
1057 if (vd.getProp<sigma>(p)*sqrt(vd.getProp<
D>(p).diagonal()[i]) > 5.0)
1059 if (smaller > 5.0/sqrt(vd.getProp<
D>(p).diagonal()[i]))
1060 {smaller = 5.0/sqrt(vd.getProp<
D>(p).diagonal()[i]);}
1063 if (smaller != std::numeric_limits<double>::max())
1064 {vd.getProp<sigma>(p) = smaller;}
1067 vd.getProp<sigma>(p) = vd.getProp<sigma>(p)*exp((cs/d_amps)*(vd.getProp<path_s>(p).norm()/chiN - 1));
1073 EMatrixXd trUp = vd.getProp<Cov_m>(p).triangularView<Eigen::Upper>();
1074 EMatrixXd trDw = vd.getProp<Cov_m>(p).triangularView<Eigen::StrictlyUpper>();
1075 vd.getProp<Cov_m>(p) = trUp + trDw.transpose();
1078 Eigen::SelfAdjointEigenSolver<EMatrixXd> eig_solver;
1080 eig_solver.compute(vd.getProp<Cov_m>(p));
1082 for (
size_t i = 0 ; i < eig_solver.eigenvalues().size() ; i++)
1083 {vd.getProp<
D>(p).diagonal()[i] = sqrt(eig_solver.eigenvalues()[i]);}
1084 vd.getProp<B>(p) = eig_solver.eigenvectors();
1087 for (
size_t i = 0 ; i < dim ; i++)
1089 if (vd.getProp<B>(p)(0,i) < 0)
1090 {vd.getProp<B>(p).col(i) = - vd.getProp<B>(p).col(i);}
1093 EMatrixXd tmp = vd.getProp<B>(p).transpose();
1097 for (
size_t i = 0 ; i < dim ; i++)
1098 {vd.getPos(p)[i] = xmean(i);}
1100 vd.getProp<sigma>(p) = adjust_sigma(vd.getProp<sigma>(p),vd.getProp<Cov_m>(p));
1103 bool stop_tol =
true;
1104 bool stop_tolX =
true;
1105 for (
size_t i = 0 ; i < dim ; i++)
1107 stop_tol &= (vd.getProp<sigma>(p)*std::max(fabs(vd.getProp<path_c>(p)(i)),sqrt(vd.getProp<Cov_m>(p)(i,i)))) < stopTolX;
1108 stop_tolX &= vd.getProp<sigma>(p)*sqrt(vd.getProp<
D>(p).diagonal()[i]) > stopTolUpX;
1111 vd.getProp<stop>(p) = stop_tol | stop_tolX;
1114 if (f_obj.get(0).f == f_obj.get(std::ceil(0.7*lambda)).f )
1116 vd.getProp<sigma>(p) = vd.getProp<sigma>(p)*exp(0.2+cs/d_amps);
1117 std::cout <<
"warning: flat fitness, consider reformulating the objective";
1120 vd.getProp<stop>(p) =
true;
1123 if (vd.getProp<stop>(p) ==
true)
1124 {std::cout <<
"Stopped" << std::endl;}
1126 if (restart_cma && vd.getProp<stop>(p) ==
true)
1128 std::cout <<
"------- Restart #" << std::endl;
1130 std::cout <<
"---------------------------------" << std::endl;
1131 std::cout <<
"Best: " << best <<
" " << fun_eval << std::endl;
1132 std::cout <<
"---------------------------------" << std::endl;
1134 vd.getProp<last_restart>(p) = step;
1135 vd.getProp<xold>(p).setZero();
1137 for (
size_t i = 0 ; i < vd.getProp<
D>(p).diagonal().size() ; i++)
1138 {vd.getProp<
D>(p).diagonal()[i] = 1.0;}
1139 vd.getProp<B>(p).resize(dim,dim);
1140 vd.getProp<B>(p).setIdentity();
1141 vd.getProp<Cov_m>(p) = vd.getProp<B>(p)*vd.getProp<
D>(p)*vd.getProp<
D>(p)*vd.getProp<B>(p);
1142 vd.getProp<path_s>(p).resize(dim);
1143 vd.getProp<path_s>(p).setZero(dim);
1144 vd.getProp<path_c>(p).resize(dim);
1145 vd.getProp<path_c>(p).setZero(dim);
1146 vd.getProp<stop>(p) =
false;
1147 vd.getProp<iniphase>(p) =
true;
1148 vd.getProp<last_restart>(p) = 0;
1149 vd.getProp<sigma>(p) = 2.0;
1152 for (
size_t i = 0 ; i < dim ; i++)
1155 vd.getPos(p)[i] = 10.0*(double)rand() / RAND_MAX - 5.0;
1160 for (
size_t i = 0 ; i < hist_size ; i++)
1161 {vd.getProp<fithist>(p)[i] = -1.0;}
1162 vd.getProp<fithist>(p)[0] = 1.0;
1163 vd.getProp<validfit>(p) = 0.0;
1169 auto & v_cl = create_vcluster();
1191 int main(
int argc,
char* argv[])
1194 openfpm_init(&argc,&argv);
1196 auto & v_cl = create_vcluster();
1201 for (
size_t i = 0 ; i < dim ; i++)
1209 for (
size_t i = 0 ; i < dim ; i++)
1210 {bc[i] = NON_PERIODIC;};
1217 particle_type vd(16,domain,bc,g);
1221 stop_fitness = 1e-10;
1222 size_t stopeval = 1e3*dim*dim;
1228 cc = 4.0 / (dim+4.0);
1229 cs = (mu_eff+2.0) / (double(dim)+mu_eff+3.0);
1230 ccov = (1.0/mu_eff) * 2.0/((dim+1.41)*(dim+1.41)) +
1231 (1.0 - 1.0/mu_eff)* std::min(1.0,(2.0*mu_eff-1.0)/((dim+2.0)*(dim+2.0) + mu_eff));
1232 d_amps = 1 + 2*std::max(0.0, sqrt((mu_eff-1.0)/(dim+1))-1) + cs;
1234 chiN = sqrt(dim)*(1.0-1.0/(4.0*dim)+1.0/(21.0*dim*dim));
1240 int seed = 24756*v_cl.rank()*v_cl.rank() + time(NULL);
1243 auto it = vd.getDomainIterator();
1249 for (
size_t i = 0 ; i < dim ; i++)
1252 vd.getPos(p)[i] = 10.0*(double)rand() / RAND_MAX - 5.0;
1255 vd.getProp<sigma>(p) = 2.0;
1259 vd.getProp<
D>(p).resize(dim);
1260 for (
size_t i = 0 ; i < vd.getProp<
D>(p).diagonal().size() ; i++)
1261 {vd.getProp<
D>(p).diagonal()[i] = 1.0;}
1262 vd.getProp<B>(p).resize(dim,dim);
1263 vd.getProp<B>(p).setIdentity();
1264 vd.getProp<Cov_m>(p) = vd.getProp<B>(p)*vd.getProp<
D>(p)*vd.getProp<
D>(p)*vd.getProp<B>(p);
1265 vd.getProp<path_s>(p).resize(dim);
1266 vd.getProp<path_s>(p).setZero(dim);
1267 vd.getProp<path_c>(p).resize(dim);
1268 vd.getProp<path_c>(p).setZero(dim);
1269 vd.getProp<stop>(p) =
false;
1270 vd.getProp<iniphase>(p) =
true;
1271 vd.getProp<last_restart>(p) = 0;
1275 for (
size_t i = 0 ; i < hist_size ; i++)
1276 {vd.getProp<fithist>(p)[i] = -1.0;}
1277 vd.getProp<fithist>(p)[0] = 1.0;
1278 vd.getProp<validfit>(p) = 0.0;
1284 if (v_cl.rank() == 0)
1285 {std::cout <<
"Starting PS-CMA-ES" << std::endl;}
1290 best = std::numeric_limits<double>::max();
1295 size_t fun_eval = 0;
1297 while (fun_eval < max_fun_eval && best > 120.000001)
1300 cma_step(vd,i+1,best,best_i,best_sol,fun_eval);
1305 if (v_cl.rank() == 0)
1307 std::cout <<
"Best solution: " << best <<
" with " << fun_eval << std::endl;
1308 std::cout <<
"at: " << std::endl;
1310 for (
size_t i = 0 ; i < best_sol.
size() ; i++)
1312 std::cout << best_sol.get(i) <<
" ";
Derivative second order on h (spacing)
__device__ __host__ void setHigh(int i, T val)
set the high interval of the box
__device__ __host__ void setLow(int i, T val)
set the low interval of the box
KeyT const ValueT ValueT OffsetIteratorT OffsetIteratorT int
[in] The number of segments that comprise the sorting data
This class represent an N-dimensional box.
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...
It model an expression expr1 + ... exprn.
Implementation of 1-D std::vector like structure.