00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_ANNEAL_MT_H 00024 #define O2SCL_ANNEAL_MT_H 00025 00026 #include <o2scl/gsl_rnga.h> 00027 #include <o2scl/ovector_tlate.h> 00028 #include <o2scl/multi_funct.h> 00029 00030 // for boost::bind() 00031 #include <boost/bind.hpp> 00032 // for boost::thread 00033 #include <boost/thread/thread.hpp> 00034 00035 #ifndef DOXYGENP 00036 namespace o2scl { 00037 #endif 00038 00039 /** 00040 \brief Multidimensional minimization by simulated annealing 00041 (Boost multi-threaded version) 00042 00043 This header-only class additionally requires the Boost 00044 libraries. It performs simulated annealing using an arbitrary 00045 number of processors using <tt>boost::thread</tt>, which is 00046 closely related to the standard Unix pthread library. It works 00047 very similarly to \ref gsl_anneal, it performs \ref ntrial 00048 evaluations over each processor, then applying the metropolis 00049 algorithm to the results from all of the processors at the end. 00050 00051 Because <tt>np</tt> function calls are happening simultaneously, 00052 where <tt>np</tt> is the number of processors, <tt>np</tt> 00053 copies of the function parameters of type <tt>param_t</tt> must 00054 also be specified. The user-defined function to minimize must 00055 also be thread-safe, allowing multiple threads to call the 00056 function at once (albeit given different parameters). The 00057 default type for these <tt>np</tt> copies of the parameters of 00058 type <tt>param_t</tt> is <tt>std::vector<param_t></tt>. 00059 00060 This works particularly well for functions which are not trivial 00061 to evaluate, i.e. functions where the execution time is more 00062 longer than the bookkeeping that \ref anneal_mt performs between 00063 trials. For functions which satisfy this requirement, this 00064 algorithm scales nearly linearly with the number of processors. 00065 00066 Verbosse I/O for this class happens only outside the theads 00067 unless the user places I/O in the streams in the function 00068 that is specified. 00069 00070 \future There may be a good way to remove the function indirection 00071 here to make this class a bit faster. 00072 */ 00073 template<class param_t, class param_vec_t=std::vector<param_t>, 00074 class func_t=multi_funct<param_t>, class vec_t=ovector_base, 00075 class alloc_vec_t=ovector, class alloc_t=ovector_alloc, 00076 class rng_t=gsl_rnga> class anneal_mt { 00077 00078 public: 00079 00080 anneal_mt() { 00081 ntrial=100; 00082 nproc=1; 00083 outs=&std::cout; 00084 ins=&std::cin; 00085 tolx=1.0e-6; 00086 T_start=1.0; 00087 T_dec=1.5; 00088 step_dec=1.5; 00089 min_step_ratio=100.0; 00090 step_vec.allocate(1); 00091 step_vec[0]=1.0; 00092 out_best=false; 00093 out_step_changes=false; 00094 } 00095 00096 virtual ~anneal_mt() { 00097 step_vec.free(); 00098 } 00099 00100 /// \name Basic usage 00101 //@{ 00102 /** \brief Calculate the minimum \c fmin of \c func w.r.t the 00103 array \c x0 of size \c nv using \c np threads. 00104 */ 00105 virtual int mmin(size_t nv, vec_t &x0, double &fmin, 00106 func_t &func, size_t np, param_vec_t &pars) { 00107 00108 if (nv==0) { 00109 O2SCL_ERR2_RET("Tried to minimize over zero variables ", 00110 " in anneal_mt::mmin().",gsl_einval); 00111 } 00112 if (np==0) { 00113 O2SCL_ERR2_RET("Tried to use 0 threads in ", 00114 "anneal_mt::mmin().",gsl_einval); 00115 } 00116 00117 allocate(nv,np); 00118 params=&pars; 00119 f=&func; 00120 fmin=0.0; 00121 00122 double E, best_E, T, old_E; 00123 int i, iter=0; 00124 size_t j; 00125 00126 for(j=0;j<nv;j++) { 00127 x[j]=x0[j]; 00128 best_x[j]=x0[j]; 00129 } 00130 00131 E=func(nv,x,pars[0]); 00132 best_E=E; 00133 00134 // Setup initial temperature and step sizes 00135 start(nv,T); 00136 00137 bool done=false; 00138 00139 boost::thread **thrd; 00140 00141 while (!done) { 00142 00143 // Copy old value of x for next() function 00144 00145 for(j=0;j<nv;j++) old_x[j]=x[j]; 00146 old_E=E; 00147 00148 size_t nmoves=0; 00149 00150 for (i=0;i<ntrial;++i) { 00151 00152 // Determine the stepsize, create and execute the new threads 00153 thrd=new boost::thread *[np]; 00154 for(size_t ip=0;ip<np;ip++) { 00155 new_x[ip]=x; 00156 thrd[ip]=new boost::thread 00157 (boost::bind(&anneal_mt::func_wrapper,this,ip)); 00158 } 00159 // Wait until all the threads are done 00160 for(size_t ip=0;ip<np;ip++) { 00161 thrd[ip]->join(); 00162 } 00163 // Delete the threads and continue 00164 for(size_t ip=0;ip<np;ip++) { 00165 delete thrd[ip]; 00166 } 00167 delete[] thrd; 00168 // Process the results from each thread 00169 for(size_t ip=0;ip<np;ip++) { 00170 00171 // Store best value obtained so far 00172 if(new_E[ip]<=best_E){ 00173 for(j=0;j<nv;j++) best_x[j]=new_x[ip][j]; 00174 best_E=new_E[ip]; 00175 if (this->verbose>0 && out_best) { 00176 std::cout << "Best: " << best_x << " " << best_E << std::endl; 00177 } 00178 } 00179 00180 // Take the crucial step: see if the new point is accepted 00181 // or not, as determined by the boltzmann probability 00182 if (new_E[ip]<E) { 00183 for(j=0;j<nv;j++) x[j]=new_x[ip][j]; 00184 E=new_E[ip]; 00185 nmoves++; 00186 } else if (this->def_rng.random() < exp(-(new_E[ip]-E)/(boltz*T)) ) { 00187 for(j=0;j<nv;j++) x[j]=new_x[ip][j]; 00188 E=new_E[ip]; 00189 nmoves++; 00190 } 00191 } 00192 00193 } 00194 00195 if (this->verbose>0) { 00196 print_iter(nv,best_x,best_E,iter,T,"anneal_mt"); 00197 iter++; 00198 } 00199 00200 // See if we're finished and proceed to the next step 00201 next(nv,old_x,old_E,x,E,T,nmoves,done); 00202 00203 } 00204 00205 for(j=0;j<nv;j++) x0[j]=best_x[j]; 00206 fmin=best_E; 00207 00208 free(np); 00209 00210 return 0; 00211 } 00212 //@} 00213 00214 /// \name Iteration control 00215 //@{ 00216 /// Determine how to change the minimization for the next iteration 00217 virtual int next(size_t nv, vec_t &x_old, double min_old, vec_t &x_new, 00218 double min_new, double &T, size_t n_moves, 00219 bool &finished) { 00220 00221 if (T/T_dec<this->tolx) { 00222 finished=true; 00223 return gsl_success; 00224 } 00225 if (n_moves==0) { 00226 bool changed=false; 00227 for(size_t i=0;i<nv;i++) { 00228 if (i<step_vec.size() && step_vec[i]>this->tolx*min_step_ratio) { 00229 step_vec[i]/=step_dec; 00230 changed=true; 00231 } 00232 } 00233 if (changed && verbose>0 && out_step_changes) { 00234 std::cout << "Step sizes changed: " << step_vec << std::endl; 00235 } 00236 } 00237 T/=T_dec; 00238 return gsl_success; 00239 } 00240 00241 /// Setup initial temperature and stepsize 00242 virtual int start(size_t nv, double &T) { 00243 T=T_start; 00244 return gsl_success; 00245 } 00246 //@} 00247 00248 /// \name Parameters 00249 //@{ 00250 /// Boltzmann factor (default 1.0). 00251 double boltz; 00252 00253 /// Number of iterations 00254 int ntrial; 00255 00256 /// Output control 00257 int verbose; 00258 00259 /// Output step size changes (default false) 00260 bool out_step_changes; 00261 00262 /// Output best point (default false) 00263 bool out_best; 00264 00265 /// The independent variable tolerance (default \f$ 10^{-6} \f$ ) 00266 double tolx; 00267 00268 /// Initial temperature (default 1.0) 00269 double T_start; 00270 00271 /// Factor to decrease temperature by (default 1.5) 00272 double T_dec; 00273 00274 /// Factor to decrease step size by (default 1.5) 00275 double step_dec; 00276 00277 /// Ratio between minimum step size and \ref tolx (default 100.0) 00278 double min_step_ratio; 00279 //@} 00280 00281 /// The default random number generator 00282 rng_t def_rng; 00283 00284 /// Return string denoting type ("anneal_mt") 00285 virtual const char *type() { return "anneal_mt"; } 00286 00287 /// Set streams for verbose I/O 00288 int set_verbose_stream(std::ostream &out, std::istream &in) { 00289 outs=&out; 00290 ins=∈ 00291 return 0; 00292 } 00293 00294 /** 00295 \brief Print out iteration information. 00296 00297 Depending on the value of the variable verbose, this prints out 00298 the iteration information. If verbose=0, then no information is 00299 printed. If verbose>0, then after each iteration, the present 00300 values of x and y are output to std::cout along with the 00301 iteration number. Also, if verbose>0, every time a new smallest 00302 function value is found, the location and the function value is 00303 output. If verbose>=2 then each iteration waits for a character 00304 between each trial. 00305 */ 00306 virtual int print_iter(size_t nv, vec_t &xx, double y, int iter, 00307 double tptr, std::string comment) 00308 { 00309 if (this->verbose<=0) return 0; 00310 00311 size_t i; 00312 char ch; 00313 00314 (*this->outs) << comment << " Iteration: " << iter << std::endl; 00315 text_out_file outsf(this->outs,79); 00316 outsf.word_out("x:"); 00317 for(i=0;i<nv;i++) outsf.double_out(xx[i]); 00318 outsf.end_line(); 00319 (*this->outs) << "y: " << y << " Tptr: " << tptr << std::endl; 00320 if (this->verbose>1) { 00321 (*this->outs) << "Press a key and type enter to continue. "; 00322 (*this->ins) >> ch; 00323 } 00324 00325 return 0; 00326 } 00327 00328 /// Set the step sizes 00329 template<class vec2_t> int set_step(size_t nv, vec2_t &stepv) { 00330 if (nv>0) { 00331 step_vec.free(); 00332 step_vec.allocate(nv); 00333 for(size_t i=0;i<nv;i++) step_vec[i]=stepv[i]; 00334 } 00335 return 0; 00336 } 00337 00338 #ifndef DOXYGEN_INTERNAL 00339 00340 protected: 00341 00342 /// The function wrapper executed by thread with index \c ip 00343 void func_wrapper(size_t ip) { 00344 step(nvar,new_x[ip]); 00345 (*f)(nvar,new_x[ip],new_E[ip],(*params)[ip]); 00346 } 00347 00348 /// Stream for verbose output 00349 std::ostream *outs; 00350 00351 /// Stream for verbose input 00352 std::istream *ins; 00353 00354 /// The number of threads to run 00355 size_t nproc; 00356 00357 /// The user-specified parameters 00358 param_vec_t *params; 00359 00360 /// The number of variables over which we minimize 00361 size_t nvar; 00362 00363 /// The function to minimize 00364 func_t *f; 00365 00366 /// \name Storage for present, next, and best vectors 00367 //@{ 00368 alloc_vec_t x, *new_x, best_x, new_E, old_x; 00369 //@} 00370 00371 /// Allocation object 00372 alloc_t ao; 00373 00374 /// Vector of step sizes 00375 ovector step_vec; 00376 00377 /** \brief Allocate memory for a minimizer over \c n dimensions 00378 with stepsize \c step 00379 */ 00380 virtual int allocate(size_t nv, size_t np) { 00381 nvar=nv; 00382 nproc=np; 00383 ao.allocate(x,nvar); 00384 new_x=new alloc_vec_t[np]; 00385 ao.allocate(old_x,nvar); 00386 for(size_t i=0;i<np;i++) { 00387 ao.allocate(new_x[i],nvar); 00388 } 00389 ao.allocate(new_E,np); 00390 ao.allocate(best_x,nvar); 00391 return 0; 00392 } 00393 00394 /// Free allocated memory 00395 virtual int free(size_t np) { 00396 ao.free(x); 00397 ao.free(old_x); 00398 ao.free(new_E); 00399 ao.free(best_x); 00400 for(size_t i=0;i<np;i++) { 00401 ao.free(new_x[i]); 00402 } 00403 delete[] new_x; 00404 return 0; 00405 } 00406 00407 /// Make a step to a new attempted minimum 00408 virtual int step(size_t nv, vec_t &sx) { 00409 size_t nstep=step_vec.size(); 00410 for(size_t i=0;i<nv;i++) { 00411 double u=this->def_rng.random(); 00412 sx[i]=(2.0*u-1.0)*step_vec[i%nstep]+sx[i]; 00413 } 00414 return 0; 00415 } 00416 00417 #endif 00418 00419 }; 00420 00421 #ifndef DOXYGENP 00422 } 00423 #endif 00424 00425 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page