Extensions
LineSearches.jl
Manopt can be used with line search algorithms implemented in LineSearches.jl. This can be illustrated by the following example of optimizing Rosenbrock function constrained to the unit sphere.
using Manopt, Manifolds, LineSearches
# define objective function and its gradient
p = [1.0, 100.0]
function rosenbrock(::AbstractManifold, x)
val = zero(eltype(x))
for i in 1:(length(x) - 1)
val += (p[1] - x[i])^2 + p[2] * (x[i + 1] - x[i]^2)^2
end
return val
end
function rosenbrock_grad!(M::AbstractManifold, storage, x)
storage .= 0.0
for i in 1:(length(x) - 1)
storage[i] += -2.0 * (p[1] - x[i]) - 4.0 * p[2] * (x[i + 1] - x[i]^2) * x[i]
storage[i + 1] += 2.0 * p[2] * (x[i + 1] - x[i]^2)
end
project!(M, storage, x, storage)
return storage
end
# define constraint
n_dims = 5
M = Manifolds.Sphere(n_dims)
# set initial point
x0 = vcat(zeros(n_dims - 1), 1.0)
# use LineSearches.jl HagerZhang method with Manopt.jl quasiNewton solver
ls_hz = Manopt.LineSearchesStepsize(M, LineSearches.HagerZhang())
x_opt = quasi_Newton(
M,
rosenbrock,
rosenbrock_grad!,
x0;
stepsize=ls_hz,
evaluation=InplaceEvaluation(),
stopping_criterion=StopAfterIteration(1000) | StopWhenGradientNormLess(1e-6),
return_state=true,
)
# Solver state for `Manopt.jl`s Quasi Newton Method
After 10 iterations
## Parameters
* direction update: limited memory InverseBFGS (size 5)and ParallelTransport() as vector transport.
* retraction method: ExponentialRetraction()
* vector transport method: ParallelTransport()
## Stepsize
LineSearchesStepsize(HagerZhang{Float64, Base.RefValue{Bool}}
delta: Float64 0.1
sigma: Float64 0.9
alphamax: Float64 Inf
rho: Float64 5.0
epsilon: Float64 1.0e-6
gamma: Float64 0.66
linesearchmax: Int64 50
psi3: Float64 0.1
display: Int64 0
mayterminate: Base.RefValue{Bool}
cache: Nothing nothing
; retraction_method=ExponentialRetraction(), vector_transport_method=ParallelTransport())
## Stopping criterion
Stop When _one_ of the following are fulfilled:
Max Iteration 1000: not reached
|grad f| < 1.0e-6: reached
Overall: reached
This indicates convergence: Yes
In general this defines the following new stepsize
Manopt.LineSearchesStepsize
β TypeLineSearchesStepsize <: Stepsize
Wrapper for line searches available in the LineSearches.jl
library.
Constructors
LineSearchesStepsize(M::AbstractManifold, linesearch; kwargs...
LineSearchesStepsize(
linesearch;
retraction_method=ExponentialRetraction(),
vector_transport_method=ParallelTransport(),
)
Wrap linesearch
(for example HagerZhang
or MoreThuente
). The initial step selection from Linesearches.jl is not yet supported and the value 1.0 is used.
Keyword Arguments
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 transports
Manifolds.jl
Loading Manifolds.jl
introduces the following additional functions
Manopt.max_stepsize
β Methodmax_stepsize(M::FixedRankMatrices, p)
Return a reasonable guess of maximum step size on FixedRankMatrices
following the choice of typical distance in Matlab Manopt, the dimension of M
. See this note
Manopt.max_stepsize
β Methodmax_stepsize(M::Hyperrectangle, p)
The default maximum stepsize for Hyperrectangle
manifold with corners is maximum of distances from p
to each boundary.
Manopt.max_stepsize
β Methodmax_stepsize(M::TangentBundle, p)
Tangent bundle has injectivity radius of either infinity (for flat manifolds) or 0 (for non-flat manifolds). This makes a guess of what a reasonable maximum stepsize on a tangent bundle might be.
ManifoldsBase.mid_point
β Functionmid_point(M, p, q, x)
mid_point!(M, y, p, q, x)
Compute the mid point between p
and q
. If there is more than one mid point of (not necessarily minimizing) geodesics (for example on the sphere), the one nearest to x
is returned (in place of y
).
Internally, Manopt.jl
provides the two additional functions to choose some Euclidean space when needed as
Manopt.Rn
β FunctionRn(args; kwargs...)
Rn(s::Symbol=:Manifolds, args; kwargs...)
A small internal helper function to choose a Euclidean space. By default, this uses the DefaultManifold
unless you load a more advanced Euclidean space like Euclidean
from Manifolds.jl
Manopt.Rn_default
β FunctionRn_default()
Specify a default value to dispatch Rn
on. This default is set to Manifolds
, indicating, that when this package is loded, it is the preferred package to ask for a vector space space.
The default within Manopt.jl
is to use the DefaultManifold
from ManifoldsBase.jl
. If you load Manifolds.jl
this switches to using Euclidan
.
JuMP.jl
Manopt can be used using the JuMP.jl interface. The manifold is provided in the @variable
macro. Note that until now, only variables (points on manifolds) are supported, that are arrays, especially structs do not yet work. The algebraic expression of the objective function is specified in the @objective
macro. The descent_state_type
attribute specifies the solver.
using JuMP, Manopt, Manifolds
model = Model(Manopt.Optimizer)
# Change the solver with this option, `GradientDescentState` is the default
set_attribute("descent_state_type", GradientDescentState)
@variable(model, U[1:2, 1:2] in Stiefel(2, 2), start = 1.0)
@objective(model, Min, sum((A - U) .^ 2))
optimize!(model)
solution_summary(model)
Interface functions
Manopt.JuMP_ArrayShape
β Typestruct ArrayShape{N} <: JuMP.AbstractShape
Shape of an Array{T,N}
of size size
.
Manopt.JuMP_VectorizedManifold
β Typestruct VectorizedManifold{M} <: MOI.AbstractVectorSet
manifold::M
end
Representation of points of manifold
as a vector of R^n
where n
is MOI.dimension(VectorizedManifold(manifold))
.
MathOptInterface.dimension
β MethodMOI.dimension(set::VectorizedManifold)
Return the representation side of points on the (vectorized in representation) manifold. As the MOI variables are real, this means if the representation_size
yields (in product) n
, this refers to the vectorized point / tangent vector from (a subset of $β^n$).
Manopt.JuMP_Optimizer
β TypeManopt.JuMP_Optimizer()
Creates a new optimizer object for the MathOptInterface (MOI). An alias Manopt.JuMP_Optimizer
is defined for convenience.
The minimization of a function f(X)
of an array X[1:n1,1:n2,...]
over a manifold M
starting at X0
, can be modeled as follows:
using JuMP
model = Model(Manopt.JuMP_Optimizer)
@variable(model, X[i1=1:n1,i2=1:n2,...] in M, start = X0[i1,i2,...])
@objective(model, Min, f(X))
The optimizer assumes that M
has a Array
shape described by ManifoldsBase.representation_size
.
MathOptInterface.empty!
β MethodMOI.empty!(model::ManoptJuMPExt.Optimizer)
Clear all model data from model
but keep the options
set.
MathOptInterface.supports
β MethodMOI.supports(::Optimizer, attr::MOI.RawOptimizerAttribute)
Return a Bool
indicating whether attr.name
is a valid option name for Manopt
.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, attr::MOI.RawOptimizerAttribute)
Return last value
set by MOI.set(model, attr, value)
.
MathOptInterface.set
β MethodMOI.get(model::Optimizer, attr::MOI.RawOptimizerAttribute)
Set the value for the keyword argument attr.name
to give for the constructor model.options[DESCENT_STATE_TYPE]
.
MathOptInterface.supports_incremental_interface
β MethodMOI.supports_incremental_interface(::JuMP_Optimizer)
Return true
indicating that Manopt.JuMP_Optimizer
implements MOI.add_constrained_variables
and MOI.set
for MOI.ObjectiveFunction
so it can be used with JuMP.direct_model
and does not require a MOI.Utilities.CachingOptimizer
. See MOI.supports_incremental_interface
.
MathOptInterface.copy_to
β MethodMOI.copy_to(dest::Optimizer, src::MOI.ModelLike)
Because supports_incremental_interface(dest)
is true
, this simply uses MOI.Utilities.default_copy_to
and copies the variables with MOI.add_constrained_variables
and the objective sense with MOI.set
.
MathOptInterface.supports_add_constrained_variables
β MethodMOI.supports_add_constrained_variables(::JuMP_Optimizer, ::Type{<:VectorizedManifold})
Return true
indicating that Manopt.JuMP_Optimizer
support optimization on variables constrained to belong in a vectorized manifold Manopt.JuMP_VectorizedManifold
.
MathOptInterface.add_constrained_variables
β MethodMOI.add_constrained_variables(model::Optimizer, set::VectorizedManifold)
Add MOI.dimension(set)
variables constrained in set
and return the list of variable indices that can be used to reference them as well a constraint index for the constraint enforcing the membership of the variables in the Manopt.JuMP_VectorizedManifold
set
.
MathOptInterface.is_valid
β MethodMOI.is_valid(model::Optimizer, vi::MOI.VariableIndex)
Return whether vi
is a valid variable index.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, ::MOI.NumberOfVariables)
Return the number of variables added in the model, this corresponds to the MOI.dimension
of the Manopt.JuMP_VectorizedManifold
.
MathOptInterface.supports
β MethodMOI.supports(::Manopt.JuMP_Optimizer, attr::MOI.RawOptimizerAttribute)
Return true
indicating that Manopt.JuMP_Optimizer
supports starting values for the variables.
MathOptInterface.set
β Methodfunction MOI.set(
model::Optimizer,
::MOI.VariablePrimalStart,
vi::MOI.VariableIndex,
value::Union{Real,Nothing},
)
Set the starting value of the variable of index vi
to value
. Note that if value
is nothing
then it essentially unset any previous starting values set and hence MOI.optimize!
unless another starting value is set.
MathOptInterface.set
β MethodMOI.set(model::Optimizer, ::MOI.ObjectiveSense, sense::MOI.OptimizationSense)
Modify the objective sense to either MOI.MAX_SENSE
, MOI.MIN_SENSE
or MOI.FEASIBILITY_SENSE
.
MathOptInterface.set
β MethodMOI.set(model::Optimizer, ::MOI.ObjectiveFunction{F}, func::F) where {F}
Set the objective function as func
for model
.
MathOptInterface.supports
β MethodMOI.supports(::Optimizer, ::Union{MOI.ObjectiveSense,MOI.ObjectiveFunction})
Return true
indicating that Optimizer
supports being set the objective sense (that is, min, max or feasibility) and the objective function.
JuMP.build_variable
β MethodJuMP.build_variable(::Function, func, m::ManifoldsBase.AbstractManifold)
Build a JuMP.VariablesConstrainedOnCreation
object containing variables and the Manopt.JuMP_VectorizedManifold
in which they should belong as well as the shape
that can be used to go from the vectorized MOI representation to the shape of the manifold, that is, Manopt.JuMP_ArrayShape
.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, ::MOI.ResultCount)
Return 0
if optimize!
hasn't been called yet and 1
otherwise indicating that one solution is available.
MathOptInterface.get
β MethodMOI.get(::Optimizer, ::MOI.SolverName)
Return the name of the Optimizer
with the value of the descent_state_type
option.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, attr::MOI.ObjectiveValue)
Return the value of the objective function evaluated at the solution.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, ::MOI.PrimalStatus)
Return MOI.NO_SOLUTION
if optimize!
hasn't been called yet and MOI.FEASIBLE_POINT
otherwise indicating that a solution is available to query with MOI.VariablePrimalStart
.
MathOptInterface.get
β MethodMOI.get(::Optimizer, ::MOI.DualStatus)
Returns MOI.NO_SOLUTION
indicating that there is no dual solution available.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, ::MOI.ResultCount)
Return MOI.OPTIMIZE_NOT_CALLED
if optimize!
hasn't been called yet and MOI.LOCALLY_SOLVED
otherwise indicating that the solver has solved the problem to local optimality the value of MOI.RawStatusString
for more details on why the solver stopped.
MathOptInterface.get
β MethodMOI.get(::Optimizer, ::MOI.SolverVersion)
Return the version of the Manopt solver, it corresponds to the version of Manopt.jl.
MathOptInterface.get
β MethodMOI.set(model::Optimizer, ::MOI.ObjectiveSense, sense::MOI.OptimizationSense)
Return the objective sense, defaults to MOI.FEASIBILITY_SENSE
if no sense has already been set.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, attr::MOI.VariablePrimal, vi::MOI.VariableIndex)
Return the value of the solution for the variable of index vi
.
MathOptInterface.get
β MethodMOI.get(model::Optimizer, ::MOI.RawStatusString)
Return a String
containing Manopt.get_reason
without the ending newline character.