All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
anneal_gsl.h
Go to the documentation of this file.
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 /* siman/siman.c
24  *
25  * Copyright (C) 2007 Brian Gough
26  * Copyright (C) 1996, 1997, 1998, 1999, 2000 Mark Galassi
27  *
28  * This program is free software; you can redistribute it and/or modify
29  * it under the terms of the GNU General Public License as published by
30  * the Free Software Foundation; either version 3 of the License, or (at
31  * your option) any later version.
32  *
33  * This program is distributed in the hope that it will be useful, but
34  * WITHOUT ANY WARRANTY; without even the implied warranty of
35  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
36  * General Public License for more details.
37  *
38  * You should have received a copy of the GNU General Public License
39  * along with this program; if not, write to the Free Software
40  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
41  * 02110-1301, USA.
42  */
43 #ifndef O2SCL_ANNEAL_GSL_H
44 #define O2SCL_ANNEAL_GSL_H
45 
46 /** \file anneal_gsl.h
47  \brief File defining \ref o2scl::anneal_gsl
48 */
49 #include <random>
50 
51 #include <boost/numeric/ublas/vector.hpp>
52 #include <boost/numeric/ublas/matrix.hpp>
53 
54 #include <o2scl/anneal.h>
55 #include <o2scl/multi_funct.h>
56 
57 #ifndef DOXYGEN_NO_O2NS
58 namespace o2scl {
59 #endif
60 
61  /** \brief Multidimensional minimization by simulated annealing (GSL)
62 
63  This class is a modification of simulated annealing as
64  implemented in GSL in the function \c gsl_siman_solve(). It acts
65  as a generic multidimensional minimizer for any function given a
66  generic temperature schedule specified by the user.
67 
68  There are a large variety of strategies for choosing the
69  temperature evolution. To offer the user the largest possible
70  flexibility, the temperature evolution is controlled by a
71  the virtual functions start() and next() which can be freely
72  changed by creating a child class which overwrites these
73  functions.
74 
75  The simulated annealing algorithm proposes a displacement of one
76  coordinate of the previous point by
77  \f[
78  x_{i,\mathrm{new}} = \mathrm{step\_size}_i
79  (2 u_i - 1) + x_{i,\mathrm{old}}
80  \f]
81  where the \f$u_i\f$ are random numbers between 0 and 1. The
82  displacement is accepted or rejected based on the Metropolis
83  method. The random number generator is set in the parent,
84  anneal.
85 
86  The default behavior is as follows: Initially, the step sizes
87  are chosen to be 1.0 (or whatever was recently specified in \ref
88  set_step() ) and the temperature to be \ref T_start (default
89  1.0). Each iteration decreases the temperature by a factor of
90  \ref T_dec (default 1.5) for each step, and the minimizer is
91  finished when the next decrease would bring the temperature
92  below multi_min::tol_abs. If none of the multi_min::ntrial steps in
93  a particular iteration changes the value of the minimum, and the
94  step sizes are greater than \ref min_step_ratio (default 100)
95  times multi_min::tol_abs, then the step sizes are decreased by a
96  factor of \ref step_dec (default 1.5) for the next iteration.
97 
98  If \ref mmin_base::verbose is greater than zero, then \ref
99  mmin() will print out information and/or request a keypress
100  after the function iterations for each temperature.
101 
102  An example demonstrating the usage of this class is given in
103  <tt>examples/ex_anneal.cpp</tt> and in the \ref ex_anneal_sect .
104 
105  The form of the user-specified function is as in \ref
106  multi_funct11 has a "function value" which is the value of the
107  function (given in the third argument as a number of type \c
108  double), and a "return value" (the integer return value). The
109  initial function evaluation which is performed at the
110  user-specified initial guess must give 0 as the return value. If
111  subsequent function evaluations have a non-zero return value,
112  then the resulting point is ignored and a new point is selected.
113 
114  This class thus can sometimes handle constrained minimization
115  problems. If the user ensures that the function's return value
116  is non-zero when the function is evaluated outside the allowed
117  region, the minimizer will not accept any steps which take the
118  minimizer outside the allowed region. Note that this should be
119  done with care, however, as this approach may cause convergence
120  problems with sufficiently difficult functions or constraints.
121 
122  See also a multi-threaded version of this class in \ref
123  anneal_mt.
124 
125  \future There's x0, old_x, new_x, best_x, and x? There's probably
126  some duplication here which could be avoided.
127 
128  \future
129  - Implement a more general simulated annealing routine
130  which would allow the solution of discrete problems like the
131  Traveling Salesman problem.
132  \comment
133  This could probably be done just by creating a parent abstract
134  base class which doesn't have any reference to step() and
135  step_vec.
136  \endcomment
137 
138  */
139  template<class func_t=multi_funct11,
141  class rng_t=int, class rng_dist_t=rng_gsl>
142  class anneal_gsl : public anneal_base<func_t,vec_t,rng_t,rng_dist_t> {
143 
144  public:
145 
147 
148  anneal_gsl() {
149  boltz=1.0;
150  step_vec.resize(1);
151  step_vec[0]=1.0;
152  T_start=1.0;
153  T_dec=1.5;
154  step_dec=1.5;
155  min_step_ratio=100.0;
156  }
157 
158  virtual ~anneal_gsl() {
159  }
160 
161  /** \brief Calculate the minimum \c fmin of \c func w.r.t the
162  array \c x0 of size \c nvar.
163  */
164  virtual int mmin(size_t nvar, vec_t &x0, double &fmin,
165  func_t &func) {
166 
167  if (nvar==0) {
168  O2SCL_ERR2("Tried to minimize over zero variables ",
169  " in anneal_gsl::mmin().",exc_einval);
170  }
171 
172  fmin=0.0;
173 
174  allocate(nvar);
175 
176  double E, new_E, best_E, T, old_E;
177  int i, iter=0;
178  size_t j;
179 
180  for(j=0;j<nvar;j++) {
181  x[j]=x0[j];
182  best_x[j]=x0[j];
183  }
184 
185  E=func(nvar,x);
186  best_E=E;
187 
188  // Setup initial temperature and step sizes
189  start(nvar,T);
190 
191  bool done=false;
192 
193  while (!done) {
194 
195  // Copy old value of x for next() function
196  for(j=0;j<nvar;j++) old_x[j]=x[j];
197  old_E=E;
198 
199  size_t nmoves=0;
200 
201  for (i=0;i<this->ntrial;++i) {
202  for (j=0;j<nvar;j++) new_x[j]=x[j];
203 
204  step(new_x,nvar);
205  new_E=func(nvar,new_x);
206 
207  // Store best value obtained so far
208  if(new_E<=best_E){
209  for(j=0;j<nvar;j++) best_x[j]=new_x[j];
210  best_E=new_E;
211  }
212 
213  // Take the crucial step: see if the new point is accepted
214  // or not, as determined by the Boltzmann probability
215  if (new_E<E) {
216  for(j=0;j<nvar;j++) x[j]=new_x[j];
217  E=new_E;
218  nmoves++;
219  } else {
220  double r=this->rng_dist(this->rng);
221  if (r < exp(-(new_E-E)/(boltz*T))) {
222  for(j=0;j<nvar;j++) x[j]=new_x[j];
223  E=new_E;
224  nmoves++;
225  }
226  }
227 
228  }
229 
230  if (this->verbose>0) {
231  this->print_iter(nvar,best_x,best_E,iter,T,"anneal_gsl");
232  iter++;
233  }
234 
235  // See if we're finished and proceed to the next step
236  next(nvar,old_x,old_E,x,E,T,nmoves,done);
237 
238  }
239 
240  for(j=0;j<nvar;j++) x0[j]=best_x[j];
241  fmin=best_E;
242 
243  return 0;
244  }
245 
246  /// Return string denoting type ("anneal_gsl")
247  virtual const char *type() { return "anneal_gsl"; }
248 
249  /** \brief Boltzmann factor (default 1.0).
250  */
251  double boltz;
252 
253  /// Set the step sizes
254  template<class vec2_t> int set_step(size_t nv, vec2_t &stepv) {
255  if (nv>0) {
256  step_vec.resize(nv);
257  for(size_t i=0;i<nv;i++) step_vec[i]=stepv[i];
258  }
259  return 0;
260  }
261 
262  /// Initial temperature (default 1.0)
263  double T_start;
264 
265  /// Factor to decrease temperature by (default 1.5)
266  double T_dec;
267 
268  /// Factor to decrease step size by (default 1.5)
269  double step_dec;
270 
271  /// Ratio between minimum step size and \ref tol_abs (default 100.0)
273 
274  /** \brief Copy constructor
275  */
279 
280  boltz=ag.boltz;
281  T_start=ag.T_start;
282  T_dec=ag.T_dec;
283  step_dec=ag.step_dec;
284  min_step_ratio=ag.min_step_ratio;
285  step_vec=ag.step_vec;
286  x=ag.x;
287  new_x=ag.new_x;
288  best_x=ag.best_x;
289  old_x=ag.old_x;
290 
291  }
292 
293  /** \brief Copy constructor from operator=
294  */
297  if (this != &ag) {
299  boltz=ag.boltz;
300  T_start=ag.T_start;
301  T_dec=ag.T_dec;
302  step_dec=ag.step_dec;
303  min_step_ratio=ag.min_step_ratio;
304  step_vec=ag.step_vec;
305  x=ag.x;
306  new_x=ag.new_x;
307  best_x=ag.best_x;
308  old_x=ag.old_x;
309  }
310  return *this;
311  }
312 
313 #ifndef DOXYGEN_INTERNAL
314 
315  protected:
316 
317  /// \name Storage for present, next, and best vectors
318  //@{
319  ubvector x, new_x, best_x, old_x;
320  //@}
321 
322  /// Vector of step sizes
324 
325  /// Determine how to change the minimization for the next iteration
326  virtual int next(size_t nvar, vec_t &x_old, double min_old,
327  vec_t &x_new, double min_new, double &T,
328  size_t n_moves, bool &finished) {
329 
330  if (T/T_dec<this->tol_abs) {
331  finished=true;
332  return success;
333  }
334  if (n_moves==0) {
335  for(size_t i=0;i<nvar;i++) {
336  if (i<step_vec.size() && step_vec[i]>this->tol_abs*min_step_ratio) {
337  step_vec[i]/=step_dec;
338  }
339  }
340  }
341  T/=T_dec;
342  return success;
343  }
344 
345  /// Setup initial temperature and stepsize
346  virtual int start(size_t nvar, double &T) {
347  T=T_start;
348  return success;
349  }
350 
351  /** \brief Allocate memory for a minimizer over \c n dimensions
352  with stepsize \c step and Boltzmann factor \c boltz_factor
353  */
354  virtual int allocate(size_t n, double boltz_factor=1.0) {
355  x.resize(n);
356  old_x.resize(n);
357  new_x.resize(n);
358  best_x.resize(n);
359  boltz=boltz_factor;
360  return 0;
361  }
362 
363  /** \brief Make a step to a new attempted minimum
364  */
365  virtual int step(vec_t &sx, int nvar) {
366  size_t nstep=step_vec.size();
367  for(int i=0;i<nvar;i++) {
368  double u=this->rng_dist(this->rng);
369  sx[i]=(2.0*u-1.0)*step_vec[i%nstep]+sx[i];
370  }
371  return 0;
372  }
373 
374 #endif
375 
376  };
377 
378 #ifndef DOXYGEN_NO_O2NS
379 }
380 #endif
381 
382 #endif
virtual int step(vec_t &sx, int nvar)
Make a step to a new attempted minimum.
Definition: anneal_gsl.h:365
double step_dec
Factor to decrease step size by (default 1.5)
Definition: anneal_gsl.h:269
invalid argument supplied by user
Definition: err_hnd.h:59
virtual int allocate(size_t n, double boltz_factor=1.0)
Allocate memory for a minimizer over n dimensions with stepsize step and Boltzmann factor boltz_facto...
Definition: anneal_gsl.h:354
virtual int next(size_t nvar, 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_gsl.h:326
virtual const char * type()
Return string denoting type ("anneal_gsl")
Definition: anneal_gsl.h:247
ubvector step_vec
Vector of step sizes.
Definition: anneal_gsl.h:323
virtual int start(size_t nvar, double &T)
Setup initial temperature and stepsize.
Definition: anneal_gsl.h:346
rng_dist_t rng_dist
The default random number distribution.
Definition: anneal.h:116
Simulated annealing base.
Definition: anneal.h:64
double T_dec
Factor to decrease temperature by (default 1.5)
Definition: anneal_gsl.h:266
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
Definition: err_hnd.h:281
rng_t rng
The default random number generator.
Definition: anneal.h:119
int set_step(size_t nv, vec2_t &stepv)
Set the step sizes.
Definition: anneal_gsl.h:254
double min_step_ratio
Ratio between minimum step size and tol_abs (default 100.0)
Definition: anneal_gsl.h:272
double T_start
Initial temperature (default 1.0)
Definition: anneal_gsl.h:263
double tol_abs
The independent variable tolerance.
Definition: mmin.h:202
Multidimensional minimization by simulated annealing (GSL)
Definition: anneal_gsl.h:142
anneal_base< func_t, vec_t, rng_t, rng_dist_t > & operator=(const anneal_base< func_t, vec_t, rng_t, rng_dist_t > &ab)
Copy constructor from operator=.
Definition: anneal.h:138
double boltz
Boltzmann factor (default 1.0).
Definition: anneal_gsl.h:251
std::function< double(size_t, const boost::numeric::ublas::vector< double > &)> multi_funct11
Multi-dimensional function typedef.
Definition: multi_funct.h:45
virtual int mmin(size_t nvar, vec_t &x0, double &fmin, func_t &func)
Calculate the minimum fmin of func w.r.t the array x0 of size nvar.
Definition: anneal_gsl.h:164
Success.
Definition: err_hnd.h:47
virtual int print_iter(size_t nv, vec_t &x, double y, int iter, double tptr, std::string comment)
Print out iteration information.
Definition: anneal.h:94
int ntrial
Maximum number of iterations.
Definition: mmin.h:196

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..