Mesh adaptive direct search (MADS)

Manopt.mesh_adaptive_direct_searchFunction
mesh_adaptive_direct_search(M, f, p=rand(M); kwargs...)
mesh_adaptive_direct_search(M, mco::AbstractManifoldCostObjective, p=rand(M); kwargs..)
mesh_adaptive_direct_search!(M, f, p; kwargs...)
mesh_adaptive_direct_search!(M, mco::AbstractManifoldCostObjective, p; kwargs..)

The Mesh Adaptive Direct Search (MADS) algorithm minimizes an objective function $f: \mathcal M → ℝ$ on the manifold M. The algorithm constructs an implicit mesh in the tangent space $T_{p}\mathcal M$ at the current candidate $p$. Each iteration consists of a search step and a poll step.

The search step selects points from the implicit mesh and attempts to find an improved candidate solution that reduces the value of $f$. If the search step fails to generate an improved candidate solution, the poll step is performed. It consists of a local exploration on the current implicit mesh in the neighbourhood of the current iterate.

Input

  • M::AbstractManifold: a Riemannian manifold $\mathcal M$
  • f: a cost function $f: \mathcal M→ ℝ$ implemented as (M, p) -> v
  • p: a point on the manifold $\mathcal M$

Keyword arguments

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.

source
Manopt.mesh_adaptive_direct_search!Function
mesh_adaptive_direct_search(M, f, p=rand(M); kwargs...)
mesh_adaptive_direct_search(M, mco::AbstractManifoldCostObjective, p=rand(M); kwargs..)
mesh_adaptive_direct_search!(M, f, p; kwargs...)
mesh_adaptive_direct_search!(M, mco::AbstractManifoldCostObjective, p; kwargs..)

The Mesh Adaptive Direct Search (MADS) algorithm minimizes an objective function $f: \mathcal M → ℝ$ on the manifold M. The algorithm constructs an implicit mesh in the tangent space $T_{p}\mathcal M$ at the current candidate $p$. Each iteration consists of a search step and a poll step.

The search step selects points from the implicit mesh and attempts to find an improved candidate solution that reduces the value of $f$. If the search step fails to generate an improved candidate solution, the poll step is performed. It consists of a local exploration on the current implicit mesh in the neighbourhood of the current iterate.

Input

  • M::AbstractManifold: a Riemannian manifold $\mathcal M$
  • f: a cost function $f: \mathcal M→ ℝ$ implemented as (M, p) -> v
  • p: a point on the manifold $\mathcal M$

Keyword arguments

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.

source

State

Manopt.MeshAdaptiveDirectSearchStateType
MeshAdaptiveDirectSearchState <: AbstractManoptSolverState

Fields

  • p::P: a point on the manifold $\mathcal M$storing the current iterate
  • mesh_size: the current (internal) mesh size
  • scale_mesh: the current scaling of the internal mesh size, yields the actual mesh size used
  • max_stepsize: an upper bound for the longest step taken in looking for a candidate in either poll or search
  • poll_size
  • stop::StoppingCriterion: a functor indicating that the stopping criterion is fulfilled
  • poll::[AbstractMeshPollFunction]: a poll step (functor) to perform
  • search::[AbstractMeshSearchFunction}(@ref) a search step (functor) to perform
source

Poll

Manopt.AbstractMeshPollFunctionType
AbstractMeshPollFunction

An abstract type for common “poll” strategies in the mesh_adaptive_direct_search solver. A subtype of this The functor has to fulfil

  • be callable as poll!(problem, mesh_size; kwargs...) and modify the state

as well as to provide functions

  • is_successful(poll!) that indicates whether the last poll was successful in finding a new candidate
  • get_basepoint(poll!) that returns the base point at which the mesh is build
  • get_candidate(poll!) that returns the last found candidate if the poll was successful. Otherwise the base point is returned
  • get_descent_direction(poll!) the the vector that points from the base point to the candidate. If the last poll was not successful, the zero vector is returned
  • update_basepoint!(M, poll!, p) that updates the base point to p and all necessary internal data to a new point to build a mesh at

The kwargs... could include

  • scale_mesh=1.0: to rescale the mesh globally
  • max_stepsize=Inf: avoid exceeding a step size beyond this value, e.g. injectivity radius. any vector longer than this should be shortened to the provided maximum step size.
source
Manopt.LowerTriangularAdaptivePollType
LowerTriangularAdaptivePoll <: AbstractMeshPollFunction

Generate a mesh (poll step) based on Section 6 and 7 of [Dre07], with two small modifications:

  • The mesh can be scaled globally so instead of $Δ_0^m=1$ a certain different scale is used
  • Any poll direction can be rescaled if it is too long. This is to not exceed the injectivity radius for example.

Functor

(p::LowerTriangularAdaptivePoll)(problem, mesh_size; scale_mesh=1.0, max_stepsize=inf)

Fields

  • base_point::P: a point on the manifold, where the mesh is build in the tangent space
  • basis: a basis of the current tangent space with respect to which the mesh is stored
  • candidate::P: a memory for a new point/candidate
  • mesh: a vector of tangent vectors storing the mesh.
  • random_vector: a $d$-dimensional random vector $b_l$`
  • random_index: a random index $ι$
  • retraction_method::AbstractRetractionMethod: a retraction $\operatorname{retr}$ to use, see the section on retractions
  • vector_transport_method::AbstractVectorTransportMethodP: a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transports
  • X::T the last successful poll direction stored as a tangent vector. initialised to the zero vector and reset to the zero vector after moving to a new tangent space.

Constructor

LowerTriangularAdaptivePoll(M, p=rand(M); kwargs...)

Keyword arguments

source

as well as the internal functions

Manopt.AbstractMeshSearchFunctionType
AbstractMeshSearchFunction

Should be callable as search!(problem, mesh_size, p, X; kwargs...)

where X is the last successful poll direction from the tangent space at p if that exists and the zero vector otherwise.

Besides that the following functions should be implemented

  • is_successful(search!) that indicates whether the last search was successful in finding a new candidate
  • get_candidate(search!) that returns the last found candidate
source
Manopt.DefaultMeshAdaptiveDirectSearchType
DefaultMeshAdaptiveDirectSearch <: AbstractMeshSearchFunction

Functor

(s::DefaultMeshAdaptiveDirectSearch)(problem, mesh_size::Real, X; scale_mesh::Real=1.0, max_stepsize::Real=inf)

Fields

  • q: a temporary memory for a point on the manifold
  • X: information to perform the search, e.g. the last direction found by poll.
  • last_search_improved::Bool indicate whether the last search was successful, i.e. improved the cost.
  • retraction_method::AbstractRetractionMethod: a retraction $\operatorname{retr}$ to use, see the section on retractions

Constructor

DefaultMeshAdaptiveDirectSearch(M::AbstractManifold, p=rand(M); kwargs...)

Keyword arguments

source

as well as the internal functions

Additional stopping criteria

Technical details

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

  • A retract!(M, q, p, X); it is recommended to set the default_retraction_method to a favourite retraction. If this default is set, a retraction_method= does not have to be specified.
  • Within the default initialization rand(M) is used to generate the initial population
  • A vector_transport_to!M, Y, p, X, q); it is recommended to set the default_vector_transport_method to a favourite retraction. If this default is set, a vector_transport_method= does not have to be specified.

Literature

[Dre07]
D. W. Dreisigmeyer. Direct Search Alogirthms over Riemannian Manifolds (Optimization Online, 2007).