libcmaes 0.10.2
A C++11 library for stochastic optimization with CMA-ES
Loading...
Searching...
No Matches
libcmaes::CMAParameters< TGenoPheno > Class Template Reference

Parameters for various flavors of the CMA-ES algorithm. More...

#include <libcmaes/cmaparameters.h>

Inheritance diagram for libcmaes::CMAParameters< TGenoPheno >:
libcmaes::Parameters< TGenoPheno >

Public Member Functions

 CMAParameters (const int &dim, const double *x0, const double &sigma, const int &lambda=-1, const uint64_t &seed=0, const TGenoPheno &gp=TGenoPheno())
 Constructor.
 
 CMAParameters (const std::vector< double > &x0, const double &sigma, const int &lambda=-1, const uint64_t &seed=0, const TGenoPheno &gp=TGenoPheno())
 Constructor.
 
 CMAParameters (const std::vector< double > &x0, const std::vector< double > &sigma, const int &lambda=-1, const std::vector< double > &lbounds=std::vector< double >(), const std::vector< double > &ubounds=std::vector< double >(), const uint64_t &seed=0)
 Constructor.
 
void initialize_parameters ()
 initialize required parameters based on dim, lambda, x0 and sigma.
 
void reset_as_fixed (const int &k)
 
void set_noisy ()
 adapt parameters for noisy objective function.
 
void set_algo (const int &algo)
 sets the optimization algorithm. Note: overrides Parameters::set_algo
 
void set_str_algo (const std::string &algo)
 sets the optimization algorithm.
 
double get_sigma_init () const
 returns initial sigma value
 
void set_gradient (const bool &gradient)
 activates the gradient injection scheme. If no gradient function is defined, injects a numerical gradient solution instead Note: overrides Parameters::set_gradient
 
void set_sep ()
 fix parameters for sep-CMA-ES, using only the diagonal of covariance matrix.
 
bool is_sep () const
 whether algorithm leverages separability.
 
void set_vd ()
 activates VD decomposition.
 
bool is_vd () const
 whether algorithm uses vd update.
 
void set_fixed_p (const int &index, const double &value)
 freezes a parameter to a given value in genotype during optimization. Adapts some generic parameters as well.
 
void unset_fixed_p (const int &index)
 unfreezes a parameter.
 
void set_restarts (const int &nrestarts)
 sets the maximum number of restarts (applies to IPOP and BIPOP).
 
int get_restarts () const
 get the number of restarts (applies to IPOP and BIPOP).
 
void set_lazy_update (const bool &lz)
 sets the lazy update (i.e. updates the eigenvalues every few steps).
 
bool get_lazy_update ()
 get lazy update status.
 
void set_elitism (const int &e)
 sets elitism: 0 -> no elitism 1 -> elitism: reinjects the best-ever seen solution 2 -> initial elitism: reinject x0 as long as it is not improved upon 3 -> initial elitism on restart: restart if best encountered solution is not the the final solution and reinjects the best solution until the population has better fitness, in its majority
 
void set_stopping_criteria (const int &criteria, const bool &active)
 all stopping criteria are active by default, this allows to control them
 
void set_tpa (const int &b)
 activates / deactivates two-point adaptation step-size mechanism. Overrides parameters::set_tpa by automatically setting dsigma value.
 
void set_tpa_dsigma (const double &d)
 sets dsigma value, use with care.
 
- Public Member Functions inherited from libcmaes::Parameters< TGenoPheno >
 Parameters ()
 empty constructor.
 
 Parameters (const int &dim, const double *x0, const int &lambda=-1, const uint64_t &seed=0, const TGenoPheno &gp=GenoPheno< NoBoundStrategy >())
 constructor
 
void set_x0 (const double &x0)
 sets initial objective function parameter values to x0 across all dimensions
 
void set_x0 (const double *x0)
 sets initial objective function parameter values to array x0
 
void set_x0 (const dVec &x0)
 sets initial objective function parameter values from Eigen vector
 
void set_x0 (const double &x0min, const double &x0max)
 sets bounds on initial objective function parameter values. Bounds are the same across all dimensions, and initial value is sampled uniformly within these bounds.
 
void set_x0 (const double *x0min, const double *x0max)
 sets bounds on initial objective function parameter values. Initial value is sampled uniformly within these bounds.
 
void set_x0 (const std::vector< double > &x0min, const std::vector< double > &x0max)
 sets bounds on initial objective function parameter values. Initial value is sampled uniformly within these bounds.
 
void set_x0 (const dVec &x0min, const dVec &x0max)
 sets bounds on initial objective function parameter values. Initial value is sampled uniformly within these bounds.
 
dVec get_x0min () const
 returns lower bound on x0 vector
 
dVec get_x0max () const
 returns upper bound on x0 vector
 
void set_fixed_p (const int &index, const double &value)
 freezes a parameter to a given value during optimization.
 
void unset_fixed_p (const int &index)
 unfreezes a parameter.
 
void set_max_iter (const int &maxiter)
 sets the maximum number of iterations allowed for the optimization.
 
int get_max_iter () const
 returns maximum number of iterations
 
void set_max_fevals (const int &fevals)
 sets the maximum budget of objective function calls allowed for the optimization.
 
int get_max_fevals () const
 returns maximum budget of objective function calls
 
void set_ftarget (const double &val)
 sets the objective function target value when known.
 
void reset_ftarget ()
 resets the objective function target value to its inactive state.
 
double get_ftarget () const
 returns objective function target value.
 
void set_seed (const int &seed)
 sets random generator's seed, 0 is special value to generate random seed.
 
int get_seed () const
 returns random generator's seed.
 
void set_ftolerance (const double &v)
 sets function tolerance as stopping criteria for TolHistFun: monitors the difference in function value over iterations and stops optimization when below tolerance.
 
double get_ftolerance () const
 returns function tolerance
 
void set_xtolerance (const double &v)
 sets parameter tolerance as stopping criteria for TolX.
 
double get_xtolerance () const
 returns parameter tolerance
 
int lambda () const
 returns lambda, number of offsprings per generation
 
int dim () const
 returns the problem's dimension
 
void set_quiet (const bool &quiet)
 sets the quiet mode (no output from the library) for the optimization at hand
 
bool quiet () const
 returns whether the quiet mode is on.
 
void set_algo (const int &algo)
 sets the optimization algorithm.
 
int get_algo () const
 returns which algorithm is set for the optimization at hand.
 
void set_gp (const TGenoPheno &gp)
 sets the genotype/phenotype transform object.
 
TGenoPheno get_gp () const
 returns the current genotype/phenotype transform object.
 
void set_fplot (const std::string &fplot)
 sets the output filename (activates the output to file).
 
void set_full_fplot (const bool &b)
 activates / deactivates the full output (for legacy plotting).
 
std::string get_fplot () const
 returns the current output filename.
 
void set_gradient (const bool &gradient)
 activates the gradient injection scheme. If no gradient function is defined, injects a numerical gradient solution instead
 
bool get_gradient () const
 returns whether the gradient injection scheme is activated.
 
void set_edm (const bool &edm)
 activates computation of expected distance to minimum when optimization has completed
 
bool get_edm () const
 returns whether edm is activated.
 
void set_mt_feval (const bool &mt)
 activate / deactivate the parallel evaluation of objective function
 
bool get_mt_feval () const
 returns whether the parallel evaluation of objective function is activated
 
void set_max_hist (const int &m)
 sets maximum history size, allows to keep memory requirements fixed.
 
void set_maximize (const bool &maximize)
 active internal maximization scheme (simply returns -f instead of f)
 
bool get_maximize () const
 returns whether the maximization mode is enabled
 
void set_initial_fvalue (const bool &b)
 whether to compute initial objective function value (i.e. at x0)
 
void set_uh (const bool &b)
 activates / deactivates uncertainty handling scheme.
 
bool get_uh () const
 get uncertainty handling status.
 
void set_tpa (const int &b)
 activates / deactivates two-point adaptation step-size mechanism
 
int get_tpa () const
 get two-point adapation step-size mechanism status.
 

Private Attributes

int _mu
 
dVec _weights
 
double _csigma
 
double _c1
 
double _cmu
 
double _cc
 
double _muw
 
double _dsigma
 
double _fact_ps
 
double _fact_pc
 
double _chi
 
double _sigma_init
 
int _nrestarts = 9
 
bool _lazy_update
 
double _lazy_value
 
double _cm
 
double _alphacov
 
double _alphaminusold
 
double _deltamaxsigma
 
double _lambdamintarget
 
double _alphaminusmin
 
bool _sep = false
 
bool _vd = false
 
bool _elitist = false
 
bool _initial_elitist = false
 
bool _initial_elitist_on_restart = false
 
std::map< int, bool_stoppingcrit
 

Friends

class CMASolutions
 
template<class U , class V >
class CMAStrategy
 
template<class U , class V , class W >
class ESOStrategy
 
template<class U >
class CMAStopCriteria
 
template<class U , class V >
class IPOPCMAStrategy
 
template<class U , class V >
class BIPOPCMAStrategy
 
class CovarianceUpdate
 
class ACovarianceUpdate
 
template<class U >
class errstats
 
class VDCMAUpdate
 

Additional Inherited Members

- Protected Attributes inherited from libcmaes::Parameters< TGenoPheno >
int _dim
 
int _lambda = -1
 
int _max_iter = -1
 
int _max_fevals = -1
 
bool _quiet = true
 
std::string _fplot = ""
 
bool _full_fplot = false
 
dVec _x0min
 
dVec _x0max
 
double _ftarget = -std::numeric_limits<double>::infinity()
 
double _ftolerance = 1e-12
 
double _xtol = 1e-12
 
uint64_t _seed = 0
 
int _algo = 0
 
bool _with_gradient =false
 
bool _with_edm =false
 
std::unordered_map< int, double_fixed_p
 
TGenoPheno _gp
 
bool _mt_feval = false
 
int _max_hist = -1
 
bool _maximize = false
 
bool _initial_fvalue = false
 
bool _uh = false
 
double _rlambda
 
double _epsuh = 1e-7
 
double _thetauh = 0.2
 
double _csuh = 1.0
 
double _alphathuh = 1.0
 
int _tpa = 1
 
double _tpa_csigma = 0.3
 
- Static Protected Attributes inherited from libcmaes::Parameters< TGenoPheno >
static std::map< std::string, int_algos = {{"cmaes",0},{"ipop",1},{"bipop",2},{"acmaes",3},{"aipop",4},{"abipop",5},{"sepcmaes",6},{"sepipop",7},{"sepbipop",8},{"sepacmaes",9},{"sepipop",10},{"sepbipop",11},{"vdcma",12},{"vdipopcma",13},{"vdbipopcma",14}}
 

Detailed Description

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
class libcmaes::CMAParameters< TGenoPheno >

Parameters for various flavors of the CMA-ES algorithm.

Constructor & Destructor Documentation

◆ CMAParameters() [1/3]

template<class TGenoPheno >
libcmaes::CMAParameters< TGenoPheno >::CMAParameters ( const int dim,
const double x0,
const double sigma,
const int lambda = -1,
const uint64_t seed = 0,
const TGenoPheno gp = TGenoPheno() 
)

Constructor.

Parameters
dimproblem dimensions
x0initial search point
sigmainitial distribution step size (positive, otherwise automatically set)
lambdanumber of offsprings sampled at each step
seedinitial random seed, useful for reproducing results (if unspecified, automatically generated from current time)
gpgenotype / phenotype object
sepwhether to use sep-CMA-ES, using diagonal covariance matrix (modifies covariance default learning rate)

◆ CMAParameters() [2/3]

template<class TGenoPheno >
libcmaes::CMAParameters< TGenoPheno >::CMAParameters ( const std::vector< double > &  x0,
const double sigma,
const int lambda = -1,
const uint64_t seed = 0,
const TGenoPheno gp = TGenoPheno() 
)

Constructor.

Parameters
x0initial search point as vector of problem dimension
sigmainitial distribution step size (positive, otherwise automatically set)
lambdanumber of offsprings sampled at each step
seedinitial random seed, useful for reproducing results (if unspecified, automatically generated from current time)
gpgenotype / phenotype object
sepwhether to use sep-CMA-ES, using diagonal covariance matrix (modifies covariance default learning rate)

◆ CMAParameters() [3/3]

template<class TGenoPheno >
libcmaes::CMAParameters< TGenoPheno >::CMAParameters ( const std::vector< double > &  x0,
const std::vector< double > &  sigma,
const int lambda = -1,
const std::vector< double > &  lbounds = std::vector<double>(),
const std::vector< double > &  ubounds = std::vector<double>(),
const uint64_t seed = 0 
)

Constructor.

Parameters
x0initial search point as vector of problem dimension
sigmavector of initial distribution step sizes (positive, otherwise automatically set)
lambdanumber of offsprings sampled at each step
seedinitial random seed, useful for reproducing results (if unspecified, automatically generated from current time)
gpgenotype / phenotype object
sepwhether to use sep-CMA-ES, using diagonal covariance matrix (modifies covariance default learning rate)

Member Function Documentation

◆ get_lazy_update()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::get_lazy_update ( )
inline

get lazy update status.

Parameters
whetherlazy update is activated

◆ get_restarts()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
int libcmaes::CMAParameters< TGenoPheno >::get_restarts ( ) const
inline

get the number of restarts (applies to IPOP and BIPOP).

Returns
number of restarts

◆ get_sigma_init()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::get_sigma_init ( ) const
inline

returns initial sigma value

Returns
initial sigma value

◆ is_sep()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::is_sep ( ) const
inline

whether algorithm leverages separability.

Returns
separability status

◆ is_vd()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::is_vd ( ) const
inline

whether algorithm uses vd update.

Returns
vd status

◆ set_algo()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
void libcmaes::CMAParameters< TGenoPheno >::set_algo ( const int algo)
inline

sets the optimization algorithm. Note: overrides Parameters::set_algo

Parameters
algofrom CMAES_DEFAULT, IPOP_CMAES, BIPOP_CMAES, aCMAES, aIPOP_CMAES, aBIPOP_CMAES, sepCMAES, sepIPOP_CMAES, sepBIPOP_CMAES, sepaCMAES, sepaIPOP_CMAES, sepaBIPOP_CMAES, VD_CMAES, VD_IPOP_CMAES, VD_BIPOP_CMAES

◆ set_fixed_p()

template<class TGenoPheno >
void libcmaes::CMAParameters< TGenoPheno >::set_fixed_p ( const int index,
const double value 
)

freezes a parameter to a given value in genotype during optimization. Adapts some generic parameters as well.

Parameters
indexdimension index of the parameter to be frozen
valuefrozen value of the parameter

◆ set_gradient()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
void libcmaes::CMAParameters< TGenoPheno >::set_gradient ( const bool gradient)
inline

activates the gradient injection scheme. If no gradient function is defined, injects a numerical gradient solution instead Note: overrides Parameters::set_gradient

Parameters
gradienttrue/false

◆ set_lazy_update()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
void libcmaes::CMAParameters< TGenoPheno >::set_lazy_update ( const bool lz)
inline

sets the lazy update (i.e. updates the eigenvalues every few steps).

Parameters
lzwhether to activate the lazy update

◆ set_restarts()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
void libcmaes::CMAParameters< TGenoPheno >::set_restarts ( const int nrestarts)
inline

sets the maximum number of restarts (applies to IPOP and BIPOP).

Parameters
nrestartsmaximum number of restarts

◆ set_stopping_criteria()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
void libcmaes::CMAParameters< TGenoPheno >::set_stopping_criteria ( const int criteria,
const bool active 
)
inline

all stopping criteria are active by default, this allows to control them

Parameters
criteriastopping criteria CMAStopCritType, see cmastopcriteria.h
activewhether to activate this criteria

◆ set_str_algo()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
void libcmaes::CMAParameters< TGenoPheno >::set_str_algo ( const std::string &  algo)
inline

sets the optimization algorithm.

Parameters
algoas string from cmaes,ipop,bipop,acmaes,aipop,abipop,sepcmaes,sepipop,sepbipop,sepacmaes,sepaipop,sepabipop,vdcma,vdipopcma,vdbipopcma

◆ set_tpa()

template<class TGenoPheno >
void libcmaes::CMAParameters< TGenoPheno >::set_tpa ( const int b)

activates / deactivates two-point adaptation step-size mechanism. Overrides parameters::set_tpa by automatically setting dsigma value.

Parameters
b0: no, 1: auto, 2: yes

◆ set_tpa_dsigma()

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
void libcmaes::CMAParameters< TGenoPheno >::set_tpa_dsigma ( const double d)
inline

sets dsigma value, use with care.

Parameters
ddsigma

◆ unset_fixed_p()

template<class TGenoPheno >
void libcmaes::CMAParameters< TGenoPheno >::unset_fixed_p ( const int index)

unfreezes a parameter.

Parameters
indexdimenion index of the parameter to unfreeze

Member Data Documentation

◆ _alphacov

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_alphacov
private

= 2 (active CMA only)

◆ _alphaminusmin

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_alphaminusmin
private

= 1 (active CMA only)

◆ _alphaminusold

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_alphaminusold
private

in [0,1] (active CMA only)

◆ _c1

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_c1
private

covariance matrix learning rate for the rank one update using pc.

◆ _cc

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_cc
private

cumulation constant for pc.

◆ _chi

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_chi
private

norm of N(0,I)

◆ _cm

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_cm
private

learning rate for the mean.

◆ _cmu

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_cmu
private

covariance matrix learning reate for the rank mu update.

◆ _csigma

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_csigma
private

cumulation constant for step size.

◆ _deltamaxsigma

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_deltamaxsigma
private

infinite (active CMA only)

◆ _dsigma

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_dsigma
private

step size damping factor.

◆ _elitist

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::_elitist = false
private

re-inject the best-ever seen solution.

◆ _initial_elitist

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::_initial_elitist = false
private

re-inject x0.

◆ _initial_elitist_on_restart

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::_initial_elitist_on_restart = false
private

activate the restart from and re-injection of the best seen solution if not the final one.

◆ _lambdamintarget

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_lambdamintarget
private

= 0.66 (active CMA only)

◆ _lazy_update

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::_lazy_update
private

covariance lazy update.

◆ _lazy_value

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_lazy_value
private

reference trigger for lazy update.

◆ _mu

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
int libcmaes::CMAParameters< TGenoPheno >::_mu
private

number of candidate solutions used to update the distribution parameters.

◆ _muw

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_muw
private

\sum^\mu _weights .

◆ _nrestarts

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
int libcmaes::CMAParameters< TGenoPheno >::_nrestarts = 9
private

maximum number of restart, when applicable.

◆ _sep

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
bool libcmaes::CMAParameters< TGenoPheno >::_sep = false
private

whether to use diagonal covariance matrix.

◆ _sigma_init

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
double libcmaes::CMAParameters< TGenoPheno >::_sigma_init
private

initial sigma value.

◆ _stoppingcrit

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
std::map<int,bool> libcmaes::CMAParameters< TGenoPheno >::_stoppingcrit
private

control list of stopping criteria.

◆ _weights

template<class TGenoPheno = GenoPheno<NoBoundStrategy>>
dVec libcmaes::CMAParameters< TGenoPheno >::_weights
private

offsprings weighting scheme.


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