How to do Constrained Optimization

Ronny Bergmann

This tutorial is a short introduction to using solvers for constraint optimisation in Manopt.jl.

Introduction

A constraint optimisation problem is given by

\[\tag{P} \begin{align*} \operatorname*{arg\,min}_{p\in\mathcal M} & f(p)\\ \text{such that} &\quad g(p) \leq 0\\ &\quad h(p) = 0,\\ \end{align*}\]

where $f\colon \mathcal M → ℝ$ is a cost function, and $g\colon \mathcal M → ℝ^m$ and $h\colon \mathcal M → ℝ^n$ are the inequality and equality constraints, respectively. The $\leq$ and $=$ in (P) are meant elementwise.

This can be seen as a balance between moving constraints into the geometry of a manifold $\mathcal M$ and keeping some, since they can be handled well in algorithms, see [BH19], [LB19] for details.

using Distributions, LinearAlgebra, Manifolds, Manopt, Random
Random.seed!(42);

In this tutorial we want to look at different ways to specify the problem and its implications. We start with specifying an example problems to illustrayte the different available forms.

We will consider the problem of a Nonnegative PCA, cf. Section 5.1.2 in [LB19]

let $v_0 ∈ ℝ^d$, $\lVert v_0 \rVert=1$ be given spike signal, that is a signal that is sparse with only $s=\lfloor δd \rfloor$ nonzero entries.

\[Z = \sqrt{σ} v_0v_0^{\mathrm{T}}+N,\]

where $\sigma$ is a signal-to-noise ratio and $N$ is a matrix with random entries, where the diagonal entries are distributed with zero mean and standard deviation $1/d$ on the off-diagonals and $2/d$ on the daigonal

d = 150; # dimension of v0
σ = 0.1^2; # SNR
δ = 0.1; s = Int(floor(δ * d)); # Sparsity
S = sample(1:d, s; replace=false);
v0 =  [i ∈ S ? 1 / sqrt(s) : 0.0 for i in 1:d];
N = rand(Normal(0, 1 / d), (d, d)); N[diagind(N, 0)] .= rand(Normal(0, 2 / d), d);
Z = Z = sqrt(σ) * v0 * transpose(v0) + N;

In order to recover $v_0$ we consider the constrained optimisation problem on the sphere $\mathcal S^{d-1}$ given by

\[\begin{align*} \operatorname*{arg\,min}_{p\in\mathcal S^{d-1}} & -p^{\mathrm{T}}Zp^{\mathrm{T}}\\ \text{such that} &\quad p \geq 0\\ \end{align*}\]

or in the previous notation $f(p) = -p^{\mathrm{T}}Zp^{\mathrm{T}}$ and $g(p) = -p$. We first initialize the manifold under consideration

M = Sphere(d - 1)
Sphere(149, ℝ)

A first Augmented Lagrangian Run

We first defined $f$ and $g$ as usual functions

f(M, p) = -transpose(p) * Z * p;
g(M, p) = -p;

since $f$ is a functions defined in the embedding $ℝ^d$ as well, we obtain its gradient by projection.

grad_f(M, p) = project(M, p, -transpose(Z) * p - Z * p);

For the constraints this is a little more involved, since each function $g_i = g(p)_i = p_i$ has to return its own gradient. These are again in the embedding just $\operatorname{grad} g_i(p) = -e_i$ the $i$ th unit vector. We can project these again onto the tangent space at $p$:

grad_g(M, p) = project.(
    Ref(M), Ref(p), [[i == j ? -1.0 : 0.0 for j in 1:d] for i in 1:d]
);

We further start in a random point:

p0 = rand(M);

Let’s check a few things for the initial point

f(M, p0)
0.005747604833124234

How much the function g is positive

maximum(g(M, p0))
0.17885478285466855

Now as a first method we can just call the Augmented Lagrangian Method with a simple call:

@time v1 = augmented_Lagrangian_method(
    M, f, grad_f, p0; g=g, grad_g=grad_g,
    debug=[:Iteration, :Cost, :Stop, " | ", (:Change, "Δp : %1.5e"), 20, "\n"],
    stopping_criterion = StopAfterIteration(300) | (
        StopWhenSmallerOrEqual(:ϵ, 1e-5) & StopWhenChangeLess(1e-8)
    )
);
Initial f(x): 0.005748 | 
# 20    f(x): -0.123842 | Δp : 9.99682e-01
# 40    f(x): -0.123842 | Δp : 8.13541e-07
# 60    f(x): -0.123842 | Δp : 7.85694e-04
The value of the variable (ϵ) is smaller than or equal to its threshold (1.0e-5).
The algorithm performed a step with a change (1.7450108123172955e-15) less than 9.77237220955808e-6.
 16.843524 seconds (43.34 M allocations: 32.293 GiB, 10.65% gc time, 37.25% compilation time)

Now we have both a lower function value and the point is nearly within the constraints, … up to numerical inaccuracies

f(M, v1)
-0.12384244779997305
maximum( g(M, v1) )
7.912675333644102e-18

A faster Augmented Lagrangian Run

Now this is a little slow, so we can modify two things, that we will directly do both – but one could also just change one of these – :

  1. Gradients should be evaluated in place, so for example
grad_f!(M, X, p) = project!(M, X, p, -transpose(Z) * p - Z * p);
  1. The constraints are currently always evaluated all together, since the function grad_g always returns a vector of gradients. We first change the constraints function into a vector of functions. We further change the gradient both into a vector of gradient functions $\operatorname{grad} g_i, i=1,\ldots,d$, as well as gradients that are computed in place.
g2 = [(M, p) -> -p[i] for i in 1:d];
grad_g2! = [
    (M, X, p) -> project!(M, X, p, [i == j ? -1.0 : 0.0 for j in 1:d]) for i in 1:d
];

We obtain

@time v2 = augmented_Lagrangian_method(
        M, f, grad_f!, p0; g=g2, grad_g=grad_g2!, evaluation=InplaceEvaluation(),
        debug=[:Iteration, :Cost, :Stop, " | ", (:Change, "Δp : %1.5e"), 20, "\n"],
        stopping_criterion = StopAfterIteration(300) | (
          StopWhenSmallerOrEqual(:ϵ, 1e-5) & StopWhenChangeLess(1e-8)
        )
    );
Initial f(x): 0.005748 | 
# 20    f(x): -0.123842 | Δp : 9.99544e-01
# 40    f(x): -0.123842 | Δp : 1.92065e-03
# 60    f(x): -0.123842 | Δp : 4.84931e-06
The value of the variable (ϵ) is smaller than or equal to its threshold (1.0e-5).
The algorithm performed a step with a change (2.7435918100802105e-17) less than 9.77237220955808e-6.
  3.547284 seconds (6.52 M allocations: 3.728 GiB, 6.70% gc time, 41.27% compilation time)

As a technical remark: Note that (by default) the change to InplaceEvaluations affects both the constrained solver as well as the inner solver of the subproblem in each iteration.

f(M, v2)
-0.12384239276300012
maximum(g(M, v2))
2.2466899389459647e-18

These are the very similar to the previous values but the solver took much less time and less memory allocations.

Exact Penalty Method

As a second solver, we have the Exact Penalty Method, which currenlty is available with two smoothing variants, which make an inner solver for smooth optimisationm, that is by default again [quasi Newton] possible: LogarithmicSumOfExponentials and LinearQuadraticHuber. We compare both here as well. The first smoothing technique is the default, so we can just call

@time v3 = exact_penalty_method(
    M, f, grad_f!, p0; g=g2, grad_g=grad_g2!, evaluation=InplaceEvaluation(),
    debug=[:Iteration, :Cost, :Stop, " | ", :Change, 50, "\n"],
);
Initial f(x): 0.005748 | 
# 50    f(x): -0.123071 | Last Change: 0.981116
# 100   f(x): -0.123840 | Last Change: 0.014124
The value of the variable (ϵ) is smaller than or equal to its threshold (1.0e-6).
The algorithm performed a step with a change (2.202641515349944e-7) less than 1.0e-6.
  2.383160 seconds (5.78 M allocations: 3.123 GiB, 7.71% gc time, 64.51% compilation time)

We obtain a similar cost value as for the Augmented Lagrangian Solver above, but here the constraint is actually fulfilled and not just numerically “on the boundary”.

f(M, v3)
-0.12384029692539944
maximum(g(M, v3))
-3.582398293370528e-6

The second smoothing technique is often beneficial, when we have a lot of constraints (in the above mentioned vectorial manner), since we can avoid several gradient evaluations for the constraint functions here. This leads to a faster iteration time.

@time v4 = exact_penalty_method(
    M, f, grad_f!, p0; g=g2, grad_g=grad_g2!,
    evaluation=InplaceEvaluation(),
    smoothing=LinearQuadraticHuber(),
    debug=[:Iteration, :Cost, :Stop, " | ", :Change, 50, "\n"],
);
Initial f(x): 0.005748 | 
# 50    f(x): -0.123845 | Last Change: 0.009235
# 100   f(x): -0.123843 | Last Change: 0.000107
The value of the variable (ϵ) is smaller than or equal to its threshold (1.0e-6).
The algorithm performed a step with a change (3.586352489111338e-7) less than 1.0e-6.
  1.557075 seconds (2.76 M allocations: 514.648 MiB, 5.08% gc time, 79.85% compilation time)

For the result we see the same behaviour as for the other smoothing.

f(M, v4)
-0.12384258173223292
maximum(g(M, v4))
2.7028045565194566e-8

Comparing to the unconstraint solver

We can compare this to the global optimum on the sphere, which is the unconstraint optimisation problem; we can just use Quasi Newton.

Note that this is much faster, since every iteration of the algorithms above does a quasi-Newton call as well.

@time w1 = quasi_Newton(
    M, f, grad_f!, p0; evaluation=InplaceEvaluation()
);
  0.706571 seconds (634.12 k allocations: 61.701 MiB, 3.18% gc time, 96.56% compilation time)
f(M, w1)
-0.14021901809807297

But for sure here the constraints here are not fulfilled and we have veru positive entries in $g(w_1)$

maximum(g(M, w1))
0.11762414497055226

Literature

[BH19]
R. Bergmann and R. Herzog. Intrinsic formulation of KKT conditions and constraint qualifications on smooth manifolds. SIAM Journal on Optimization 29, 2423–2444 (2019), arXiv: [1804.06214](https://arxiv.org/abs/1804.06214).
[LB19]