A manifold objective

The Objective describes that actual cost function and all its properties.

Manopt.AbstractManifoldObjectiveType
AbstractManifoldObjective{E<:AbstractEvaluationType}

Describe the collection of the optimization function $f: \mathcal M → ℝ$ (or even a vectorial range) and its corresponding elements, which might for example be a gradient or (one or more) proximal maps.

All these elements should usually be implemented as functions (M, p) -> ..., or (M, X, p) -> ... that is

  • the first argument of these functions should be the manifold M they are defined on
  • the argument X is present, if the computation is performed in-place of X (see InplaceEvaluation)
  • the argument p is the place the function ($f$ or one of its elements) is evaluated at.

the type T indicates the global AbstractEvaluationType.

source

Which has two main different possibilities for its containing functions concerning the evaluation mode, not necessarily the cost, but for example gradient in an AbstractManifoldGradientObjective.

Decorators for objectives

An objective can be decorated using the following trait and function to initialize

Manopt.dispatch_objective_decoratorFunction
dispatch_objective_decorator(o::AbstractManoptSolverState)

Indicate internally, whether an AbstractManifoldObjective o to be of decorating type, it stores (encapsulates) an object in itself, by default in the field o.objective.

Decorators indicate this by returning Val{true} for further dispatch.

The default is Val{false}, so by default an state is not decorated.

source
Manopt.decorate_objective!Function
decorate_objective!(M, o::AbstractManifoldObjective)

decorate the AbstractManifoldObjectiveo with specific decorators.

Optional arguments

optional arguments provide necessary details on the decorators. A specific one is used to activate certain decorators.

  • cache=missing: specify a cache. Currently :Simple is supported and :LRU if you load LRUCache.jl. For this case a tuple specifying what to cache and how many can be provided, has to be specified. For example (:LRU, [:Cost, :Gradient], 10) states that the last 10 used cost function evaluations and gradient evaluations should be stored. See objective_cache_factory for details.
  • count=missing: specify calls to the objective to be called, see ManifoldCountObjective for the full list
  • objective_type=:Riemannian: specify that an objective is :Riemannian or :Euclidean. The :Euclidean symbol is equivalent to specifying it as :Embedded, since in the end, both refer to converting an objective from the embedding (whether its Euclidean or not) to the Riemannian one.

See also

objective_cache_factory

source

Embedded objectives

Manopt.EmbeddedManifoldObjectiveType
EmbeddedManifoldObjective{P, T, E, O2, O1<:AbstractManifoldObjective{E}} <:
   AbstractDecoratedManifoldObjective{E,O2}

Declare an objective to be defined in the embedding. This also declares the gradient to be defined in the embedding, and especially being the Riesz representer with respect to the metric in the embedding. The types can be used to still dispatch on also the undecorated objective type O2.

Fields

  • objective: the objective that is defined in the embedding
  • p=nothing: a point in the embedding.
  • X=nothing: a tangent vector in the embedding

When a point in the embedding p is provided, embed! is used in place of this point to reduce memory allocations. Similarly X is used when embedding tangent vectors

source

Cache objective

Since single function calls, for example to the cost or the gradient, might be expensive, a simple cache objective exists as a decorator, that caches one cost value or gradient.

It can be activated/used with the cache= keyword argument available for every solver.

Manopt.reset_counters!Function
reset_counters(co::ManifoldCountObjective, value::Integer=0)

Reset all values in the count objective to value.

source
Manopt.objective_cache_factoryFunction
objective_cache_factory(M::AbstractManifold, o::AbstractManifoldObjective, cache::Symbol)

Generate a cached variant of the AbstractManifoldObjective o on the AbstractManifold M based on the symbol cache.

The following caches are available

  • :Simple generates a SimpleManifoldCachedObjective
  • :LRU generates a ManifoldCachedObjective where you should use the form (:LRU, [:Cost, :Gradient]) to specify what should be cached or (:LRU, [:Cost, :Gradient], 100) to specify the cache size. Here this variant defaults to (:LRU, [:Cost, :Gradient], 100), caching up to 100 cost and gradient values.[1]
source
objective_cache_factory(M::AbstractManifold, o::AbstractManifoldObjective, cache::Tuple{Symbol, Array, Array})
objective_cache_factory(M::AbstractManifold, o::AbstractManifoldObjective, cache::Tuple{Symbol, Array})

Generate a cached variant of the AbstractManifoldObjective o on the AbstractManifold M based on the symbol cache[1], where the second element cache[2] are further arguments to the cache and the optional third is passed down as keyword arguments.

For all available caches see the simpler variant with symbols.

source

A simple cache

A first generic cache is always available, but it only caches one gradient and one cost function evaluation (for the same point).

Manopt.SimpleManifoldCachedObjectiveType
 SimpleManifoldCachedObjective{O<:AbstractManifoldGradientObjective{E,TC,TG}, P, T,C} <: AbstractManifoldGradientObjective{E,TC,TG}

Provide a simple cache for an AbstractManifoldGradientObjective that is for a given point p this cache stores a point p and a gradient $\operatorname{grad} f(p)$ in X as well as a cost value $f(p)$ in c.

Both X and c are accompanied by booleans to keep track of their validity.

Constructor

SimpleManifoldCachedObjective(M::AbstractManifold, obj::AbstractManifoldGradientObjective; kwargs...)

Keyword arguments

  • p=rand(M): a point on the manifold to initialize the cache with
  • X=get_gradient(M, obj, p) or zero_vector(M,p): a tangent vector to store the gradient in, see also initialize=
  • c=[get_cost](@ref)(M, obj, p)or0.0: a value to store the cost function ininitialize`
  • initialized=true: whether to initialize the cached X and c or not.
source

A generic cache

For the more advanced cache, you need to implement some type of cache yourself, that provides a get! and implement init_caches. This is for example provided if you load LRUCache.jl. Then you obtain

Manopt.ManifoldCachedObjectiveType
ManifoldCachedObjective{E,P,O<:AbstractManifoldObjective{<:E},C<:NamedTuple{}} <: AbstractDecoratedManifoldObjective{E,P}

Create a cache for an objective, based on a NamedTuple that stores some kind of cache.

Constructor

ManifoldCachedObjective(M, o::AbstractManifoldObjective, caches::Vector{Symbol}; kwargs...)

Create a cache for the AbstractManifoldObjective where the Symbols in caches indicate, which function evaluations to cache.

Supported symbols

SymbolCaches calls to (incl. ! variants)Comment
:Costget_cost
:EqualityConstraintget_equality_constraint(M, p, i)
:EqualityConstraintsget_equality_constraint(M, p, :)
:GradEqualityConstraintget_grad_equality_constrainttangent vector per (p,i)
:GradInequalityConstraintget_inequality_constrainttangent vector per (p,i)
:Gradientget_gradient(M,p)tangent vectors
:Hessianget_hessiantangent vectors
:InequalityConstraintget_inequality_constraint(M, p, j)
:InequalityConstraintsget_inequality_constraint(M, p, :)
:Preconditionerget_preconditionertangent vectors
:ProximalMapget_proximal_mappoint per (p,λ,i)
:StochasticGradientsget_gradientsvector of tangent vectors
:StochasticGradientget_gradient(M, p, i)tangent vector per (p,i)
:SubGradientget_subgradienttangent vectors
:SubtrahendGradientget_subtrahend_gradienttangent vectors

Keyword arguments

  • p=rand(M): the type of the keys to be used in the caches. Defaults to the default representation on M.
  • value=get_cost(M, objective, p): the type of values for numeric values in the cache
  • X=zero_vector(M,p): the type of values to be cached for gradient and Hessian calls.
  • cache=[:Cost]: a vector of symbols indicating which function calls should be cached.
  • cache_size=10: number of (least recently used) calls to cache
  • cache_sizes=Dict{Symbol,Int}(): a named tuple or dictionary specifying the sizes individually for each cache.
source
Manopt.init_cachesFunction
init_caches(caches, T::Type{LRU}; kwargs...)

Given a vector of symbols caches, this function sets up the NamedTuple of caches, where T is the type of cache to use.

Keyword arguments

  • p=rand(M): a point on a manifold, to both infer its type for keys and initialize caches
  • value=0.0: a value both typing and initialising number-caches, the default is for (Float) values like the cost.
  • X=zero_vector(M, p): a tangent vector at p to both type and initialize tangent vector caches
  • cache_size=10: a default cache size to use
  • cache_sizes=Dict{Symbol,Int}(): a dictionary of sizes for the caches to specify different (non-default) sizes
source
init_caches(M::AbstractManifold, caches, T; kwargs...)

Given a vector of symbols caches, this function sets up the NamedTuple of caches for points/vectors on M, where T is the type of cache to use.

source

Count objective

Manopt.ManifoldCountObjectiveType
ManifoldCountObjective{E,P,O<:AbstractManifoldObjective,I<:Integer} <: AbstractDecoratedManifoldObjective{E,P}

A wrapper for any AbstractManifoldObjective of type O to count different calls to parts of the objective.

Fields

  • counts a dictionary of symbols mapping to integers keeping the counted values
  • objective the wrapped objective

Supported symbols

SymbolCounts calls to (incl. ! variants)Comment
:Costget_cost
:EqualityConstraintget_equality_constraintrequires vector of counters
:EqualityConstraintsget_equality_constraintwhen evaluating all of them with :
:GradEqualityConstraintget_grad_equality_constraintrequires vector of counters
:GradEqualityConstraintsget_grad_equality_constraintwhen evaluating all of them with :
:GradInequalityConstraintget_inequality_constraintrequires vector of counters
:GradInequalityConstraintsget_inequality_constraintwhen evaluating all of them with :
:Gradientget_gradient(M,p)
:Hessianget_hessian
:InequalityConstraintget_inequality_constraintrequires vector of counters
:InequalityConstraintsget_inequality_constraintwhen evaluating all of them with :
:Preconditionerget_preconditioner
:ProximalMapget_proximal_map
:StochasticGradientsget_gradients
:StochasticGradientget_gradient(M, p, i)
:SubGradientget_subgradient
:SubtrahendGradientget_subtrahend_gradient

Constructors

ManifoldCountObjective(objective::AbstractManifoldObjective, counts::Dict{Symbol, <:Integer})

Initialise the ManifoldCountObjective to wrap objective initializing the set of counts

ManifoldCountObjective(M::AbstractManifold, objective::AbstractManifoldObjective, count::AbstractVecor{Symbol}, init=0)

Count function calls on objective using the symbols in count initialising all entries to init.

source

Internal decorators

Manopt.ReturnManifoldObjectiveType
ReturnManifoldObjective{E,O2,O1<:AbstractManifoldObjective{E}} <:
   AbstractDecoratedManifoldObjective{E,O2}

A wrapper to indicate that get_solver_result should return the inner objective.

The types are such that one can still dispatch on the undecorated type O2 of the original objective as well.

source

Specific Objective typed and their access functions

Cost objective

Manopt.ManifoldCostObjectiveType
ManifoldCostObjective{T, TC} <: AbstractManifoldCostObjective{T, TC}

specify an AbstractManifoldObjective that does only have information about the cost function $f: \mathbb M → ℝ$ implemented as a function (M, p) -> c to compute the cost value c at p on the manifold M.

  • cost: a function $f: \mathcal M → ℝ$ to minimize

Constructors

ManifoldCostObjective(f)

Generate a problem. While this Problem does not have any allocating functions, the type T can be set for consistency reasons with other problems.

Used with

NelderMead, particle_swarm

source

Access functions

Manopt.get_costFunction
get_cost(amp::AbstractManoptProblem, p)

evaluate the cost function f stored within the AbstractManifoldObjective of an AbstractManoptProblem amp at the point p.

source
get_cost(M::AbstractManifold, obj::AbstractManifoldObjective, p)

evaluate the cost function f defined on M stored within the AbstractManifoldObjective at the point p.

source
get_cost(M::AbstractManifold, mco::AbstractManifoldCostObjective, p)

Evaluate the cost function from within the AbstractManifoldCostObjective on M at p.

By default this implementation assumed that the cost is stored within mco.cost.

source
get_cost(TpM, trmo::TrustRegionModelObjective, X)

Evaluate the tangent space TrustRegionModelObjective

\[m(X) = f(p) + ⟨\operatorname{grad} f(p), X ⟩_p + \frac{1}{2} ⟨\operatorname{Hess} f(p)[X], X⟩_p.\]

source
get_cost(TpM, trmo::AdaptiveRagularizationWithCubicsModelObjective, X)

Evaluate the tangent space AdaptiveRagularizationWithCubicsModelObjective

\[m(X) = f(p) + ⟨\operatorname{grad} f(p), X ⟩_p + \frac{1}{2} ⟨\operatorname{Hess} f(p)[X], X⟩_p + \frac{σ}{3} \lVert X \rVert^3,\]

at X, cf. Eq. (33) in [ABBC20].

source
get_cost(TpM::TangentSpace, slso::SymmetricLinearSystemObjective, X)

evaluate the cost

\[f(X) = \frac{1}{2} \lVert \mathcal A[X] + b \rVert_{p}^2,\qquad X ∈ T_{p}\mathcal M,\]

at X.

source
get_cost(M::AbstractManifold, sgo::ManifoldStochasticGradientObjective, p, i)

Evaluate the ith summand of the cost.

If you use a single function for the stochastic cost, then only the index ì=1` is available to evaluate the whole cost.

source
get_cost(M::AbstractManifold,emo::EmbeddedManifoldObjective, p)

Evaluate the cost function of an objective defined in the embedding by first embedding p before calling the cost function stored in the EmbeddedManifoldObjective.

source

and internally

Manopt.get_cost_functionFunction
get_cost_function(amco::AbstractManifoldCostObjective)

return the function to evaluate (just) the cost $f(p)=c$ as a function (M,p) -> c.

source

Gradient objectives

Manopt.ManifoldGradientObjectiveType
ManifoldGradientObjective{T<:AbstractEvaluationType} <: AbstractManifoldGradientObjective{T}

specify an objective containing a cost and its gradient

Fields

  • cost: a function $f: \mathcal M → ℝ$
  • gradient!!: the gradient $\operatorname{grad}f: \mathcal M → \mathcal T\mathcal M$ of the cost function $f$.

Depending on the AbstractEvaluationType T the gradient can have to forms

Constructors

ManifoldGradientObjective(cost, gradient; evaluation=AllocatingEvaluation())

Used with

gradient_descent, conjugate_gradient_descent, quasi_Newton

source
Manopt.ManifoldAlternatingGradientObjectiveType
ManifoldAlternatingGradientObjective{E<:AbstractEvaluationType,TCost,TGradient} <: AbstractManifoldGradientObjective{E}

An alternating gradient objective consists of

  • a cost function $F(x)$
  • a gradient $\operatorname{grad}F$ that is either
    • given as one function $\operatorname{grad}F$ returning a tangent vector X on M or
    • an array of gradient functions $\operatorname{grad}F_i$, ì=1,…,n s each returning a component of the gradient
    which might be allocating or mutating variants, but not a mix of both.
Note

This Objective is usually defined using the ProductManifold from Manifolds.jl, so Manifolds.jl to be loaded.

Constructors

ManifoldAlternatingGradientObjective(F, gradF::Function;
    evaluation=AllocatingEvaluation()
)
ManifoldAlternatingGradientObjective(F, gradF::AbstractVector{<:Function};
    evaluation=AllocatingEvaluation()
)

Create a alternating gradient problem with an optional cost and the gradient either as one function (returning an array) or a vector of functions.

source
Manopt.ManifoldStochasticGradientObjectiveType
ManifoldStochasticGradientObjective{T<:AbstractEvaluationType} <: AbstractManifoldGradientObjective{T}

A stochastic gradient objective consists of

  • a(n optional) cost function $f(p) = \displaystyle\sum_{i=1}^n f_i(p)$
  • an array of gradients, $\operatorname{grad}f_i(p), i=1,\ldots,n$ which can be given in two forms
    • as one single function $(\mathcal M, p) ↦ (X_1,…,X_n) ∈ (T_p\mathcal M)^n$
    • as a vector of functions $\bigl( (\mathcal M, p) ↦ X_1, …, (\mathcal M, p) ↦ X_n\bigr)$.

Where both variants can also be provided as InplaceEvaluation functions (M, X, p) -> X, where X is the vector of X1,...,Xn and (M, X1, p) -> X1, ..., (M, Xn, p) -> Xn, respectively.

Constructors

ManifoldStochasticGradientObjective(
    grad_f::Function;
    cost=Missing(),
    evaluation=AllocatingEvaluation()
)
ManifoldStochasticGradientObjective(
    grad_f::AbstractVector{<:Function};
    cost=Missing(), evaluation=AllocatingEvaluation()
)

Create a Stochastic gradient problem with the gradient either as one function (returning an array of tangent vectors) or a vector of functions (each returning one tangent vector).

The optional cost can also be given as either a single function (returning a number) pr a vector of functions, each returning a value.

Used with

stochastic_gradient_descent

Note that this can also be used with a gradient_descent, since the (complete) gradient is just the sums of the single gradients.

source
Manopt.NonlinearLeastSquaresObjectiveType
NonlinearLeastSquaresObjective{T<:AbstractEvaluationType} <: AbstractManifoldObjective{T}

A type for nonlinear least squares problems. T is a AbstractEvaluationType for the F and Jacobian functions.

Specify a nonlinear least squares problem

Fields

  • f a function $f: \mathcal M → ℝ^d$ to minimize
  • jacobian!! Jacobian of the function $f$
  • jacobian_tangent_basis the basis of tangent space used for computing the Jacobian.
  • num_components number of values returned by f (equal to d).

Depending on the AbstractEvaluationType T the function $F$ has to be provided:

  • as a functions (M::AbstractManifold, p) -> v that allocates memory for v itself for an AllocatingEvaluation,
  • as a function (M::AbstractManifold, v, p) -> v that works in place of v for a InplaceEvaluation.

Also the Jacobian $jacF!!$ is required:

  • as a functions (M::AbstractManifold, p; basis_domain::AbstractBasis) -> v that allocates memory for v itself for an AllocatingEvaluation,
  • as a function (M::AbstractManifold, v, p; basis_domain::AbstractBasis) -> v that works in place of v for an InplaceEvaluation.

Constructors

NonlinearLeastSquaresProblem(M, F, jacF, num_components; evaluation=AllocatingEvaluation(), jacobian_tangent_basis=DefaultOrthonormalBasis())

See also

LevenbergMarquardt, LevenbergMarquardtState

source

There is also a second variant, if just one function is responsible for computing the cost and the gradient

Manopt.ManifoldCostGradientObjectiveType
ManifoldCostGradientObjective{T} <: AbstractManifoldObjective{T}

specify an objective containing one function to perform a combined computation of cost and its gradient

Fields

  • costgrad!!: a function that computes both the cost $f: \mathcal M → ℝ$ and its gradient $\operatorname{grad}f: \mathcal M → \mathcal T\mathcal M$

Depending on the AbstractEvaluationType T the gradient can have to forms

Constructors

ManifoldCostGradientObjective(costgrad; evaluation=AllocatingEvaluation())

Used with

gradient_descent, conjugate_gradient_descent, quasi_Newton

source

Access functions

Manopt.get_gradientFunction
get_gradient(s::AbstractManoptSolverState)

return the (last stored) gradient within AbstractManoptSolverStates`. By default also undecorates the state beforehand

source
get_gradient(amp::AbstractManoptProblem, p)
get_gradient!(amp::AbstractManoptProblem, X, p)

evaluate the gradient of an AbstractManoptProblem amp at the point p.

The evaluation is done in place of X for the !-variant.

source
get_gradient(M::AbstractManifold, mgo::AbstractManifoldGradientObjective{T}, p)
get_gradient!(M::AbstractManifold, X, mgo::AbstractManifoldGradientObjective{T}, p)

evaluate the gradient of a AbstractManifoldGradientObjective{T} mgo at p.

The evaluation is done in place of X for the !-variant. The T=AllocatingEvaluation problem might still allocate memory within. When the non-mutating variant is called with a T=InplaceEvaluation memory for the result is allocated.

Note that the order of parameters follows the philosophy of Manifolds.jl, namely that even for the mutating variant, the manifold is the first parameter and the (in-place) tangent vector X comes second.

source
get_gradient(agst::AbstractGradientSolverState)

return the gradient stored within gradient options. THe default returns agst.X.

source
get_gradient(M::AbstractManifold, vgf::VectorGradientFunction, p, i)
get_gradient(M::AbstractManifold, vgf::VectorGradientFunction, p, i, range)
get_gradient!(M::AbstractManifold, X, vgf::VectorGradientFunction, p, i)
get_gradient!(M::AbstractManifold, X, vgf::VectorGradientFunction, p, i, range)

Evaluate the gradients of the vector function vgf on the manifold M at p and the values given in range, specifying the representation of the gradients.

Since i is assumed to be a linear index, you can provide

  • a single integer
  • a UnitRange to specify a range to be returned like 1:3
  • a BitVector specifying a selection
  • a AbstractVector{<:Integer} to specify indices
  • : to return the vector of all gradients
source
get_gradient(TpM, trmo::TrustRegionModelObjective, X)

Evaluate the gradient of the TrustRegionModelObjective

\[\operatorname{grad} m(X) = \operatorname{grad} f(p) + \operatorname{Hess} f(p)[X].\]

source
get_gradient(TpM, trmo::AdaptiveRagularizationWithCubicsModelObjective, X)

Evaluate the gradient of the AdaptiveRagularizationWithCubicsModelObjective

\[\operatorname{grad} m(X) = \operatorname{grad} f(p) + \operatorname{Hess} f(p)[X] + σ\lVert X \rVert X,\]

at X, cf. Eq. (37) in [ABBC20].

source
get_gradient(TpM::TangentSpace, slso::SymmetricLinearSystemObjective, X)
get_gradient!(TpM::TangentSpace, Y, slso::SymmetricLinearSystemObjective, X)

evaluate the gradient of

\[f(X) = \frac{1}{2} \lVert \mathcal A[X] + b \rVert_{p}^2,\qquad X ∈ T_{p}\mathcal M,\]

Which is $\operatorname{grad} f(X) = \mathcal A[X]+b$. This can be computed in-place of Y.

source
get_gradient(M::AbstractManifold, sgo::ManifoldStochasticGradientObjective, p, k)
get_gradient!(M::AbstractManifold, sgo::ManifoldStochasticGradientObjective, Y, p, k)

Evaluate one of the summands gradients $\operatorname{grad}f_k$, $k∈\{1,…,n\}$, at x (in place of Y).

If you use a single function for the stochastic gradient, that works in-place, then get_gradient is not available, since the length (or number of elements of the gradient required for allocation) can not be determined.

source
get_gradient(M::AbstractManifold, sgo::ManifoldStochasticGradientObjective, p)
get_gradient!(M::AbstractManifold, sgo::ManifoldStochasticGradientObjective, X, p)

Evaluate the complete gradient $\operatorname{grad} f = \displaystyle\sum_{i=1}^n \operatorname{grad} f_i(p)$ at p (in place of X).

If you use a single function for the stochastic gradient, that works in-place, then get_gradient is not available, since the length (or number of elements of the gradient required for allocation) can not be determined.

source
get_gradient(M::AbstractManifold, emo::EmbeddedManifoldObjective, p)
get_gradient!(M::AbstractManifold, X, emo::EmbeddedManifoldObjective, p)

Evaluate the gradient function of an objective defined in the embedding, that is embed p before calling the gradient function stored in the EmbeddedManifoldObjective.

The returned gradient is then converted to a Riemannian gradient calling riemannian_gradient.

source
Manopt.get_gradientsFunction
get_gradients(M::AbstractManifold, sgo::ManifoldStochasticGradientObjective, p)
get_gradients!(M::AbstractManifold, X, sgo::ManifoldStochasticGradientObjective, p)

Evaluate all summands gradients $\{\operatorname{grad}f_i\}_{i=1}^n$ at p (in place of X).

If you use a single function for the stochastic gradient, that works in-place, then get_gradient is not available, since the length (or number of elements of the gradient) can not be determined.

source

and internally

Manopt.get_gradient_functionFunction
get_gradient_function(amgo::AbstractManifoldGradientObjective, recursive=false)

return the function to evaluate (just) the gradient $\operatorname{grad} f(p)$, where either the gradient function using the decorator or without the decorator is used.

By default recursive is set to false, since usually to just pass the gradient function somewhere, one still wants for example the cached one or the one that still counts calls.

Depending on the AbstractEvaluationType E this is a function

source

Internal helpers

Manopt.get_gradient_from_Jacobian!Function
get_gradient_from_Jacobian!(
    M::AbstractManifold,
    X,
    nlso::NonlinearLeastSquaresObjective{InplaceEvaluation},
    p,
    Jval=zeros(nlso.num_components, manifold_dimension(M)),
)

Compute gradient of NonlinearLeastSquaresObjective nlso at point p in place of X, with temporary Jacobian stored in the optional argument Jval.

source

Subgradient objective

Manopt.ManifoldSubgradientObjectiveType
ManifoldSubgradientObjective{T<:AbstractEvaluationType,C,S} <:AbstractManifoldCostObjective{T, C}

A structure to store information about a objective for a subgradient based optimization problem

Fields

  • cost: the function $f$ to be minimized
  • subgradient: a function returning a subgradient $∂f$ of $f$

Constructor

ManifoldSubgradientObjective(f, ∂f)

Generate the ManifoldSubgradientObjective for a subgradient objective, consisting of a (cost) function f(M, p) and a function ∂f(M, p) that returns a not necessarily deterministic element from the subdifferential at p on a manifold M.

source

Access functions

Manopt.get_subgradientFunction
X = get_subgradient(M::AbstractManifold, sgo::AbstractManifoldGradientObjective, p)
get_subgradient!(M::AbstractManifold, X, sgo::AbstractManifoldGradientObjective, p)

Evaluate the subgradient, which for the case of a objective having a gradient, means evaluating the gradient itself.

While in general, the result might not be deterministic, for this case it is.

source
get_subgradient(amp::AbstractManoptProblem, p)
get_subgradient!(amp::AbstractManoptProblem, X, p)

evaluate the subgradient of an AbstractManoptProblem amp at point p.

The evaluation is done in place of X for the !-variant. The result might not be deterministic, one element of the subdifferential is returned.

source
X = get_subgradient(M;;AbstractManifold, sgo::ManifoldSubgradientObjective, p)
get_subgradient!(M;;AbstractManifold, X, sgo::ManifoldSubgradientObjective, p)

Evaluate the (sub)gradient of a ManifoldSubgradientObjective sgo at the point p.

The evaluation is done in place of X for the !-variant. The result might not be deterministic, one element of the subdifferential is returned.

source

Proximal map objective

Manopt.ManifoldProximalMapObjectiveType
ManifoldProximalMapObjective{E<:AbstractEvaluationType, TC, TP, V <: Vector{<:Integer}} <: AbstractManifoldCostObjective{E, TC}

specify a problem for solvers based on the evaluation of proximal maps, which represents proximal maps $\operatorname{prox}_{λf_i}$ for summands $f = f_1 + f_2+ … + f_N$ of the cost function $f$.

Fields

  • cost: a function $f:\mathcal M→ℝ$ to minimize
  • proxes: proximal maps $\operatorname{prox}_{λf_i}:\mathcal M → \mathcal M$ as functions (M, λ, p) -> q or in-place (M, q, λ, p).
  • number_of_proxes: number of proximal maps per function, to specify when one of the maps is a combined one such that the proximal maps functions return more than one entry per function, you have to adapt this value. if not specified, it is set to one prox per function.

Constructor

ManifoldProximalMapObjective(f, proxes_f::Union{Tuple,AbstractVector}, numer_of_proxes=onex(length(proxes));
   evaluation=Allocating)

Generate a proximal problem with a tuple or vector of funtions, where by default every function computes a single prox of one component of $f$.

ManifoldProximalMapObjective(f, prox_f); evaluation=Allocating)

Generate a proximal objective for $f$ and its proxial map $\operatorname{prox}_{λf}$

See also

cyclic_proximal_point, get_cost, get_proximal_map

source

Access functions

Manopt.get_proximal_mapFunction
q = get_proximal_map(M::AbstractManifold, mpo::ManifoldProximalMapObjective, λ, p)
get_proximal_map!(M::AbstractManifold, q, mpo::ManifoldProximalMapObjective, λ, p)
q = get_proximal_map(M::AbstractManifold, mpo::ManifoldProximalMapObjective, λ, p, i)
get_proximal_map!(M::AbstractManifold, q, mpo::ManifoldProximalMapObjective, λ, p, i)

evaluate the (ith) proximal map of ManifoldProximalMapObjective p at the point p of p.M with parameter $λ>0$.

source

Hessian objective

Manopt.ManifoldHessianObjectiveType
ManifoldHessianObjective{T<:AbstractEvaluationType,C,G,H,Pre} <: AbstractManifoldHessianObjective{T,C,G,H}

specify a problem for Hessian based algorithms.

Fields

  • cost: a function $f:\mathcal M→ℝ$ to minimize
  • gradient: the gradient $\operatorname{grad}f:\mathcal M → \mathcal T\mathcal M$ of the cost function $f$
  • hessian: the Hessian $\operatorname{Hess}f(x)[⋅]: \mathcal T_{x} \mathcal M → \mathcal T_{x} \mathcal M$ of the cost function $f$
  • preconditioner: the symmetric, positive definite preconditioner as an approximation of the inverse of the Hessian of $f$, a map with the same input variables as the hessian to numerically stabilize iterations when the Hessian is ill-conditioned

Depending on the AbstractEvaluationType T the gradient and can have to forms

Constructor

ManifoldHessianObjective(f, grad_f, Hess_f, preconditioner = (M, p, X) -> X;
    evaluation=AllocatingEvaluation())

See also

truncated_conjugate_gradient_descent, trust_regions

source

Access functions

Manopt.get_hessianFunction
Y = get_hessian(amp::AbstractManoptProblem{T}, p, X)
get_hessian!(amp::AbstractManoptProblem{T}, Y, p, X)

evaluate the Hessian of an AbstractManoptProblem amp at p applied to a tangent vector X, computing $\operatorname{Hess}f(q)[X]$, which can also happen in-place of Y.

source
get_hessian(M::AbstractManifold, vgf::VectorHessianFunction, p, X, i)
get_hessian(M::AbstractManifold, vgf::VectorHessianFunction, p, X, i, range)
get_hessian!(M::AbstractManifold, X, vgf::VectorHessianFunction, p, X, i)
get_hessian!(M::AbstractManifold, X, vgf::VectorHessianFunction, p, X, i, range)

Evaluate the Hessians of the vector function vgf on the manifold M at p in direction X and the values given in range, specifying the representation of the gradients.

Since i is assumed to be a linear index, you can provide

  • a single integer
  • a UnitRange to specify a range to be returned like 1:3
  • a BitVector specifying a selection
  • a AbstractVector{<:Integer} to specify indices
  • : to return the vector of all gradients
source
get_hessian(TpM, trmo::TrustRegionModelObjective, X)

Evaluate the Hessian of the TrustRegionModelObjective

\[\operatorname{Hess} m(X)[Y] = \operatorname{Hess} f(p)[Y].\]

source
get_Hessian(TpM::TangentSpace, slso::SymmetricLinearSystemObjective, X, V)
get_Hessian!(TpM::TangentSpace, W, slso::SymmetricLinearSystemObjective, X, V)

evaluate the Hessian of

\[f(X) = \frac{1}{2} \lVert \mathcal A[X] + b \rVert_{p}^2,\qquad X ∈ T_{p}\mathcal M,\]

Which is $\operatorname{Hess} f(X)[Y] = \mathcal A[V]$. This can be computed in-place of W.

source
get_hessian(M::AbstractManifold, emo::EmbeddedManifoldObjective, p, X)
get_hessian!(M::AbstractManifold, Y, emo::EmbeddedManifoldObjective, p, X)

Evaluate the Hessian of an objective defined in the embedding, that is embed p and X before calling the Hessian function stored in the EmbeddedManifoldObjective.

The returned Hessian is then converted to a Riemannian Hessian calling riemannian_Hessian.

source
Manopt.get_preconditionerFunction
get_preconditioner(amp::AbstractManoptProblem, p, X)

evaluate the symmetric, positive definite preconditioner (approximation of the inverse of the Hessian of the cost function f) of a AbstractManoptProblem amps objective at the point p applied to a tangent vector X.

source
get_preconditioner(M::AbstractManifold, mho::ManifoldHessianObjective, p, X)

evaluate the symmetric, positive definite preconditioner (approximation of the inverse of the Hessian of the cost function F) of a ManifoldHessianObjective mho at the point p applied to a tangent vector X.

source

and internally

Primal-dual based objectives

Manopt.AbstractPrimalDualManifoldObjectiveType
AbstractPrimalDualManifoldObjective{E<:AbstractEvaluationType,C,P} <: AbstractManifoldCostObjective{E,C}

A common abstract super type for objectives that consider primal-dual problems.

source
Manopt.PrimalDualManifoldObjectiveType
PrimalDualManifoldObjective{T<:AbstractEvaluationType} <: AbstractPrimalDualManifoldObjective{T}

Describes an Objective linearized or exact Chambolle-Pock algorithm, cf. [BHS+21], [CP11]

Fields

All fields with !! can either be in-place or allocating functions, which should be set depending on the evaluation= keyword in the constructor and stored in T <: AbstractEvaluationType.

  • cost: $F + G(Λ(⋅))$ to evaluate interim cost function values
  • linearized_forward_operator!!: linearized operator for the forward operation in the algorithm $DΛ$
  • linearized_adjoint_operator!!: the adjoint differential $(DΛ)^* : \mathcal N → T\mathcal M$
  • prox_f!!: the proximal map belonging to $f$
  • prox_G_dual!!: the proximal map belonging to $g_n^*$
  • Λ!!: the forward operator (if given) $Λ: \mathcal M → \mathcal N$

Either the linearized operator $DΛ$ or $Λ$ are required usually.

Constructor

PrimalDualManifoldObjective(cost, prox_f, prox_G_dual, adjoint_linearized_operator;
    linearized_forward_operator::Union{Function,Missing}=missing,
    Λ::Union{Function,Missing}=missing,
    evaluation::AbstractEvaluationType=AllocatingEvaluation()
)

The last optional argument can be used to provide the 4 or 5 functions as allocating or mutating (in place computation) ones. Note that the first argument is always the manifold under consideration, the mutated one is the second.

source
Manopt.PrimalDualManifoldSemismoothNewtonObjectiveType
PrimalDualManifoldSemismoothNewtonObjective{E<:AbstractEvaluationType, TC, LO, ALO, PF, DPF, PG, DPG, L} <: AbstractPrimalDualManifoldObjective{E, TC, PF}

Describes a Problem for the Primal-dual Riemannian semismooth Newton algorithm. [DL21]

Fields

  • cost: $F + G(Λ(⋅))$ to evaluate interim cost function values
  • linearized_operator: the linearization $DΛ(⋅)[⋅]$ of the operator $Λ(⋅)$.
  • linearized_adjoint_operator: the adjoint differential $(DΛ)^* : \mathcal N → T\mathcal M$
  • prox_F: the proximal map belonging to $F$
  • diff_prox_F: the (Clarke Generalized) differential of the proximal maps of $F$
  • prox_G_dual: the proximal map belonging to G^\ast_n`
  • diff_prox_dual_G: the (Clarke Generalized) differential of the proximal maps of $G^\ast_n$
  • Λ: the exact forward operator. This operator is required if Λ(m)=n does not hold.

Constructor

PrimalDualManifoldSemismoothNewtonObjective(cost, prox_F, prox_G_dual, forward_operator, adjoint_linearized_operator,Λ)
source

Access functions

Manopt.adjoint_linearized_operatorFunction
X = adjoint_linearized_operator(N::AbstractManifold, apdmo::AbstractPrimalDualManifoldObjective, m, n, Y)
adjoint_linearized_operator(N::AbstractManifold, X, apdmo::AbstractPrimalDualManifoldObjective, m, n, Y)

Evaluate the adjoint of the linearized forward operator of $(DΛ(m))^*[Y]$ stored within the AbstractPrimalDualManifoldObjective (in place of X). Since $Y∈T_n\mathcal N$, both $m$ and $n=Λ(m)$ are necessary arguments, mainly because the forward operator $Λ$ might be missing in p.

source
Manopt.forward_operatorFunction
q = forward_operator(M::AbstractManifold, N::AbstractManifold, apdmo::AbstractPrimalDualManifoldObjective, p)
forward_operator!(M::AbstractManifold, N::AbstractManifold, q, apdmo::AbstractPrimalDualManifoldObjective, p)

Evaluate the forward operator of $Λ(x)$ stored within the TwoManifoldProblem (in place of q).

source
Manopt.get_differential_dual_proxFunction
η = get_differential_dual_prox(N::AbstractManifold, pdsno::PrimalDualManifoldSemismoothNewtonObjective, n, τ, X, ξ)
get_differential_dual_prox!(N::AbstractManifold, pdsno::PrimalDualManifoldSemismoothNewtonObjective, η, n, τ, X, ξ)

Evaluate the differential proximal map of $G_n^*$ stored within PrimalDualManifoldSemismoothNewtonObjective

\[D\operatorname{prox}_{τG_n^*}(X)[ξ]\]

which can also be computed in place of η.

source
Manopt.get_differential_primal_proxFunction
y = get_differential_primal_prox(M::AbstractManifold, pdsno::PrimalDualManifoldSemismoothNewtonObjective σ, x)
get_differential_primal_prox!(p::TwoManifoldProblem, y, σ, x)

Evaluate the differential proximal map of $F$ stored within AbstractPrimalDualManifoldObjective

\[D\operatorname{prox}_{σF}(x)[X]\]

which can also be computed in place of y.

source
Manopt.get_dual_proxFunction
Y = get_dual_prox(N::AbstractManifold, apdmo::AbstractPrimalDualManifoldObjective, n, τ, X)
get_dual_prox!(N::AbstractManifold, apdmo::AbstractPrimalDualManifoldObjective, Y, n, τ, X)

Evaluate the proximal map of $g_n^*$ stored within AbstractPrimalDualManifoldObjective

\[ Y = \operatorname{prox}}_{τG_n^*}(X)\]

which can also be computed in place of Y.

source
Manopt.get_primal_proxFunction
q = get_primal_prox(M::AbstractManifold, p::AbstractPrimalDualManifoldObjective, σ, p)
get_primal_prox!(M::AbstractManifold, p::AbstractPrimalDualManifoldObjective, q, σ, p)

Evaluate the proximal map of $F$ stored within AbstractPrimalDualManifoldObjective

\[\operatorname{prox}_{σF}(x)\]

which can also be computed in place of y.

source
Manopt.linearized_forward_operatorFunction
Y = linearized_forward_operator(M::AbstractManifold, N::AbstractManifold, apdmo::AbstractPrimalDualManifoldObjective, m, X, n)
linearized_forward_operator!(M::AbstractManifold, N::AbstractManifold, Y, apdmo::AbstractPrimalDualManifoldObjective, m, X, n)

Evaluate the linearized operator (differential) $DΛ(m)[X]$ stored within the AbstractPrimalDualManifoldObjective (in place of Y), where n = Λ(m).

source

Constrained objective

Manopt.ConstrainedManifoldObjectiveType
ConstrainedManifoldObjective{T<:AbstractEvaluationType, C<:ConstraintType} <: AbstractManifoldObjective{T}

Describes the constrained objective

\[\begin{aligned} \operatorname*{arg\,min}_{p ∈\mathcal{M}} & f(p)\\ \text{subject to } &g_i(p)\leq0 \quad \text{ for all } i=1,…,m,\\ \quad &h_j(p)=0 \quad \text{ for all } j=1,…,n. \end{aligned}\]

Fields

  • objective: an AbstractManifoldObjective representing the unconstrained objective, that is containing cost $f$, the gradient of the cost $f$ and maybe the Hessian.
  • equality_constraints: an AbstractManifoldObjective representing the equality constraints

$h: \mathcal M → \mathbb R^n$ also possibly containing its gradient and/or Hessian

$h: \mathcal M → \mathbb R^n$ also possibly containing its gradient and/or Hessian

Constructors

ConstrainedManifoldObjective(M::AbstractManifold, f, grad_f;
    g=nothing,
    grad_g=nothing,
    h=nothing,
    grad_h=nothing;
    hess_f=nothing,
    hess_g=nothing,
    hess_h=nothing,
    equality_constraints=nothing,
    inequality_constraints=nothing,
    evaluation=AllocatingEvaluation(),
    M = nothing,
    p = isnothing(M) ? nothing : rand(M),
)

Generate the constrained objective based on all involved single functions f, grad_f, g, grad_g, h, grad_h, and optionally a Hessian for each of these. With equality_constraints and inequality_constraints you have to provide the dimension of the ranges of h and g, respectively. You can also provide a manifold M and a point p to use one evaluation of the constraints to automatically try to determine these sizes.

ConstrainedManifoldObjective(M::AbstractManifold, mho::AbstractManifoldObjective;
    equality_constraints = nothing,
    inequality_constraints = nothing
)

Generate the constrained objective either with explicit constraints $g$ and $h$, and their gradients, or in the form where these are already encapsulated in VectorGradientFunctions.

Both variants require that at least one of the constraints (and its gradient) is provided. If any of the three parts provides a Hessian, the corresponding object, that is a ManifoldHessianObjective for f or a VectorHessianFunction for g or h, respectively, is created.

source

It might be beneficial to use the adapted problem to specify different ranges for the gradients of the constraints

Manopt.ConstrainedManoptProblemType
ConstrainedProblem{
    TM <: AbstractManifold,
    O <: AbstractManifoldObjective
    HR<:Union{AbstractPowerRepresentation,Nothing},
    GR<:Union{AbstractPowerRepresentation,Nothing},
    HHR<:Union{AbstractPowerRepresentation,Nothing},
    GHR<:Union{AbstractPowerRepresentation,Nothing},
} <: AbstractManoptProblem{TM}

A constrained problem might feature different ranges for the (vectors of) gradients of the equality and inequality constraints.

The ranges are required in a few places to allocate memory and access elements correctly, they work as follows:

Assume the objective is

\[\begin{aligned} \operatorname*{arg\,min}_{p ∈\mathcal{M}} & f(p)\\ \text{subject to } &g_i(p)\leq0 \quad \text{ for all } i=1,…,m,\\ \quad &h_j(p)=0 \quad \text{ for all } j=1,…,n. \end{aligned}\]

then the gradients can (classically) be considered as vectors of the components gradients, for example $\bigl(\operatorname{grad} g_1(p), \operatorname{grad} g_2(p), …, \operatorname{grad} g_m(p) \bigr)$.

In another interpretation, this can be considered a point on the tangent space at $P = (p,…,p) \in \mathcal M^m$, so in the tangent space to the PowerManifold $\mathcal M^m$. The case where this is a NestedPowerRepresentation this agrees with the interpretation from before, but on power manifolds, more efficient representations exist.

To then access the elements, the range has to be specified. That is what this problem is for.

Constructor

ConstrainedManoptProblem(
    M::AbstractManifold,
    co::ConstrainedManifoldObjective;
    range=NestedPowerRepresentation(),
    gradient_equality_range=range,
    gradient_inequality_range=range
    hessian_equality_range=range,
    hessian_inequality_range=range
)

Creates a constrained Manopt problem specifying an AbstractPowerRepresentation for both the gradient_equality_range and the gradient_inequality_range, respectively.

source

as well as the helper functions

Manopt.AbstractConstrainedFunctorType
AbstractConstrainedFunctor{T}

A common supertype for fucntors that model constraint functions.

This supertype provides access for the fields $λ$ and $μ$, the dual variables of constraintsnof type T.

source
Manopt.AbstractConstrainedSlackFunctorType
AbstractConstrainedSlackFunctor{T,R}

A common supertype for fucntors that model constraint functions with slack.

This supertype additionally provides access for the fields

  • μ::T the dual for the inequality constraints
  • s::T the slack parametyer, and
  • β::R the the barrier parameter

which is also of typee T.

source
Manopt.LagrangianCostType
LagrangianCost{CO,T} <: AbstractConstrainedFunctor{T}

Implement the Lagrangian of a ConstrainedManifoldObjective co.

\[\mathcal L(p; μ, λ) = f(p) + \sum_{i=1}^m μ_ig_i(p) + \sum_{j=1}^n λ_jh_j(p)\]

Fields

  • co::CO, μ::T, λ::T as mentioned, where T represents a vector type.

Constructor

LagrangianCost(co, μ, λ)

Create a functor for the Lagrangian with fixed dual variables.

Example

When you directly want to evaluate the Lagrangian $\mathcal L$ you can also call

LagrangianCost(co, μ, λ)(M,p)
source
Manopt.LagrangianGradientType
LagrangianGradient{CO,T}

The gradient of the Lagrangian of a ConstrainedManifoldObjective co with respect to the variable $p$. The formula reads

\[\operatorname{grad}_p \mathcal L(p; μ, λ) = \operatorname{grad} f(p) + \sum_{i=1}^m μ_i \operatorname{grad} g_i(p) + \sum_{j=1}^n λ_j \operatorname{grad} h_j(p)\]

Fields

  • co::CO, μ::T, λ::T as mentioned, where T represents a vector type.

Constructor

LagrangianGradient(co, μ, λ)

Create a functor for the Lagrangian with fixed dual variables.

Example

When you directly want to evaluate the gradient of the Lagrangian $\operatorname{grad}_p \mathcal L$ you can also call LagrangianGradient(co, μ, λ)(M,p) or LagrangianGradient(co, μ, λ)(M,X,p) for the in-place variant.

source
Manopt.LagrangianHessianType
LagrangianHessian{CO, V, T}

The Hesian of the Lagrangian of a ConstrainedManifoldObjective co with respect to the variable $p$. The formula reads

\[\operatorname{Hess}_p \mathcal L(p; μ, λ)[X] = \operatorname{Hess} f(p) + \sum_{i=1}^m μ_i \operatorname{Hess} g_i(p)[X] + \sum_{j=1}^n λ_j \operatorname{Hess} h_j(p)[X]\]

Fields

  • co::CO, μ::T, λ::T as mentioned, where T represents a vector type.

Constructor

LagrangianHessian(co, μ, λ)

Create a functor for the Lagrangian with fixed dual variables.

Example

When you directly want to evaluate the Hessian of the Lagrangian $\operatorname{Hess}_p \mathcal L$ you can also call LagrangianHessian(co, μ, λ)(M, p, X) or LagrangianHessian(co, μ, λ)(M, Y, p, X) for the in-place variant.

source

Access functions

Manopt.get_inequality_constraintFunction
get_inequality_constraint(amp::AbstractManoptProblem, p, j=:)
get_inequality_constraint(M::AbstractManifold, co::ConstrainedManifoldObjective, p, j=:, range=NestedPowerRepresentation())

Evaluate inequality constraints of a ConstrainedManifoldObjective objective at point p and indices j (by default : which corresponds to all indices).

source
Manopt.get_grad_equality_constraintFunction
get_grad_equality_constraint(amp::AbstractManoptProblem, p, j)
get_grad_equality_constraint(M::AbstractManifold, co::ConstrainedManifoldObjective, p, j, range=NestedPowerRepresentation())
get_grad_equality_constraint!(amp::AbstractManoptProblem, X, p, j)
get_grad_equality_constraint!(M::AbstractManifold, X, co::ConstrainedManifoldObjective, p, j, range=NestedPowerRepresentation())

Evaluate the gradient or gradients of the equality constraint $(\operatorname{grad} h(p))_j$ or $\operatorname{grad} h_j(p)$,

See also the ConstrainedManoptProblem to specify the range of the gradient.

source
Manopt.get_grad_inequality_constraintFunction
get_grad_inequality_constraint(amp::AbstractManoptProblem, p, j=:)
get_grad_inequality_constraint(M::AbstractManifold, co::ConstrainedManifoldObjective, p, j=:, range=NestedPowerRepresentation())
get_grad_inequality_constraint!(amp::AbstractManoptProblem, X, p, j=:)
get_grad_inequality_constraint!(M::AbstractManifold, X, co::ConstrainedManifoldObjective, p, j=:, range=NestedPowerRepresentation())

Evaluate the gradient or gradients of the inequality constraint $(\operatorname{grad} g(p))_j$ or $\operatorname{grad} g_j(p)$,

See also the ConstrainedManoptProblem to specify the range of the gradient.

source
Manopt.get_hess_equality_constraintFunction
get_hess_equality_constraint(amp::AbstractManoptProblem, p, j=:)
get_hess_equality_constraint(M::AbstractManifold, co::ConstrainedManifoldObjective, p, j, range=NestedPowerRepresentation())
get_hess_equality_constraint!(amp::AbstractManoptProblem, X, p, j=:)
get_hess_equality_constraint!(M::AbstractManifold, X, co::ConstrainedManifoldObjective, p, j, range=NestedPowerRepresentation())

Evaluate the Hessian or Hessians of the equality constraint $(\operatorname{Hess} h(p))_j$ or $\operatorname{Hess} h_j(p)$,

See also the ConstrainedManoptProblem to specify the range of the Hessian.

source
Manopt.get_hess_inequality_constraintFunction
get_hess_inequality_constraint(amp::AbstractManoptProblem, p, X, j=:)
get_hess_inequality_constraint(M::AbstractManifold, co::ConstrainedManifoldObjective, p, j=:, range=NestedPowerRepresentation())
get_hess_inequality_constraint!(amp::AbstractManoptProblem, Y, p, j=:)
get_hess_inequality_constraint!(M::AbstractManifold, Y, co::ConstrainedManifoldObjective, p, X, j=:, range=NestedPowerRepresentation())

Evaluate the Hessian or Hessians of the inequality constraint $(\operatorname{Hess} g(p)[X])_j$ or $\operatorname{Hess} g_j(p)[X]$,

See also the ConstrainedManoptProblem to specify the range of the Hessian.

source
Manopt.is_feasibleFunction
is_feasible(M::AbstractManifold, cmo::ConstrainedManifoldObjective, p, kwargs...)

Evaluate whether a boint p on M is feasible with respect to the ConstrainedManifoldObjective cmo. That is for the provided inequality constaints $g: \mathcal M → ℝ^m$ and equality constaints $h: \mathcal M \to ℝ^m$ from within cmo, the point $p ∈ \mathcal M$ is feasible if

\[g_i(p) ≤ 0, \text{ for all } i=1,…,m\quad\text{ and }\quad h_j(p) = 0, \text{ for all } j=1,…,n.\]

Keyword arguments

  • check_point::Bool=true: whether to also verify that `p∈\mathcal M holds, using is_point
  • error::Symbol=:none: if the point is not feasible, this symbol determines how to report the error.
    • :error: throws an error
    • :info: displays the error message as an @info
    • :none: (default) the function just returns true/false
    • :warn: displays the error message as a @warning.

The keyword error= and all other kwargs... are passed on to is_point if the point is verfied (see check_point).

All other keywords are passed on to is_poi

source

Internal functions

Manopt.get_feasibility_statusFunction
get_feasibility_status(
    M::AbstractManifold,
    cmo::ConstrainedManifoldObjective,
    g = get_inequality_constraints(M, cmo, p),
    h = get_equality_constraints(M, cmo, p),
)

Generate a message about the feasibiliy of p with respect to the ConstrainedManifoldObjective. You can also provide the evaluated vectors for the values of g and h as keyword arguments, in case you had them evaluated before.

source

Vectorial objectives

Manopt.AbstractVectorFunctionType
AbstractVectorFunction{E, FT} <: Function

Represent an abstract vectorial function $f:\mathcal M → ℝ^n$ with an AbstractEvaluationType E and an AbstractVectorialType to specify the format $f$ is implemented as.

Representations of $f$

There are three different representations of $f$, which might be beneficial in one or the other situation:

For the ComponentVectorialType imagine that $f$ could also be written using its component functions,

\[f(p) = \bigl( f_1(p), f_2(p), \ldots, f_n(p) \bigr)^{\mathrm{T}}\]

In this representation f is given as a vector [f1(M,p), f2(M,p), ..., fn(M,p)] of its component functions. An advantage is that the single components can be evaluated and from this representation one even can directly read of the number n. A disadvantage might be, that one has to implement a lot of individual (component) functions.

For the FunctionVectorialType $f$ is implemented as a single function f(M, p), that returns an AbstractArray. And advantage here is, that this is a single function. A disadvantage might be, that if this is expensive even to compute a single component, all of f has to be evaluated

For the ComponentVectorialType of f, each of the component functions is a (classical) objective.

source
Manopt.VectorGradientFunctionType
VectorGradientFunction{E, FT, JT, F, J, I} <: AbstractVectorGradientFunction{E, FT, JT}

Represent a function $f:\mathcal M → ℝ^n$ including it first derivative, either as a vector of gradients of a Jacobian

And hence has a gradient $\oepratorname{grad} f_i(p) ∈ T_p\mathcal M$. Putting these gradients into a vector the same way as the functions, yields a ComponentVectorialType

\[\operatorname{grad} f(p) = \Bigl( \operatorname{grad} f_1(p), \operatorname{grad} f_2(p), …, \operatorname{grad} f_n(p) \Bigr)^{\mathrm{T}} ∈ (T_p\mathcal M)^n\]

And advantage here is, that again the single components can be evaluated individually

Fields

  • value!!: the cost function $f$, which can take different formats
  • cost_type: indicating / string data for the type of f
  • jacobian!!: the Jacobian of $f$
  • jacobian_type: indicating / storing data for the type of $J_f$
  • parameters: the number n from, the size of the vector $f$ returns.

Constructor

VectorGradientFunction(f, Jf, range_dimension;
    evaluation::AbstractEvaluationType=AllocatingEvaluation(),
    function_type::AbstractVectorialType=FunctionVectorialType(),
    jacobian_type::AbstractVectorialType=FunctionVectorialType(),
)

Create a VectorGradientFunction of f and its Jacobian (vector of gradients) Jf, where f maps into the Euclidean space of dimension range_dimension. Their types are specified by the function_type, and jacobian_type, respectively. The Jacobian can further be given as an allocating variant or an in-place variant, specified by the evaluation= keyword.

source
Manopt.VectorHessianFunctionType
VectorHessianFunction{E, FT, JT, HT, F, J, H, I} <: AbstractVectorGradientFunction{E, FT, JT}

Represent a function $f:\mathcal M → ℝ^n$ including it first derivative, either as a vector of gradients of a Jacobian, and the Hessian, as a vector of Hessians of the component functions.

Both the Jacobian and the Hessian can map into either a sequence of tangent spaces or a single tangent space of the power manifold of lenth n.

Fields

  • value!!: the cost function $f$, which can take different formats
  • cost_type: indicating / string data for the type of f
  • jacobian!!: the Jacobian of $f$
  • jacobian_type: indicating / storing data for the type of $J_f$
  • hessians!!: the Hessians of $f$ (in a component wise sense)
  • hessian_type: indicating / storing data for the type of $H_f$
  • parameters: the number n from, the size of the vector $f$ returns.

Constructor

VectorGradientFunction(f, Jf, Hess_f, range_dimension;
    evaluation::AbstractEvaluationType=AllocatingEvaluation(),
    function_type::AbstractVectorialType=FunctionVectorialType(),
    jacobian_type::AbstractVectorialType=FunctionVectorialType(),
    hessian_type::AbstractVectorialType=FunctionVectorialType(),
)

Create a VectorGradientFunction of f and its Jacobian (vector of gradients) Jf and (vector of) Hessians, where f maps into the Euclidean space of dimension range_dimension. Their types are specified by the function_type, and jacobian_type, and hessian_type, respectively. The Jacobian and Hessian can further be given as an allocating variant or an inplace-variant, specified by the evaluation= keyword.

source
Manopt.AbstractVectorialTypeType
AbstractVectorialType

An abstract type for different representations of a vectorial function $f: \mathcal M → \mathbb R^m$ and its (component-wise) gradient/Jacobian

source
Manopt.CoordinateVectorialTypeType
CoordinateVectorialType{B<:AbstractBasis} <: AbstractVectorialType

A type to indicate that gradient of the constraints is implemented as a Jacobian matrix with respect to a certain basis, that is if the constraints are given as $g: \mathcal M → ℝ^m$ with respect to a basis $\mathcal B$ of $T_p\mathcal M$, at $p∈ \mathcal M$ This can be written as $J_g(p) = (c_1^{\mathrm{T}},…,c_m^{\mathrm{T}})^{\mathrm{T}} \in ℝ^{m,d}$, that is, every row $c_i$ of this matrix is a set of coefficients such that get_coefficients(M, p, c, B) is the tangent vector $\oepratorname{grad} g_i(p)$

for example $g_i(p) ∈ ℝ^m$ or $\operatorname{grad} g_i(p) ∈ T_p\mathcal M$, $i=1,…,m$.

Fields

source
Manopt.ComponentVectorialTypeType
ComponentVectorialType <: AbstractVectorialType

A type to indicate that constraints are implemented as component functions, for example $g_i(p) ∈ ℝ^m$ or $\operatorname{grad} g_i(p) ∈ T_p\mathcal M$, $i=1,…,m$.

source
Manopt.FunctionVectorialTypeType
FunctionVectorialType <: AbstractVectorialType

A type to indicate that constraints are implemented one whole functions, for example $g(p) ∈ ℝ^m$ or $\operatorname{grad} g(p) ∈ (T_p\mathcal M)^m$.

source

Access functions

Manopt.get_valueFunction
get_value(M::AbstractManifold, vgf::AbstractVectorFunction, p[, i=:])

Evaluate the vector function VectorGradientFunction vgf at p. The range can be used to specify a potential range, but is currently only present for consistency.

The i can be a linear index, you can provide

  • a single integer
  • a UnitRange to specify a range to be returned like 1:3
  • a BitVector specifying a selection
  • a AbstractVector{<:Integer} to specify indices
  • : to return the vector of all gradients, which is also the default
source
Base.lengthMethod
length(vgf::AbstractVectorFunction)

Return the length of the vector the function $f: \mathcal M → ℝ^n$ maps into, that is the number n.

source

Internal functions

Manopt._to_iterable_indicesFunction
_to_iterable_indices(A::AbstractVector, i)

Convert index i (integer, colon, vector of indices, etc.) for array A into an iterable structure of indices.

source

Subproblem objective

This objective can be use when the objective of a sub problem solver still needs access to the (outer/main) objective.

Manopt.AbstractManifoldSubObjectiveType
AbstractManifoldSubObjective{O<:AbstractManifoldObjective} <: AbstractManifoldObjective

An abstract type for objectives of sub problems within a solver but still store the original objective internally to generate generic objectives for sub solvers.

source

Access functions

Manopt.get_objective_costFunction
get_objective_cost(M, amso::AbstractManifoldSubObjective, p)

Evaluate the cost of the (original) objective stored within the sub objective.

source
Manopt.get_objective_gradientFunction
X = get_objective_gradient(M, amso::AbstractManifoldSubObjective, p)
get_objective_gradient!(M, X, amso::AbstractManifoldSubObjective, p)

Evaluate the gradient of the (original) objective stored within the sub objective amso.

source
Manopt.get_objective_hessianFunction
Y = get_objective_Hessian(M, amso::AbstractManifoldSubObjective, p, X)
get_objective_Hessian!(M, Y, amso::AbstractManifoldSubObjective, p, X)

Evaluate the Hessian of the (original) objective stored within the sub objective amso.

source
Manopt.get_objective_preconditionerFunction
Y = get_objective_preconditioner(M, amso::AbstractManifoldSubObjective, p, X)
get_objective_preconditioner(M, Y, amso::AbstractManifoldSubObjective, p, X)

Evaluate the Hessian of the (original) objective stored within the sub objective amso.

source