Utilities and data structures
Active set
FrankWolfe.AbstractActiveSet — Type
AbstractActiveSet{AT, R, IT}Abstract type for an active set of atoms of type AT with weights of type R and iterate of type IT. An active set is typically expected to have a field weights, a field atoms, and a field x. Otherwise, all active set methods from src/active_set.jl can be overwritten.
FrankWolfe.ActiveSet — Type
ActiveSet{AT, R, IT}Represents an active set of extreme vertices collected in a FW algorithm, along with their coefficients (λ_i, a_i). R is the type of the λ_i, AT is the type of the atoms a_i. The iterate x = ∑λ_i a_i is stored in x with type IT.
FrankWolfe.DeletedVertexStorage — Type
Vertex storage to store dropped vertices or find a suitable direction in lazy settings. The algorithm will look for at most return_kth suitable atoms before returning the best. See Extra-lazification with a vertex storage for usage.
A vertex storage can be any type that implements two operations:
Base.push!(storage, atom)to add an atom to the storage.
Note that it is the storage type responsibility to ensure uniqueness of the atoms present.
storage_find_argmin_vertex(storage, direction, lazy_threshold) -> (found, vertex)
returning whether a vertex with sufficient progress was found and the vertex. It is up to the storage to remove vertices (or not) when they have been picked up.
FrankWolfe._unsafe_equal — Method
_unsafe_equal(a, b)Like isequal on arrays but without the checks. Assumes a and b have the same axes.
FrankWolfe.active_set_add_weight! — Method
active_set_add_weight!(active_set, lambda, i)Adds lambda to the weight of the ith atom in active_set.
FrankWolfe.active_set_argmin — Method
active_set_argmin(active_set::AbstractActiveSet, direction)Computes the linear minimizer in the direction on the active set. Returns (λ_i, a_i, i)
FrankWolfe.active_set_argminmax — Method
active_set_argminmax(active_set::AbstractActiveSet, direction)Computes the linear minimizer in the direction on the active set. Returns (λ_min, a_min, i_min, val_min, λ_max, a_max, i_max, val_max, val_max-val_min ≥ Φ)
FrankWolfe.active_set_initialize! — Method
active_set_initialize!(as, v)Resets the active set structure to a single vertex v with unit weight.
FrankWolfe.active_set_mul_weights! — Method
active_set_mul_weights!(active_set, lambda)Multiplies all weights in active_set by lambda.
FrankWolfe.active_set_update! — Method
active_set_update!(active_set::AbstractActiveSet, lambda, atom)Adds the atom to the active set with weight lambda or adds lambda to existing atom.
FrankWolfe.active_set_update_iterate_pairwise! — Method
active_set_update_iterate_pairwise!(active_set, x, lambda, fw_atom, away_atom)Operates x ← x + λ a_fw - λ a_aw.
FrankWolfe.active_set_update_pairwise! — Method
active_set_update_pairwise!(active_set, gamma, gamma_max, v_local_loc, a_loc, v_local, a, add_dropped_vertices, extra_vertex_storage)Updates the active set for a pairwise step with step size gamma.
FrankWolfe.active_set_update_scale! — Method
active_set_update_scale!(x, lambda, atom)Operates x ← (1-λ) x + λ a.
FrankWolfe.compute_active_set_iterate! — Method
compute_active_set_iterate!(active_set::AbstractActiveSet) -> xRecomputes from scratch the iterate x from the current weights and vertices of the active set. Returns the iterate x.
FrankWolfe.get_active_set_iterate — Method
get_active_set_iterate(active_set)Return the current iterate corresponding. Does not recompute it.
FrankWolfe.pre_computed_set_argminmax — Method
Computes the linear minimizer in the direction on the precomputedset. Precomputedset stores the vertices computed as extreme points v in each iteration.
FrankWolfe.storage_find_argmin_vertex — Method
Give the vertex v in the storage that minimizes s = direction ⋅ v and whether s achieves s ≤ lazy_threshold.
FrankWolfe.ActiveSetQuadraticPartialCaching — Type
ActiveSetQuadraticPartialCaching{AT, R, IT}Represents an active set of extreme vertices collected in a FW algorithm, along with their coefficients (λ_i, a_i). R is the type of the λ_i, AT is the type of the atoms a_i. The iterate x = ∑λ_i a_i is stored in x with type IT. The objective function is assumed to be of the form f(x)=½⟨x,Ax⟩+⟨b,x⟩+c so that the gradient is simply ∇f(x)=Ax+b.
In contrast to ActiveSetQuadraticProductCaching, this active set assumes that the Hessian matrix A is fixed up to a scaling factor λ and that the linear term b can change during the iterations. Therefore, we cache only the dot products ⟨A * x, a_i⟩ and ⟨A * a_i, a_j⟩. This active set might be used for experiments with Alternating Linear Minimization (ALM) or Split FW where we usually have chaning linear terms and Hessian matrices which only change in scale.
FrankWolfe.ActiveSetQuadraticProductCaching — Type
ActiveSetQuadraticProductCaching{AT, R, IT}Represents an active set of extreme vertices collected in a FW algorithm, along with their coefficients (λ_i, a_i). R is the type of the λ_i, AT is the type of the atoms a_i. The iterate x = ∑λ_i a_i is stored in x with type IT. The objective function is assumed to be of the form f(x)=½⟨x,Ax⟩+⟨b,x⟩+c so that the gradient is simply ∇f(x)=Ax+b.
FrankWolfe.ActiveSetQuadraticLinearSolve — Type
ActiveSetQuadraticLinearSolve{AT, R, IT}Represents an active set of extreme vertices collected in a FW algorithm, along with their coefficients (λ_i, a_i). R is the type of the λ_i, AT is the type of the atoms a_i. The iterate x = ∑λ_i a_i is stored in x with type IT. The objective function is assumed to be of the form f(x)=½⟨x,Ax⟩+⟨b,x⟩+c so that the gradient is ∇f(x)=Ax+b.
This active set stores an inner active_set that keeps track of the current set of vertices and convex decomposition. It therefore delegates all update, deletion, and addition operations to this inner active set. The weight, atoms, and x fields should only be accessed to read and are effectively the same objects as those in the inner active set. The flag wolfe_step determines whether to use a Wolfe step from the min-norm point algorithm or the normal direct solve. The Wolfe step solves the auxiliary subproblem over the affine hull of the current active set (instead of the convex hull).
The structure also contains a scheduler struct which is called with the should_solve_lp function. To define a new frequency at which the LP should be solved, one can define another scheduler struct and implement the corresponding method.
FrankWolfe.ActiveSetQuadraticLinearSolve — Method
ActiveSetQuadraticLinearSolve(tuple_values::Vector{Tuple{R,AT}}, A, b, lp_optimizer)Creates an ActiveSetQuadraticLinearSolve from the given Hessian A, linear term b and lp_optimizer by creating an inner ActiveSetQuadraticProductCaching active set.
FrankWolfe.ActiveSetQuadraticLinearSolve — Method
ActiveSetQuadraticLinearSolve(tuple_values::Vector{Tuple{R,AT}}, grad!::Function, lp_optimizer)Creates an ActiveSetQuadraticLinearSolve by computing the Hessian and linear term from grad!.
FrankWolfe.solve_quadratic_activeset_lp! — Method
solve_quadratic_activeset_lp!(as::ActiveSetQuadraticLinearSolve{AT, R, IT, H}))Solves the auxiliary LP over the current active set. The method is specialized by type H of the Hessian matrix A.
Functions and gradients
FrankWolfe.AbstractStochasticObjective — Type
AbstractStochasticObjectiveRepresents an objective function accessed through a zeroth and first-order oracle. It must implement the following functions:
compute_value(f::AbstractStochasticObjective, θ; batch_size, rng)grad!(storage, θ, x)compute_gradient(f::AbstractStochasticObjective, θ; batch_size, rng)
FrankWolfe.ObjectiveFunction — Type
ObjectiveFunctionRepresents an objective function optimized by algorithms. Subtypes of ObjectiveFunction must implement at least
compute_value(::ObjectiveFunction, x)for primal value evaluationcompute_gradient(::ObjectiveFunction, x)for gradient evaluation.
and optionally compute_value_gradient(::ObjectiveFunction, x) returning the (primal, gradient) pair. compute_gradient may always use the same storage and return a reference to it.
FrankWolfe.SimpleFunctionObjective — Type
SimpleFunctionObjective{F,G,S}An objective function built from separate primal objective f(x) and in-place gradient function grad!(storage, x). It keeps an internal storage of type s used to evaluate the gradient in-place.
FrankWolfe.StochasticObjective — Type
StochasticObjective{F, G, XT, S}(f::F, grad!::G, xs::XT, storage::S)Represents a finite sum function F(θ) = ∑_i f(θ, x_i) evaluated with stochastic gradient. f(θ, x) evaluates the loss for a single data point x and parameter θ. grad!(storage, θ, x) adds to storage the partial gradient with respect to data point x at parameter θ. xs must be an indexable iterable (Vector{Vector{Float64}} for instance). Functions using a StochasticObjective have optional keyword arguments rng, batch_size and full_evaluation controlling whether the function should be evaluated over all data points.
Note: grad! must not reset the storage to 0 before adding to it.
FrankWolfe.compute_gradient — Function
compute_gradient(f::ObjectiveFunction, x; [kwargs...])Computes the gradient of f at x. May return a reference to an internal storage.
FrankWolfe.compute_value — Function
compute_value(f::ObjectiveFunction, x; [kwargs...])Computes the objective f at x.
FrankWolfe.compute_value_gradient — Method
compute_value_gradient(f::ObjectiveFunction, x; [kwargs...])Computes in one call the pair (value, gradient) evaluated at x. By default, calls compute_value and compute_gradient with keywords kwargs passed down to both.
Callbacks
FrankWolfe.CallbackState — Type
CallbackStateMain structure created before and passed to the callback in first position.
Fields
tprimaldualdual_gaptimexvdgammafgrad!lmogradientstep_type
Custom vertex storage
Custom extreme point types
For some feasible sets, the extreme points of the feasible set returned by the LMO possess a specific structure that can be represented in an efficient manner both for storage and for common operations like scaling and addition with an iterate. They are presented below:
FrankWolfe.ScaledHotVector — Type
ScaledHotVector{T}Represents a vector of at most one value different from 0.
FrankWolfe.RankOneMatrix — Type
RankOneMatrix{T, UT, VT}Represents a rank-one matrix R = u * vt'. Composes like a charm.
FrankWolfe.NegatingArray — Type
Given an array array, NegatingArray represents -1 * array lazily.
FrankWolfe.SubspaceVector — Type
SubspaceVector{HasMultiplicities, T}Companion structure of SubspaceLMO containing three fields:
datais the full structure to be deflated,vecis a vector in the reduced subspace in which computations are performed,mulis only used to compute scalar products whenHasMultiplicities = true,
which should be avoided (when possible) for performance reasons.
Utils
FrankWolfe.muladd_memory_mode — Method
muladd_memory_mode(memory_mode::MemoryEmphasis, d, x, v)Performs d = x - v in-place or not depending on MemoryEmphasis
FrankWolfe.muladd_memory_mode — Method
(memory_mode::MemoryEmphasis, storage, x, gamma::Real, d)Performs storage = x - gamma * d in-place or not depending on MemoryEmphasis
FrankWolfe.muladd_memory_mode — Method
(memory_mode::MemoryEmphasis, x, gamma::Real, d)Performs x = x - gamma * d in-place or not depending on MemoryEmphasis
FrankWolfe.trajectory_callback — Method
trajectory_callback(storage)Callback pushing the state at each iteration to the passed storage. The state data is only the 5 first fields, usually: (t,primal,dual,dual_gap,time)
Oracle counting trackers
The following structures are wrapping given oracles to behave similarly but additionally track the number of calls.
FrankWolfe.TrackingObjective — Type
A function acting like the normal objective f but tracking the number of calls.
FrankWolfe.TrackingGradient — Type
A function acting like the normal grad! but tracking the number of calls.
FrankWolfe.TrackingLMO — Type
TrackingLMO{LMO}(lmo)An LMO wrapping another one and tracking the number of calls.
Also see the example Tracking, counters and custom callbacks for Frank Wolfe.
Update order for block-coordinate methods
Block-coordinate methods can be run with different update orders. All update orders are subtypes of FrankWolfe.BlockCoordinateUpdateOrder. They have to implement the method FrankWolfe.select_update_indices which selects which blocks to update in what order.
FrankWolfe.BlockCoordinateUpdateOrder — Type
Update order for a block-coordinate method. A BlockCoordinateUpdateOrder must implement
select_update_indices(::BlockCoordinateUpdateOrder, s::CallbackState, dual_gaps)FrankWolfe.select_update_indices — Function
select_update_indices(::BlockCoordinateUpdateOrder, s::CallbackState, dual_gaps)Returns a list of lists of the block indices. Each sublist represents one round of updates in an iteration. The indices in a list show which blocks should be updated parallely in one round. For example, a full update is given by [1:l] and a blockwise update by [[i] for i=1:l], where l is the number of blocks.
FrankWolfe.FullUpdate — Type
The full update initiates a parallel update of all blocks in one single round.
FrankWolfe.CyclicUpdate — Type
The cyclic update initiates a sequence of update rounds. In each round only one block is updated. The order of the blocks is determined by the given order of the LMOs.
FrankWolfe.StochasticUpdate — Type
The stochastic update initiates a sequence of update rounds. In each round only one block is updated. The order of the blocks is a random.
FrankWolfe.LazyUpdate — Type
The Lazy update order is discussed in "Flexible block-iterative analysis for the Frank-Wolfe algorithm," by Braun, Pokutta, & Woodstock (2024). 'lazyblock' is an index of a computationally expensive block to update; 'refreshrate' describes the frequency at which we perform a full activation; and 'blocksize' describes the number of "faster" blocks (i.e., those excluding 'lazyblock') activated (chosen uniformly at random) during each of the "faster" iterations; for more detail, see the article. If 'block_size' is unspecified, this defaults to
Note: This methodology is currently only proven to work with 'FrankWolfe.Shortstep' linesearches and a (not-yet implemented) adaptive method; see the article for details.
Update step for block-coordinate Frank-Wolfe
Block-coordinate Frank-Wolfe (BCFW) can run different FW algorithms on different blocks. All update steps are subtypes of FrankWolfe.UpdateStep and implement FrankWolfe.update_block_iterate which defines one iteration of the corresponding method.
FrankWolfe.UpdateStep — Type
Update step for block-coordinate Frank-Wolfe. These are implementations of different FW-algorithms to be used in a blockwise manner. Each update step must implement FrankWolfe.update_iterate.
FrankWolfe.update_block_iterate — Function
update_block_iterate(
step::UpdateStep,
x,
lmo,
f,
gradient,
grad!,
dual_gap,
t,
line_search,
linesearch_workspace,
memory_mode,
epsilon,
d,
)Executes one iteration of the defined FrankWolfe.UpdateStep and updates the iterate x implicitly. The function returns a tuple (dual_gap, v, d, gamma, step_type):
dual_gapis the updated FrankWolfe gapvis the used vertexdis the update directiongammais the applied step-sizestep_typeis the applied step-type
The d passed as argument is used as a container to avoid allocating d inside the function
FrankWolfe.FrankWolfeStep — Type
Implementation of the vanilla Frank-Wolfe algorithm as an update step for block-coordinate Frank-Wolfe.
FrankWolfe.BPCGStep — Type
Implementation of the blended pairwise conditional gradient (BPCG) method as an update step for block-coordinate Frank-Wolfe.
Block vector
FrankWolfe.BlockVector — Type
BlockVector{T, MT <: AbstractArray{T}, ST <: Tuple} <: AbstractVector{T}Represents a vector consisting of blocks. T is the element type of the vector, MT is the type of the underlying data array, and ST is the type of the tuple representing the sizes of each block. Each block can be accessed with the blocks field, and the sizes of the blocks are stored in the block_sizes field.
Index
FrankWolfe.AbstractActiveSetFrankWolfe.AbstractStochasticObjectiveFrankWolfe.ActiveSetFrankWolfe.ActiveSetQuadraticLinearSolveFrankWolfe.ActiveSetQuadraticLinearSolveFrankWolfe.ActiveSetQuadraticLinearSolveFrankWolfe.ActiveSetQuadraticPartialCachingFrankWolfe.ActiveSetQuadraticProductCachingFrankWolfe.BPCGStepFrankWolfe.BlockCoordinateUpdateOrderFrankWolfe.BlockVectorFrankWolfe.CallbackStateFrankWolfe.CyclicUpdateFrankWolfe.DeletedVertexStorageFrankWolfe.FrankWolfeStepFrankWolfe.FullUpdateFrankWolfe.LazyUpdateFrankWolfe.NegatingArrayFrankWolfe.ObjectiveFunctionFrankWolfe.RankOneMatrixFrankWolfe.ScaledHotVectorFrankWolfe.SimpleFunctionObjectiveFrankWolfe.StochasticObjectiveFrankWolfe.StochasticUpdateFrankWolfe.SubspaceVectorFrankWolfe.TrackingGradientFrankWolfe.TrackingLMOFrankWolfe.TrackingObjectiveFrankWolfe.UpdateStepBase.copyFrankWolfe._unsafe_equalFrankWolfe.active_set_add_weight!FrankWolfe.active_set_argminFrankWolfe.active_set_argminmaxFrankWolfe.active_set_initialize!FrankWolfe.active_set_mul_weights!FrankWolfe.active_set_update!FrankWolfe.active_set_update_iterate_pairwise!FrankWolfe.active_set_update_pairwise!FrankWolfe.active_set_update_scale!FrankWolfe.compute_active_set_iterate!FrankWolfe.compute_gradientFrankWolfe.compute_valueFrankWolfe.compute_value_gradientFrankWolfe.get_active_set_iterateFrankWolfe.muladd_memory_modeFrankWolfe.muladd_memory_modeFrankWolfe.muladd_memory_modeFrankWolfe.pre_computed_set_argminmaxFrankWolfe.select_update_indicesFrankWolfe.solve_quadratic_activeset_lp!FrankWolfe.storage_find_argmin_vertexFrankWolfe.trajectory_callbackFrankWolfe.update_block_iterate