00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef _edoEstimatorNormalAdaptive_h
00030 #define _edoEstimatorNormalAdaptive_h
00031
00032 #ifdef WITH_EIGEN
00033
00034 #include <algorithm>
00035
00036 #include<Eigen/Dense>
00037
00038 #include "edoNormalAdaptive.h"
00039 #include "edoEstimatorAdaptive.h"
00040
00047 template< typename EOT, typename D = edoNormalAdaptive<EOT> >
00048 class edoEstimatorNormalAdaptive : public edoEstimatorAdaptive< D >
00049 {
00050 public:
00051 typedef typename EOT::AtomType AtomType;
00052 typedef typename D::Vector Vector;
00053 typedef typename D::Matrix Matrix;
00054
00055 edoEstimatorNormalAdaptive( D& distrib ) :
00056 edoEstimatorAdaptive<D>( distrib ),
00057 _calls(0),
00058 _eigeneval(0)
00059 {}
00060
00061 private:
00062 Eigen::VectorXd edoCMAESweights( unsigned int pop_size )
00063 {
00064
00065 Eigen::VectorXd weights( pop_size );
00066 double sum_w = 0;
00067 for( unsigned int i = 0; i < pop_size; ++i ) {
00068 double w_i = log( pop_size + 0.5 ) - log( i + 1 );
00069 weights(i) = w_i;
00070 sum_w += w_i;
00071 }
00072
00073 weights /= sum_w;
00074
00075 assert( weights.size() == pop_size);
00076 return weights;
00077 }
00078
00079 public:
00080 void resetCalls()
00081 {
00082 _calls = 0;
00083 }
00084
00085
00086 edoNormalAdaptive<EOT> operator()( eoPop<EOT>& pop )
00087 {
00088
00089
00090
00091
00092
00093 unsigned int N = pop[0].size();
00094 unsigned int lambda = pop.size();
00095
00096
00097 _calls++;
00098
00099 unsigned int counteval = _calls * lambda;
00100
00101
00102
00103
00104 Matrix arx( N, lambda );
00105
00106
00107 for( unsigned int d = 0; d < N; ++d ) {
00108 for( unsigned int i = 0; i < lambda; ++i ) {
00109 arx(d,i) = pop[i][d];
00110 }
00111 }
00112
00113
00114 Eigen::VectorXd weights = edoCMAESweights( lambda );
00115 assert( weights.size() == lambda );
00116
00117
00118
00119
00120 double mueff = pow(weights.sum(), 2) / (weights.array().square()).sum();
00121
00122
00123 double cc = (4+mueff/N) / (N+4 + 2*mueff/N);
00124
00125
00126 double cs = (mueff+2) / (N+mueff+5);
00127
00128
00129 double c1 = 2 / (pow(N+1.3,2)+mueff);
00130
00131
00132 double cmu = 2 * (mueff-2+1/mueff) / ( pow(N+2,2)+mueff);
00133
00134
00135 double damps = 1 + 2*std::max(0.0, sqrt((mueff-1)/(N+1))-1) + cs;
00136
00137
00138
00139 D& d = this->distribution();
00140
00141
00142 Matrix invsqrtC =
00143 d.coord_sys() * d.scaling().asDiagonal().inverse()
00144 * d.coord_sys().transpose();
00145 assert( invsqrtC.innerSize() == d.coord_sys().innerSize() );
00146 assert( invsqrtC.outerSize() == d.coord_sys().outerSize() );
00147
00148
00149 double chiN = sqrt(N)*(1-1/(4*N)+1/(21*pow(N,2)));
00150
00151
00152
00153
00154
00155
00156
00157 Vector xold = d.mean();
00158 assert( xold.size() == N );
00159
00160 Vector xmean = arx * weights;
00161 assert( xmean.size() == N );
00162 d.mean( xmean );
00163
00164
00165
00166
00167
00168
00169
00170 d.path_sigma(
00171 (1.0-cs)*d.path_sigma() + sqrt(cs*(2.0-cs)*mueff)*invsqrtC*(xmean-xold)/d.sigma()
00172 );
00173
00174
00175 double hsig;
00176 if( d.path_sigma().norm()/sqrt(1.0-pow((1.0-cs),(2.0*counteval/lambda)))/chiN
00177 < 1.4 + 2.0/(N+1.0)
00178 ) {
00179 hsig = 1.0;
00180 } else {
00181 hsig = 0.0;
00182 }
00183
00184
00185 d.path_covar(
00186 (1.0-cc)*d.path_covar() + hsig*sqrt(cc*(2.0-cc)*mueff)*(xmean-xold) / d.sigma()
00187 );
00188
00189 Matrix xmu( N, lambda);
00190 xmu = xold.rowwise().replicate(lambda);
00191 assert( xmu.innerSize() == N );
00192 assert( xmu.outerSize() == lambda );
00193 Matrix artmp = (1.0/d.sigma()) * (arx - xmu);
00194
00195 assert( artmp.innerSize() == N && artmp.outerSize() == lambda );
00196
00197
00198
00199
00200
00201
00202 d.covar(
00203 (1-c1-cmu) * d.covar()
00204 + c1 * (d.path_covar()*d.path_covar().transpose()
00205 + (1-hsig) * cc*(2-cc) * d.covar())
00206 + cmu * artmp * weights.asDiagonal() * artmp.transpose()
00207 );
00208
00209
00210 d.sigma( d.sigma() * exp((cs/damps)*(d.path_sigma().norm()/chiN - 1)) );
00211
00212
00213
00214
00215
00216
00217
00218
00219 if( counteval - _eigeneval > lambda/(c1+cmu)/N/10 ) {
00220 _eigeneval = counteval;
00221
00222
00223 Matrix C = d.covar();
00224
00225
00226
00227
00228 d.covar( C );
00229
00230 Eigen::SelfAdjointEigenSolver<Matrix> eigensolver( d.covar() );
00231 d.coord_sys( eigensolver.eigenvectors() );
00232 Matrix mD = eigensolver.eigenvalues().asDiagonal();
00233 assert( mD.innerSize() == N && mD.outerSize() == N );
00234
00235
00236 mD.cwiseSqrt();
00237 d.scaling( mD.diagonal() );
00238 }
00239
00240 return d;
00241 }
00242
00243 protected:
00244
00245 unsigned int _calls;
00246 unsigned int _eigeneval;
00247
00248
00249
00250 };
00251 #endif // WITH_EIGEN
00252
00253 #endif // !_edoEstimatorNormalAdaptive_h