Mesh adaptive direct search (MADS)
Manopt.mesh_adaptive_direct_search
— Functionmesh_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
mesh_basis=
DefaultOrthonormalBasis
: a basis to generate the mesh in. The mesh is generated in coordinates of this basis in every tangent spacemax_stepsize=
injectivity_radius
(M)
: a maximum step size to take. any vector generated on the mesh is shortened to this length to avoid leaving the injectivity radius,poll::
AbstractMeshPollFunction
=
LowerTriangularAdaptivePoll
(M, copy(M,p))
: the poll function to use. Themesh_basis
(asbasis
),retraction_method
, andvector_transport_method
are passed to this default as well.retraction_method=
default_retraction_method
(M, typeof(p))
: a retraction $\operatorname{retr}$ to use, see the section on retractionsscale_mesh=
injectivity_radius
(M) / 4
: initial scaling of the meshsearch::
AbstractMeshSearchFunction
=
DefaultMeshAdaptiveDirectSearch
(M, copy(M,p))
: the search function to use. Theretraction_method
is passed to this default as well.stopping_criterion=
StopAfterIteration
(500)
|
StopWhenPollSizeLess
(1e-10)
: a functor indicating that the stopping criterion is fulfilledvector_transport_method=
default_vector_transport_method
(M, typeof(p))
: a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transportsX=
zero_vector
(M, p)
: a tangent vector at the point $p$ on the manifold $\mathcal 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.
Manopt.mesh_adaptive_direct_search!
— Functionmesh_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
mesh_basis=
DefaultOrthonormalBasis
: a basis to generate the mesh in. The mesh is generated in coordinates of this basis in every tangent spacemax_stepsize=
injectivity_radius
(M)
: a maximum step size to take. any vector generated on the mesh is shortened to this length to avoid leaving the injectivity radius,poll::
AbstractMeshPollFunction
=
LowerTriangularAdaptivePoll
(M, copy(M,p))
: the poll function to use. Themesh_basis
(asbasis
),retraction_method
, andvector_transport_method
are passed to this default as well.retraction_method=
default_retraction_method
(M, typeof(p))
: a retraction $\operatorname{retr}$ to use, see the section on retractionsscale_mesh=
injectivity_radius
(M) / 4
: initial scaling of the meshsearch::
AbstractMeshSearchFunction
=
DefaultMeshAdaptiveDirectSearch
(M, copy(M,p))
: the search function to use. Theretraction_method
is passed to this default as well.stopping_criterion=
StopAfterIteration
(500)
|
StopWhenPollSizeLess
(1e-10)
: a functor indicating that the stopping criterion is fulfilledvector_transport_method=
default_vector_transport_method
(M, typeof(p))
: a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transportsX=
zero_vector
(M, p)
: a tangent vector at the point $p$ on the manifold $\mathcal 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.MeshAdaptiveDirectSearchState
— TypeMeshAdaptiveDirectSearchState <: AbstractManoptSolverState
Fields
p::P
: a point on the manifold $\mathcal M$storing the current iteratemesh_size
: the current (internal) mesh sizescale_mesh
: the current scaling of the internal mesh size, yields the actual mesh size usedmax_stepsize
: an upper bound for the longest step taken in looking for a candidate in either poll or searchpoll_size
stop::StoppingCriterion
: a functor indicating that the stopping criterion is fulfilledpoll::
[AbstractMeshPollFunction
]: a poll step (functor) to performsearch::
[AbstractMeshSearchFunction
}(@ref) a search step (functor) to perform
Poll
Manopt.AbstractMeshPollFunction
— TypeAbstractMeshPollFunction
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 candidateget_basepoint(poll!)
that returns the base point at which the mesh is buildget_candidate(poll!)
that returns the last found candidate if the poll was successful. Otherwise the base point is returnedget_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 returnedupdate_basepoint!(M, poll!, p)
that updates the base point top
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 globallymax_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.
Manopt.LowerTriangularAdaptivePoll
— TypeLowerTriangularAdaptivePoll <: 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 spacebasis
: a basis of the current tangent space with respect to which the mesh is storedcandidate::P
: a memory for a new point/candidatemesh
: 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 retractionsvector_transport_method::AbstractVectorTransportMethodP
: a vector transport $\mathcal T_{⋅←⋅}$ to use, see the section on vector transportsX::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
basis=
DefaultOrthonormalBasis
retraction_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 transportsX=
zero_vector
(M, p)
: a tangent vector at the point $p$ on the manifold $\mathcal M$
as well as the internal functions
Manopt.get_descent_direction
— Methodget_descent_direction(ltap::LowerTriangularAdaptivePoll)
Return the direction of the last LowerTriangularAdaptivePoll
that yields a descent of the cost. If the poll was not successful, the zero vector is returned
Manopt.is_successful
— Methodis_successful(ltap::LowerTriangularAdaptivePoll)
Return whether the last LowerTriangularAdaptivePoll
step was successful
Manopt.get_candidate
— Methodget_candidate(ltap::LowerTriangularAdaptivePoll)
Return the candidate of the last successful LowerTriangularAdaptivePoll
. If the poll was unsuccessful, the base point is returned.
Manopt.get_basepoint
— Methodget_basepoint(ltap::LowerTriangularAdaptivePoll)
Return the base point of the tangent space, where the mash for the LowerTriangularAdaptivePoll
is build in.
Manopt.update_basepoint!
— Methodupdate_basepoint!(M, ltap::LowerTriangularAdaptivePoll, p)
Update the base point of the LowerTriangularAdaptivePoll
. This especially also updates the basis, that is used to build a (new) mesh.
Search
Manopt.AbstractMeshSearchFunction
— TypeAbstractMeshSearchFunction
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 candidateget_candidate(search!)
that returns the last found candidate
Manopt.DefaultMeshAdaptiveDirectSearch
— TypeDefaultMeshAdaptiveDirectSearch <: 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 manifoldX
: 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
- `X::T=zero_vector(M, p)
retraction_method=
default_retraction_method
(M, typeof(p))
: a retraction $\operatorname{retr}$ to use, see the section on retractions
as well as the internal functions
Manopt.is_successful
— Methodis_successful(dmads::DefaultMeshAdaptiveDirectSearch)
Return whether the last DefaultMeshAdaptiveDirectSearch
was successful.
Manopt.get_candidate
— Methodget_candidate(dmads::DefaultMeshAdaptiveDirectSearch)
Return the last candidate a DefaultMeshAdaptiveDirectSearch
found
Additional stopping criteria
Manopt.StopWhenPollSizeLess
— TypeStopWhenPollSizeLess <: StoppingCriterion
stores a threshold when to stop looking at the poll mesh size of an MeshAdaptiveDirectSearchState
.
Constructor
StopWhenPollSizeLess(ε)
initialize the stopping criterion to a threshold ε
.
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 thedefault_retraction_method
to a favourite retraction. If this default is set, aretraction_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 thedefault_vector_transport_method
to a favourite retraction. If this default is set, avector_transport_method=
does not have to be specified.
Literature
- [Dre07]
- D. W. Dreisigmeyer. Direct Search Alogirthms over Riemannian Manifolds (Optimization Online, 2007).