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), projections, 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}
; 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.LineSearchesStepsizeType
LineSearchesStepsize <: Stepsize

Wrapper for line searches available in the LineSearches.jl library.

Constructors

LineSearchesStepsize(
    M::AbstractManifold,
    linesearch;
    retraction_method::AbstractRetractionMethod=default_retraction_method(M),
    vector_transport_method::AbstractVectorTransportMethod=default_vector_transport_method(M),
)
LineSearchesStepsize(
    linesearch;
    retraction_method::AbstractRetractionMethod=ExponentialRetraction(),
    vector_transport_method::AbstractVectorTransportMethod=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. The retraction used for determining the line along which the search is performed can be provided as retraction_method. Gradient vectors are transported between points using vector_transport_method.

source

Manifolds.jl

Loading Manifolds.jl introduces the following additional functions

ManifoldsBase.mid_pointFunction
mid_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).

source
Manopt.max_stepsizeMethod
max_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.

source
Manopt.max_stepsizeMethod
max_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

source

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_VectorizedManifoldType
struct 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)).

source
MathOptInterface.dimensionMethod
MOI.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$).

source
Manopt.JuMP_OptimizerType
Manopt.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.

source
MathOptInterface.supportsMethod
MOI.supports(::Optimizer, attr::MOI.RawOptimizerAttribute)

Return a Bool indicating whether attr.name is a valid option name for Manopt.

source
MathOptInterface.getMethod
MOI.get(model::Optimizer, attr::MOI.RawOptimizerAttribute)

Return last value set by MOI.set(model, attr, value).

source
MathOptInterface.setMethod
MOI.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].

source
MathOptInterface.copy_toMethod
MOI.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.

source
MathOptInterface.supportsMethod
MOI.supports(::Manopt.JuMP_Optimizer, attr::MOI.RawOptimizerAttribute)

Return true indicating that Manopt.JuMP_Optimizer supports starting values for the variables.

source
MathOptInterface.setMethod
function 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.

source
MathOptInterface.setMethod
MOI.set(model::Optimizer, ::MOI.ObjectiveSense, sense::MOI.OptimizationSense)

Modify the objective sense to either MOI.MAX_SENSE, MOI.MIN_SENSE or MOI.FEASIBILITY_SENSE.

source
MathOptInterface.setMethod
MOI.set(model::Optimizer, ::MOI.ObjectiveFunction{F}, func::F) where {F}

Set the objective function as func for model.

source
MathOptInterface.supportsMethod
MOI.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.

source
JuMP.build_variableMethod
JuMP.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.

source
MathOptInterface.getMethod
MOI.get(model::Optimizer, ::MOI.ResultCount)

Return 0 if optimize! hasn't been called yet and 1 otherwise indicating that one solution is available.

source
MathOptInterface.getMethod
MOI.get(::Optimizer, ::MOI.SolverName)

Return the name of the Optimizer with the value of the descent_state_type option.

source
MathOptInterface.getMethod
MOI.get(model::Optimizer, attr::MOI.ObjectiveValue)

Return the value of the objective function evaluated at the solution.

source
MathOptInterface.getMethod
MOI.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.

source
MathOptInterface.getMethod
MOI.get(::Optimizer, ::MOI.DualStatus)

Returns MOI.NO_SOLUTION indicating that there is no dual solution available.

source
MathOptInterface.getMethod
MOI.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.

source
MathOptInterface.getMethod
MOI.get(::Optimizer, ::MOI.SolverVersion)

Return the version of the Manopt solver, it corresponds to the version of Manopt.jl.

source
MathOptInterface.getMethod
MOI.set(model::Optimizer, ::MOI.ObjectiveSense, sense::MOI.OptimizationSense)

Return the objective sense, defaults to MOI.FEASIBILITY_SENSE if no sense has already been set.

source
MathOptInterface.getMethod
MOI.get(model::Optimizer, attr::MOI.VariablePrimal, vi::MOI.VariableIndex)

Return the value of the solution for the variable of index vi.

source
MathOptInterface.getMethod
MOI.get(model::Optimizer, ::MOI.RawStatusString)

Return a String containing Manopt.get_reason without the ending newline character.

source