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 pointp
σ=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 theStopWhenPopulationCostConcentrated
, similar to absolute difference between function values at subsequent pointstol_x=1e-12
: tolerance for theStopWhenPopulationStronglyConcentrated
, 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 fulfilledretraction_method=
default_retraction_method
(M, typeof(p))
: a retraction $\operatorname{retr}$ to use, see the section on retractionsvector_transport_method=
default_vector_transport_method
(M, typeof(p))
: a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transportsbasis
(DefaultOrthonormalBasis()
) basis used to represent covariance inrng=default_rng()
: random number generator for generating new points onM
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} <: AbstractManoptSolverState
State of covariance matrix adaptation evolution strategy.
Fields
p::P
: a point on the manifold $\mathcal M$ storing the best point found so farp_obj
objective value atp
μ
parent numberλ
population sizeμ_eff
variance effective selection mass for the meanc_1
learning rate for the rank-one updatec_c
decay rate for cumulation path for the rank-one updatec_μ
learning rate for the rank-μ updatec_σ
decay rate for the cumulation path for the step-size controlc_m
learning rate for the meand_σ
damping parameter for step-size updatepopulation
population of the current generationys_c
coordinates of random vectors for the current generationcovariance_matrix
coordinates of the covariance matrixcovariance_matrix_eigen
eigen decomposition ofcovariance_matrix
covariance_matrix_cond
condition number ofcovariance_matrix
, updated after eigen decompositionbest_fitness_current_gen
best fitness value of individuals in the current generationmedian_fitness_current_gen
median fitness value of individuals in the current generationworst_fitness_current_gen
worst fitness value of individuals in the current generationp_m
point around which the search for new candidates is doneσ
step sizep_σ
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 RNGbuffer
buffer for random number generation andwmean_y_c
of lengthn_coords
e_mv_norm
expected value of norm of then_coords
-variable standard normal distributionrecombination_weights
recombination weights used for updating covariance matrixretraction_method::AbstractRetractionMethod
: a retraction $\operatorname{retr}$ to use, see the section on retractionsstop::StoppingCriterion
: a functor indicating that the stopping criterion is fulfilledvector_transport_method::AbstractVectorTransportMethodP
: a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transportsbasis
a real coefficient basis for covariance matrixrng
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
Stopping criteria
Manopt.StopWhenBestCostInGenerationConstant
— TypeStopWhenBestCostInGenerationConstant <: 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
.
Manopt.StopWhenCovarianceIllConditioned
— TypeStopWhenCovarianceIllConditioned <: StoppingCriterion
Stop CMA-ES if condition number of covariance matrix exceeds threshold
. This corresponds to ConditionCov
condition from [Han23].
Manopt.StopWhenEvolutionStagnates
— TypeStopWhenEvolutionStagnates{TParam<:Real} <: StoppingCriterion
The 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} <: 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)
Manopt.StopWhenPopulationDiverges
— TypeStopWhenPopulationDiverges{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].
Manopt.StopWhenPopulationStronglyConcentrated
— TypeStopWhenPopulationStronglyConcentrated{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].
Fields
tol
the tolerance to verify againstat_iteration
an 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_method
to 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_method
to 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)
andget_vector!
(M, X, p, c, b)
with respect to theAbstractBasis
b
provided, which isDefaultOrthonormalBasis
by default from thebasis=
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).