Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
o2scl::anneal_gsl< func_t, vec_t, rng_t > Class Template Reference

Multidimensional minimization by simulated annealing (GSL) More...

#include <anneal_gsl.h>

Inheritance diagram for o2scl::anneal_gsl< func_t, vec_t, rng_t >:
o2scl::anneal_base< func_t, vec_t, rng_t > o2scl::mmin_base< func_t, func_t, vec_t >

Public Types

typedef boost::numeric::ublas::vector< double > ubvector
 

Public Member Functions

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.
 
virtual const char * type ()
 Return string denoting type ("anneal_gsl")
 
template<class vec2_t >
int set_step (size_t nv, vec2_t &stepv)
 Set the step sizes.
 
 anneal_gsl (const anneal_gsl< func_t, vec_t, rng_t > &ag)
 Copy constructor.
 
anneal_gsl< func_t, vec_t, rng_t > & operator= (const anneal_gsl< func_t, vec_t, rng_t > &ag)
 Copy constructor from operator=.
 
- Public Member Functions inherited from o2scl::anneal_base< func_t, vec_t, rng_t >
virtual int print_iter (size_t nv, vec_t &x, double y, int iter, double tptr, std::string comment)
 Print out iteration information. More...
 
 anneal_base (const anneal_base< func_t, vec_t, rng_t > &ab)
 Copy constructor.
 
anneal_base< func_t, vec_t, rng_t > & operator= (const anneal_base< func_t, vec_t, rng_t > &ab)
 Copy constructor from operator=.
 
- Public Member Functions inherited from o2scl::mmin_base< func_t, func_t, vec_t >
 mmin_base (const mmin_base< func_t, func_t, vec_t > &mb)
 Copy constructor.
 
int set_verbose_stream (std::ostream &out, std::istream &in)
 Set streams for verbose I/O. More...
 
virtual int mmin_de (size_t nvar, vec_t &x, double &fmin, func_t &func, func_t &dfunc)
 Calculate the minimum min of func w.r.t. the array x of size nvar with gradient dfunc.
 
int print_iter (size_t nv, vec2_t &x, double y, int iter, double value, double limit, std::string comment)
 Print out iteration information. More...
 
const char * type ()
 Return string denoting type ("mmin_base")
 
mmin_base< func_t, func_t, vec_t > & operator= (const mmin_base< func_t, func_t, vec_t > &mb)
 Copy constructor from operator=.
 

Public Attributes

double boltz
 Boltzmann factor (default 1.0).
 
double T_start
 Initial temperature (default 1.0)
 
double T_dec
 Factor to decrease temperature by (default 1.5)
 
double step_dec
 Factor to decrease step size by (default 1.5)
 
double min_step_ratio
 Ratio between minimum step size and tol_abs (default 100.0)
 
- Public Attributes inherited from o2scl::anneal_base< func_t, vec_t, rng_t >
rng_t rng
 The default random number generator.
 
o2scl::prob_dens_uniform dist
 The random distribution object.
 
- Public Attributes inherited from o2scl::mmin_base< func_t, func_t, vec_t >
int verbose
 Output control.
 
int ntrial
 Maximum number of iterations.
 
double tol_rel
 Function value tolerance.
 
double tol_abs
 The independent variable tolerance.
 
int last_ntrial
 The number of iterations for in the most recent minimization.
 
bool err_nonconv
 If true, call the error handler if the routine does not "converge".
 

Protected Member Functions

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, double best_E, bool &finished)
 Determine how to change the minimization for the next iteration.
 
virtual int start (size_t nvar, double &T)
 Setup initial temperature and stepsize.
 
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_factor.
 
virtual int step (vec_t &sx, int nvar)
 Make a step to a new attempted minimum.
 

Protected Attributes

double step_norm
 Normalization for step.
 
ubvector step_vec
 Vector of step sizes.
 
Storage for points in parameter space
ubvector x
 Current point.
 
ubvector new_x
 Proposed next point.
 
ubvector best_x
 Optimum point over all iterations.
 
ubvector old_x
 Last point from previous iteration.
 
- Protected Attributes inherited from o2scl::mmin_base< func_t, func_t, vec_t >
std::ostream * outs
 Stream for verbose output.
 
std::istream * ins
 Stream for verbose input.
 

Detailed Description

template<class func_t = multi_funct11, class vec_t = boost::numeric::ublas::vector<double>, class rng_t = rng_gsl>
class o2scl::anneal_gsl< func_t, vec_t, rng_t >

This class is a modification of simulated annealing as implemented in GSL in the function gsl_siman_solve(). It acts as a generic multidimensional minimizer for any function given a generic temperature schedule specified by the user.

There are a large variety of strategies for choosing the temperature evolution. To offer the user the largest possible flexibility, the temperature evolution is controlled by a the virtual functions start() and next() which can be freely changed by creating a child class which overwrites these functions.

The simulated annealing algorithm proposes a displacement of one coordinate of the previous point by

\[ x_{i,\mathrm{new}} = \mathrm{step\_size}_i (2 u_i - 1) + x_{i,\mathrm{old}} \]

where the $u_i$ are random numbers between 0 and 1. The displacement is accepted or rejected based on the Metropolis method. The random number generator is set in the parent, anneal.

The default behavior is as follows: Initially, the step sizes are chosen to be 1.0 (or whatever was recently specified in set_step() ) and the temperature to be T_start (default 1.0). Each iteration decreases the temperature by a factor of T_dec (default 1.5) for each step, and the minimizer is finished when the next decrease would bring the temperature below o2scl::mmin_base::tol_abs. If none of the mmin_base::ntrial steps in a particular iteration changes the value of the minimum, and the step sizes are greater than min_step_ratio (default 100) times o2scl::mmin_base::tol_abs, then the step sizes are decreased by a factor of step_dec (default 1.5) for the next iteration.

If o2scl::mmin_base::verbose is greater than zero, then mmin() will print out information and/or request a keypress after the function iterations for each temperature.

An example demonstrating the usage of this class is given in examples/ex_anneal.cpp and in the Simulated annealing example .

See also a multi-threaded version of this class in anneal_mt.

Idea for Future:
Implement a more general simulated annealing routine which would allow the solution of discrete problems like the Traveling Salesman problem.

Definition at line 147 of file anneal_gsl.h.


The documentation for this class was generated from the following file:

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).