cma (version 1.1.06 $Revision: 4129 $ $Date: 2015-01-23 20:13:51 +0100 (Fri, 23 Jan 2015) $)
index
cma.py

Module cma implements the CMA-ES (Covariance Matrix Adaptation
Evolution Strategy).
 
CMA-ES is a stochastic optimizer for robust non-linear non-convex
derivative- and function-value-free numerical optimization.
 
This implementation can be used with Python versions >= 2.6, namely
2.6, 2.7, 3.3, 3.4.
 
CMA-ES searches for a minimizer (a solution x in :math:`R^n`) of an
objective function f (cost function), such that f(x) is minimal.
Regarding f, only a passably reliable ranking of the candidate
solutions in each iteration is necessary. Neither the function values
itself, nor the gradient of f need to be available or do matter (like
in the downhill simplex Nelder-Mead algorithm). Some termination
criteria however depend on actual f-values.
 
Two interfaces are provided:
 
  - function `fmin(func, x0, sigma0,...)`
        runs a complete minimization
        of the objective function func with CMA-ES.
 
  - class `CMAEvolutionStrategy`
      allows for minimization such that the control of the iteration
      loop remains with the user.
 
 
Used packages:
 
    - unavoidable: `numpy` (see `barecmaes2.py` if `numpy` is not
      available),
    - avoidable with small changes: `time`, `sys`
    - optional: `matplotlib.pyplot` (for `plot` etc., highly
      recommended), `pprint` (pretty print), `pickle` (in class
      `Sections`), `doctest`, `inspect`, `pygsl` (never by default)
 
Install
-------
The file ``cma.py`` only needs to be visible in the python path (e.g. in
the current working directory).
 
The preferred way of (system-wide) installation is calling
 
    pip install cma
 
from the command line.
 
The ``cma.py`` file can also be installed from the
system shell terminal command line by::
 
    python cma.py --install
 
which solely calls the ``setup`` function from the standard
``distutils.core`` package for installation. If the ``setup.py``
file is been provided with ``cma.py``, the standard call is
 
    python setup.py cma
 
Both calls need to see ``cma.py`` in the current working directory and
might need to be preceded with ``sudo``.
 
To upgrade the currently installed version from the Python Package Index,
and also for first time installation, type in the system shell::
 
    pip install --upgrade cma
 
Testing
-------
From the system shell::
 
    python cma.py --test
 
or from the Python shell ``ipython``::
 
    run cma.py --test
 
or from any python shell
 
    import cma
    cma.main('--test')
 
runs ``doctest.testmod(cma)`` showing only exceptions (and not the
tests that fail due to small differences in the output) and should
run without complaints in about between 20 and 100 seconds.
 
Example
-------
From a python shell::
 
    import cma
    help(cma)  # "this" help message, use cma? in ipython
    help(cma.fmin)
    help(cma.CMAEvolutionStrategy)
    help(cma.CMAOptions)
    cma.CMAOptions('tol')  # display 'tolerance' termination options
    cma.CMAOptions('verb') # display verbosity options
    res = cma.fmin(cma.Fcts.tablet, 15 * [1], 1)
    res[0]  # best evaluated solution
    res[5]  # mean solution, presumably better with noise
 
:See: `fmin()`, `CMAOptions`, `CMAEvolutionStrategy`
 
:Author: Nikolaus Hansen, 2008-2015
:Contributor: Petr Baudis, 2014
 
:License: BSD 3-Clause, see below.

 
Modules
       
collections
numpy
matplotlib.pyplot
sys
time

 
Classes
       
__builtin__.dict(__builtin__.object)
CMAOptions
__builtin__.object
BaseDataLogger
CMADataLogger
BestSolution
CMAAdaptSigmaBase
CMAAdaptSigmaCSA
CMAAdaptSigmaDistanceProportional
CMAAdaptSigmaMedianImprovement
CMAAdaptSigmaNone
CMAAdaptSigmaTPA
ElapsedTime
GenoPheno
MathHelperFunctions
Misc
NoiseHandler
OOOptimizer
CMAEvolutionStrategy
Rotation
Sections
_abcoll.MutableMapping(_abcoll.Mapping)
DerivedDictBase
SolutionDict
CMASolutionDict
BoundaryHandlerBase(__builtin__.object)
BoundNone
BoundPenalty
BoundTransform

 
class BaseDataLogger(__builtin__.object)
    "abstract" base class for a data logger that can be used with an `OOOptimizer`
 
Details: attribute `modulo` is used in ``OOOptimizer.optimize``
 
  Methods defined here:
add(self, optim=None, more_data=[])
abstract method, add a "data point" from the state of `optim` into the
logger, the argument `optim` can be omitted if it was `register()`-ed before,
acts like an event handler
data(self)
return logged data in a dictionary (not implemented)
disp(self)
display some data trace (not implemented)
plot(self)
plot data (not implemented)
register(self, optim)
abstract method, register an optimizer `optim`, only needed if `add()` is
called without a value for the `optim` argument

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class BestSolution(__builtin__.object)
    container to keep track of the best solution seen
 
  Methods defined here:
__init__(self, x=None, f=inf, evals=None)
initialize the best solution with `x`, `f`, and `evals`.
Better solutions have smaller `f`-values.
get(self)
return ``(x, f, evals)``
update(self, arx, xarchive=None, arf=None, evals=None)
checks for better solutions in list `arx`.
 
Based on the smallest corresponding value in `arf`,
alternatively, `update` may be called with a `BestSolution`
instance like ``update(another_best_solution)`` in which case
the better solution becomes the current best.
 
`xarchive` is used to retrieve the genotype of a solution.

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class BoundNone(BoundaryHandlerBase)
    # ____________________________________________________________
# ____________________________________________________________
 
 
Method resolution order:
BoundNone
BoundaryHandlerBase
__builtin__.object

Methods defined here:
__init__(self, bounds=None)
is_in_bounds(self, x)

Methods inherited from BoundaryHandlerBase:
__call__(self, solutions, *args, **kwargs)
return penalty or list of penalties, by default zero(s).
 
This interface seems too specifically tailored to the derived
BoundPenalty class, it should maybe change.
get_bounds(self, which, dimension)
``get_bounds('lower', 8)`` returns the lower bounds in 8-D
has_bounds(self)
return True, if any variable is bounded
inverse(self, y, copy_if_changed=True, copy_always=False)
repair(self, x, copy_if_changed=True, copy_always=False)
projects infeasible values on the domain bound, might be
overwritten by derived class
to_dim_times_two(self, bounds)
return boundaries in format ``[[lb0, ub0], [lb1, ub1], ...]``,
as used by ``BoxConstraints...`` class.
update(self, *args, **kwargs)

Data descriptors inherited from BoundaryHandlerBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class BoundPenalty(BoundaryHandlerBase)
    Computes the boundary penalty. Must be updated each iteration,
using the `update` method.
 
Details
-------
The penalty computes like ``sum(w[i] * (x[i]-xfeas[i])**2)``,
where `xfeas` is the closest feasible (in-bounds) solution from `x`.
The weight `w[i]` should be updated during each iteration using
the update method.
 
Example:
 
>>> import cma
>>> cma.fmin(cma.felli, 6 * [1], 1,
...          {
...              'boundary_handling': 'BoundPenalty',
...              'bounds': [-1, 1],
...              'fixed_variables': {0: 0.012, 2:0.234}
...          })
 
Reference: Hansen et al 2009, A Method for Handling Uncertainty...
IEEE TEC, with addendum, see
http://www.lri.fr/~hansen/TEC2009online.pdf
 
 
Method resolution order:
BoundPenalty
BoundaryHandlerBase
__builtin__.object

Methods defined here:
__call__(self, x, archive, gp)
returns the boundary violation penalty for `x` ,where `x` is a
single solution or a list or array of solutions.
__init__(self, bounds=None)
Argument bounds can be `None` or ``bounds[0]`` and ``bounds[1]``
are lower  and upper domain boundaries, each is either `None` or
a scalar or a list or array of appropriate size.
feasible_ratio(self, solutions)
counts for each coordinate the number of feasible values in
``solutions`` and returns an array of length ``len(solutions[0])``
with the ratios.
 
`solutions` is a list or array of repaired ``Solution``
instances,
repair(self, x, copy_if_changed=True, copy_always=False)
sets out-of-bounds components of ``x`` on the bounds.
update(self, function_values, es)
updates the weights for computing a boundary penalty.
 
Arguments
---------
`function_values`
    all function values of recent population of solutions
`es`
    `CMAEvolutionStrategyobject instance, in particular
    mean and variances and the methods from the attribute
    `gp` of type `GenoPheno` are used.

Methods inherited from BoundaryHandlerBase:
get_bounds(self, which, dimension)
``get_bounds('lower', 8)`` returns the lower bounds in 8-D
has_bounds(self)
return True, if any variable is bounded
inverse(self, y, copy_if_changed=True, copy_always=False)
is_in_bounds(self, x)
not yet tested
to_dim_times_two(self, bounds)
return boundaries in format ``[[lb0, ub0], [lb1, ub1], ...]``,
as used by ``BoxConstraints...`` class.

Data descriptors inherited from BoundaryHandlerBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class BoundTransform(BoundaryHandlerBase)
    Handles boundary by a smooth, piecewise linear and quadratic
transformation into the feasible domain.
 
>>> import cma
>>> veq = cma.Mh.vequals_approximately
>>> b = cma.BoundTransform([None, 1])
>>> assert b.bounds == [[None], [1]]
>>> assert veq(b.repair([0, 1, 1.2]), array([ 0., 0.975, 0.975]))
>>> assert b.is_in_bounds([0, 0.5, 1])
>>> assert veq(b.transform([0, 1, 2]), [ 0.   ,  0.975,  0.2  ])
>>> o=cma.fmin(cma.fcts.sphere, 6 * [-2], 0.5, options={
...    'boundary_handling': 'BoundTransform ',
...    'bounds': [[], 5 * [-1] + [inf]] })
>>> assert o[1] < 5 + 1e-8
>>> import numpy as np
>>> b = cma.BoundTransform([-np.random.rand(120), np.random.rand(120)])
>>> for i in range(100):
...     x = (-i-1) * np.random.rand(120) + i * np.random.randn(120)
...     x_to_b = b.repair(x)
...     x2 = b.inverse(x_to_b)
...     x2_to_b = b.repair(x2)
...     x3 = b.inverse(x2_to_b)
...     x3_to_b = b.repair(x3)
...     assert veq(x_to_b, x2_to_b)
...     assert veq(x2, x3)
...     assert veq(x2_to_b, x3_to_b)
 
Details: this class uses ``class BoxConstraintsLinQuadTransformation``
 
 
Method resolution order:
BoundTransform
BoundaryHandlerBase
__builtin__.object

Methods defined here:
__init__(self, bounds=None)
Argument bounds can be `None` or ``bounds[0]`` and ``bounds[1]``
are lower and upper domain boundaries, each is either `None` or
a scalar or a list or array of appropriate size.
inverse(self, x, copy_if_changed=True, copy_always=False)
inverse transform of ``x`` from the bounded domain.
repair(self, x, copy_if_changed=True, copy_always=False)
transforms ``x`` into the bounded domain.
 
``copy_always`` option might disappear.
transform(self, x)

Methods inherited from BoundaryHandlerBase:
__call__(self, solutions, *args, **kwargs)
return penalty or list of penalties, by default zero(s).
 
This interface seems too specifically tailored to the derived
BoundPenalty class, it should maybe change.
get_bounds(self, which, dimension)
``get_bounds('lower', 8)`` returns the lower bounds in 8-D
has_bounds(self)
return True, if any variable is bounded
is_in_bounds(self, x)
not yet tested
to_dim_times_two(self, bounds)
return boundaries in format ``[[lb0, ub0], [lb1, ub1], ...]``,
as used by ``BoxConstraints...`` class.
update(self, *args, **kwargs)

Data descriptors inherited from BoundaryHandlerBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAAdaptSigmaBase(__builtin__.object)
    step-size adaptation base class, implementing hsig functionality
via an isotropic evolution path.
 
  Methods defined here:
__init__(self, *args, **kwargs)
hsig(self, es)
return "OK-signal" for rank-one update, `True` (OK) or `False`
(stall rank-one update), based on the length of an evolution path
initialize_base(self, es)
set parameters and state variable based on dimension,
mueff and possibly further options.
update(self, es, **kwargs)
update ``es.sigma``

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAAdaptSigmaCSA(CMAAdaptSigmaBase)
    
Method resolution order:
CMAAdaptSigmaCSA
CMAAdaptSigmaBase
__builtin__.object

Methods defined here:
__init__(self)
postpone initialization to a method call where dimension and mueff should be known.
initialize(self, es)
set parameters and state variable based on dimension,
mueff and possibly further options.
update(self, es, **kwargs)

Methods inherited from CMAAdaptSigmaBase:
hsig(self, es)
return "OK-signal" for rank-one update, `True` (OK) or `False`
(stall rank-one update), based on the length of an evolution path
initialize_base(self, es)
set parameters and state variable based on dimension,
mueff and possibly further options.

Data descriptors inherited from CMAAdaptSigmaBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAAdaptSigmaDistanceProportional(CMAAdaptSigmaBase)
    artificial setting of ``sigma`` for test purposes, e.g.
to simulate optimal progress rates.
 
 
Method resolution order:
CMAAdaptSigmaDistanceProportional
CMAAdaptSigmaBase
__builtin__.object

Methods defined here:
__init__(self, coefficient=1.2)
update(self, es, **kwargs)

Methods inherited from CMAAdaptSigmaBase:
hsig(self, es)
return "OK-signal" for rank-one update, `True` (OK) or `False`
(stall rank-one update), based on the length of an evolution path
initialize_base(self, es)
set parameters and state variable based on dimension,
mueff and possibly further options.

Data descriptors inherited from CMAAdaptSigmaBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAAdaptSigmaMedianImprovement(CMAAdaptSigmaBase)
    Compares median fitness against a fitness percentile of the previous iteration,
see Ait ElHara et al, GECCO 2013.
 
 
Method resolution order:
CMAAdaptSigmaMedianImprovement
CMAAdaptSigmaBase
__builtin__.object

Methods defined here:
__init__(self)
initialize(self, es)
update(self, es, **kwargs)

Methods inherited from CMAAdaptSigmaBase:
hsig(self, es)
return "OK-signal" for rank-one update, `True` (OK) or `False`
(stall rank-one update), based on the length of an evolution path
initialize_base(self, es)
set parameters and state variable based on dimension,
mueff and possibly further options.

Data descriptors inherited from CMAAdaptSigmaBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAAdaptSigmaNone(CMAAdaptSigmaBase)
    
Method resolution order:
CMAAdaptSigmaNone
CMAAdaptSigmaBase
__builtin__.object

Methods defined here:
update(self, es, **kwargs)
no update, ``es.sigma`` remains constant.
 
:param es: ``CMAEvolutionStrategy`` class instance
:param kwargs: whatever else is needed to update ``es.sigma``

Methods inherited from CMAAdaptSigmaBase:
__init__(self, *args, **kwargs)
hsig(self, es)
return "OK-signal" for rank-one update, `True` (OK) or `False`
(stall rank-one update), based on the length of an evolution path
initialize_base(self, es)
set parameters and state variable based on dimension,
mueff and possibly further options.

Data descriptors inherited from CMAAdaptSigmaBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAAdaptSigmaTPA(CMAAdaptSigmaBase)
    two point adaptation for step-size sigma. Relies on a specific
sampling of the first two offspring, whose objective function
value ranks are used to decide on the step-size change.
 
Example
=======
 
>>> import cma
>>> cma.CMAOptions('adapt').pprint()
>>> es = cma.CMAEvolutionStrategy(10 * [0.2], 0.1, {'AdaptSigma': cma.CMAAdaptSigmaTPA, 'ftarget': 1e-8})
>>> es.optimize(cma.fcts.rosen)
>>> assert 'ftarget' in es.stop()
>>> assert es.result()[1] <= 1e-8
>>> assert es.result()[2] < 6500  # typically < 5500
 
References: loosely based on Hansen 2008, CMA-ES with Two-Point
Step-Size Adaptation, more tightly based on an upcoming paper by
Hansen et al.
 
 
Method resolution order:
CMAAdaptSigmaTPA
CMAAdaptSigmaBase
__builtin__.object

Methods defined here:
__init__(self, dimension=None, opts=None)
initialize(self, N=None, opts=None)
update(self, es, function_values, **kwargs)
the first and second value in ``function_values``
must reflect two mirrored solutions sampled
in direction / in opposite direction of
the previous mean shift, respectively.

Methods inherited from CMAAdaptSigmaBase:
hsig(self, es)
return "OK-signal" for rank-one update, `True` (OK) or `False`
(stall rank-one update), based on the length of an evolution path
initialize_base(self, es)
set parameters and state variable based on dimension,
mueff and possibly further options.

Data descriptors inherited from CMAAdaptSigmaBase:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMADataLogger(BaseDataLogger)
    data logger for class `CMAEvolutionStrategy`. The logger is
identified by its name prefix and (over-)writes or reads according
data files. Therefore, the logger must be considered as *global* variable
with unpredictable side effects, if two loggers with the same name
and on the same working folder are used at the same time.
 
Examples
========
::
 
    import cma
    es = cma.CMAEvolutionStrategy(...)
    logger = cma.CMADataLogger().register(es)
    while not es.stop():
        ...
        logger.add()  # add can also take an argument
 
    logger.plot() # or a short cut can be used:
    cma.plot()  # plot data from logger with default name
 
 
    logger2 = cma.CMADataLogger('just_another_filename_prefix').load()
    logger2.plot()
    logger2.disp()
 
::
 
    import cma
    from matplotlib.pylab import *
    res = cma.fmin(cma.Fcts.sphere, rand(10), 1e-0)
    logger = res[-1]  # the CMADataLogger
    logger.load()  # by "default" data are on disk
    semilogy(logger.f[:,0], logger.f[:,5])  # plot f versus iteration, see file header
    show()
 
Details
=======
After loading data, the logger has the attributes `xmean`, `xrecent`,
`std`, `f`, `D` and `corrspec` corresponding to ``xmean``,
``xrecentbest``, ``stddev``, ``fit``, ``axlen`` and ``axlencorr``
filename trails.
 
:See: `disp()`, `plot()`
 
 
Method resolution order:
CMADataLogger
BaseDataLogger
__builtin__.object

Methods defined here:
__init__(self, name_prefix=u'outcmaes', modulo=1, append=False)
initialize logging of data from a `CMAEvolutionStrategy`
instance, default ``modulo=1`` means logging with each call
add(self, es=None, more_data=[], modulo=None)
append some logging data from `CMAEvolutionStrategy` class instance `es`,
if ``number_of_times_called % modulo`` equals to zero, never if ``modulo==0``.
 
The sequence ``more_data`` must always have the same length.
 
When used for a different optimizer class, this function can be
(easily?) adapted by changing the assignments under INTERFACE
in the implemention.
closefig(self)
data(self)
return dictionary with data.
 
If data entries are None or incomplete, consider calling
``.load().data()`` to (re-)load the data from files first.
disp(self, idx=100)
displays selected data from (files written by) the class `CMADataLogger`.
 
Arguments
---------
   `idx`
       indices corresponding to rows in the data file;
       if idx is a scalar (int), the first two, then every idx-th,
       and the last three rows are displayed. Too large index values are removed.
 
Example
-------
>>> import cma, numpy as np
>>> res = cma.fmin(cma.fcts.elli, 7 * [0.1], 1, {'verb_disp':1e9})  # generate data
>>> assert res[1] < 1e-9
>>> assert res[2] < 4400
>>> l = cma.CMADataLogger()  # == res[-1], logger with default name, "points to" above data
>>> l.disp([0,-1])  # first and last
>>> l.disp(20)  # some first/last and every 20-th line
>>> l.disp(np.r_[0:999999:100, -1]) # every 100-th and last
>>> l.disp(np.r_[0, -10:0]) # first and ten last
>>> cma.disp(l.name_prefix, np.r_[0::100, -10:])  # the same as l.disp(...)
 
Details
-------
The data line with the best f-value is displayed as last line.
 
:See: `disp()`
disp_header(self)
downsampling(self, factor=10, first=3, switch=True, verbose=True)
rude downsampling of a `CMADataLogger` data file by `factor`,
keeping also the first `first` entries. This function is a
stump and subject to future changes. Return self.
 
Arguments
---------
   - `factor` -- downsampling factor
   - `first` -- keep first `first` entries
   - `switch` -- switch the new logger to the downsampled logger
        original_name+'down'
 
Details
-------
``self.name_prefix+'down'`` files are written
 
Example
-------
::
 
    import cma
    cma.downsampling()  # takes outcmaes* files
    cma.plot('outcmaesdown')
initialize(self, modulo=None)
reset logger, overwrite original files, `modulo`: log only every modulo call
load(self, filenameprefix=None)
load (or reload) data from output files, `load()` is called in
`plot()` and `disp()`.
 
Argument `filenameprefix` is the filename prefix of data to be
loaded (six files), by default ``'outcmaes'``.
 
Return self with (added) attributes `xrecent`, `xmean`,
`f`, `D`, `std`, 'corrspec'
plot(self, fig=None, iabscissa=1, iteridx=None, plot_mean=False, foffset=1e-19, x_opt=None, fontsize=9)
plot data from a `CMADataLogger` (using the files written 
by the logger).
 
Arguments
---------
    `fig`
        figure number, by default 325
    `iabscissa`
        ``0==plot`` versus iteration count,
        ``1==plot`` versus function evaluation number
    `iteridx`
        iteration indices to plot
 
Return `CMADataLogger` itself.
 
Examples
--------
::
 
    import cma
    logger = cma.CMADataLogger()  # with default name
    # try to plot the "default logging" data (e.g.
    #   from previous fmin calls, which is essentially what
    #   also cma.plot() does)
    logger.plot()
    cma.savefig('fig325.png')  # save current figure
    logger.closefig()
 
Dependencies: matlabplotlib/pyplot.
plot_all(self, fig=None, iabscissa=1, iteridx=None, foffset=1e-19, x_opt=None, fontsize=9)
plot data from a `CMADataLogger` (using the files written by the logger).
 
Arguments
---------
    `fig`
        figure number, by default 425
    `iabscissa`
        ``0==plot`` versus iteration count,
        ``1==plot`` versus function evaluation number
    `iteridx`
        iteration indices to plot
 
Return `CMADataLogger` itself.
 
Examples
--------
::
 
    import cma
    logger = cma.CMADataLogger()  # with default name
    # try to plot the "default logging" data (e.g.
    #   from previous fmin calls, which is essentially what
    #   also cma.plot() does)
    logger.plot_all()
    cma.savefig('fig425.png')  # save current figure
    logger.closefig()
 
Dependencies: matlabplotlib/pyplot.
plot_axes_scaling(self, iabscissa=1)
plot_correlations(self, iabscissa=1)
spectrum of correlation matrix and largest correlation
plot_divers(self, iabscissa=1, foffset=1e-19)
plot fitness, sigma, axis ratio...
 
:param iabscissa: 0 means vs evaluations, 1 means vs iterations
:param foffset: added to f-value
 
:See: `plot()`
plot_mean(self, iabscissa=1, x_opt=None, annotations=None)
plot_stds(self, iabscissa=1)
plot_xrecent(self, iabscissa=1, x_opt=None, annotations=None)
register(self, es, append=None, modulo=None)
register a `CMAEvolutionStrategy` instance for logging,
``append=True`` appends to previous data logged under the same name,
by default previous data are overwritten.
save_to(self, nameprefix, switch=False)
saves logger data to a different set of files, for
``switch=True`` also the loggers name prefix is switched to
the new value
select_data(self, iteration_indices)
keep only data of `iteration_indices`

Data and other attributes defined here:
default_prefix = u'outcmaes'

Data descriptors inherited from BaseDataLogger:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAEvolutionStrategy(OOOptimizer)
    CMA-ES stochastic optimizer class with ask-and-tell interface.
 
Calling Sequences
=================
 
    es = CMAEvolutionStrategy(x0, sigma0)
 
    es = CMAEvolutionStrategy(x0, sigma0, opts)
 
    es = CMAEvolutionStrategy(x0, sigma0).optimize(objective_fct)
 
    res = CMAEvolutionStrategy(x0, sigma0,
                            opts).optimize(objective_fct).result()
 
Arguments
=========
    `x0`
        initial solution, starting point. `x0` is given as "phenotype"
        which means, if::
 
            opts = {'transformation': [transform, inverse]}
 
        is given and ``inverse is None``, the initial mean is not
        consistent with `x0` in that ``transform(mean)`` does not
        equal to `x0` unless ``transform(mean)`` equals ``mean``.
    `sigma0`
        initial standard deviation.  The problem variables should
        have been scaled, such that a single standard deviation
        on all variables is useful and the optimum is expected to
        lie within about `x0` +- ``3*sigma0``. See also options
        `scaling_of_variables`. Often one wants to check for
        solutions close to the initial point. This allows,
        for example, for an easier check of consistency of the
        objective function and its interfacing with the optimizer.
        In this case, a much smaller `sigma0` is advisable.
    `opts`
        options, a dictionary with optional settings,
        see class `CMAOptions`.
 
Main interface / usage
======================
The interface is inherited from the generic `OOOptimizer`
class (see also there). An object instance is generated from
 
    es = cma.CMAEvolutionStrategy(8 * [0.5], 0.2)
 
The least verbose interface is via the optimize method::
 
    es.optimize(objective_func)
    res = es.result()
 
More verbosely, the optimization is done using the
methods ``stop``, ``ask``, and ``tell``::
 
    while not es.stop():
        solutions = es.ask()
        es.tell(solutions, [cma.fcts.rosen(s) for s in solutions])
        es.disp()
    es.result_pretty()
 
 
where ``ask`` delivers new candidate solutions and ``tell`` updates
the ``optim`` instance by passing the respective function values
(the objective function ``cma.fcts.rosen`` can be replaced by any
properly defined objective function, see ``cma.fcts`` for more
examples).
 
To change an option, for example a termination condition to
continue the optimization, call
 
    es.opts.set({'tolfacupx': 1e4})
 
The class `CMAEvolutionStrategy` also provides::
 
    (solutions, func_values) = es.ask_and_eval(objective_func)
 
and an entire optimization can also be written like::
 
    while not es.stop():
        es.tell(*es.ask_and_eval(objective_func))
 
Besides for termination criteria, in CMA-ES only the ranks of the
`func_values` are relevant.
 
Attributes and Properties
=========================
    - `inputargs` -- passed input arguments
    - `inopts` -- passed options
    - `opts` -- actually used options, some of them can be changed any
      time via ``opts.set``, see class `CMAOptions`
    - `popsize` -- population size lambda, number of candidate
       solutions returned by `ask()`
    - `logger` -- a `CMADataLogger` instance utilized by `optimize`
 
Examples
========
Super-short example, with output shown:
 
>>> import cma
>>> # construct an object instance in 4-D, sigma0=1:
>>> es = cma.CMAEvolutionStrategy(4 * [1], 1, {'seed':234})
(4_w,8)-CMA-ES (mu_w=2.6,w_1=52%) in dimension 4 (seed=234)
>>>
>>> # optimize the ellipsoid function
>>> es.optimize(cma.fcts.elli, verb_disp=1)
Iterat #Fevals   function value     axis ratio  sigma   minstd maxstd min:sec
    1       8 2.093015112685775e+04 1.0e+00 9.27e-01  9e-01  9e-01 0:0.0
    2      16 4.964814235917688e+04 1.1e+00 9.54e-01  9e-01  1e+00 0:0.0
    3      24 2.876682459926845e+05 1.2e+00 1.02e+00  9e-01  1e+00 0:0.0
  100     800 6.809045875281943e-01 1.3e+02 1.41e-02  1e-04  1e-02 0:0.2
  200    1600 2.473662150861846e-10 8.0e+02 3.08e-05  1e-08  8e-06 0:0.5
  233    1864 2.766344961865341e-14 8.6e+02 7.99e-07  8e-11  7e-08 0:0.6
>>>
>>> cma.pprint(es.result())
(array([ -1.98546755e-09,  -1.10214235e-09,   6.43822409e-11,
        -1.68621326e-11]),
 4.5119610261406537e-16,
 1666,
 1672,
 209,
 array([ -9.13545269e-09,  -1.45520541e-09,  -6.47755631e-11,
        -1.00643523e-11]),
 array([  3.20258681e-08,   3.15614974e-09,   2.75282215e-10,
         3.27482983e-11]))
>>> assert es.result()[1] < 1e-9
>>> help(es.result)
Help on method result in module cma:
 
result(self) method of cma.CMAEvolutionStrategy instance
    return ``(xbest, f(xbest), evaluations_xbest, evaluations, iterations, pheno(xmean), effective_stds)``
 
 
The optimization loop can also be written explicitly.
 
>>> import cma
>>> es = cma.CMAEvolutionStrategy(4 * [1], 1)
>>> while not es.stop():
...    X = es.ask()
...    es.tell(X, [cma.fcts.elli(x) for x in X])
...    es.disp()
<output omitted>
 
achieving the same result as above.
 
An example with lower bounds (at zero) and handling infeasible
solutions:
 
>>> import cma
>>> import numpy as np
>>> es = cma.CMAEvolutionStrategy(10 * [0.2], 0.5, {'bounds': [0, np.inf]})
>>> while not es.stop():
...     fit, X = [], []
...     while len(X) < es.popsize:
...         curr_fit = None
...         while curr_fit in (None, np.NaN):
...             x = es.ask(1)[0]
...             curr_fit = cma.fcts.somenan(x, cma.fcts.elli) # might return np.NaN
...         X.append(x)
...         fit.append(curr_fit)
...     es.tell(X, fit)
...     es.logger.add()
...     es.disp()
<output omitted>
>>>
>>> assert es.result()[1] < 1e-9
>>> assert es.result()[2] < 9000  # by internal termination
>>> # es.logger.plot()  # will plot data
>>> # cma.show()  # display plot window
 
An example with user-defined transformation, in this case to realize
a lower bound of 2.
 
>>> es = cma.CMAEvolutionStrategy(5 * [3], 1,
...                 {"transformation": [lambda x: x**2+2, None]})
>>> es.optimize(cma.fcts.rosen)
<output omitted>
>>> assert cma.fcts.rosen(es.result()[0]) < 1e-6 + 5.530760944396627e+02
>>> assert es.result()[2] < 3300
 
The inverse transformation is (only) necessary if the `BoundPenalty`
boundary handler is used at the same time.
 
The ``CMAEvolutionStrategy`` class also provides a default logger
(cave: files are overwritten when the logger is used with the same
filename prefix):
 
>>> import cma
>>> es = cma.CMAEvolutionStrategy(4 * [0.2], 0.5, {'verb_disp': 0})
>>> es.logger.disp_header()  # to understand the print of disp
Iterat Nfevals  function value    axis ratio maxstd   minstd
>>> while not es.stop():
...     X = es.ask()
...     es.tell(X, [cma.fcts.sphere(x) for x in X])
...     es.logger.add()  # log current iteration
...     es.logger.disp([-1])  # display info for last iteration
1      8 2.72769793021748e+03 1.0e+00 4.05e-01 3.99e-01
2     16 6.58755537926063e+03 1.1e+00 4.00e-01 3.39e-01
<output ommitted>
193   1544 3.15195320957214e-15 1.2e+03 3.70e-08 3.45e-11
>>> es.logger.disp_header()
Iterat Nfevals  function value    axis ratio maxstd   minstd
>>> # es.logger.plot() # will make a plot
 
Example implementing restarts with increasing popsize (IPOP), output
is not displayed:
 
>>> import cma, numpy as np
>>>
>>> # restart with increasing population size (IPOP)
>>> bestever = cma.BestSolution()
>>> for lam in 10 * 2**np.arange(8):  # 10, 20, 40, 80, ..., 10 * 2**7
...     es = cma.CMAEvolutionStrategy('6 - 8 * np.random.rand(9)',  # 9-D
...                                   5,  # initial std sigma0
...                                   {'popsize': lam,  # options
...                                    'verb_append': bestever.evalsall})
...     logger = cma.CMADataLogger().register(es, append=bestever.evalsall)
...     while not es.stop():
...         X = es.ask()    # get list of new solutions
...         fit = [cma.fcts.rastrigin(x) for x in X]  # evaluate each solution
...         es.tell(X, fit) # besides for termination only the ranking in fit is used
...
...         # display some output
...         logger.add()  # add a "data point" to the log, writing in files
...         es.disp()  # uses option verb_disp with default 100
...
...     print('termination:', es.stop())
...     cma.pprint(es.best.__dict__)
...
...     bestever.update(es.best)
...
...     # show a plot
...     # logger.plot();
...     if bestever.f < 1e-8:  # global optimum was hit
...         break
<output omitted>
>>> assert es.result()[1] < 1e-8
 
On the Rastrigin function, usually after five restarts the global
optimum is located.
 
Using the ``multiprocessing`` module, we can evaluate the function in
parallel with a simple modification of the example (however
multiprocessing seems not always reliable)::
 
    try:
        import multiprocessing as mp
        import cma
        es = cma.CMAEvolutionStrategy(22 * [0.0], 1.0, {'maxiter':10})
        pool = mp.Pool(es.popsize)
        while not es.stop():
            X = es.ask()
            f_values = pool.map_async(cma.felli, X).get()
            # use chunksize parameter as es.popsize/len(pool)?
            es.tell(X, f_values)
            es.disp()
            es.logger.add()
    except ImportError:
        pass
 
The final example shows how to resume:
 
>>> import cma, pickle
>>>
>>> es = cma.CMAEvolutionStrategy(12 * [0.1],  # a new instance, 12-D
...                               0.5)         # initial std sigma0
>>> es.optimize(cma.fcts.rosen, iterations=100)
>>> pickle.dump(es, open('saved-cma-object.pkl', 'wb'))
>>> print('saved')
>>> del es  # let's start fresh
>>>
>>> es = pickle.load(open('saved-cma-object.pkl', 'rb'))
>>> print('resumed')
>>> es.optimize(cma.fcts.rosen, verb_disp=200)
>>> assert es.result()[2] < 15000
>>> cma.pprint(es.result())
 
Details
=======
The following two enhancements are implemented, the latter is turned
on by default only for very small population size.
 
*Active CMA* is implemented with option ``CMA_active`` and
conducts an update of the covariance matrix with negative weights.
The negative update is implemented, such that positive definiteness
is guarantied. The update is applied after the default update and
only before the covariance matrix is decomposed, which limits the
additional computational burden to be at most a factor of three
(typically smaller). A typical speed up factor (number of
f-evaluations) is between 1.1 and two.
 
References: Jastrebski and Arnold, CEC 2006, Glasmachers et al, GECCO 2010.
 
*Selective mirroring* is implemented with option ``CMA_mirrors``
in the method ``get_mirror()``. Only the method `ask_and_eval()`
(used by `fmin`) will then sample selectively mirrored vectors. In
selective mirroring, only the worst solutions are mirrored. With
the default small number of mirrors, *pairwise selection* (where at
most one of the two mirrors contribute to the update of the
distribution mean) is implicitly guarantied under selective
mirroring and therefore not explicitly implemented.
 
References: Brockhoff et al, PPSN 2010, Auger et al, GECCO 2011.
 
:See: `fmin()`, `OOOptimizer`, `CMAOptions`, `plot()`, `ask()`,
    `tell()`, `ask_and_eval()`
 
 
Method resolution order:
CMAEvolutionStrategy
OOOptimizer
__builtin__.object

Methods defined here:
__init__(self, x0, sigma0, inopts={})
see class `CMAEvolutionStrategy`
ask(self, number=None, xmean=None, sigma_fac=1, gradf=None, args=())
get new candidate solutions, sampled from a multi-variate
normal distribution and transformed to f-representation
(phenotype) to be evaluated.
 
Arguments
---------
    `number`
        number of returned solutions, by default the
        population size ``popsize`` (AKA ``lambda``).
    `xmean`
        distribution mean, phenotyp?
    `sigma_fac`
        multiplier for internal sample width (standard
        deviation)
    `gradf`
        gradient, ``len(gradf(x)) == len(x)``, if
        ``gradf is not None`` the third solution in the
        returned list is "sampled" in supposedly Newton
        direction ``dot(C, gradf(xmean, *args))``.
    `args`
        additional arguments passed to gradf
 
Return
------
A list of N-dimensional candidate solutions to be evaluated
 
Example
-------
>>> import cma
>>> es = cma.CMAEvolutionStrategy([0,0,0,0], 0.3)
>>> while not es.stop() and es.best.f > 1e-6:  # my_desired_target_f_value
...     X = es.ask()  # get list of new solutions
...     fit = [cma.fcts.rosen(x) for x in X]  # call function rosen with each solution
...     es.tell(X, fit)  # feed values
 
:See: `ask_and_eval`, `ask_geno`, `tell`
ask_and_eval(self, func, args=(), gradf=None, number=None, xmean=None, sigma_fac=1, evaluations=1, aggregation=<function median>, kappa=1)
samples `number` solutions and evaluates them on `func`, where
each solution `s` is resampled until ``self.is_feasible(s, func(s)) is True``.
 
Arguments
---------
    `func`
        objective function, ``func(x)`` returns a scalar
    `args`
        additional parameters for `func`
    `gradf`
        gradient of objective function, ``g = gradf(x, *args)``
        must satisfy ``len(g) == len(x)``
    `number`
        number of solutions to be sampled, by default
        population size ``popsize`` (AKA lambda)
    `xmean`
        mean for sampling the solutions, by default ``self.mean``.
    `sigma_fac`
        multiplier for sampling width, standard deviation, for example
        to get a small perturbation of solution `xmean`
    `evaluations`
        number of evaluations for each sampled solution
    `aggregation`
        function that aggregates `evaluations` values to
        as single value.
    `kappa`
        multiplier used for the evaluation of the solutions, in
        that ``func(m + kappa*(x - m))`` is the f-value for x.
 
Return
------
``(X, fit)``, where
    X -- list of solutions
    fit -- list of respective function values
 
Details
-------
While ``not self.is_feasible(x, func(x))``new solutions are sampled. By
default ``self.is_feasible == cma.feasible == lambda x, f: f not in (None, np.NaN)``.
The argument to `func` can be freely modified within `func`.
 
Depending on the ``CMA_mirrors`` option, some solutions are not sampled
independently but as mirrors of other bad solutions. This is a simple
derandomization that can save 10-30% of the evaluations in particular
with small populations, for example on the cigar function.
 
Example
-------
>>> import cma
>>> x0, sigma0 = 8*[10], 1  # 8-D
>>> es = cma.CMAEvolutionStrategy(x0, sigma0)
>>> while not es.stop():
...     X, fit = es.ask_and_eval(cma.fcts.elli)  # handles NaN with resampling
...     es.tell(X, fit)  # pass on fitness values
...     es.disp(20) # print every 20-th iteration
>>> print('terminated on ' + str(es.stop()))
<output omitted>
 
A single iteration step can be expressed in one line, such that
an entire optimization after initialization becomes
::
 
    while not es.stop():
        es.tell(*es.ask_and_eval(cma.fcts.elli))
ask_geno(self, number=None, xmean=None, sigma_fac=1)
get new candidate solutions in genotyp, sampled from a
multi-variate normal distribution.
 
Arguments are
    `number`
        number of returned solutions, by default the
        population size `popsize` (AKA lambda).
    `xmean`
        distribution mean
    `sigma_fac`
        multiplier for internal sample width (standard
        deviation)
 
`ask_geno` returns a list of N-dimensional candidate solutions
in genotyp representation and is called by `ask`.
 
Details: updates the sample distribution and might change
the geno-pheno transformation during this update.
 
:See: `ask`, `ask_and_eval`
clip_or_fit_solutions(self, pop, idx)
make sure that solutions fit to sample distribution, this interface will probably change.
 
In particular the frequency of long vectors appearing in pop[idx] - self.mean is limited.
copy_constructor(self, es)
correlation_matrix(self)
# ____________________________________________________________
# ____________________________________________________________
decompose_C(self)
eigen-decompose self.C and update self.dC, self.C, self.B.
 
Known bugs: this might give a runtime error with
CMA_diagonal / separable option on.
disp(self, modulo=None)
prints some single-line infos according to `disp_annotation()`,
if ``iteration_counter % modulo == 0``
disp_annotation(self)
print annotation for `disp()`
eval_mean(self, func, args=())
evaluate the distribution mean, this is not (yet) effective
in terms of termination or display
feedForResume(self, X, function_values)
Given all "previous" candidate solutions and their respective
function values, the state of a `CMAEvolutionStrategyobject
can be reconstructed from this history. This is the purpose of
function `feedForResume`.
 
Arguments
---------
    `X`
      (all) solution points in chronological order, phenotypic
      representation. The number of points must be a multiple
      of popsize.
    `function_values`
      respective objective function values
 
Details
-------
`feedForResume` can be called repeatedly with only parts of
the history. The part must have the length of a multiple
of the population size.
`feedForResume` feeds the history in popsize-chunks into `tell`.
The state of the random number generator might not be
reconstructed, but this would be only relevant for the future.
 
Example
-------
::
 
    import cma
 
    # prepare
    (x0, sigma0) = ... # initial values from previous trial
    X = ... # list of generated solutions from a previous trial
    f = ... # respective list of f-values
 
    # resume
    es = cma.CMAEvolutionStrategy(x0, sigma0)
    es.feedForResume(X, f)
 
    # continue with func as objective function
    while not es.stop():
       X = es.ask()
       es.tell(X, [func(x) for x in X])
 
Credits to Dirk Bueche and Fabrice Marchal for the feeding idea.
 
:See: class `CMAEvolutionStrategy` for a simple dump/load to resume
get_mirror(self, x, preserve_length=False)
return ``pheno(self.mean - (geno(x) - self.mean))``.
 
>>> import cma
>>> es = cma.CMAEvolutionStrategy(cma.np.random.randn(3), 1)
>>> x = cma.np.random.randn(3)
>>> assert cma.Mh.vequals_approximately(es.mean - (x - es.mean), es.get_mirror(x, preserve_length=True))
>>> x = es.ask(1)[0]
>>> vals = (es.get_mirror(x) - es.mean) / (x - es.mean)
>>> assert cma.Mh.equals_approximately(sum(vals), len(vals) * vals[0])
 
TODO: this implementation is yet experimental.
 
TODO: this implementation includes geno-pheno transformation,
however in general GP-transformation should be separated from
specific code.
 
Selectively mirrored sampling improves to a moderate extend but
overadditively with active CMA for quite understandable reasons.
 
Optimal number of mirrors are suprisingly small: 1,2,3 for
maxlam=7,13,20 where 3,6,10 are the respective maximal possible
mirrors that must be clearly suboptimal.
get_selective_mirrors(self, number=None, pop_sorted=None)
get mirror genotypic directions of the `number` worst
solution, based on ``pop_sorted`` attribute (from last
iteration).
 
Details:
Takes the last ``number=sp.lam_mirr`` entries in
``pop_sorted=self.pop_sorted`` as solutions to be mirrored.
inject(self, solutions)
inject a genotypic solution. The solution is used as direction
relative to the distribution mean to compute a new candidate
solution returned in method `ask_geno` which in turn is used in
method `ask`.
 
>>> import cma
>>> es = cma.CMAEvolutionStrategy(4 * [1], 2)
>>> while not es.stop():
...     es.inject([4 * [0.0]])
...     X = es.ask()
...     break
>>> assert X[0][0] == X[0][1]
mahalanobis_norm(self, dx)
compute the Mahalanobis norm that is induced by the adapted
sample distribution, covariance matrix ``C`` times ``sigma**2``,
including ``sigma_vec``. The expected Mahalanobis distance to
the sample mean is about ``sqrt(dimension)``.
 
Argument
--------
A *genotype* difference `dx`.
 
Example
-------
>>> import cma, numpy
>>> es = cma.CMAEvolutionStrategy(numpy.ones(10), 1)
>>> xx = numpy.random.randn(2, 10)
>>> d = es.mahalanobis_norm(es.gp.geno(xx[0]-xx[1]))
 
`d` is the distance "in" the true sample distribution,
sampled points have a typical distance of ``sqrt(2*es.N)``,
where ``es.N`` is the dimension, and an expected distance of
close to ``sqrt(N)`` to the sample mean. In the example,
`d` is the Euclidean distance, because C = I and sigma = 1.
multiplyC(self, alpha)
multiply C with a scalar and update all related internal variables (dC, D,...)
plot(self)
prepare_injection_directions(self)
provide genotypic directions for TPA and selective mirroring,
with no specific length normalization, to be used in the
coming iteration.
 
Details:
This method is called in the end of `tell`. The result is
assigned to ``self.pop_injection_directions`` and used in
`ask_geno`.
 
TODO: should be rather appended?
random_rescale_to_mahalanobis(self, x)
change `x` like for injection, all on genotypic level
random_rescaling_factor_to_mahalanobis_size(self, y)
``self.mean + self.random_rescaling_factor_to_mahalanobis_size(y)``
is guarantied to appear like from the sample distribution.
readProperties(self)
reads dynamic parameters from property file (not implemented)
repair_genotype(self, x, copy_if_changed=False)
make sure that solutions fit to the sample distribution, this interface will probably change.
 
In particular the frequency of x - self.mean being long is limited.
result(self)
return::
 
(xbest, f(xbest), evaluations_xbest, evaluations, iterations,
    pheno(xmean), effective_stds)
result_pretty(self, number_of_runs=0, time_str=None, fbestever=None)
pretty print result.
 
Returns ``self.result()``
stop(self, check=True)
return a dictionary with the termination status.
With ``check==False``, the termination conditions are not checked
and the status might not reflect the current situation.
tell(self, solutions, function_values, check_points=None, copy=False)
pass objective function values to prepare for next
iteration. This core procedure of the CMA-ES algorithm updates
all state variables, in particular the two evolution paths, the
distribution mean, the covariance matrix and a step-size.
 
Arguments
---------
    `solutions`
        list or array of candidate solution points (of
        type `numpy.ndarray`), most presumably before
        delivered by method `ask()` or `ask_and_eval()`.
    `function_values`
        list or array of objective function values
        corresponding to the respective points. Beside for termination
        decisions, only the ranking of values in `function_values`
        is used.
    `check_points`
        If ``check_points is None``, only solutions that are not generated
        by `ask()` are possibly clipped (recommended). ``False`` does not clip
        any solution (not recommended).
        If ``True``, clips solutions that realize long steps (i.e. also
        those that are unlikely to be generated with `ask()`). `check_points`
        can be a list of indices to be checked in solutions.
    `copy`
        ``solutions`` can be modified in this routine, if ``copy is False``
 
Details
-------
`tell()` updates the parameters of the multivariate
normal search distribution, namely covariance matrix and
step-size and updates also the attributes ``countiter`` and
``countevals``. To check the points for consistency is quadratic
in the dimension (like sampling points).
 
Bugs
----
The effect of changing the solutions delivered by `ask()`
depends on whether boundary handling is applied. With boundary
handling, modifications are disregarded. This is necessary to
apply the default boundary handling that uses unrepaired
solutions but might change in future.
 
Example
-------
::
 
    import cma
    func = cma.fcts.elli  # choose objective function
    es = cma.CMAEvolutionStrategy(cma.np.random.rand(10), 1)
    while not es.stop():
       X = es.ask()
       es.tell(X, [func(x) for x in X])
    es.result()  # where the result can be found
 
:See: class `CMAEvolutionStrategy`, `ask()`, `ask_and_eval()`, `fmin()`
updateBD(self)
update internal variables for sampling the distribution with the
current covariance matrix C. This method is O(N^3), if C is not diagonal.
update_exponential(self, Z, eta, BDpair=None)
exponential update of C that guarantees positive definiteness, that is,
instead of the assignment ``C = C + eta * Z``,
we have ``C = C**.5 * exp(eta * C**-.5 * Z * C**-.5) * C**.5``.
 
Parameter `Z` should have expectation zero, e.g. sum(w[i] * z[i] * z[i].T) - C
if E z z.T = C.
 
Parameter `eta` is the learning rate, for ``eta == 0`` nothing is updated.
 
This function conducts two eigendecompositions, assuming that
B and D are not up to date, unless `BDpair` is given. Given BDpair,
B is the eigensystem and D is the vector of sqrt(eigenvalues), one
eigendecomposition is omitted.
 
Reference: Glasmachers et al 2010, Exponential Natural Evolution Strategies

Data descriptors defined here:
popsize
number of samples by default returned by` ask()`

Methods inherited from OOOptimizer:
initialize(self)
(re-)set to the initial state
optimize(self, objective_fct, iterations=None, min_iterations=1, args=(), verb_disp=None, logger=None, call_back=None)
find minimizer of `objective_fct`.
 
CAVEAT: the return value for `optimize` has changed to ``self``.
 
Arguments
---------
 
    `objective_fct`
        function be to minimized
    `iterations`
        number of (maximal) iterations, while ``not self.stop()``
    `min_iterations`
        minimal number of iterations, even if ``not self.stop()``
    `args`
        arguments passed to `objective_fct`
    `verb_disp`
        print to screen every `verb_disp` iteration, if ``None``
        the value from ``self.logger`` is "inherited", if
        available.
    ``logger``
        a `BaseDataLogger` instance, which must be compatible
        with the type of ``self``.
    ``call_back``
        call back function called like ``call_back(self)`` or
        a list of call back functions.
 
``return self``, that is, the `OOOptimizer` instance.
 
Example
-------
>>> import cma
>>> es = cma.CMAEvolutionStrategy(7 * [0.1], 0.5
...              ).optimize(cma.fcts.rosen, verb_disp=100)
(4_w,9)-CMA-ES (mu_w=2.8,w_1=49%) in dimension 7 (seed=630721393)
Iterat #Fevals   function value    axis ratio  sigma   minstd maxstd min:sec
    1       9 3.163954777181882e+01 1.0e+00 4.12e-01  4e-01  4e-01 0:0.0
    2      18 3.299006223906629e+01 1.0e+00 3.60e-01  3e-01  4e-01 0:0.0
    3      27 1.389129389866704e+01 1.1e+00 3.18e-01  3e-01  3e-01 0:0.0
  100     900 2.494847340045985e+00 8.6e+00 5.03e-02  2e-02  5e-02 0:0.3
  200    1800 3.428234862999135e-01 1.7e+01 3.77e-02  6e-03  3e-02 0:0.5
  300    2700 3.216640032470860e-04 5.6e+01 6.62e-03  4e-04  9e-03 0:0.8
  400    3600 6.155215286199821e-12 6.6e+01 7.44e-06  1e-07  4e-06 0:1.1
  438    3942 1.187372505161762e-14 6.0e+01 3.27e-07  4e-09  9e-08 0:1.2
  438    3942 1.187372505161762e-14 6.0e+01 3.27e-07  4e-09  9e-08 0:1.2
('termination by', {'tolfun': 1e-11})
('best f-value =', 1.1189867885201275e-14)
('solution =', array([ 1.        ,  1.        ,  1.        ,  0.99999999,  0.99999998,
        0.99999996,  0.99999992]))
>>> print(es.result()[0])
array([ 1.          1.          1.          0.99999999  0.99999998  0.99999996
  0.99999992])

Data descriptors inherited from OOOptimizer:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class CMAOptions(__builtin__.dict)
    ``CMAOptions()`` returns a dictionary with the available options
and their default values for class ``CMAEvolutionStrategy``.
 
``CMAOptions('pop')`` returns a subset of recognized options that
contain 'pop' in there keyword name or (default) value or description.
 
``CMAOptions(opts)`` returns the subset of recognized options in
``dict(opts)``.
 
Option values can be "written" in a string and, when passed to fmin
or CMAEvolutionStrategy, are evaluated using "N" and "popsize" as
known values for dimension and population size (sample size, number
of new solutions per iteration). All default option values are such
a string.
 
Details
-------
``CMAOptions`` entries starting with ``tol`` are termination
"tolerances".
 
For `tolstagnation`, the median over the first and the second half
of at least `tolstagnation` iterations are compared for both, the
per-iteration best and per-iteration median function value.
 
Example
-------
::
 
    import cma
    cma.CMAOptions('tol')
 
is a shortcut for cma.CMAOptions().match('tol') that returns all options
that contain 'tol' in their name or description.
 
To set an option::
 
    import cma
    opts = cma.CMAOptions()
    opts.set('tolfun', 1e-12)
    opts['tolx'] = 1e-11
 
:See: `fmin`(), `CMAEvolutionStrategy`, `_CMAParameters`
 
 
Method resolution order:
CMAOptions
__builtin__.dict
__builtin__.object

Methods defined here:
__call__(self, key, default=None, loc=None)
evaluate and return the value of option `key` on the fly, or
returns those options whose name or description contains `key`,
case disregarded.
 
Details
-------
Keys that contain `filename` are not evaluated.
For ``loc==None``, `self` is used as environment
but this does not define ``N``.
 
:See: `eval()`, `evalall()`
__init__(self, s=None, unchecked=False)
return an `CMAOptions` instance, either with the default
options, if ``s is None``, or with all options whose name or
description contains `s`, if `s` is a string (case is
disregarded), or with entries from dictionary `s` as options,
not complemented with default options or settings
 
Returns: see above.
check(self, options=None)
check for ambiguous keys and move attributes into dict
check_attributes(self, opts=None)
check for attributes and moves them into the dictionary
check_values(self, options=None)
complement(self)
add all missing options with their default values
corrected_key(self, key)
return the matching valid key, if ``key.lower()`` is a unique
starting sequence to identify the valid key, ``else None``
eval(self, key, default=None, loc=None, correct_key=True)
Evaluates and sets the specified option value in
environment `loc`. Many options need ``N`` to be defined in
`loc`, some need `popsize`.
 
Details
-------
Keys that contain 'filename' are not evaluated.
For `loc` is None, the self-dict is used as environment
 
:See: `evalall()`, `__call__`
evalall(self, loc=None, defaults=None)
Evaluates all option values in environment `loc`.
 
:See: `eval()`
init(self, dict_or_str, val=None, warn=True)
initialize one or several options.
 
Arguments
---------
    `dict_or_str`
        a dictionary if ``val is None``, otherwise a key.
        If `val` is provided `dict_or_str` must be a valid key.
    `val`
        value for key
 
Details
-------
Only known keys are accepted. Known keys are in `CMAOptions.defaults()`
match(self, s=u'')
return all options that match, in the name or the description,
with string `s`, case is disregarded.
 
Example: ``cma.CMAOptions().match('verb')`` returns the verbosity
options.
pp(self)
pprint(self, linebreak=80)
print_ = pprint(self, linebreak=80)
printme = pprint(self, linebreak=80)
set(self, dic, val=None, force=False)
set can assign versatile options from
`CMAOptions.versatile_options()` with a new value, use `init()`
for the others.
 
Arguments
---------
    `dic`
        either a dictionary or a key. In the latter
        case, `val` must be provided
    `val`
        value for `key`, approximate match is sufficient
    `force`
        force setting of non-versatile options, use with caution
 
This method will be most probably used with the ``opts`` attribute of
a `CMAEvolutionStrategy` instance.
settable(self)
return the subset of those options that are settable at any
time.
 
Settable options are in `versatile_options()`, but the
list might be incomplete.

Static methods defined here:
defaults()
return a dictionary with default option values and description
merge(self, dict_=None)
not is use so far, see check()
versatile_options()
return list of options that can be changed at any time (not
only be initialized), however the list might not be entirely up
to date.
 
The string ' #v ' in the default value indicates a 'versatile'
option that can be changed any time.

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

Methods inherited from __builtin__.dict:
__cmp__(...)
x.__cmp__(y) <==> cmp(x,y)
__contains__(...)
D.__contains__(k) -> True if D has a key k, else False
__delitem__(...)
x.__delitem__(y) <==> del x[y]
__eq__(...)
x.__eq__(y) <==> x==y
__ge__(...)
x.__ge__(y) <==> x>=y
__getattribute__(...)
x.__getattribute__('name') <==> x.name
__getitem__(...)
x.__getitem__(y) <==> x[y]
__gt__(...)
x.__gt__(y) <==> x>y
__iter__(...)
x.__iter__() <==> iter(x)
__le__(...)
x.__le__(y) <==> x<=y
__len__(...)
x.__len__() <==> len(x)
__lt__(...)
x.__lt__(y) <==> x<y
__ne__(...)
x.__ne__(y) <==> x!=y
__repr__(...)
x.__repr__() <==> repr(x)
__setitem__(...)
x.__setitem__(i, y) <==> x[i]=y
__sizeof__(...)
D.__sizeof__() -> size of D in memory, in bytes
clear(...)
D.clear() -> None.  Remove all items from D.
copy(...)
D.copy() -> a shallow copy of D
fromkeys(...)
dict.fromkeys(S[,v]) -> New dict with keys from S and values equal to v.
v defaults to None.
get(...)
D.get(k[,d]) -> D[k] if k in D, else d.  d defaults to None.
has_key(...)
D.has_key(k) -> True if D has a key k, else False
items(...)
D.items() -> list of D's (key, value) pairs, as 2-tuples
iteritems(...)
D.iteritems() -> an iterator over the (key, value) items of D
iterkeys(...)
D.iterkeys() -> an iterator over the keys of D
itervalues(...)
D.itervalues() -> an iterator over the values of D
keys(...)
D.keys() -> list of D's keys
pop(...)
D.pop(k[,d]) -> v, remove specified key and return the corresponding value.
If key is not found, d is returned if given, otherwise KeyError is raised
popitem(...)
D.popitem() -> (k, v), remove and return some (key, value) pair as a
2-tuple; but raise KeyError if D is empty.
setdefault(...)
D.setdefault(k[,d]) -> D.get(k,d), also set D[k]=d if k not in D
update(...)
D.update([E, ]**F) -> None.  Update D from dict/iterable E and F.
If E present and has a .keys() method, does:     for k in E: D[k] = E[k]
If E present and lacks .keys() method, does:     for (k, v) in E: D[k] = v
In either case, this is followed by: for k in F: D[k] = F[k]
values(...)
D.values() -> list of D's values
viewitems(...)
D.viewitems() -> a set-like object providing a view on D's items
viewkeys(...)
D.viewkeys() -> a set-like object providing a view on D's keys
viewvalues(...)
D.viewvalues() -> an object providing a view on D's values

Data and other attributes inherited from __builtin__.dict:
__hash__ = None
__new__ = <built-in method __new__ of type object>
T.__new__(S, ...) -> a new object with type S, a subtype of T

 
class CMASolutionDict(SolutionDict)
    
Method resolution order:
CMASolutionDict
SolutionDict
DerivedDictBase
_abcoll.MutableMapping
_abcoll.Mapping
_abcoll.Sized
_abcoll.Iterable
_abcoll.Container
__builtin__.object

Methods defined here:
__init__(self, *args, **kwargs)
insert(self, key, geno=None, iteration=None, fitness=None, value=None)
insert an entry with key ``key`` and value
``value if value is not None else {'geno':key}`` and
``self[key]['kwarg'] = kwarg if kwarg is not None`` for the further kwargs.

Data and other attributes defined here:
__abstractmethods__ = frozenset([])

Methods inherited from SolutionDict:
__delitem__(self, key)
remove only most current key-entry
__getitem__(self, key)
defines self[key]
__setitem__(self, key, value)
defines self[key] = value
key(self, x)
truncate(self, max_len, min_iter)

Methods inherited from DerivedDictBase:
__contains__(self, key)
__iter__(self)
__len__(self)

Methods inherited from _abcoll.MutableMapping:
clear(self)
D.clear() -> None.  Remove all items from D.
pop(self, key, default=<object object>)
D.pop(k[,d]) -> v, remove specified key and return the corresponding value.
If key is not found, d is returned if given, otherwise KeyError is raised.
popitem(self)
D.popitem() -> (k, v), remove and return some (key, value) pair
as a 2-tuple; but raise KeyError if D is empty.
setdefault(self, key, default=None)
D.setdefault(k[,d]) -> D.get(k,d), also set D[k]=d if k not in D
update(*args, **kwds)
D.update([E, ]**F) -> None.  Update D from mapping/iterable E and F.
If E present and has a .keys() method, does:     for k in E: D[k] = E[k]
If E present and lacks .keys() method, does:     for (k, v) in E: D[k] = v
In either case, this is followed by: for k, v in F.items(): D[k] = v

Methods inherited from _abcoll.Mapping:
__eq__(self, other)
__ne__(self, other)
get(self, key, default=None)
D.get(k[,d]) -> D[k] if k in D, else d.  d defaults to None.
items(self)
D.items() -> list of D's (key, value) pairs, as 2-tuples
iteritems(self)
D.iteritems() -> an iterator over the (key, value) items of D
iterkeys(self)
D.iterkeys() -> an iterator over the keys of D
itervalues(self)
D.itervalues() -> an iterator over the values of D
keys(self)
D.keys() -> list of D's keys
values(self)
D.values() -> list of D's values

Data and other attributes inherited from _abcoll.Mapping:
__hash__ = None

Class methods inherited from _abcoll.Sized:
__subclasshook__(cls, C) from abc.ABCMeta

Data descriptors inherited from _abcoll.Sized:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

Data and other attributes inherited from _abcoll.Sized:
__metaclass__ = <class 'abc.ABCMeta'>
Metaclass for defining Abstract Base Classes (ABCs).
 
Use this metaclass to create an ABC.  An ABC can be subclassed
directly, and then acts as a mix-in class.  You can also register
unrelated concrete classes (even built-in classes) and unrelated
ABCs as 'virtual subclasses' -- these and their descendants will
be considered subclasses of the registering ABC by the built-in
issubclass() function, but the registering ABC won't show up in
their MRO (Method Resolution Order) nor will method
implementations defined by the registering ABC be callable (not
even via super()).

 
class DerivedDictBase(_abcoll.MutableMapping)
    for conveniently adding "features" to a dictionary. The actual
dictionary is in ``self.data``. Copy-paste
and modify setitem, getitem, and delitem, if necessary.
 
Details: This is the clean way to subclass build-in dict.
 
 
Method resolution order:
DerivedDictBase
_abcoll.MutableMapping
_abcoll.Mapping
_abcoll.Sized
_abcoll.Iterable
_abcoll.Container
__builtin__.object

Methods defined here:
__contains__(self, key)
__delitem__(self, key)
__getitem__(self, key)
defines self[key]
__init__(self, *args, **kwargs)
__iter__(self)
__len__(self)
__setitem__(self, key, value)
defines self[key] = value

Data and other attributes defined here:
__abstractmethods__ = frozenset([])

Methods inherited from _abcoll.MutableMapping:
clear(self)
D.clear() -> None.  Remove all items from D.
pop(self, key, default=<object object>)
D.pop(k[,d]) -> v, remove specified key and return the corresponding value.
If key is not found, d is returned if given, otherwise KeyError is raised.
popitem(self)
D.popitem() -> (k, v), remove and return some (key, value) pair
as a 2-tuple; but raise KeyError if D is empty.
setdefault(self, key, default=None)
D.setdefault(k[,d]) -> D.get(k,d), also set D[k]=d if k not in D
update(*args, **kwds)
D.update([E, ]**F) -> None.  Update D from mapping/iterable E and F.
If E present and has a .keys() method, does:     for k in E: D[k] = E[k]
If E present and lacks .keys() method, does:     for (k, v) in E: D[k] = v
In either case, this is followed by: for k, v in F.items(): D[k] = v

Methods inherited from _abcoll.Mapping:
__eq__(self, other)
__ne__(self, other)
get(self, key, default=None)
D.get(k[,d]) -> D[k] if k in D, else d.  d defaults to None.
items(self)
D.items() -> list of D's (key, value) pairs, as 2-tuples
iteritems(self)
D.iteritems() -> an iterator over the (key, value) items of D
iterkeys(self)
D.iterkeys() -> an iterator over the keys of D
itervalues(self)
D.itervalues() -> an iterator over the values of D
keys(self)
D.keys() -> list of D's keys
values(self)
D.values() -> list of D's values

Data and other attributes inherited from _abcoll.Mapping:
__hash__ = None

Class methods inherited from _abcoll.Sized:
__subclasshook__(cls, C) from abc.ABCMeta

Data descriptors inherited from _abcoll.Sized:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

Data and other attributes inherited from _abcoll.Sized:
__metaclass__ = <class 'abc.ABCMeta'>
Metaclass for defining Abstract Base Classes (ABCs).
 
Use this metaclass to create an ABC.  An ABC can be subclassed
directly, and then acts as a mix-in class.  You can also register
unrelated concrete classes (even built-in classes) and unrelated
ABCs as 'virtual subclasses' -- these and their descendants will
be considered subclasses of the registering ABC by the built-in
issubclass() function, but the registering ABC won't show up in
their MRO (Method Resolution Order) nor will method
implementations defined by the registering ABC be callable (not
even via super()).

 
class ElapsedTime(__builtin__.object)
    using ``time.clock`` with overflow handling to measure CPU time.
 
Example:
 
>>> clock = ElapsedTime()  # clock starts here
>>> t1 = clock()  # get elapsed CPU time
 
Details: 32-bit C overflows after int(2**32/1e6) == 4294s about 72 min
 
  Methods defined here:
__call__(self)
__init__(self)
reset = __init__(self)

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class GenoPheno(__builtin__.object)
    Genotype-phenotype transformation.
 
Method `pheno` provides the transformation from geno- to phenotype,
that is from the internal representation to the representation used
in the objective function. Method `geno` provides the "inverse" pheno-
to genotype transformation. The geno-phenotype transformation comprises,
in this order:
 
   - insert fixed variables (with the phenotypic and therefore quite
     possibly "wrong" values)
   - affine linear transformation (first scaling then shift)
   - user-defined transformation
   - repair (e.g. into feasible domain due to boundaries)
   - assign fixed variables their original phenotypic value
 
By default all transformations are the identity. The repair is only applied,
if the transformation is given as argument to the method `pheno`.
 
``geno`` is only necessary, if solutions have been injected.
 
  Methods defined here:
__init__(self, dim, scaling=None, typical_x=None, fixed_values=None, tf=None)
return `GenoPheno` instance with phenotypic dimension `dim`.
 
Keyword Arguments
-----------------
    `scaling`
        the diagonal of a scaling transformation matrix, multipliers
        in the genotyp-phenotyp transformation, see `typical_x`
    `typical_x`
        ``pheno = scaling*geno + typical_x``
    `fixed_values`
        a dictionary of variable indices and values, like ``{0:2.0, 2:1.1}``,
        that are not subject to change, negative indices are ignored
        (they act like incommenting the index), values are phenotypic
        values.
    `tf`
        list of two user-defined transformation functions, or `None`.
 
        ``tf[0]`` is a function that transforms the internal representation
        as used by the optimizer into a solution as used by the
        objective function. ``tf[1]`` does the back-transformation.
        For example::
 
            tf_0 = lambda x: [xi**2 for xi in x]
            tf_1 = lambda x: [abs(xi)**0.5 fox xi in x]
 
        or "equivalently" without the `lambda` construct::
 
            def tf_0(x):
                return [xi**2 for xi in x]
            def tf_1(x):
                return [abs(xi)**0.5 fox xi in x]
 
        ``tf=[tf_0, tf_1]`` is a reasonable way to guaranty that only positive
        values are used in the objective function.
 
Details
-------
If ``tf_0`` is not the identity and ``tf_1`` is ommitted,
the genotype of ``x0`` cannot be computed consistently and
"injection" of phenotypic solutions is likely to lead to
unexpected results.
geno(self, y, from_bounds=None, copy_if_changed=True, copy_always=False, repair=None, archive=None)
maps the phenotypic input argument into the genotypic space,
that is, computes essentially the inverse of ``pheno``.
 
By default a copy is made only to prevent to modify ``y``.
 
The inverse of the user-defined transformation (if any)
is only needed if external solutions are injected, it is not
applied to the initial solution x0.
 
Details
=======
``geno`` searches first in ``archive`` for the genotype of
``y`` and returns the found value, typically unrepaired.
Otherwise, first ``from_bounds`` is applied, to revert a
projection into the bound domain (if necessary) and ``pheno``
is reverted. ``repair`` is applied last, and is usually the
method ``CMAEvolutionStrategy.repair_genotype`` that limits the
Mahalanobis norm of ``geno(y) - mean``.
pheno(self, x, into_bounds=None, copy=True, copy_always=False, archive=None, iteration=None)
maps the genotypic input argument into the phenotypic space,
see help for class `GenoPheno`
 
Details
-------
If ``copy``, values from ``x`` are copied if changed under the transformation.

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
Mh = class MathHelperFunctions(__builtin__.object)
    static convenience math helper functions, if the function name
is preceded with an "a", a numpy array is returned
 
  Static methods defined here:
aclamp(x, upper)
amax(vec, vec_or_scalar)
amin(vec_or_scalar, vec_or_scalar2)
aminmax(val, min_val, max_val)
apos(x, lower=0)
clips argument (scalar or array) from below at lower
cauchy_with_variance_one()
equals_approximately(a, b, eps=1e-12)
expms(A, eig=<function eigh>)
matrix exponential for a symmetric matrix
max(vec, vec_or_scalar)
min(a, b)
minmax(val, min_val, max_val)
norm(vec, expo=2)
prctile(data, p_vals=[0, 25, 50, 75, 100], sorted_=False)
``prctile(data, 50)`` returns the median, but p_vals can
also be a sequence.
 
Provides for small samples better values than matplotlib.mlab.prctile,
however also slower.
sround(nb)
return stochastic round: floor(nb) + (rand()<remainder(nb))
standard_finite_cauchy(size=1)
vequals_approximately(a, b, eps=1e-12)

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class Misc(__builtin__.object)
     Static methods defined here:
eig(C)
eigendecomposition of a symmetric matrix, much slower than
`numpy.linalg.eigh`, return ``(EVals, Basis)``, the eigenvalues
and an orthonormal basis of the corresponding eigenvectors, where
 
    ``Basis[i]``
        the i-th row of ``Basis``
    columns of ``Basis``, ``[Basis[j][i] for j in range(len(Basis))]``
        the i-th eigenvector with eigenvalue ``EVals[i]``
likelihood(x, m=None, Cinv=None, sigma=1, detC=None)
return likelihood of x for the normal density N(m, sigma**2 * Cinv**-1)
loglikelihood(self, x, previous=False)
return log-likelihood of `x` regarding the current sample distribution

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

Data and other attributes defined here:
MathHelperFunctions = <class 'cma.MathHelperFunctions'>
static convenience math helper functions, if the function name
is preceded with an "a", a numpy array is returned

 
class NoiseHandler(__builtin__.object)
    Noise handling according to [Hansen et al 2009, A Method for
Handling Uncertainty in Evolutionary Optimization...]
 
The interface of this class is yet versatile and subject to changes.
 
The noise handling follows closely [Hansen et al 2009] in the
measurement part, but the implemented treatment is slightly
different: for ``noiseS > 0``, ``evaluations`` (time) and sigma are
increased by ``alpha``. For ``noiseS < 0``, ``evaluations`` (time)
is decreased by ``alpha**(1/4)``.
 
The (second) parameter ``evaluations`` defines the maximal number
of evaluations for a single fitness computation. If it is a list,
the smallest element defines the minimal number and if the list has
three elements, the median value is the start value for
``evaluations``.
 
``NoiseHandler`` serves to control the noise via steps-size
increase and number of re-evaluations, for example via ``fmin`` or
with ``ask_and_eval()``.
 
Examples
--------
Minimal example together with `fmin` on a non-noisy function:
 
>>> import cma
>>> cma.fmin(cma.felli, 7 * [1], 1, noise_handler=cma.NoiseHandler(7))
 
in dimension 7 (which needs to be given tice). More verbose example
in the optimization loop with a noisy function defined in ``func``:
 
>>> import cma, numpy as np
>>> func = lambda x: cma.fcts.sphere(x) * (1 + 4 * np.random.randn() / len(x))  # cma.Fcts.noisysphere
>>> es = cma.CMAEvolutionStrategy(np.ones(10), 1)
>>> nh = cma.NoiseHandler(es.N, maxevals=[1, 1, 30])
>>> while not es.stop():
...     X, fit_vals = es.ask_and_eval(func, evaluations=nh.evaluations)
...     es.tell(X, fit_vals)  # prepare for next iteration
...     es.sigma *= nh(X, fit_vals, func, es.ask)  # see method __call__
...     es.countevals += nh.evaluations_just_done  # this is a hack, not important though
...     es.logger.add(more_data = [nh.evaluations, nh.noiseS])  # add a data point
...     es.disp()
...     # nh.maxevals = ...  it might be useful to start with smaller values and then increase
>>> print(es.stop())
>>> print(es.result()[-2])  # take mean value, the best solution is totally off
>>> assert sum(es.result()[-2]**2) < 1e-9
>>> print(X[np.argmin(fit_vals)])  # not bad, but probably worse than the mean
>>> # es.logger.plot()
 
 
The command ``logger.plot()`` will plot the logged data.
 
The noise options of `fmin()` control a `NoiseHandler` instance
similar to this example. The command ``cma.CMAOptions('noise')``
lists in effect the parameters of `__init__` apart from
``aggregate``.
 
Details
-------
The parameters reevals, theta, c_s, and alpha_t are set differently
than in the original publication, see method `__init__()`. For a
very small population size, say popsize <= 5, the measurement
technique based on rank changes is likely to fail.
 
Missing Features
----------------
In case no noise is found, ``self.lam_reeval`` should be adaptive
and get at least as low as 1 (however the possible savings from this
are rather limited). Another option might be to decide during the
first call by a quantitative analysis of fitness values whether
``lam_reeval`` is set to zero. More generally, an automatic noise
mode detection might also set the covariance matrix learning rates
to smaller values.
 
:See: `fmin()`, `CMAEvolutionStrategy.ask_and_eval()`
 
  Methods defined here:
__call__(self, X, fit, func, ask=None, args=())
proceed with noise measurement, set anew attributes ``evaluations``
(proposed number of evaluations to "treat" noise) and ``evaluations_just_done``
and return a factor for increasing sigma.
 
Parameters
----------
    `X`
        a list/sequence/vector of solutions
    `fit`
        the respective list of function values
    `func`
        the objective function, ``fit[i]`` corresponds to ``func(X[i], *args)``
    `ask`
        a method to generate a new, slightly disturbed solution. The argument
        is (only) mandatory if ``epsilon`` is not zero, see `__init__()`.
    `args`
        optional additional arguments to `func`
 
Details
-------
Calls the methods ``reeval()``, ``update_measure()`` and ``treat()`` in this order.
``self.evaluations`` is adapted within the method `treat()`.
__init__(self, N, maxevals=[1, 1, 1], aggregate=<function median>, reevals=None, epsilon=1e-07, parallel=False)
parameters are
 
`N`
    dimension, (only) necessary to adjust the internal
    "alpha"-parameters
`maxevals`
    maximal value for ``self.evaluations``, where
    ``self.evaluations`` function calls are aggregated for
    noise treatment. With ``maxevals == 0`` the noise
    handler is (temporarily) "switched off". If `maxevals`
    is a list, min value and (for >2 elements) median are
    used to define minimal and initial value of
    ``self.evaluations``. Choosing ``maxevals > 1`` is only
    reasonable, if also the original ``fit`` values (that
    are passed to `__call__`) are computed by aggregation of
    ``self.evaluations`` values (otherwise the values are
    not comparable), as it is done within `fmin()`.
`aggregate`
    function to aggregate single f-values to a 'fitness', e.g.
    ``np.median``.
`reevals`
    number of solutions to be reevaluated for noise
    measurement, can be a float, by default set to ``2 +
    popsize/20``, where ``popsize = len(fit)`` in
    ``__call__``. zero switches noise handling off.
`epsilon`
    multiplier for perturbation of the reevaluated solutions
`parallel`
    a single f-call with all resampled solutions
 
:See: `fmin()`, `CMAOptions`, `CMAEvolutionStrategy.ask_and_eval()`
get_evaluations(self)
return ``self.evaluations``, the number of evalutions to get a single fitness measurement
indices(self, fit)
return the set of indices to be reevaluated for noise
measurement.
 
Given the first values are the earliest, this is a useful policy also
with a time changing objective.
reeval(self, X, fit, func, ask, args=())
store two fitness lists, `fit` and ``fitre`` reevaluating some
solutions in `X`.
``self.evaluations`` evaluations are done for each reevaluated
fitness value.
See `__call__()`, where `reeval()` is called.
treat(self)
adapt self.evaluations depending on the current measurement value
and return ``sigma_fac in (1.0, self.alphasigma)``
update_measure(self)
updated noise level measure using two fitness lists ``self.fit`` and
``self.fitre``, return ``self.noiseS, all_individual_measures``.
 
Assumes that `self.idx` contains the indices where the fitness
lists differ

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class OOOptimizer(__builtin__.object)
    "abstract" base class for an Object Oriented Optimizer interface.
 
 Relevant methods are `__init__`, `ask`, `tell`, `stop`, `result`,
 and `optimize`. Only `optimize` is fully implemented in this base
 class.
 
Examples
--------
All examples minimize the function `elli`, the output is not shown.
(A preferred environment to execute all examples is ``ipython`` in
``%pylab`` mode.)
 
First we need::
 
    from cma import CMAEvolutionStrategy
    # CMAEvolutionStrategy derives from the OOOptimizer class
    felli = lambda x: sum(1e3**((i-1.)/(len(x)-1.)*x[i])**2 for i in range(len(x)))
 
The shortest example uses the inherited method
`OOOptimizer.optimize()`::
 
    es = CMAEvolutionStrategy(8 * [0.1], 0.5).optimize(felli)
 
The input parameters to `CMAEvolutionStrategy` are specific to this
inherited class. The remaining functionality is based on interface
defined by `OOOptimizer`. We might have a look at the result::
 
    print(es.result()[0])  # best solution and
    print(es.result()[1])  # its function value
 
In order to display more exciting output we do::
 
    es.logger.plot()  # if matplotlib is available
 
Virtually the same example can be written with an explicit loop
instead of using `optimize()`. This gives the necessary insight into
the `OOOptimizer` class interface and entire control over the
iteration loop::
 
    optim = CMAEvolutionStrategy(9 * [0.5], 0.3)
    # a new CMAEvolutionStrategy instance
 
    # this loop resembles optimize()
    while not optim.stop():  # iterate
        X = optim.ask()      # get candidate solutions
        f = [felli(x) for x in X]  # evaluate solutions
        #  in case do something else that needs to be done
        optim.tell(X, f)     # do all the real "update" work
        optim.disp(20)       # display info every 20th iteration
        optim.logger.add()   # log another "data line"
 
    # final output
    print('termination by', optim.stop())
    print('best f-value =', optim.result()[1])
    print('best solution =', optim.result()[0])
    optim.logger.plot()  # if matplotlib is available
 
Details
-------
Most of the work is done in the method `tell(...)`. The method
`result()` returns more useful output.
 
  Methods defined here:
__init__(self, xstart, **more_args)
``xstart`` is a mandatory argument
ask(self, gradf=None, **more_args)
abstract method, AKA "get" or "sample_distribution", deliver
new candidate solution(s), a list of "vectors"
disp(self, modulo=None)
abstract method, display some iteration infos if
``self.iteration_counter % modulo == 0``
initialize(self)
(re-)set to the initial state
optimize(self, objective_fct, iterations=None, min_iterations=1, args=(), verb_disp=None, logger=None, call_back=None)
find minimizer of `objective_fct`.
 
CAVEAT: the return value for `optimize` has changed to ``self``.
 
Arguments
---------
 
    `objective_fct`
        function be to minimized
    `iterations`
        number of (maximal) iterations, while ``not self.stop()``
    `min_iterations`
        minimal number of iterations, even if ``not self.stop()``
    `args`
        arguments passed to `objective_fct`
    `verb_disp`
        print to screen every `verb_disp` iteration, if ``None``
        the value from ``self.logger`` is "inherited", if
        available.
    ``logger``
        a `BaseDataLogger` instance, which must be compatible
        with the type of ``self``.
    ``call_back``
        call back function called like ``call_back(self)`` or
        a list of call back functions.
 
``return self``, that is, the `OOOptimizer` instance.
 
Example
-------
>>> import cma
>>> es = cma.CMAEvolutionStrategy(7 * [0.1], 0.5
...              ).optimize(cma.fcts.rosen, verb_disp=100)
(4_w,9)-CMA-ES (mu_w=2.8,w_1=49%) in dimension 7 (seed=630721393)
Iterat #Fevals   function value    axis ratio  sigma   minstd maxstd min:sec
    1       9 3.163954777181882e+01 1.0e+00 4.12e-01  4e-01  4e-01 0:0.0
    2      18 3.299006223906629e+01 1.0e+00 3.60e-01  3e-01  4e-01 0:0.0
    3      27 1.389129389866704e+01 1.1e+00 3.18e-01  3e-01  3e-01 0:0.0
  100     900 2.494847340045985e+00 8.6e+00 5.03e-02  2e-02  5e-02 0:0.3
  200    1800 3.428234862999135e-01 1.7e+01 3.77e-02  6e-03  3e-02 0:0.5
  300    2700 3.216640032470860e-04 5.6e+01 6.62e-03  4e-04  9e-03 0:0.8
  400    3600 6.155215286199821e-12 6.6e+01 7.44e-06  1e-07  4e-06 0:1.1
  438    3942 1.187372505161762e-14 6.0e+01 3.27e-07  4e-09  9e-08 0:1.2
  438    3942 1.187372505161762e-14 6.0e+01 3.27e-07  4e-09  9e-08 0:1.2
('termination by', {'tolfun': 1e-11})
('best f-value =', 1.1189867885201275e-14)
('solution =', array([ 1.        ,  1.        ,  1.        ,  0.99999999,  0.99999998,
        0.99999996,  0.99999992]))
>>> print(es.result()[0])
array([ 1.          1.          1.          0.99999999  0.99999998  0.99999996
  0.99999992])
result(self)
abstract method, return ``(x, f(x), ...)``, that is, the
minimizer, its function value, ...
stop(self)
abstract method, return satisfied termination conditions in
a dictionary like ``{'termination reason': value, ...}``,
for example ``{'tolfun': 1e-12}``, or the empty dictionary ``{}``.
The implementation of `stop()` should prevent an infinite
loop.
tell(self, solutions, function_values)
abstract method, AKA "update", pass f-values and prepare for
next iteration

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class Rotation(__builtin__.object)
    Rotation class that implements an orthogonal linear transformation,
one for each dimension.
 
By default reach ``Rotation`` instance provides a different "random"
but fixed rotation. This class is used to implement non-separable
test functions, most conveniently via `FFWrapper.RotatedFitness`.
 
Example:
 
>>> import cma, numpy as np
>>> R = cma.Rotation()
>>> R2 = cma.Rotation() # another rotation
>>> x = np.array((1,2,3))
>>> print(R(R(x), inverse=1))
[ 1.  2.  3.]
 
See: `FFWrapper.RotatedFitness`
 
  Methods defined here:
__call__(self, x, inverse=False)
Rotates the input array `x` with a fixed rotation matrix
(``self.dicMatrices['str(len(x))']``)
__init__(self, seed=None)
by default a random but fixed rotation, different for each instance

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

Data and other attributes defined here:
dicMatrices = {}

 
class Sections(__builtin__.object)
    plot sections through an objective function.
 
A first rational thing to do, when facing an (expensive)
application. By default 6 points in each coordinate are evaluated.
This class is still experimental.
 
Examples
--------
 
>>> import cma, numpy as np
>>> s = cma.Sections(cma.Fcts.rosen, np.zeros(3)).do(plot=False)
>>> s.do(plot=False)  # evaluate the same points again, i.e. check for noise
>> try:
...     s.plot()
... except:
...     print('plotting failed: matplotlib.pyplot package missing?')
 
Details
-------
Data are saved after each function call during `do()`. The filename
is attribute ``name`` and by default ``str(func)``, see `__init__()`.
 
A random (orthogonal) basis can be generated with
``cma.Rotation()(np.eye(3))``.
 
CAVEAT: The default name is unique in the function name, but it
should be unique in all parameters of `__init__()` but `plot_cmd`
and `load`. If, for example, a different basis is chosen, either
the name must be changed or the ``.pkl`` file containing the
previous data must first be renamed or deleted.
 
``s.res`` is a dictionary with an entry for each "coordinate" ``i``
and with an entry ``'x'``, the middle point. Each entry ``i`` is
again a dictionary with keys being different dx values and the
value being a sequence of f-values. For example ``s.res[2][0.1] ==
[0.01, 0.01]``, which is generated using the difference vector ``s
.basis[2]`` like
 
``s.res[2][dx] += func(s.res['x'] + dx * s.basis[2])``.
 
:See: `__init__()`
 
  Methods defined here:
__init__(self, func, x, args=(), basis=None, name=None, plot_cmd=<function plot>, load=True)
Parameters
----------
    `func`
        objective function
    `x`
        point in search space, middle point of the sections
    `args`
        arguments passed to `func`
    `basis`
        evaluated points are ``func(x + locations[j] * basis[i])
        for i in len(basis) for j in len(locations)``,
        see `do()`
    `name`
        filename where to save the result
    `plot_cmd`
        command used to plot the data, typically matplotlib pyplots `plot` or `semilogy`
    `load`
        load previous data from file ``str(func) + '.pkl'``
do(self, repetitions=1, locations=array([-0.5, -0.3, -0.1, 0.1, 0.3, 0.5]), plot=True)
generates, plots and saves function values ``func(y)``,
where ``y`` is 'close' to `x` (see `__init__()`). The data are stored in
the ``res`` attribute and the class instance is saved in a file
with (the weired) name ``str(func)``.
 
Parameters
----------
    `repetitions`
        for each point, only for noisy functions is >1 useful. For
        ``repetitions==0`` only already generated data are plotted.
    `locations`
        coordinated wise deviations from the middle point given in `__init__`
flattened(self)
return flattened data ``(x, f)`` such that for the sweep through
coordinate ``i`` we have for data point ``j`` that ``f[i][j] == func(x[i][j])``
load(self, name=None)
load from file
plot(self, plot_cmd=None, tf=<function <lambda>>)
plot the data we have, return ``self``
save(self, name=None)
save to file

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class SolutionDict(DerivedDictBase)
    dictionary with computation of an hash key.
 
The hash key is generated from the inserted solution and a stack of
previously inserted same solutions is provided. Each entry is meant
to store additional information related to the solution.
 
    >>> import cma, numpy as np
    >>> d = cma.SolutionDict()
    >>> x = np.array([1,2,4])
    >>> d[x] = {'f': sum(x**2), 'iteration': 1}
    >>> assert d[x]['iteration'] == 1
    >>> assert d.get(x) == (d[x] if d.key(x) in d.keys() else None)
 
TODO: data_with_same_key behaves like a stack (see setitem and
delitem), but rather should behave like a queue?! A queue is less
consistent with the operation self[key] = ..., if
self.data_with_same_key[key] is not empty.
 
TODO: iteration key is used to clean up without error management
 
 
Method resolution order:
SolutionDict
DerivedDictBase
_abcoll.MutableMapping
_abcoll.Mapping
_abcoll.Sized
_abcoll.Iterable
_abcoll.Container
__builtin__.object

Methods defined here:
__delitem__(self, key)
remove only most current key-entry
__getitem__(self, key)
defines self[key]
__init__(self, *args, **kwargs)
__setitem__(self, key, value)
defines self[key] = value
key(self, x)
truncate(self, max_len, min_iter)

Data and other attributes defined here:
__abstractmethods__ = frozenset([])

Methods inherited from DerivedDictBase:
__contains__(self, key)
__iter__(self)
__len__(self)

Methods inherited from _abcoll.MutableMapping:
clear(self)
D.clear() -> None.  Remove all items from D.
pop(self, key, default=<object object>)
D.pop(k[,d]) -> v, remove specified key and return the corresponding value.
If key is not found, d is returned if given, otherwise KeyError is raised.
popitem(self)
D.popitem() -> (k, v), remove and return some (key, value) pair
as a 2-tuple; but raise KeyError if D is empty.
setdefault(self, key, default=None)
D.setdefault(k[,d]) -> D.get(k,d), also set D[k]=d if k not in D
update(*args, **kwds)
D.update([E, ]**F) -> None.  Update D from mapping/iterable E and F.
If E present and has a .keys() method, does:     for k in E: D[k] = E[k]
If E present and lacks .keys() method, does:     for (k, v) in E: D[k] = v
In either case, this is followed by: for k, v in F.items(): D[k] = v

Methods inherited from _abcoll.Mapping:
__eq__(self, other)
__ne__(self, other)
get(self, key, default=None)
D.get(k[,d]) -> D[k] if k in D, else d.  d defaults to None.
items(self)
D.items() -> list of D's (key, value) pairs, as 2-tuples
iteritems(self)
D.iteritems() -> an iterator over the (key, value) items of D
iterkeys(self)
D.iterkeys() -> an iterator over the keys of D
itervalues(self)
D.itervalues() -> an iterator over the values of D
keys(self)
D.keys() -> list of D's keys
values(self)
D.values() -> list of D's values

Data and other attributes inherited from _abcoll.Mapping:
__hash__ = None

Class methods inherited from _abcoll.Sized:
__subclasshook__(cls, C) from abc.ABCMeta

Data descriptors inherited from _abcoll.Sized:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

Data and other attributes inherited from _abcoll.Sized:
__metaclass__ = <class 'abc.ABCMeta'>
Metaclass for defining Abstract Base Classes (ABCs).
 
Use this metaclass to create an ABC.  An ABC can be subclassed
directly, and then acts as a mix-in class.  You can also register
unrelated concrete classes (even built-in classes) and unrelated
ABCs as 'virtual subclasses' -- these and their descendants will
be considered subclasses of the registering ABC by the built-in
issubclass() function, but the registering ABC won't show up in
their MRO (Method Resolution Order) nor will method
implementations defined by the registering ABC be callable (not
even via super()).

 
Functions
       
closefig = close(*args)
Close a figure window.
 
``close()`` by itself closes the current figure
 
``close(h)`` where *h* is a :class:`Figure` instance, closes that figure
 
``close(num)`` closes figure number *num*
 
``close(name)`` where *name* is a string, closes figure with that label
 
``close('all')`` closes all the figure windows
disp(name=None, idx=None)
displays selected data from (files written by) the class `CMADataLogger`.
 
The call ``cma.disp(name, idx)`` is a shortcut for ``cma.CMADataLogger(name).disp(idx)``.
 
Arguments
---------
    `name`
        name of the logger, filename prefix, `None` evaluates to
        the default ``'outcmaes'``
    `idx`
        indices corresponding to rows in the data file; by
        default the first five, then every 100-th, and the last
        10 rows. Too large index values are removed.
 
Examples
--------
::
 
   import cma, numpy
   # assume some data are available from previous runs
   cma.disp(None,numpy.r_[0,-1])  # first and last
   cma.disp(None,numpy.r_[0:1e9:100,-1]) # every 100-th and last
   cma.disp(idx=numpy.r_[0,-10:0]) # first and ten last
   cma.disp(idx=numpy.r_[0:1e9:1e3,-10:0])
 
:See: `CMADataLogger.disp()`
felli(x)
unbound test function, needed to test multiprocessor
fmin(objective_function, x0, sigma0, options=None, args=(), gradf=None, restarts=0, restart_from_best=u'False', incpopsize=2, eval_initial_x=False, noise_handler=None, noise_change_sigma_exponent=1, noise_kappa_exponent=0, bipop=False)
functional interface to the stochastic optimizer CMA-ES
for non-convex function minimization.
 
Calling Sequences
=================
    ``fmin(objective_function, x0, sigma0)``
        minimizes `objective_function` starting at `x0` and with standard deviation
        `sigma0` (step-size)
    ``fmin(objective_function, x0, sigma0, options={'ftarget': 1e-5})``
        minimizes `objective_function` up to target function value 1e-5, which
        is typically useful for benchmarking.
    ``fmin(objective_function, x0, sigma0, args=('f',))``
        minimizes `objective_function` called with an additional argument ``'f'``.
    ``fmin(objective_function, x0, sigma0, options={'ftarget':1e-5, 'popsize':40})``
        uses additional options ``ftarget`` and ``popsize``
    ``fmin(objective_function, esobj, None, options={'maxfevals': 1e5})``
        uses the `CMAEvolutionStrategyobject instance `esobj` to optimize
        `objective_function`, similar to `esobj.optimize()`.
 
Arguments
=========
    `objective_function`
        function to be minimized. Called as ``objective_function(x,
        *args)``. `x` is a one-dimensional `numpy.ndarray`.
        `objective_function` can return `numpy.NaN`,
        which is interpreted as outright rejection of solution `x`
        and invokes an immediate resampling and (re-)evaluation
        of a new solution not counting as function evaluation.
    `x0`
        list or `numpy.ndarray`, initial guess of minimum solution
        before the application of the geno-phenotype transformation
        according to the ``transformation`` option.  It can also be
        a string holding a Python expression that is evaluated
        to yield the initial guess - this is important in case
        restarts are performed so that they start from different
        places.  Otherwise `x0` can also be a `cma.CMAEvolutionStrategy`
        object instance, in that case `sigma0` can be ``None``.
    `sigma0`
        scalar, initial standard deviation in each coordinate.
        `sigma0` should be about 1/4th of the search domain width
        (where the optimum is to be expected). The variables in
        `objective_function` should be scaled such that they
        presumably have similar sensitivity.
        See also option `scaling_of_variables`.
    `options`
        a dictionary with additional options passed to the constructor
        of class ``CMAEvolutionStrategy``, see ``cma.CMAOptions()``
        for a list of available options.
    ``args=()``
        arguments to be used to call the `objective_function`
    ``gradf``
        gradient of f, where ``len(gradf(x, *args)) == len(x)``.
        `gradf` is called once in each iteration if
        ``gradf is not None``.
    ``restarts=0``
        number of restarts with increasing population size, see also
        parameter `incpopsize`, implementing the IPOP-CMA-ES restart
        strategy, see also parameter `bipop`; to restart from
        different points (recommended), pass `x0` as a string.
    ``restart_from_best=False``
        which point to restart from
    ``incpopsize=2``
        multiplier for increasing the population size `popsize` before
        each restart
    ``eval_initial_x=None``
        evaluate initial solution, for `None` only with elitist option
    ``noise_handler=None``
        a ``NoiseHandler`` instance or ``None``, a simple usecase is
        ``cma.fmin(f, 6 * [1], 1, noise_handler=cma.NoiseHandler(6))``
        see ``help(cma.NoiseHandler)``.
    ``noise_change_sigma_exponent=1``
        exponent for sigma increment for additional noise treatment
    ``noise_evaluations_as_kappa``
        instead of applying reevaluations, the "number of evaluations"
        is (ab)used as scaling factor kappa (experimental).
    ``bipop``
        if True, run as BIPOP-CMA-ES; BIPOP is a special restart
        strategy switching between two population sizings - small
        (like the default CMA, but with more focused search) and
        large (progressively increased as in IPOP). This makes the
        algorithm perform well both on functions with many regularly
        or irregularly arranged local optima (the latter by frequently
        restarting with small populations).  For the `bipop` parameter
        to actually take effect, also select non-zero number of
        (IPOP) restarts; the recommended setting is ``restarts<=9``
        and `x0` passed as a string.  Note that small-population
        restarts do not count into the total restart count.
 
Optional Arguments
==================
All values in the `options` dictionary are evaluated if they are of
type `str`, besides `verb_filenameprefix`, see class `CMAOptions` for
details. The full list is available via ``cma.CMAOptions()``.
 
>>> import cma
>>> cma.CMAOptions()
 
Subsets of options can be displayed, for example like
``cma.CMAOptions('tol')``, or ``cma.CMAOptions('bound')``,
see also class `CMAOptions`.
 
Return
======
Return the list provided by `CMAEvolutionStrategy.result()` appended
with termination conditions, an `OOOptimizer` and a `BaseDataLogger`::
 
    res = es.result() + (es.stop(), es, logger)
 
where
    - ``res[0]`` (``xopt``) -- best evaluated solution
    - ``res[1]`` (``fopt``) -- respective function value
    - ``res[2]`` (``evalsopt``) -- respective number of function evaluations
    - ``res[3]`` (``evals``) -- number of overall conducted objective function evaluations
    - ``res[4]`` (``iterations``) -- number of overall conducted iterations
    - ``res[5]`` (``xmean``) -- mean of the final sample distribution
    - ``res[6]`` (``stds``) -- effective stds of the final sample distribution
    - ``res[-3]`` (``stop``) -- termination condition(s) in a dictionary
    - ``res[-2]`` (``cmaes``) -- class `CMAEvolutionStrategy` instance
    - ``res[-1]`` (``logger``) -- class `CMADataLogger` instance
 
Details
=======
This function is an interface to the class `CMAEvolutionStrategy`. The
latter class should be used when full control over the iteration loop
of the optimizer is desired.
 
Examples
========
The following example calls `fmin` optimizing the Rosenbrock function
in 10-D with initial solution 0.1 and initial step-size 0.5. The
options are specified for the usage with the `doctest` module.
 
>>> import cma
>>> # cma.CMAOptions()  # returns all possible options
>>> options = {'CMA_diagonal':100, 'seed':1234, 'verb_time':0}
>>>
>>> res = cma.fmin(cma.fcts.rosen, [0.1] * 10, 0.5, options)
(5_w,10)-CMA-ES (mu_w=3.2,w_1=45%) in dimension 10 (seed=1234)
   Covariance matrix is diagonal for 10 iterations (1/ccov=29.0)
Iterat #Fevals   function value     axis ratio  sigma   minstd maxstd min:sec
    1      10 1.264232686260072e+02 1.1e+00 4.40e-01  4e-01  4e-01
    2      20 1.023929748193649e+02 1.1e+00 4.00e-01  4e-01  4e-01
    3      30 1.214724267489674e+02 1.2e+00 3.70e-01  3e-01  4e-01
  100    1000 6.366683525319511e+00 6.2e+00 2.49e-02  9e-03  3e-02
  200    2000 3.347312410388666e+00 1.2e+01 4.52e-02  8e-03  4e-02
  300    3000 1.027509686232270e+00 1.3e+01 2.85e-02  5e-03  2e-02
  400    4000 1.279649321170636e-01 2.3e+01 3.53e-02  3e-03  3e-02
  500    5000 4.302636076186532e-04 4.6e+01 4.78e-03  3e-04  5e-03
  600    6000 6.943669235595049e-11 5.1e+01 5.41e-06  1e-07  4e-06
  650    6500 5.557961334063003e-14 5.4e+01 1.88e-07  4e-09  1e-07
termination on tolfun : 1e-11
final/bestever f-value = 5.55796133406e-14 2.62435631419e-14
mean solution:  [ 1.          1.00000001  1.          1.
    1.          1.00000001  1.00000002  1.00000003 ...]
std deviation: [ 3.9193387e-09  3.7792732e-09  4.0062285e-09  4.6605925e-09
    5.4966188e-09   7.4377745e-09   1.3797207e-08   2.6020765e-08 ...]
>>>
>>> print('best solutions fitness = %f' % (res[1]))
best solutions fitness = 2.62435631419e-14
>>> assert res[1] < 1e-12
 
The above call is pretty much equivalent with the slightly more
verbose call::
 
    es = cma.CMAEvolutionStrategy([0.1] * 10, 0.5,
                options=options).optimize(cma.fcts.rosen)
 
The following example calls `fmin` optimizing the Rastrigin function
in 3-D with random initial solution in [-2,2], initial step-size 0.5
and the BIPOP restart strategy (that progressively increases population).
The options are specified for the usage with the `doctest` module.
 
>>> import cma
>>> # cma.CMAOptions()  # returns all possible options
>>> options = {'seed':12345, 'verb_time':0, 'ftarget': 1e-8}
>>>
>>> res = cma.fmin(cma.fcts.rastrigin, '2. * np.random.rand(3) - 1', 0.5,
...                options, restarts=9, bipop=True)
(3_w,7)-aCMA-ES (mu_w=2.3,w_1=58%) in dimension 3 (seed=12345)
Iterat #Fevals   function value    axis ratio  sigma  minstd maxstd min:sec
    1       7 1.633489455763566e+01 1.0e+00 4.35e-01  4e-01  4e-01
    2      14 9.762462950258016e+00 1.2e+00 4.12e-01  4e-01  4e-01
    3      21 2.461107851413725e+01 1.4e+00 3.78e-01  3e-01  4e-01
  100     700 9.949590571272680e-01 1.7e+00 5.07e-05  3e-07  5e-07
  123     861 9.949590570932969e-01 1.3e+00 3.93e-06  9e-09  1e-08
termination on tolfun=1e-11
final/bestever f-value = 9.949591e-01 9.949591e-01
mean solution: [  9.94958638e-01  -7.19265205e-10   2.09294450e-10]
std deviation: [  8.71497860e-09   8.58994807e-09   9.85585654e-09]
[...]
(4_w,9)-aCMA-ES (mu_w=2.8,w_1=49%) in dimension 3 (seed=12349)
Iterat #Fevals   function value    axis ratio  sigma  minstd maxstd min:sec
    1  5342.0 2.114883315350800e+01 1.0e+00 3.42e-02  3e-02  4e-02
    2  5351.0 1.810102940125502e+01 1.4e+00 3.79e-02  3e-02  4e-02
    3  5360.0 1.340222457448063e+01 1.4e+00 4.58e-02  4e-02  6e-02
   50  5783.0 8.631491965616078e-09 1.6e+00 2.01e-04  8e-06  1e-05
termination on ftarget=1e-08 after 4 restarts
final/bestever f-value = 8.316963e-09 8.316963e-09
mean solution: [ -3.10652459e-06   2.77935436e-06  -4.95444519e-06]
std deviation: [  1.02825265e-05   8.08348144e-06   8.47256408e-06]
 
In either case, the method::
 
    cma.plot();
 
(based on `matplotlib.pyplot`) produces a plot of the run and, if
necessary::
 
    cma.show()
 
shows the plot in a window. Finally::
 
    cma.savefig('myfirstrun')  # savefig from matplotlib.pyplot
 
will save the figure in a png.
 
We can use the gradient like
 
>>> import cma
>>> res = cma.fmin(cma.fcts.rosen, np.zeros(10), 0.1,
...             options = {'ftarget':1e-8,},
...             gradf=cma.fcts.grad_rosen,
...         )
>>> assert cma.fcts.rosen(res[0]) < 1e-8
>>> assert res[2] < 3600  # 1% are > 3300
>>> assert res[3] < 3600  # 1% are > 3300
 
:See: `CMAEvolutionStrategy`, `OOOptimizer.optimize(), `plot()`,
    `CMAOptions`, `scipy.optimize.fmin()`
is_feasible(x, f)
default to check feasibility, see also ``cma_default_options``
main(argv=None)
to install and/or test from the command line use::
 
    python cma.py [options | func dim sig0 [optkey optval][optkey optval]...]
 
with options being
 
``--test`` (or ``-t``) to run the doctest, ``--test -v`` to get (much) verbosity.
 
``install`` to install cma.py (uses setup from distutils.core).
 
``--doc`` for more infos.
 
Or start Python or (even better) ``ipython`` and::
 
    import cma
    cma.main('--test')
    help(cma)
    help(cma.fmin)
    res = fmin(cma.fcts.rosen, 10 * [0], 1)
    cma.plot()
 
Examples
========
Testing with the local python distribution from a command line
in a folder where ``cma.py`` can be found::
 
    python cma.py --test
 
And a single run on the Rosenbrock function::
 
    python cma.py rosen 10 1  # dimension initial_sigma
    python cma.py plot
 
In the python shell::
 
    import cma
    cma.main('--test')
plot(name=None, fig=None, abscissa=1, iteridx=None, plot_mean=False, foffset=1e-19, x_opt=None, fontsize=9)
plot data from files written by a `CMADataLogger`,
the call ``cma.plot(name, **argsdict)`` is a shortcut for
``cma.CMADataLogger(name).plot(**argsdict)``
 
Arguments
---------
    `name`
        name of the logger, filename prefix, None evaluates to
        the default 'outcmaes'
    `fig`
        filename or figure number, or both as a tuple (any order)
    `abscissa`
        0==plot versus iteration count,
        1==plot versus function evaluation number
    `iteridx`
        iteration indices to plot
 
Return `None`
 
Examples
--------
::
 
   cma.plot();  # the optimization might be still
                # running in a different shell
   cma.savefig('fig325.png')
   cma.closefig()
 
   cdl = cma.CMADataLogger().downsampling().plot()
   # in case the file sizes are large
 
Details
-------
Data from codes in other languages (C, Java, Matlab, Scilab) have the same
format and can be plotted just the same.
 
:See: `CMADataLogger`, `CMADataLogger.plot()`
pprint(to_be_printed)
nicely formated print
savefig(*args, **kwargs)
Save the current figure.
 
Call signature::
 
  savefig(fname, dpi=None, facecolor='w', edgecolor='w',
          orientation='portrait', papertype=None, format=None,
          transparent=False, bbox_inches=None, pad_inches=0.1,
          frameon=None)
 
The output formats available depend on the backend being used.
 
Arguments:
 
  *fname*:
    A string containing a path to a filename, or a Python
    file-like object, or possibly some backend-dependent object
    such as :class:`~matplotlib.backends.backend_pdf.PdfPages`.
 
    If *format* is *None* and *fname* is a string, the output
    format is deduced from the extension of the filename. If
    the filename has no extension, the value of the rc parameter
    ``savefig.format`` is used.
 
    If *fname* is not a string, remember to specify *format* to
    ensure that the correct backend is used.
 
Keyword arguments:
 
  *dpi*: [ *None* | ``scalar > 0`` ]
    The resolution in dots per inch.  If *None* it will default to
    the value ``savefig.dpi`` in the matplotlibrc file.
 
  *facecolor*, *edgecolor*:
    the colors of the figure rectangle
 
  *orientation*: [ 'landscape' | 'portrait' ]
    not supported on all backends; currently only on postscript output
 
  *papertype*:
    One of 'letter', 'legal', 'executive', 'ledger', 'a0' through
    'a10', 'b0' through 'b10'. Only supported for postscript
    output.
 
  *format*:
    One of the file extensions supported by the active
    backend.  Most backends support png, pdf, ps, eps and svg.
 
  *transparent*:
    If *True*, the axes patches will all be transparent; the
    figure patch will also be transparent unless facecolor
    and/or edgecolor are specified via kwargs.
    This is useful, for example, for displaying
    a plot on top of a colored background on a web page.  The
    transparency of these patches will be restored to their
    original values upon exit of this function.
 
  *frameon*:
    If *True*, the figure patch will be colored, if *False*, the
    figure background will be transparent.  If not provided, the
    rcParam 'savefig.frameon' will be used.
 
  *bbox_inches*:
    Bbox in inches. Only the given portion of the figure is
    saved. If 'tight', try to figure out the tight bbox of
    the figure.
 
  *pad_inches*:
    Amount of padding around the figure when bbox_inches is
    'tight'.
 
  *bbox_extra_artists*:
    A list of extra artists that will be considered when the
    tight bbox is calculated.
show()
unitdoctest()
is used to describe test cases and might in future become helpful
as an experimental tutorial as well. The main testing feature at the
moment is by doctest with ``cma._test()`` or conveniently by
``python cma.py --test``. With the ``--verbose`` option added, the
results will always slightly differ and many "failed" test cases
might be reported.
 
A simple first overall test:
    >>> import cma
    >>> res = cma.fmin(cma.fcts.elli, 3*[1], 1,
    ...                {'CMA_diagonal':2, 'seed':1, 'verb_time':0})
    (3_w,7)-CMA-ES (mu_w=2.3,w_1=58%) in dimension 3 (seed=1)
       Covariance matrix is diagonal for 2 iterations (1/ccov=7.0)
    Iterat #Fevals   function value     axis ratio  sigma   minstd maxstd min:sec
        1       7 1.453161670768570e+04 1.2e+00 1.08e+00  1e+00  1e+00
        2      14 3.281197961927601e+04 1.3e+00 1.22e+00  1e+00  2e+00
        3      21 1.082851071704020e+04 1.3e+00 1.24e+00  1e+00  2e+00
      100     700 8.544042012075362e+00 1.4e+02 3.18e-01  1e-03  2e-01
      200    1400 5.691152415221861e-12 1.0e+03 3.82e-05  1e-09  1e-06
      220    1540 3.890107746209078e-15 9.5e+02 4.56e-06  8e-11  7e-08
    termination on tolfun : 1e-11
    final/bestever f-value = 3.89010774621e-15 2.52273602735e-15
    mean solution:  [ -4.63614606e-08  -3.42761465e-10   1.59957987e-11]
    std deviation: [  6.96066282e-08   2.28704425e-09   7.63875911e-11]
 
Test on the Rosenbrock function with 3 restarts. The first trial only
finds the local optimum, which happens in about 20% of the cases.
 
    >>> import cma
    >>> res = cma.fmin(cma.fcts.rosen, 4*[-1], 1,
    ...                options={'ftarget':1e-6, 'verb_time':0,
    ...                    'verb_disp':500, 'seed':3},
    ...                restarts=3)
    (4_w,8)-CMA-ES (mu_w=2.6,w_1=52%) in dimension 4 (seed=3)
    Iterat #Fevals   function value     axis ratio  sigma   minstd maxstd min:sec
        1       8 4.875315645656848e+01 1.0e+00 8.43e-01  8e-01  8e-01
        2      16 1.662319948123120e+02 1.1e+00 7.67e-01  7e-01  8e-01
        3      24 6.747063604799602e+01 1.2e+00 7.08e-01  6e-01  7e-01
      184    1472 3.701428610430019e+00 4.3e+01 9.41e-07  3e-08  5e-08
    termination on tolfun : 1e-11
    final/bestever f-value = 3.70142861043 3.70142861043
    mean solution:  [-0.77565922  0.61309336  0.38206284  0.14597202]
    std deviation: [  2.54211502e-08   3.88803698e-08   4.74481641e-08   3.64398108e-08]
    (8_w,16)-CMA-ES (mu_w=4.8,w_1=32%) in dimension 4 (seed=4)
    Iterat #Fevals   function value     axis ratio  sigma   minstd maxstd min:sec
        1    1489 2.011376859371495e+02 1.0e+00 8.90e-01  8e-01  9e-01
        2    1505 4.157106647905128e+01 1.1e+00 8.02e-01  7e-01  7e-01
        3    1521 3.548184889359060e+01 1.1e+00 1.02e+00  8e-01  1e+00
      111    3249 6.831867555502181e-07 5.1e+01 2.62e-02  2e-04  2e-03
    termination on ftarget : 1e-06
    final/bestever f-value = 6.8318675555e-07 1.18576673231e-07
    mean solution:  [ 0.99997004  0.99993938  0.99984868  0.99969505]
    std deviation: [ 0.00018973  0.00038006  0.00076479  0.00151402]
    >>> assert res[1] <= 1e-6
 
Notice the different termination conditions. Termination on the target
function value ftarget prevents further restarts.
 
Test of scaling_of_variables option
 
    >>> import cma
    >>> opts = cma.CMAOptions()
    >>> opts['seed'] = 456
    >>> opts['verb_disp'] = 0
    >>> opts['CMA_active'] = 1
    >>> # rescaling of third variable: for searching in  roughly
    >>> #   x0 plus/minus 1e3*sigma0 (instead of plus/minus sigma0)
    >>> opts['scaling_of_variables'] = [1, 1, 1e3, 1]
    >>> res = cma.fmin(cma.fcts.rosen, 4 * [0.1], 0.1, opts)
    termination on tolfun : 1e-11
    final/bestever f-value = 2.68096173031e-14 1.09714829146e-14
    mean solution:  [ 1.00000001  1.00000002  1.00000004  1.00000007]
    std deviation: [  3.00466854e-08   5.88400826e-08   1.18482371e-07   2.34837383e-07]
 
The printed std deviations reflect the actual value in the parameters
of the function (not the one in the internal representation which 
can be different).
 
Test of CMA_stds scaling option.
 
    >>> import cma
    >>> opts = cma.CMAOptions()
    >>> s = 5 * [1]
    >>> s[0] = 1e3
    >>> opts.set('CMA_stds', s)
    >>> opts.set('verb_disp', 0)
    >>> res = cma.fmin(cma.fcts.cigar, 5 * [0.1], 0.1, opts)
    >>> assert res[1] < 1800
 
:See: cma.main(), cma._test()

 
Data
        Fcts = <cma.FitnessFunctions object>
__all__ = (u'main', u'fmin', u'fcts', u'Fcts', u'felli', u'rotate', u'pprint', u'plot', u'disp', u'show', u'savefig', u'closefig', u'use_archives', u'is_feasible', u'unitdoctest', u'DerivedDictBase', u'SolutionDict', u'CMASolutionDict', u'BestSolution', u'BoundNone', ...)
__author__ = u'Nikolaus Hansen'
__docformat__ = u'reStructuredText'
__version__ = u'1.1.06 $Revision: 4129 $ $Date: 2015-01-23 20:13:51 +0100 (Fri, 23 Jan 2015) $'
fcts = <cma.FitnessFunctions object>
rotate = <cma.Rotation object>
use_archives = True

 
Author
        Nikolaus Hansen