All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
anneal_mt.h
1 /*
2  -------------------------------------------------------------------
3 
4  Copyright (C) 2006-2014, Andrew W. Steiner
5 
6  This file is part of O2scl.
7 
8  O2scl is free software; you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version.
12 
13  O2scl is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with O2scl. If not, see <http://www.gnu.org/licenses/>.
20 
21  -------------------------------------------------------------------
22 */
23 #ifndef O2SCL_ANNEAL_MT_H
24 #define O2SCL_ANNEAL_MT_H
25 
26 #include <o2scl/rng_gsl.h>
27 #include <o2scl/multi_funct.h>
28 
29 // for boost::bind()
30 #include <boost/bind.hpp>
31 // for boost::thread
32 #include <boost/thread/thread.hpp>
33 
34 #ifndef DOXYGEN_NO_O2NS
35 namespace o2scl {
36 #endif
37 
38  /** \brief Multidimensional minimization by simulated annealing
39  (Boost multi-threaded version)
40 
41  This header-only class additionally requires the Boost
42  libraries. It performs simulated annealing using an arbitrary
43  number of processors using <tt>boost::thread</tt>, which is
44  closely related to the standard Unix pthread library. It works
45  very similarly to \ref anneal_gsl, it performs \ref ntrial
46  evaluations over each processor, then applying the metropolis
47  algorithm to the results from all of the processors at the end.
48 
49  Because <tt>np</tt> function calls are happening simultaneously,
50  where <tt>np</tt> is the number of processors, <tt>np</tt>
51  copies of the function parameters of type <tt>param_t</tt> must
52  also be specified. The user-defined function to minimize must
53  also be thread-safe, allowing multiple threads to call the
54  function at once (albeit given different parameters). The
55  default type for these <tt>np</tt> copies of the parameters of
56  type <tt>param_t</tt> is <tt>std::vector<param_t></tt>.
57 
58  This works particularly well for functions which are not trivial
59  to evaluate, i.e. functions where the execution time is more
60  longer than the bookkeeping that \ref anneal_mt performs between
61  trials. For functions which satisfy this requirement, this
62  algorithm scales nearly linearly with the number of processors.
63 
64  Verbose I/O for this class happens only outside the theads
65  unless the user places I/O in the streams in the function that
66  is specified.
67 
68  \future There may be a good way to remove the function indirection
69  here to make this class a bit faster.
70  */
71  template<class param_t, class param_vec_t=std::vector<param_t>,
72  class func_t=multi_funct<param_t>, class vec_t=ubvector,
73  class alloc_vec_t=ubvector, class alloc_t=ubvector_alloc,
74  class rng_t=rng_gsl> class anneal_mt {
75 
76  public:
77 
78  anneal_mt() {
79  ntrial=100;
80  nproc=1;
81  outs=&std::cout;
82  ins=&std::cin;
83  tolx=1.0e-6;
84  T_start=1.0;
85  T_dec=1.5;
86  step_dec=1.5;
87  min_step_ratio=100.0;
88  step_vec.allocate(1);
89  step_vec[0]=1.0;
90  out_best=false;
91  out_step_changes=false;
92  }
93 
94  virtual ~anneal_mt() {
95  step_vec.free();
96  }
97 
98  /// \name Basic usage
99  //@{
100  /** \brief Calculate the minimum \c fmin of \c func w.r.t the
101  array \c x0 of size \c nv using \c np threads.
102  */
103  virtual int mmin(size_t nv, vec_t &x0, double &fmin,
104  func_t &func, size_t np, param_vec_t &pars) {
105 
106  if (nv==0) {
107  O2SCL_ERR2("Tried to minimize over zero variables ",
108  " in anneal_mt::mmin().",gsl_einval);
109  }
110  if (np==0) {
111  O2SCL_ERR2("Tried to use zero threads in ",
112  "anneal_mt::mmin().",gsl_einval);
113  }
114 
115  allocate(nv,np);
116  params=&pars;
117  f=&func;
118  fmin=0.0;
119 
120  double E, best_E, T, old_E;
121  int i, iter=0;
122  size_t j;
123 
124  for(j=0;j<nv;j++) {
125  x[j]=x0[j];
126  best_x[j]=x0[j];
127  }
128 
129  E=func(nv,x,pars[0]);
130  best_E=E;
131 
132  // Setup initial temperature and step sizes
133  start(nv,T);
134 
135  bool done=false;
136 
137  boost::thread **thrd;
138 
139  while (!done) {
140 
141  // Copy old value of x for next() function
142 
143  for(j=0;j<nv;j++) old_x[j]=x[j];
144  old_E=E;
145 
146  size_t nmoves=0;
147 
148  for (i=0;i<ntrial;++i) {
149 
150  // Determine the stepsize, create and execute the new threads
151  thrd=new boost::thread *[np];
152  for(size_t ip=0;ip<np;ip++) {
153  new_x[ip]=x;
154  thrd[ip]=new boost::thread
155  (boost::bind(&anneal_mt::func_wrapper,this,ip));
156  }
157  // Wait until all the threads are done
158  for(size_t ip=0;ip<np;ip++) {
159  thrd[ip]->join();
160  }
161  // Delete the threads and continue
162  for(size_t ip=0;ip<np;ip++) {
163  delete thrd[ip];
164  }
165  delete[] thrd;
166  // Process the results from each thread
167  for(size_t ip=0;ip<np;ip++) {
168 
169  // Store best value obtained so far
170  if(new_E[ip]<=best_E){
171  for(j=0;j<nv;j++) best_x[j]=new_x[ip][j];
172  best_E=new_E[ip];
173  if (this->verbose>0 && out_best) {
174  std::cout << "Best: " << best_x << " " << best_E << std::endl;
175  }
176  }
177 
178  // Take the crucial step: see if the new point is accepted
179  // or not, as determined by the boltzmann probability
180  if (new_E[ip]<E) {
181  for(j=0;j<nv;j++) x[j]=new_x[ip][j];
182  E=new_E[ip];
183  nmoves++;
184  } else if (this->def_rng.random() < exp(-(new_E[ip]-E)/(boltz*T)) ) {
185  for(j=0;j<nv;j++) x[j]=new_x[ip][j];
186  E=new_E[ip];
187  nmoves++;
188  }
189  }
190 
191  }
192 
193  if (this->verbose>0) {
194  print_iter(nv,best_x,best_E,iter,T,"anneal_mt");
195  iter++;
196  }
197 
198  // See if we're finished and proceed to the next step
199  next(nv,old_x,old_E,x,E,T,nmoves,done);
200 
201  }
202 
203  for(j=0;j<nv;j++) x0[j]=best_x[j];
204  fmin=best_E;
205 
206  free(np);
207 
208  return 0;
209  }
210  //@}
211 
212  /// \name Iteration control
213  //@{
214  /// Determine how to change the minimization for the next iteration
215  virtual int next(size_t nv, vec_t &x_old, double min_old, vec_t &x_new,
216  double min_new, double &T, size_t n_moves,
217  bool &finished) {
218 
219  if (T/T_dec<this->tolx) {
220  finished=true;
221  return gsl_success;
222  }
223  if (n_moves==0) {
224  bool changed=false;
225  for(size_t i=0;i<nv;i++) {
226  if (i<step_vec.size() && step_vec[i]>this->tolx*min_step_ratio) {
227  step_vec[i]/=step_dec;
228  changed=true;
229  }
230  }
231  if (changed && verbose>0 && out_step_changes) {
232  std::cout << "Step sizes changed: " << step_vec << std::endl;
233  }
234  }
235  T/=T_dec;
236  return gsl_success;
237  }
238 
239  /// Setup initial temperature and stepsize
240  virtual int start(size_t nv, double &T) {
241  T=T_start;
242  return gsl_success;
243  }
244  //@}
245 
246  /// \name Parameters
247  //@{
248  /// Boltzmann factor (default 1.0).
249  double boltz;
250 
251  /// Number of iterations
252  int ntrial;
253 
254  /// Output control
255  int verbose;
256 
257  /// Output step size changes (default false)
259 
260  /// Output best point (default false)
261  bool out_best;
262 
263  /// The independent variable tolerance (default \f$ 10^{-6} \f$ )
264  double tolx;
265 
266  /// Initial temperature (default 1.0)
267  double T_start;
268 
269  /// Factor to decrease temperature by (default 1.5)
270  double T_dec;
271 
272  /// Factor to decrease step size by (default 1.5)
273  double step_dec;
274 
275  /// Ratio between minimum step size and \ref tolx (default 100.0)
277  //@}
278 
279  /// The default random number generator
280  rng_t def_rng;
281 
282  /// Return string denoting type ("anneal_mt")
283  virtual const char *type() { return "anneal_mt"; }
284 
285  /// Set streams for verbose I/O
286  int set_verbose_stream(std::ostream &out, std::istream &in) {
287  outs=&out;
288  ins=&in;
289  return 0;
290  }
291 
292  /** \brief Print out iteration information.
293 
294  Depending on the value of the variable verbose, this prints out
295  the iteration information. If verbose=0, then no information is
296  printed. If verbose>0, then after each iteration, the present
297  values of x and y are output to std::cout along with the
298  iteration number. Also, if verbose>0, every time a new smallest
299  function value is found, the location and the function value is
300  output. If verbose>=2 then each iteration waits for a character
301  between each trial.
302  */
303  virtual int print_iter(size_t nv, vec_t &xx, double y, int iter,
304  double tptr, std::string comment)
305  {
306  if (this->verbose<=0) return 0;
307 
308  size_t i;
309  char ch;
310 
311  (*this->outs) << comment << " Iteration: " << iter << std::endl;
312  text_out_file outsf(this->outs,79);
313  outsf.word_out("x:");
314  for(i=0;i<nv;i++) outsf.double_out(xx[i]);
315  outsf.end_line();
316  (*this->outs) << "y: " << y << " Tptr: " << tptr << std::endl;
317  if (this->verbose>1) {
318  (*this->outs) << "Press a key and type enter to continue. ";
319  (*this->ins) >> ch;
320  }
321 
322  return 0;
323  }
324 
325  /// Set the step sizes
326  template<class vec2_t> int set_step(size_t nv, vec2_t &stepv) {
327  if (nv>0) {
328  step_vec.free();
329  step_vec.allocate(nv);
330  for(size_t i=0;i<nv;i++) step_vec[i]=stepv[i];
331  }
332  return 0;
333  }
334 
335 #ifndef DOXYGEN_INTERNAL
336 
337  protected:
338 
339  /// The function wrapper executed by thread with index \c ip
340  void func_wrapper(size_t ip) {
341  step(nvar,new_x[ip]);
342  (*f)(nvar,new_x[ip],new_E[ip],(*params)[ip]);
343  }
344 
345  /// Stream for verbose output
346  std::ostream *outs;
347 
348  /// Stream for verbose input
349  std::istream *ins;
350 
351  /// The number of threads to run
352  size_t nproc;
353 
354  /// The user-specified parameters
355  param_vec_t *params;
356 
357  /// The number of variables over which we minimize
358  size_t nvar;
359 
360  /// The function to minimize
361  func_t *f;
362 
363  /// \name Storage for present, next, and best vectors
364  //@{
365  alloc_vec_t x, *new_x, best_x, new_E, old_x;
366  //@}
367 
368  /// Allocation object
369  alloc_t ao;
370 
371  /// Vector of step sizes
372  ubvector step_vec;
373 
374  /** \brief Allocate memory for a minimizer over \c n dimensions
375  with stepsize \c step
376  */
377  virtual int allocate(size_t nv, size_t np) {
378  nvar=nv;
379  nproc=np;
380  ao.allocate(x,nvar);
381  new_x=new alloc_vec_t[np];
382  ao.allocate(old_x,nvar);
383  for(size_t i=0;i<np;i++) {
384  ao.allocate(new_x[i],nvar);
385  }
386  ao.allocate(new_E,np);
387  ao.allocate(best_x,nvar);
388  return 0;
389  }
390 
391  /// Free allocated memory
392  virtual int free(size_t np) {
393  ao.free(x);
394  ao.free(old_x);
395  ao.free(new_E);
396  ao.free(best_x);
397  for(size_t i=0;i<np;i++) {
398  ao.free(new_x[i]);
399  }
400  delete[] new_x;
401  return 0;
402  }
403 
404  /// Make a step to a new attempted minimum
405  virtual int step(size_t nv, vec_t &sx) {
406  size_t nstep=step_vec.size();
407  for(size_t i=0;i<nv;i++) {
408  double u=this->def_rng.random();
409  sx[i]=(2.0*u-1.0)*step_vec[i%nstep]+sx[i];
410  }
411  return 0;
412  }
413 
414 #endif
415 
416  };
417 
418 #ifndef DOXYGEN_NO_O2NS
419 }
420 #endif
421 
422 #endif
virtual int free(size_t np)
Free allocated memory.
Definition: anneal_mt.h:392
virtual const char * type()
Return string denoting type ("anneal_mt")
Definition: anneal_mt.h:283
virtual int step(size_t nv, vec_t &sx)
Make a step to a new attempted minimum.
Definition: anneal_mt.h:405
alloc_t ao
Allocation object.
Definition: anneal_mt.h:369
param_vec_t * params
The user-specified parameters.
Definition: anneal_mt.h:355
ubvector step_vec
Vector of step sizes.
Definition: anneal_mt.h:372
double min_step_ratio
Ratio between minimum step size and tolx (default 100.0)
Definition: anneal_mt.h:276
virtual int print_iter(size_t nv, vec_t &xx, double y, int iter, double tptr, std::string comment)
Print out iteration information.
Definition: anneal_mt.h:303
double tolx
The independent variable tolerance (default )
Definition: anneal_mt.h:264
virtual int next(size_t nv, vec_t &x_old, double min_old, vec_t &x_new, double min_new, double &T, size_t n_moves, bool &finished)
Determine how to change the minimization for the next iteration.
Definition: anneal_mt.h:215
func_t * f
The function to minimize.
Definition: anneal_mt.h:361
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
Definition: err_hnd.h:281
size_t nproc
The number of threads to run.
Definition: anneal_mt.h:352
int ntrial
Number of iterations.
Definition: anneal_mt.h:252
size_t nvar
The number of variables over which we minimize.
Definition: anneal_mt.h:358
std::istream * ins
Stream for verbose input.
Definition: anneal_mt.h:349
bool out_best
Output best point (default false)
Definition: anneal_mt.h:261
std::ostream * outs
Stream for verbose output.
Definition: anneal_mt.h:346
int set_verbose_stream(std::ostream &out, std::istream &in)
Set streams for verbose I/O.
Definition: anneal_mt.h:286
double T_dec
Factor to decrease temperature by (default 1.5)
Definition: anneal_mt.h:270
int verbose
Output control.
Definition: anneal_mt.h:255
Multidimensional minimization by simulated annealing (Boost multi-threaded version) ...
Definition: anneal_mt.h:74
int set_step(size_t nv, vec2_t &stepv)
Set the step sizes.
Definition: anneal_mt.h:326
bool out_step_changes
Output step size changes (default false)
Definition: anneal_mt.h:258
virtual int mmin(size_t nv, vec_t &x0, double &fmin, func_t &func, size_t np, param_vec_t &pars)
Calculate the minimum fmin of func w.r.t the array x0 of size nv using np threads.
Definition: anneal_mt.h:103
double T_start
Initial temperature (default 1.0)
Definition: anneal_mt.h:267
double step_dec
Factor to decrease step size by (default 1.5)
Definition: anneal_mt.h:273
double boltz
Boltzmann factor (default 1.0).
Definition: anneal_mt.h:249
void func_wrapper(size_t ip)
The function wrapper executed by thread with index ip.
Definition: anneal_mt.h:340
rng_t def_rng
The default random number generator.
Definition: anneal_mt.h:280
virtual int start(size_t nv, double &T)
Setup initial temperature and stepsize.
Definition: anneal_mt.h:240
virtual int allocate(size_t nv, size_t np)
Allocate memory for a minimizer over n dimensions with stepsize step.
Definition: anneal_mt.h:377

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..