Main.check_gradients
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 either a point in the intersection, $x \in P \cap Q$, or a pair of points, $(x_P, x_Q) \in 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")
Main.check_gradients
Setting up objective, gradient and linear minimization oracles
Since we only consider the feasibility problem the objective function as well as the gradient are 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
We run Alternating Linear Minimization (ALM) with FrankWolfe.block_coordinate_frank_wolfe
. This method allows three different update orders, 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).
LAMBDA: 1.0
Block coordinate Frank-Wolfe (BCFW).
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: DataType[FrankWolfe.Adaptive{Float64, Int64}, FrankWolfe.Adaptive{Float64, Int64}] EPSILON: 1.0e-6 MAXITERATION: 10000 TYPE: Float64
MOMENTUM: nothing GRADIENTTYPE: FrankWolfe.BlockVector{Float64, Vector{Float64}, Tuple{Int64}} UPDATE_ORDER: FrankWolfe.FullUpdate() UPDATE_STEP: DataType[FrankWolfe.FrankWolfeStep, FrankWolfe.FrankWolfeStep]
[ Info: In memory_mode memory iterates are written back into x0!
-------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec
-------------------------------------------------------------------------------------------------
I 1 3.778464e+00 -4.022154e+01 4.400000e+01 0.000000e+00 Inf
----------------
Infeas
----------------
3.778464e+00
FW 1000 5.000047e-02 4.988536e-02 1.151053e-04 7.017723e-01 1.424964e+03
5.000047e-02
FW 2000 5.000000e-02 4.999558e-02 4.423306e-06 7.113388e-01 2.811600e+03
5.000000e-02
Last 2445 5.000000e-02 4.999898e-02 1.022986e-06 8.821796e-01 2.771545e+03
-------------------------------------------------------------------------------------------------
5.000000e-02
----------------
Alternating Linear Minimization (ALM).
LAMBDA: 1.0
Block coordinate Frank-Wolfe (BCFW).
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: DataType[FrankWolfe.Adaptive{Float64, Int64}, FrankWolfe.Adaptive{Float64, Int64}] EPSILON: 1.0e-6 MAXITERATION: 10000 TYPE: Float64
MOMENTUM: nothing GRADIENTTYPE: FrankWolfe.BlockVector{Float64, Vector{Float64}, Tuple{Int64}} UPDATE_ORDER: FrankWolfe.CyclicUpdate() UPDATE_STEP: DataType[FrankWolfe.FrankWolfeStep, FrankWolfe.FrankWolfeStep]
[ Info: In memory_mode memory iterates are written back into x0!
-------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec
-------------------------------------------------------------------------------------------------
I 1 3.778464e+00 -4.022154e+01 4.400000e+01 0.000000e+00 Inf
----------------
Infeas
----------------
3.778464e+00
FW 1000 5.000047e-02 4.988833e-02 1.121389e-04 8.125300e-03 1.230724e+05
5.000047e-02
FW 2000 5.000000e-02 4.999569e-02 4.308009e-06 1.640463e-02 1.219168e+05
5.000000e-02
Last 2446 5.000000e-02 4.999898e-02 1.018029e-06 2.227696e-02 1.097996e+05
-------------------------------------------------------------------------------------------------
5.000000e-02
----------------
Alternating Linear Minimization (ALM).
LAMBDA: 1.0
Block coordinate Frank-Wolfe (BCFW).
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: DataType[FrankWolfe.Adaptive{Float64, Int64}, FrankWolfe.Adaptive{Float64, Int64}] EPSILON: 1.0e-6 MAXITERATION: 10000 TYPE: Float64
MOMENTUM: nothing GRADIENTTYPE: FrankWolfe.BlockVector{Float64, Vector{Float64}, Tuple{Int64}} UPDATE_ORDER: FrankWolfe.StochasticUpdate() UPDATE_STEP: DataType[FrankWolfe.FrankWolfeStep, FrankWolfe.FrankWolfeStep]
[ Info: In memory_mode memory iterates are written back into x0!
-------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec
-------------------------------------------------------------------------------------------------
I 1 3.778464e+00 -4.022154e+01 4.400000e+01 0.000000e+00 Inf
----------------
Infeas
----------------
3.778464e+00
FW 1000 5.000043e-02 4.988980e-02 1.106340e-04 9.367196e-02 1.067555e+04
5.000043e-02
FW 2000 5.000000e-02 4.999606e-02 3.941984e-06 1.073940e-01 1.862302e+04
5.000000e-02
Last 2416 5.000000e-02 4.999900e-02 1.002753e-06 1.133112e-01 2.132181e+04
-------------------------------------------------------------------------------------------------
5.000000e-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).
LAMBDA: 1.0
Away-step Frank-Wolfe Algorithm.
MEMORY_MODE: FrankWolfe.InplaceEmphasis() STEPSIZE: Adaptive{Float64, Int64} EPSILON: 1.0e-6 MAXITERATION: 10000 TYPE: Float64
GRADIENTTYPE: FrankWolfe.BlockVector{Float64, Vector{Float64}, Tuple{Int64}} LAZY: false lazy_tolerance: 2.0 MOMENTUM: nothing AWAYSTEPS: true
LMO: FrankWolfe.ProductLMO{2, Tuple{FrankWolfe.ProbabilitySimplexOracle{Float64}, FrankWolfe.ScaledBoundLInfNormBall{Float64, 1, Vector{Float64}, Vector{Float64}}}}
----------------------------------------------------------------------------------------------------------------
Type Iteration Primal Dual Dual Gap Time It/sec #ActiveSet
----------------------------------------------------------------------------------------------------------------
I 1 2.300000e+01 -Inf Inf 0.000000e+00 Inf 2
----------------
Infeas
----------------
2.010582e+00
Last 171 5.000000e-02 4.999908e-02 9.177388e-07 9.940259e-01 1.720277e+02 84
----------------------------------------------------------------------------------------------------------------
5.000000e-02
----------------
PP 171 5.000000e-02 4.999908e-02 9.177388e-07 1.079442e+00 1.584152e+02 84
----------------------------------------------------------------------------------------------------------------
5.000000e-02
----------------
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: Vector{Vector{Float64}}
[ Info: In memory_mode memory iterates are written back into x0!
----------------------------------------------------------------------------------
Type Iteration Dual Gap Infeas Time It/sec
----------------------------------------------------------------------------------
I 1 4.081717e-01 2.040858e-01 0.000000e+00 Inf
FW 100 1.045308e-04 5.000040e-02 1.233165e+00 8.109217e+01
FW 200 2.524997e-05 5.000002e-02 1.788905e+00 1.118003e+02
FW 300 1.122811e-05 5.000000e-02 2.457983e+00 1.220513e+02
FW 400 6.334101e-06 5.000000e-02 3.213475e+00 1.244758e+02
FW 500 4.028960e-06 5.000000e-02 4.031201e+00 1.240325e+02
FW 600 2.823896e-06 5.000000e-02 4.907775e+00 1.222550e+02
FW 700 2.088682e-06 5.000000e-02 5.843668e+00 1.197878e+02
FW 800 1.569986e-06 5.000000e-02 6.780057e+00 1.179931e+02
FW 900 1.258972e-06 5.000000e-02 7.765538e+00 1.158967e+02
FW 1000 9.944106e-07 5.000000e-02 8.794483e+00 1.137077e+02
Last 1000 9.944106e-07 5.000000e-02 8.800321e+00 1.136322e+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.