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 optimzation 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_esFunction
cma_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

Optional

  • 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))
  • retraction_method: (default_retraction_method(M, typeof(p_m)))
  • vector_transport_method: (default_vector_transport_method(M, typeof(p_m)))
  • basis (DefaultOrthonormalBasis()) basis used to represent covariance in
  • rng: (default_rng()) random number generator for generating new points on M

Output

the obtained (approximate) minimizer $p^*$. To obtain the whole final state of the solver, see get_solver_return for details.

source

State

Manopt.CMAESStateType
CMAESState{P,T} <: AbstractManoptSolverState

State of covariance matrix adaptation evolution strategy.

Fields

  • p the best point found so far
  • p_obj objective value at p
  • μ parent number
  • λ population size
  • μ_eff variance effective selection mass for the mean
  • c_1 learning rate for the rank-one update
  • c_c decay 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_m learning rate for the mean
  • d_σ damping parameter for step-size update
  • stop stopping criteria, StoppingCriterion
  • population population of the current generation
  • ys_c coordinates of random vectors for the current generation
  • covariance_matrix coordinates of the covariance matrix
  • covariance_matrix_eigen eigendecomposition of covariance_matrix
  • covariance_matrix_cond condition number of covariance_matrix, updated after eigendecomposition
  • best_fitness_current_gen best fitness value of individuals in the current generation
  • median_fitness_current_gen median fitness value of individuals in the current generation
  • worst_fitness_current_gen worst fitness value of individuals in the current generation
  • p_m point around which we search for new candidates
  • σ step size
  • p_σ coordinates of a vector in $T_{p_m} \mathcal M$
  • p_c coordinates of a vector in $T_{p_m} \mathcal M$
  • deviations standard deviations of coordinate RNG
  • buffer buffer for random number generation and wmean_y_c of length n_coords
  • e_mv_norm expected value of norm of the n_coords-variable standard normal distribution
  • recombination_weights recombination weights used for updating covariance matrix
  • retraction_method an AbstractRetractionMethod
  • vector_transport_method a vector transport to use
  • basis a real coefficient basis for covariance matrix
  • rng RNG 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=DefaultOrthonormalBasis(),
    rng::TRng=default_rng(),
) where {
    P,
    TParams<:Real,
    TStopping<:StoppingCriterion,
    TRetraction<:AbstractRetractionMethod,
    TVTM<:AbstractVectorTransportMethod,
    TB<:AbstractBasis,
    TRng<:AbstractRNG,
}

See also

cma_es

source

Stopping Criteria

Manopt.StopWhenBestCostInGenerationConstantType
StopWhenBestCostInGenerationConstant <: StoppingCriterion

Stop 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.

source
Manopt.StopWhenEvolutionStagnatesType
StopWhenEvolutionStagnates{TParam<:Real} <: StoppingCriterion

The best and median fitness in each iteraion 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.

source
Manopt.StopWhenPopulationCostConcentratedType
StopWhenPopulationCostConcentrated{TParam<:Real} <: StoppingCriterion

Stop 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)
source
Manopt.StopWhenPopulationDivergesType
StopWhenPopulationDiverges{TParam<:Real} <: StoppingCriterion

Stop 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].

source
Manopt.StopWhenPopulationStronglyConcentratedType
StopWhenPopulationStronglyConcentrated{TParam<:Real} <: StoppingCriterion

Stop 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].

source

Technical details

The cma_es solver requires the following functions of a manifold to be available

Internal helpers

You may add new methods to eigenvector_transport! if you know a more optimized implementation for your manifold.

Manopt.eigenvector_transport!Function
eigenvector_transport!(
    M::AbstractManifold,
    matrix_eigen::Eigen,
    p,
    q,
    basis::AbstractBasis,
    vtm::AbstractVectorTransportMethod,
)

Transport the matrix with matrix_eig eigendecomposition 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.

source

Literature