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 See also the unthreaded version of this class in \ref gsl_anneal. 00044 00045 \future There may be a good way to remove the function indirection 00046 here to make this class a bit faster. 00047 */ 00048 template<class param_t, class param_vec_t=std::vector<param_t>, 00049 class func_t=multi_funct<param_t>, class vec_t=ovector_base, 00050 class alloc_vec_t=ovector, class alloc_t=ovector_alloc, 00051 class rng_t=gsl_rnga> class anneal_mt { 00052 00053 public: 00054 00055 anneal_mt() { 00056 ntrial=100; 00057 nstep=0; 00058 nproc=1; 00059 outs=&std::cout; 00060 ins=&std::cin; 00061 } 00062 00063 virtual ~anneal_mt() { 00064 if (nstep>0) { 00065 delete[] step_sizes; 00066 } 00067 } 00068 00069 /// \name Basic usage 00070 //@{ 00071 /** \brief Calculate the minimum \c fmin of \c func w.r.t the 00072 array \c x0 of size \c nv using \c np threads. 00073 */ 00074 virtual int mmin(size_t nv, vec_t &x0, double &fmin, 00075 func_t &func, size_t np, param_vec_t &pars) { 00076 00077 if (nv==0) { 00078 O2SCL_ERR2_RET("Tried to minimize over zero variables ", 00079 " in anneal_mt::mmin().",gsl_einval); 00080 } 00081 if (np==0) { 00082 O2SCL_ERR2_RET("Tried to use 0 threads in ", 00083 "anneal_mt::mmin().",gsl_einval); 00084 } 00085 00086 allocate(nv,np); 00087 params=&pars; 00088 f=&func; 00089 fmin=0.0; 00090 00091 double E, best_E, T, old_E; 00092 int i, iter=0; 00093 size_t j; 00094 00095 for(j=0;j<nv;j++) { 00096 x[j]=x0[j]; 00097 best_x[j]=x0[j]; 00098 } 00099 00100 E=func(nv,x,pars[0]); 00101 best_E=E; 00102 00103 // Setup initial temperature and step sizes 00104 start(nv,T); 00105 00106 bool done=false; 00107 00108 boost::thread **thrd; 00109 00110 while (!done) { 00111 00112 // Copy old value of x for next() function 00113 00114 for(j=0;j<nv;j++) old_x[j]=x[j]; 00115 old_E=E; 00116 00117 size_t nmoves=0; 00118 00119 for (i=0;i<this->ntrial;++i) { 00120 00121 // Determine the stepsize and call the function 00122 thrd=new boost::thread *[np]; 00123 for(size_t ip=0;ip<np;ip++) { 00124 new_x[ip]=x; 00125 thrd[ip]=new boost::thread 00126 (boost::bind(&anneal_mt::func_wrapper,this,ip)); 00127 } 00128 // Wait until all the threads are done 00129 for(size_t ip=0;ip<np;ip++) { 00130 thrd[ip]->join(); 00131 } 00132 // Delete the threads and continue 00133 for(size_t ip=0;ip<np;ip++) { 00134 delete thrd[ip]; 00135 } 00136 delete[] thrd; 00137 // Process the results from each thread 00138 for(size_t ip=0;ip<np;ip++) { 00139 00140 // Store best value obtained so far 00141 if(new_E[ip]<=best_E){ 00142 for(j=0;j<nv;j++) best_x[j]=new_x[ip][j]; 00143 best_E=new_E[ip]; 00144 std::cout << "Best: " << best_x << " " << best_E << std::endl; 00145 } 00146 00147 // Take the crucial step: see if the new point is accepted 00148 // or not, as determined by the boltzmann probability 00149 if (new_E[ip]<E) { 00150 for(j=0;j<nv;j++) x[j]=new_x[ip][j]; 00151 E=new_E[ip]; 00152 nmoves++; 00153 } else if (this->def_rng.random() < exp(-(new_E[ip]-E)/(boltz*T)) ) { 00154 for(j=0;j<nv;j++) x[j]=new_x[ip][j]; 00155 E=new_E[ip]; 00156 nmoves++; 00157 } 00158 } 00159 00160 } 00161 00162 if (this->verbose>0) { 00163 print_iter(nv,best_x,best_E,iter,T,"gsl_anneal"); 00164 iter++; 00165 } 00166 00167 // See if we're finished and proceed to the next step 00168 next(nv,old_x,old_E,x,E,T,nmoves,done); 00169 00170 } 00171 00172 for(j=0;j<nv;j++) x0[j]=best_x[j]; 00173 fmin=best_E; 00174 00175 free(np); 00176 00177 return 0; 00178 } 00179 //@} 00180 00181 /// \name Iteration control 00182 //@{ 00183 /// Determine how to change the minimization for the next iteration 00184 virtual int next(size_t nv, vec_t &x_old, double min_old, vec_t &x_new, 00185 double min_new, double &T, size_t n_moves, 00186 bool &finished) { 00187 00188 if (T/2.0<this->tolx) { 00189 finished=true; 00190 return gsl_success; 00191 } 00192 if (n_moves==0) { 00193 for(size_t i=0;i<nv;i++) { 00194 if (step_sizes[i]>this->tolx*1.0e2) { 00195 step_sizes[i]/=1.5; 00196 } 00197 } 00198 } 00199 T/=1.5; 00200 return gsl_success; 00201 } 00202 00203 /// Setup initial temperature and stepsize 00204 virtual int start(size_t nv, double &T) { 00205 T=1.0; 00206 nstep=nv; 00207 step_sizes=new double[nv]; 00208 for(size_t i=0;i<nv;i++) { 00209 step_sizes[i]=10.0; 00210 } 00211 00212 return gsl_success; 00213 } 00214 //@} 00215 00216 /// \name Parameters 00217 //@{ 00218 /// Boltzmann factor (default 1.0). 00219 double boltz; 00220 00221 /// Number of iterations 00222 int ntrial; 00223 00224 /// Output control 00225 int verbose; 00226 00227 /// The independent variable tolerance 00228 double tolx; 00229 //@} 00230 00231 /// The default random number generator 00232 rng_t def_rng; 00233 00234 /// Return string denoting type ("anneal_mt") 00235 virtual const char *type() { return "anneal_mt"; } 00236 00237 /** 00238 \brief Print out iteration information. 00239 00240 Depending on the value of the variable verbose, this prints out 00241 the iteration information. If verbose=0, then no information is 00242 printed, while if verbose>1, then after each iteration, the 00243 present values of x and y are output to std::cout along with the 00244 iteration number. If verbose>=2 then each iteration waits for a 00245 character. 00246 */ 00247 virtual int print_iter(size_t nv, vec_t &xx, double y, int iter, 00248 double tptr, std::string comment) 00249 { 00250 if (this->verbose<=0) return 0; 00251 00252 size_t i; 00253 char ch; 00254 00255 (*this->outs) << comment << " Iteration: " << iter << std::endl; 00256 text_out_file outsf(this->outs,79); 00257 outsf.word_out("x:"); 00258 for(i=0;i<nv;i++) outsf.double_out(xx[i]); 00259 outsf.end_line(); 00260 (*this->outs) << "y: " << y << " Tptr: " << tptr << std::endl; 00261 if (this->verbose>1) { 00262 (*this->outs) << "Press a key and type enter to continue. "; 00263 (*this->ins) >> ch; 00264 } 00265 00266 return 0; 00267 } 00268 00269 #ifndef DOXYGEN_INTERNAL 00270 00271 protected: 00272 00273 /// The function wrapper executed by thread with index \c ip 00274 void func_wrapper(size_t ip) { 00275 step(new_x[ip],nvar); 00276 (*f)(nvar,new_x[ip],new_E[ip],(*params)[ip]); 00277 } 00278 00279 /// Stream for verbose output 00280 std::ostream *outs; 00281 00282 /// Stream for verbose input 00283 std::istream *ins; 00284 00285 /// The number of threads to run 00286 size_t nproc; 00287 00288 /// The user-specified parameters 00289 param_vec_t *params; 00290 00291 /// The number of variables over which we minimize 00292 size_t nvar; 00293 00294 /// The function to minimize 00295 func_t *f; 00296 00297 /// \name Storage for present, next, and best vectors 00298 //@{ 00299 alloc_vec_t x, *new_x, best_x, new_E, old_x; 00300 //@} 00301 00302 /// Allocation object 00303 alloc_t ao; 00304 00305 /// Number of step sizes 00306 size_t nstep; 00307 00308 /// Step sizes 00309 double *step_sizes; 00310 00311 /** \brief Allocate memory for a minimizer over \c n dimensions 00312 with stepsize \c step 00313 */ 00314 virtual int allocate(size_t nv, size_t np) { 00315 nvar=nv; 00316 nproc=np; 00317 ao.allocate(x,nvar); 00318 new_x=new alloc_vec_t[np]; 00319 ao.allocate(old_x,nvar); 00320 for(size_t i=0;i<np;i++) { 00321 ao.allocate(new_x[i],nvar); 00322 } 00323 ao.allocate(new_E,np); 00324 ao.allocate(best_x,nvar); 00325 return 0; 00326 } 00327 00328 /// Free allocated memory 00329 virtual int free(size_t np) { 00330 ao.free(x); 00331 ao.free(old_x); 00332 ao.free(new_E); 00333 ao.free(best_x); 00334 for(size_t i=0;i<np;i++) { 00335 ao.free(new_x[i]); 00336 } 00337 delete[] new_x; 00338 return 0; 00339 } 00340 00341 /// Make a step to a new attempted minimum 00342 virtual int step(vec_t &sx, int nv) { 00343 for(int i=0;i<nv;i++) { 00344 double u=this->def_rng.random(); 00345 sx[i]=(2.0*u-1.0)*step_sizes[i%nstep]+sx[i]; 00346 } 00347 return 0; 00348 } 00349 00350 #endif 00351 00352 }; 00353 00354 #ifndef DOXYGENP 00355 } 00356 #endif 00357 00358 #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