Convex bundle method
Manopt.convex_bundle_method — Functionconvex_bundle_method(M, f, ∂f, p)perform a convex bundle method $p_{j+1} = \mathrm{retr}(p_k, -g_k)$, where $\mathrm{retr}$ is a retraction and
\[g_k = \sum_{j\in J_k} λ_j^k \mathrm{P}_{p_k←q_j}X_{q_j},\]
$p_k$ is the last serious iterate, $X_{q_j} ∈ ∂f(q_j)$, and the $λ_j^k$ are solutions to the quadratic subproblem provided by the convex_bundle_method_subsolver.
Though the subdifferential might be set valued, the argument ∂f should always return one element from the subdifferential, but not necessarily deterministic.
For more details, see [BHJ24].
Input
M: a manifold $\mathcal M$f: a cost function $f:\mathcal M→ℝ$ to minimize∂f: the subgradient $∂f: \mathcal M → T\mathcal M$ of f restricted to always only returning one value/element from the subdifferential. This function can be passed as an allocation function(M, p) -> Xor a mutating function(M, X, p) -> X, seeevaluation.p: (rand(M)) an initial value $p_0 ∈ \mathcal M$
Optional
atol_λ: (eps()) tolerance parameter for the convex coefficients in λ.atol_errors: (eps()) tolerance parameter for the linearization errors.m: (1e-3) the parameter to test the decrease of the cost: $f(q_{k+1}) \le f(p_k) + m \xi$.diameter: (50.0) estimate for the diameter of the level set of the objective function at the starting point.domain: ((M, p) -> isfinite(f(M, p))) a function to that evaluates to true when the current candidate is in the domain of the objectivef, and false otherwise, for example domain = (M, p) -> p ∈ dom f(M, p) ? true : false.k_max: upper bound on the sectional curvature of the manifold.evaluation: (AllocatingEvaluation) specify whether the subgradient works by allocation (default) form∂f(M, q)orInplaceEvaluationin place, that is of the form∂f!(M, X, p).inverse_retraction_method: (default_inverse_retraction_method(M, typeof(p))) an inverse retraction method to useretraction_method: (default_retraction_method(M, typeof(p))) aretraction(M, p, X)to use.stopping_criterion: (StopWhenLagrangeMultiplierLess(1e-8; names=["-ξ"])) a functor, seeStoppingCriterion, indicating when to stopvector_transport_method: (default_vector_transport_method(M, typeof(p))) a vector transport method to usesub_problem: a function evaluating with new allocations that solves the sub problem onMgiven the last serious iteratep_last_serious, the linearization errorslinearization_errors, and the transported subgradientstransported_subgradients
Output
the obtained (approximate) minimizer $p^*$, see get_solver_return for details
Manopt.convex_bundle_method! — Functionconvex_bundle_method!(M, f, ∂f, p)perform a bundle method $p_{j+1} = \mathrm{retr}(p_k, -g_k)$ in place of p.
Input
M: a manifold $\mathcal M$f: a cost function $f:\mathcal M→ℝ$ to minimize∂f: the (sub)gradient $∂f:\mathcal M→ T\mathcal M$ of F restricted to always only returning one value/element from the subdifferential. This function can be passed as an allocation function(M, p) -> Xor a mutating function(M, X, p) -> X, seeevaluation.p: an initial value $p_0=p ∈ \mathcal M$
for more details and all optional parameters, see convex_bundle_method.
State
Manopt.ConvexBundleMethodState — TypeConvexBundleMethodState <: AbstractManoptSolverStateStores option values for a convex_bundle_method solver.
Fields
atol_λ: (eps()) tolerance parameter for the convex coefficients in λatol_errors: (eps()) tolerance parameter for the linearization errorsbundle: bundle that collects each iterate with the computed subgradient at the iteratebundle_cap: (25) the maximal number of elements the bundle is allowed to rememberdiameter: (50.0) estimate for the diameter of the level set of the objective function at the starting pointdomain: ((M, p) -> isfinite(f(M, p))) a function to that evaluates to true when the current candidate is in the domain of the objectivef, and false otherwise, for exampledomain = (M, p) -> p ∈ dom f(M, p) ? true : falseg: descent directioninverse_retraction_method: the inverse retraction to use withink_max: upper bound on the sectional curvature of the manifoldlinearization_errors: linearization errors at the last serious stepm: (1e-3) the parameter to test the decrease of the cost: $f(q_{k+1}) \le f(p_k) + m \xi$.p: current candidate pointp_last_serious: last serious iterateretraction_method: the retraction to use withinstop: aStoppingCriteriontransported_subgradients: subgradients of the bundle that are transported top_last_seriousvector_transport_method: the vector transport method to use withinX: (zero_vector(M, p)) the current element from the possible subgradients atpthat was last evaluated.stepsize: (ConstantStepsize(M)) aStepsizeε: convex combination of the linearization errorsλ: convex coefficients that solve the subproblemξ: the stopping parameter given by $ξ = -\lvert g\rvert^2 – ε$sub_problem: ([convex_bundle_method_subsolver]) a function that solves the sub problem onMgiven the last serious iteratep_last_serious, the linearization errorslinearization_errors, and the transported subgradientstransported_subgradients,sub_state: anAbstractManoptSolverStatefor the subsolver
Constructor
ConvexBundleMethodState(M::AbstractManifold, p; kwargs...)with keywords for all fields with defaults besides p_last_serious which obtains the same type as p. You can use for example X= to specify the type of tangent vector to use
Stopping criteria
Manopt.StopWhenLagrangeMultiplierLess — TypeStopWhenLagrangeMultiplierLess <: StoppingCriterionStopping Criteria for Lagrange multipliers.
Currently these are meant for the convex_bundle_method and proximal_bundle_method, where based on the Lagrange multipliers an approximate (sub)gradient $g$ and an error estimate $ε$ is computed.
The mode=:both requires that both $ε$ and $\lvert g \rvert$ are smaller than their tolerances for the convex_bundle_method, and that $c$ and $\lvert d \rvert$ are smaller than their tolerances for the proximal_bundle_method.
The mode=:estimate requires that, for the convex_bundle_method $-ξ = \lvert g \rvert^2 + ε$ is less than a given tolerance. For the proximal_bundle_method, the equation reads $-ν = μ \lvert d \rvert^2 + c$.
Constructors
StopWhenLagrangeMultiplierLess(tolerance=1e-6; mode::Symbol=:estimate, names=nothing)Create the stopping criterion for one of the modes mentioned. Note that tolerance can be a single number for the :estimate case, but a vector of two values is required for the :both mode. Here the first entry specifies the tolerance for $ε$ ($c$), the second the tolerance for $\lvert g \rvert$ ($\lvert d \rvert$), respectively.
Debug functions
Manopt.DebugWarnIfLagrangeMultiplierIncreases — TypeDebugWarnIfLagrangeMultiplierIncreases <: DebugActionprint a warning if the Lagrange parameter based value $-ξ$ of the bundle method increases.
Constructor
DebugWarnIfLagrangeMultiplierIncreases(warn=:Once; tol=1e2)Initialize the warning to warning level (:Once) and introduce a tolerance for the test of 1e2.
The warn level can be set to :Once to only warn the first time the cost increases, to :Always to report an increase every time it happens, and it can be set to :No to deactivate the warning, then this DebugAction is inactive. All other symbols are handled as if they were :Always:
Helpers and internal functions
Manopt.convex_bundle_method_subsolver — Functionλ = convex_bundle_method_subsolver(M, p_last_serious, linearization_errors, transported_subgradients)
convex_bundle_method_subsolver!(M, λ, p_last_serious, linearization_errors, transported_subgradients)solver for the subproblem of the convex bundle method at the last serious iterate $p_k$ given the current linearization errors $c_j^k$, and transported subgradients $\mathrm{P}_{p_k←q_j} X_{q_j}$.
The computation can also be done in-place of λ.
The subproblem for the convex bundle method is
\[\begin{align*} \operatorname*{arg\,min}_{λ ∈ ℝ^{\lvert J_k\rvert}}& \frac{1}{2} \Bigl\lVert \sum_{j ∈ J_k} λ_j \mathrm{P}_{p_k←q_j} X_{q_j} \Bigr\rVert^2 + \sum_{j ∈ J_k} λ_j \, c_j^k \\ \text{s. t.}\quad & \sum_{j ∈ J_k} λ_j = 1, \quad λ_j ≥ 0 \quad \text{for all } j ∈ J_k, \end{align*}\]
where $J_k = \{j ∈ J_{k-1} \ | \ λ_j > 0\} \cup \{k\}$. See [BHJ24] for more details
A default subsolver based on RipQP.jl and QuadraticModels is available if these two packages are loaded.
Manopt.DomainBackTrackingStepsize — TypeDomainBackTrackingStepsize <: StepsizeImplement a backtrack as long as we are $q = \operatorname{retr}_p(X)$ yields a point closer to $p$ than $\lVert X \rVert_p$ or $q$ is not on the domain. For the domain this step size requires a ConvexBundleMethodState
Literature
- [BHJ24]
- R. Bergmann, R. Herzog and H. Jasa. The Riemannian Convex Bundle Method, preprint (2024), arXiv:2402.13670.