edoEstimatorNormalAdaptive.h
00001 /*
00002 The Evolving Distribution Objects framework (EDO) is a template-based,
00003 ANSI-C++ evolutionary computation library which helps you to write your
00004 own estimation of distribution algorithms.
00005 
00006 This library is free software; you can redistribute it and/or
00007 modify it under the terms of the GNU Lesser General Public
00008 License as published by the Free Software Foundation; either
00009 version 2.1 of the License, or (at your option) any later version.
00010 
00011 This library is distributed in the hope that it will be useful,
00012 but WITHOUT ANY WARRANTY; without even the implied warranty of
00013 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00014 Lesser General Public License for more details.
00015 
00016 You should have received a copy of the GNU Lesser General Public
00017 License along with this library; if not, write to the Free Software
00018 Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
00019 
00020 Copyright (C) 2010 Thales group
00021 */
00022 /*
00023 Authors:
00024     Johann Dréo <johann.dreo@thalesgroup.com>
00025     Pierre Savéant <pierre.saveant@thalesgroup.com>
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; // column vectors @see edoNormalAdaptive
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         // compute recombination weights
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         // normalization of weights
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     // update the distribution reference this->distribution()
00086     edoNormalAdaptive<EOT> operator()( eoPop<EOT>& pop )
00087     {
00088 
00089         /**********************************************************************
00090          * INITIALIZATION
00091          *********************************************************************/
00092 
00093         unsigned int N = pop[0].size(); // FIXME expliciter la dimension du pb ?
00094         unsigned int lambda = pop.size();
00095 
00096         // number of calls to the operator == number of generations
00097         _calls++;
00098         // number of "evaluations" until now
00099         unsigned int counteval = _calls * lambda;
00100 
00101         // Here, if we are in canonical CMA-ES,
00102         // pop is supposed to be the mu ranked better solutions,
00103         // as the rank mu selection is supposed to have occured.
00104         Matrix arx( N, lambda );
00105 
00106         // copy the pop (most probably a vector of vectors) in a Eigen3 matrix
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]; // NOTE: pop = arx.transpose()
00110             } // dimensions
00111         } // individuals
00112 
00113         // muXone array for weighted recombination
00114         Eigen::VectorXd weights = edoCMAESweights( lambda );
00115         assert( weights.size() == lambda );
00116 
00117         // FIXME exposer les constantes dans l'interface
00118 
00119         // variance-effectiveness of sum w_i x_i
00120         double mueff = pow(weights.sum(), 2) / (weights.array().square()).sum();
00121 
00122         // time constant for cumulation for C
00123         double cc = (4+mueff/N) / (N+4 + 2*mueff/N);
00124 
00125         // t-const for cumulation for sigma control
00126         double cs = (mueff+2) / (N+mueff+5);
00127 
00128         // learning rate for rank-one update of C
00129         double c1 = 2 / (pow(N+1.3,2)+mueff);
00130 
00131         // and for rank-mu update
00132         double cmu = 2 * (mueff-2+1/mueff) / ( pow(N+2,2)+mueff);
00133 
00134         // damping for sigma
00135         double damps = 1 + 2*std::max(0.0, sqrt((mueff-1)/(N+1))-1) + cs;
00136 
00137 
00138         // shortcut to the referenced distribution
00139         D& d = this->distribution();
00140 
00141         // C^-1/2
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         // expectation of ||N(0,I)|| == norm(randn(N,1))
00149         double chiN = sqrt(N)*(1-1/(4*N)+1/(21*pow(N,2)));
00150 
00151 
00152         /**********************************************************************
00153          * WEIGHTED MEAN
00154          *********************************************************************/
00155 
00156         // compute weighted mean into xmean
00157         Vector xold = d.mean();
00158         assert( xold.size() == N );
00159         //  xmean ( N, 1 ) = arx( N, lambda ) * weights( lambda, 1 )
00160         Vector xmean = arx * weights;
00161         assert( xmean.size() == N );
00162         d.mean( xmean );
00163 
00164 
00165         /**********************************************************************
00166          * CUMULATION: UPDATE EVOLUTION PATHS
00167          *********************************************************************/
00168 
00169         // cumulation for sigma
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         // sign of h
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         // cumulation for the covariance matrix
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         // Matrix artmp = (1.0/d.sigma()) * arx - xold.colwise().replicate(lambda);
00195         assert( artmp.innerSize() == N && artmp.outerSize() == lambda );
00196 
00197 
00198         /**********************************************************************
00199          * COVARIANCE MATRIX ADAPTATION
00200          *********************************************************************/
00201 
00202         d.covar(
00203                 (1-c1-cmu) * d.covar()                            // regard old matrix
00204                 + c1 * (d.path_covar()*d.path_covar().transpose() // plus rank one update
00205                     + (1-hsig) * cc*(2-cc) * d.covar())   // minor correction if hsig==0
00206                 + cmu * artmp * weights.asDiagonal() * artmp.transpose() // plus rank mu update
00207                );
00208 
00209         // Adapt step size sigma
00210         d.sigma( d.sigma() * exp((cs/damps)*(d.path_sigma().norm()/chiN - 1)) );
00211 
00212 
00213 
00214         /**********************************************************************
00215          * DECOMPOSITION OF THE COVARIANCE MATRIX
00216          *********************************************************************/
00217 
00218         // Decomposition of C into B*diag(D.^2)*B' (diagonalization)
00219         if( counteval - _eigeneval > lambda/(c1+cmu)/N/10 ) {  // to achieve O(N^2)
00220             _eigeneval = counteval;
00221 
00222             // enforce symmetry of the covariance matrix
00223             Matrix C = d.covar();
00224             // FIXME edoEstimatorNormalAdaptive.h:213:44: erreur: expected primary-expression before ‘)’ token
00225             // copy the upper part in the lower one
00226             //C.triangularView<Eigen::Lower>() = C.adjoint();
00227             // Matrix CS = C.triangularView<Eigen::Upper>() + C.triangularView<Eigen::StrictlyUpper>().transpose();
00228             d.covar( C );
00229 
00230             Eigen::SelfAdjointEigenSolver<Matrix> eigensolver( d.covar() ); // FIXME use JacobiSVD?
00231             d.coord_sys( eigensolver.eigenvectors() );
00232             Matrix mD = eigensolver.eigenvalues().asDiagonal();
00233             assert( mD.innerSize() == N && mD.outerSize() == N );
00234 
00235             // from variance to standard deviations
00236             mD.cwiseSqrt();
00237             d.scaling( mD.diagonal() );
00238         }
00239 
00240         return d;
00241     } // operator()
00242 
00243 protected:
00244 
00245     unsigned int _calls;
00246     unsigned int _eigeneval;
00247 
00248 
00249     // D & distribution() inherited from edoEstimatorAdaptive
00250 };
00251 #endif // WITH_EIGEN
00252 
00253 #endif // !_edoEstimatorNormalAdaptive_h
 All Classes Functions Variables Typedefs