plot_sparsity (generic function with 1 method)
Alternating methods
In this example we will compare FrankWolfe.alternating_linear_minimization
and FrankWolfe.alternating_projections
for a very simple feasibility problem.
We consider the probability simplex
\[P = \{ x \in \mathbb{R}^n \colon \sum_{i=1}^n x_i = 1, x_i \geq 0 ~~ i=1,\dots,n\} ~.\]
and a scaled, shifted $\ell^{\infty}$ norm ball
\[Q = [-1,0]^n ~.\]
The goal is to find a point that lies both in $P$ and $Q$. We do this by reformulating the problem first. Instead of a finding a point in the intersection $P \cap Q$, we search for a pair of points, $(x_P, x_Q)$ in the cartesian product $P \times Q$, which attains minimal distance between $P$ and $Q$,
\[\|x_P - x_Q\|_2 = \min_{(x,y) \in P \times Q} \|x - y \|_2 ~.\]
using FrankWolfe
include("../examples/plot_utils.jl")
plot_sparsity (generic function with 1 method)
Setting up objective, gradient and linear minimization oracles
Alternating Linear Minimization (ALM) allows for an additional objective such that one can optimize over an intersection of sets instead of finding only feasible points. Since this example only considers the feasibility, we set the objective function as well as the gradient to zero.
n = 20
f(x) = 0
function grad!(storage, x)
@. storage = zero(x)
end
lmo1 = FrankWolfe.ProbabilitySimplexOracle(1.0)
lmo2 = FrankWolfe.ScaledBoundLInfNormBall(-ones(n), zeros(n))
lmos = (lmo1, lmo2)
x0 = rand(n)
target_tolerance = 1e-6
trajectories = [];
Running Alternating Linear Minimization
The method FrankWolfe.alternating_linear_minimization
is not a FrankWolfe method itself. It is a wrapper translating a problem over the intersection of multiple sets to a problem over the product space. ALM can be called with any FW method. The default choice though is FrankWolfe.block_coordinate_frank_wolfe
as it allows to update the blocks separately. There are three different update orders implemented, FullUpdate
, CyclicUpdate
and Stochasticupdate
. Accordingly both blocks are updated either simulatenously, sequentially or in random order.
for order in [FrankWolfe.FullUpdate(), FrankWolfe.CyclicUpdate(), FrankWolfe.StochasticUpdate()]
_, _, _, _, _, alm_trajectory = FrankWolfe.alternating_linear_minimization(
FrankWolfe.block_coordinate_frank_wolfe,
f,
grad!,
lmos,
x0,
update_order=order,
verbose=true,
trajectory=true,
epsilon=target_tolerance,
)
push!(trajectories, alm_trajectory)
end
Alternating Linear Minimization (ALM).
FW METHOD: block_coordinate_frank_wolfe
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: FrankWolfe.Adaptive{Float64, Int64} EPSILON: 1.0e-6 MAXITERATION: 10000
TYPE: Float64 GRADIENTTYPE: Float64
LAMBDA: 1.0
[ Info: In memory_mode memory iterates are written back into x0!
----------------------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec Dist2
----------------------------------------------------------------------------------------------------------------
I 1 1.889232e+00 -2.011077e+01 2.200000e+01 0.000000e+00 Inf 1.889232e+00
FW 1000 2.500023e-02 2.494268e-02 5.755263e-05 8.898701e-01 1.123760e+03 2.500023e-02
FW 2000 2.500000e-02 2.499779e-02 2.211653e-06 9.052215e-01 2.209404e+03 2.500000e-02
Last 2236 2.500000e-02 2.499900e-02 1.003006e-06 1.068085e+00 2.093467e+03 2.500000e-02
----------------------------------------------------------------------------------------------------------------
Alternating Linear Minimization (ALM).
FW METHOD: block_coordinate_frank_wolfe
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: FrankWolfe.Adaptive{Float64, Int64} EPSILON: 1.0e-6 MAXITERATION: 10000
TYPE: Float64 GRADIENTTYPE: Float64
LAMBDA: 1.0
[ Info: In memory_mode memory iterates are written back into x0!
----------------------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec Dist2
----------------------------------------------------------------------------------------------------------------
I 1 1.889232e+00 -2.011077e+01 2.200000e+01 0.000000e+00 Inf 1.889232e+00
FW 1000 2.500023e-02 2.494416e-02 5.606946e-05 1.823902e-01 5.482750e+03 2.500023e-02
FW 2000 2.500000e-02 2.499785e-02 2.154004e-06 2.044633e-01 9.781705e+03 2.500000e-02
Last 2237 2.500000e-02 2.499900e-02 9.981445e-07 2.100389e-01 1.065041e+04 2.500000e-02
----------------------------------------------------------------------------------------------------------------
Alternating Linear Minimization (ALM).
FW METHOD: block_coordinate_frank_wolfe
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: FrankWolfe.Adaptive{Float64, Int64} EPSILON: 1.0e-6 MAXITERATION: 10000
TYPE: Float64 GRADIENTTYPE: Float64
LAMBDA: 1.0
[ Info: In memory_mode memory iterates are written back into x0!
----------------------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec Dist2
----------------------------------------------------------------------------------------------------------------
I 1 2.866179e-01 -2.171338e+01 2.200000e+01 0.000000e+00 Inf 2.866179e-01
FW 1000 2.500018e-02 2.495033e-02 4.985431e-05 8.633746e-02 1.158246e+04 2.500018e-02
FW 2000 2.500000e-02 2.499809e-02 1.906507e-06 1.078695e-01 1.854092e+04 2.500000e-02
Last 2194 2.500000e-02 2.499900e-02 9.989094e-07 1.122315e-01 1.954888e+04 2.500000e-02
----------------------------------------------------------------------------------------------------------------
As an alternative to Block-Coordiante Frank-Wolfe (BCFW), one can also run alternating linear minimization with standard Frank-Wolfe algorithm. These methods perform then the full (simulatenous) update at each iteration. In this example we also use FrankWolfe.away_frank_wolfe
.
_, _, _, _, _, afw_trajectory = FrankWolfe.alternating_linear_minimization(
FrankWolfe.away_frank_wolfe,
f,
grad!,
lmos,
x0,
verbose=true,
trajectory=true,
epsilon=target_tolerance,
)
push!(trajectories, afw_trajectory);
Alternating Linear Minimization (ALM).
FW METHOD: away_frank_wolfe
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: FrankWolfe.Adaptive{Float64, Int64} EPSILON: 1.0e-6 MAXITERATION: 10000
TYPE: Float64 GRADIENTTYPE: Float64
LAMBDA: 1.0
[ Info: In memory_mode memory iterates are written back into x0!
-------------------------------------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec Dist2 #ActiveSet
-------------------------------------------------------------------------------------------------------------------------------
I 1 1.150000e+01 -Inf Inf 0.000000e+00 Inf 1.005291e+00 2
Last 161 2.500000e-02 2.499904e-02 9.621137e-07 4.042893e-01 3.982297e+02 2.500000e-02 78
-------------------------------------------------------------------------------------------------------------------------------
PP 161 2.500000e-02 2.499904e-02 9.621137e-07 4.846365e-01 3.322078e+02 2.500000e-02 78
-------------------------------------------------------------------------------------------------------------------------------
Running Alternating Projections
Unlike ALM, Alternating Projections (AP) is only suitable for feasibility problems. One omits the objective and gradient as parameters.
_, _, _, _, ap_trajectory = FrankWolfe.alternating_projections(
lmos,
x0,
trajectory=true,
verbose=true,
print_iter=100,
epsilon=target_tolerance,
)
push!(trajectories, ap_trajectory);
Alternating Projections.
MEMORY_MODE: FrankWolfe.InplaceEmphasis() EPSILON: 1.0e-6 MAXITERATION: 10000 TYPE: Float64
GRADIENTTYPE: FrankWolfe.BlockVector{Float64, Vector{Float64}, Tuple{Int64}}
[ Info: In memory_mode memory iterates are written back into x0!
----------------------------------------------------------------------------------
Type Iteration Dual Gap dist2 Time It/sec
----------------------------------------------------------------------------------
I 1 6.021733e-01 1.505433e-01 0.000000e+00 Inf
FW 100 9.968329e-05 2.500020e-02 1.406556e+00 7.109565e+01
FW 200 2.316528e-05 2.500001e-02 1.411063e+00 1.417371e+02
FW 300 1.122389e-05 2.500000e-02 1.414872e+00 2.120333e+02
FW 400 6.105757e-06 2.500000e-02 1.418152e+00 2.820572e+02
FW 500 3.822934e-06 2.500000e-02 1.421399e+00 3.517661e+02
FW 600 2.624313e-06 2.500000e-02 1.424755e+00 4.211250e+02
FW 700 1.881535e-06 2.500000e-02 1.428158e+00 4.901417e+02
FW 800 1.505897e-06 2.500000e-02 1.431304e+00 5.589307e+02
FW 900 1.216129e-06 2.500000e-02 1.434457e+00 6.274152e+02
Last 992 9.990101e-07 2.500000e-02 1.841607e+00 5.386599e+02
----------------------------------------------------------------------------------
Plotting the resulting trajectories
labels = ["BCFW - Full", "BCFW - Cyclic", "BCFW - Stochastic", "AFW", "AP"]
plot_trajectories(trajectories, labels, xscalelog=true)
This page was generated using Literate.jl.