Multidimensional minimization by simulated annealing (OpenMP/MPI version) More...
#include <anneal_para.h>
This class works particularly well for functions which are not trivial to evaluate, i.e. functions where the execution time is more longer than the bookkeeping that anneal_para performs between trials. For functions which satisfy this requirement, this algorithm scales nearly linearly with the number of processors.
Verbose I/O for this class happens only outside the parallel regions unless the user places I/O in the streams in the function that is specified.
Definition at line 53 of file anneal_para.h.
Public Types | |
| typedef boost::numeric::ublas::vector< double > | ubvector |
Public Types inherited from o2scl::anneal_gsl< multi_funct, boost::numeric::ublas::vector< double > > | |
| typedef boost::numeric::ublas::vector< double > | ubvector |
Public Member Functions | |
| virtual int | step (vec_t &x, vec_t &sx, int nvar, size_t ithread) |
| Make a step to a new attempted minimum. | |
Basic usage | |
| virtual int | mmin (size_t nv, vec_t &x0, double &fmin, std::vector< func_t > &func) |
| Desc. More... | |
| 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, vec_t &best_x, double best_E, int &finished) |
| Determine how to change the minimization for the next iteration. | |
| virtual int | mmin (size_t nv, vec_t &x0, double &fmin, func_t &func) |
| Desc. | |
| virtual const char * | type () |
| Return string denoting type ("anneal_para") | |
Public Member Functions inherited from o2scl::anneal_gsl< multi_funct, boost::numeric::ublas::vector< double > > | |
| anneal_gsl (const anneal_gsl< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > &ag) | |
| Copy constructor. | |
| virtual int | mmin (size_t nvar, boost::numeric::ublas::vector< double > &x0, double &fmin, multi_funct &func) |
Calculate the minimum fmin of func w.r.t the array x0 of size nvar. | |
| int | set_step (size_t nv, vec2_t &stepv) |
| Set the step sizes. | |
| anneal_gsl< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > & | operator= (const anneal_gsl< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > &ag) |
| Copy constructor from operator=. | |
Public Member Functions inherited from o2scl::anneal_base< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > | |
| anneal_base (const anneal_base< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > &ab) | |
| Copy constructor. | |
| virtual int | print_iter (size_t nv, boost::numeric::ublas::vector< double > &x, double y, int iter, double tptr, std::string comment) |
| Print out iteration information. More... | |
| anneal_base< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > & | operator= (const anneal_base< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > &ab) |
| Copy constructor from operator=. | |
Public Member Functions inherited from o2scl::mmin_base< multi_funct, multi_funct, boost::numeric::ublas::vector< double > > | |
| mmin_base (const mmin_base< multi_funct, multi_funct, boost::numeric::ublas::vector< double > > &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, boost::numeric::ublas::vector< double > &x, double &fmin, multi_funct &func, multi_funct &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< multi_funct, multi_funct, boost::numeric::ublas::vector< double > > & | operator= (const mmin_base< multi_funct, multi_funct, boost::numeric::ublas::vector< double > > &mb) |
| Copy constructor from operator=. | |
Public Attributes | |
| size_t | n_threads |
| The number of OpenMP threads. | |
| double | start_time |
| The MPI starting time. | |
| double | max_time |
| The maximum time. | |
| bool | collect_all_ranks |
| If true, obtain the global minimum over all MPI ranks (default true) | |
| std::vector< rng_gsl > | vrng |
Public Attributes inherited from o2scl::anneal_gsl< multi_funct, boost::numeric::ublas::vector< double > > | |
| 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< multi_funct, boost::numeric::ublas::vector< double >, rng_gsl > | |
| rng_gsl | rng |
| The default random number generator. | |
| o2scl::prob_dens_uniform | dist |
| The random distribution object. | |
Public Attributes inherited from o2scl::mmin_base< multi_funct, multi_funct, boost::numeric::ublas::vector< double > > | |
| 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". | |
Additional Inherited Members | |
Protected Member Functions inherited from o2scl::anneal_gsl< multi_funct, boost::numeric::ublas::vector< double > > | |
| virtual int | next (size_t nvar, boost::numeric::ublas::vector< double > &x_old, double min_old, boost::numeric::ublas::vector< double > &x_new, double min_new, double &T, size_t n_moves, boost::numeric::ublas::vector< double > &best_x, 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 | step (boost::numeric::ublas::vector< double > &sx, int nvar) |
| Make a step to a new attempted minimum. | |
Protected Attributes inherited from o2scl::anneal_gsl< multi_funct, boost::numeric::ublas::vector< double > > | |
| double | step_norm |
| Normalization for step. | |
| ubvector | step_vec |
| Vector of step sizes. | |
Protected Attributes inherited from o2scl::mmin_base< multi_funct, multi_funct, boost::numeric::ublas::vector< double > > | |
| std::ostream * | outs |
| Stream for verbose output. | |
| std::istream * | ins |
| Stream for verbose input. | |
|
inlinevirtual |
A different random number generator for each OpenMP thread
Definition at line 108 of file anneal_para.h.
Documentation generated with Doxygen. Provided under the
GNU Free Documentation License (see License Information).