Covariance matrix adaptation evolutionary strategy
The CMA-ES algorithm has been implemented based on [Han23] with basic Riemannian adaptations, related to transport of covariance matrix and its update vectors. Other attempts at adapting CMA-ES to Riemannian optimization include [CFFS10]. The algorithm is suitable for global optimization.
Covariance matrix transport between consecutive mean points is handled by eigenvector_transport! function which is based on the idea of transport of matrix eigenvectors.
Manopt.cma_es — Functioncma_es(M, f, p_m=rand(M); σ::Real=1.0, kwargs...)Perform covariance matrix adaptation evolutionary strategy search for global gradient-free randomized optimization. It is suitable for complicated non-convex functions. It can be reasonably expected to find global minimum within 3σ distance from p_m.
Implementation is based on [Han23] with basic adaptations to the Riemannian setting.
Input
- M: a manifold $\mathcal M$
- f: a cost function $f: \mathcal M→ℝ$ to find a minimizer $p^*$ for
Keyword arguments
- p_m=- rand- (M): an initial point- p
- σ=1.0: initial standard deviation
- λ: (- 4 + Int(floor(3 * log(manifold_dimension(M))))population size (can be increased for a more thorough global search but decreasing is not recommended)
- tol_fun=1e-12: tolerance for the- StopWhenPopulationCostConcentrated, similar to absolute difference between function values at subsequent points
- tol_x=1e-12: tolerance for the- StopWhenPopulationStronglyConcentrated, similar to absolute difference between subsequent point but actually computed from distribution parameters.
- stopping_criterion=- default_cma_es_stopping_criterion(M, λ; tol_fun=tol_fun, tol_x=tol_x): a functor indicating that the stopping criterion is fulfilled
- retraction_method=- default_retraction_method- (M, typeof(p)): a retraction $\operatorname{retr}$ to use, see the section on retractions
- vector_transport_method=- default_vector_transport_method- (M, typeof(p)): a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transports
- basis(- DefaultOrthonormalBasis()) basis used to represent covariance in
- rng=default_rng(): random number generator for generating new points on- M
All other keyword arguments are passed to decorate_state! for state decorators or decorate_objective! for objective decorators, respectively.
Output
The obtained approximate minimizer $p^*$. To obtain the whole final state of the solver, see get_solver_return for details, especially the return_state= keyword.
State
Manopt.CMAESState — TypeCMAESState{P,T} <: AbstractManoptSolverStateState of covariance matrix adaptation evolution strategy.
Fields
- p::P: a point on the manifold $\mathcal M$ storing the best point found so far
- p_objobjective value at- p
- μparent number
- λpopulation size
- μ_effvariance effective selection mass for the mean
- c_1learning rate for the rank-one update
- c_cdecay rate for cumulation path for the rank-one update
- c_μlearning rate for the rank-μ update
- c_σdecay rate for the cumulation path for the step-size control
- c_mlearning rate for the mean
- d_σdamping parameter for step-size update
- populationpopulation of the current generation
- ys_ccoordinates of random vectors for the current generation
- covariance_matrixcoordinates of the covariance matrix
- covariance_matrix_eigeneigen decomposition of- covariance_matrix
- covariance_matrix_condcondition number of- covariance_matrix, updated after eigen decomposition
- best_fitness_current_genbest fitness value of individuals in the current generation
- median_fitness_current_genmedian fitness value of individuals in the current generation
- worst_fitness_current_genworst fitness value of individuals in the current generation
- p_mpoint around which the search for new candidates is done
- σstep size
- p_σcoordinates of a vector in $T_{p_m}\mathcal M$
- p_ccoordinates of a vector in $T_{p_m}\mathcal M$
- deviationsstandard deviations of coordinate RNG
- bufferbuffer for random number generation and- wmean_y_cof length- n_coords
- e_mv_normexpected value of norm of the- n_coords-variable standard normal distribution
- recombination_weightsrecombination weights used for updating covariance matrix
- retraction_method::AbstractRetractionMethod: a retraction $\operatorname{retr}$ to use, see the section on retractions
- stop::StoppingCriterion: a functor indicating that the stopping criterion is fulfilled
- vector_transport_method::AbstractVectorTransportMethodP: a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transports
- basisa real coefficient basis for covariance matrix
- rngRNG for generating new points
Constructor
CMAESState(
    M::AbstractManifold,
    p_m::P,
    μ::Int,
    λ::Int,
    μ_eff::TParams,
    c_1::TParams,
    c_c::TParams,
    c_μ::TParams,
    c_σ::TParams,
    c_m::TParams,
    d_σ::TParams,
    stop::TStopping,
    covariance_matrix::Matrix{TParams},
    σ::TParams,
    recombination_weights::Vector{TParams};
    retraction_method::TRetraction=default_retraction_method(M, typeof(p_m)),
    vector_transport_method::TVTM=default_vector_transport_method(M, typeof(p_m)),
    basis::TB=default_basis(M, typeof(p_m)),
    rng::TRng=default_rng(),
) where {
    P,
    TParams<:Real,
    TStopping<:StoppingCriterion,
    TRetraction<:AbstractRetractionMethod,
    TVTM<:AbstractVectorTransportMethod,
    TB<:AbstractBasis,
    TRng<:AbstractRNG,
}See also
Stopping criteria
Manopt.StopWhenBestCostInGenerationConstant — TypeStopWhenBestCostInGenerationConstant <: StoppingCriterionStop if the range of the best objective function values of the last iteration_range generations is zero. This corresponds to EqualFUnValues condition from [Han23].
See also StopWhenPopulationCostConcentrated.
Manopt.StopWhenCovarianceIllConditioned — TypeStopWhenCovarianceIllConditioned <: StoppingCriterionStop CMA-ES if condition number of covariance matrix exceeds threshold. This corresponds to ConditionCov condition from [Han23].
Manopt.StopWhenEvolutionStagnates — TypeStopWhenEvolutionStagnates{TParam<:Real} <: StoppingCriterionThe best and median fitness in each iteration is tracked over the last 20% but at least min_size and no more than max_size iterations. Solver is stopped if in both histories the median of the most recent fraction of values is not better than the median of the oldest fraction.
Manopt.StopWhenPopulationCostConcentrated — TypeStopWhenPopulationCostConcentrated{TParam<:Real} <: StoppingCriterionStop if the range of the best objective function value in the last max_size generations and all function values in the current generation is below tol. This corresponds to TolFun condition from [Han23].
Constructor
StopWhenPopulationCostConcentrated(tol::Real, max_size::Int)Manopt.StopWhenPopulationDiverges — TypeStopWhenPopulationDiverges{TParam<:Real} <: StoppingCriterionStop if σ times maximum deviation increased by more than tol. This usually indicates a far too small σ, or divergent behavior. This corresponds to TolXUp condition from [Han23].
Manopt.StopWhenPopulationStronglyConcentrated — TypeStopWhenPopulationStronglyConcentrated{TParam<:Real} <: StoppingCriterionStop if the standard deviation in all coordinates is smaller than tol and norm of σ * p_c is smaller than tol. This corresponds to TolX condition from [Han23].
Fields
- tolthe tolerance to verify against
- at_iterationan internal field to indicate at with iteration $i \geq 0$ the tolerance was met.
Constructor
StopWhenPopulationStronglyConcentrated(tol::Real)Technical details
The cma_es solver requires the following functions of a manifold to be available
- A retract!(M, q, p, X); it is recommended to set thedefault_retraction_methodto a favourite retraction. If this default is set, aretraction_method=does not have to be specified.
- A vector_transport_to!M, Y, p, X, q); it is recommended to set thedefault_vector_transport_methodto a favourite retraction. If this default is set, avector_transport_method=does not have to be specified.
- A copyto!(M, q, p)andcopy(M,p)for points and similarlycopy(M, p, X)for tangent vectors.
- get_coordinates!- (M, Y, p, X, b)and- get_vector!- (M, X, p, c, b)with respect to the- AbstractBasis- bprovided, which is- DefaultOrthonormalBasisby default from the- basis=keyword.
- An is_flat(M).
Internal helpers
You may add new methods to eigenvector_transport! if you know a more optimized implementation for your manifold.
Manopt.eigenvector_transport! — Functioneigenvector_transport!(
    M::AbstractManifold,
    matrix_eigen::Eigen,
    p,
    q,
    basis::AbstractBasis,
    vtm::AbstractVectorTransportMethod,
)Transport the matrix with matrix_eig eigen decomposition when expanded in basis from point p to point q on M. Update matrix_eigen in-place.
(p, matrix_eig) belongs to the fiber bundle of $B = \mathcal M × SPD(n)$, where n is the (real) dimension of M. The function corresponds to the Ehresmann connection defined by vector transport vtm of eigenvectors of matrix_eigen.
Literature
- [CFFS10]
- S. Colutto, F. Fruhauf, M. Fuchs and O. Scherzer. The CMA-ES on Riemannian Manifolds to Reconstruct Shapes in 3-D Voxel Images. IEEE Transactions on Evolutionary Computation 14, 227–245 (2010).
- [Han23]
- N. Hansen. The CMA Evolution Strategy: A Tutorial. ArXiv Preprint (2023).