COSMOAccelerators.AbstractAccelerator
— TypeAbstractAccelerator
Abstract supertype for acceleration objects that can be used to speed up a fixed-point iterations g = g(x) of a nonexpansive operator g
. They must implement the following methods to communicate with the fixed-point algorithm:
- update!(aa::AbstractAccelerator, g, x, num_iter) #stores the fixed-point iterates
- accelerate!(g::AbstractVector, x, aa::AbstractAccelerator, num_iter) #recombines past iterates to determine an accelerated point and overwrites
g
- restart!(aa::AbstractAccelerator, args...; kwargs...) # algorithm tells the accelerator to restart
The algorithm has to be able to query the following information:
- was_successful(aa::AbstractAccelerator) –> Bool #indicate whether accelerate! was succesful at the last iteration
The following method is optional:
- log!(aa::AbstractAccelerator, args...; kwargs...) # algorithm tells accelerator to log certain information for debugging
COSMOAccelerators.AbstractBroydenType
— TypeAbstractBroydenType
Abstract supertype for the Broyden type of the accelerator.
COSMOAccelerators.AbstractLeastSquaresMethod
— TypeAbstractLeastSquaresMethod
Abstract supertype for least squares solution method when Type-II acceleration is used
COSMOAccelerators.AbstractMemory
— TypeAbstractMemory
Abstract supertype for the memory management of the accelerator.
COSMOAccelerators.AbstractRegularizer
— TypeAbstractRegularizer
Abstract supertype for Anderson Acceleration Type-II regularisation schemes.
COSMOAccelerators.AndersonAccelerator
— TypeAndersonAccelerator{T, BT, MT, RE} <: AbstractAccelerator
Accelerator object implementing Anderson Acceleration. Parameterized by:
- T: AbstractFloat, floating-point type
- BT: Broyden-type, i.e. Type1 or Type2
- MT: AbstractMemory, how full memory buffers are handled
- RE: AbstractRegularizer
COSMOAccelerators.EmptyAccelerator
— TypeEmptyAccelerator <: AbstractAccelerator
An accelerator that does nothing.
COSMOAccelerators.RestartedMemory
— TypeRestartedMemory
The accelerator will delete the history once the memory buffers are full.
COSMOAccelerators.RollingMemory
— TypeRollingMemory
The accelerator will append new iterates to the history and discard the oldest iterate.
COSMOAccelerators.accelerate!
— MethodRecombine past iterates to compute an accelerated point. Overwrite g with accelerated point.
COSMOAccelerators.accelerate!
— MethodRecombine past iterates to compute an accelerated point. Overwrite g with accelerated point. Uses QR-decomposition.
COSMOAccelerators.apply_memory_approach!
— MethodDepending on the AndersonAccelerator parameter AbstractMemory, dispatch on the correct method that handles the case when the memory buffers are full.
COSMOAccelerators.assemble_inv_matrix!
— MethodAnderson type-1: eta = M^{-1} X' where M = (X'F).
COSMOAccelerators.assemble_inv_matrix!
— MethodAnderson type-2: eta = M^{-1} X' where M = (F'F).
COSMOAccelerators.compute_Δf!
— MethodCompute Δf = f - flast and store result in flast.
COSMOAccelerators.empty_caches!
— MethodReset the pointer that tracks the number of valid cols in the cached history of past vectors.
COSMOAccelerators.fill_caches!
— MethodUpdate the history matrices X = [Δxi, Δxi+1, ...], G = [Δgi, Δgi+1, ...] and F = [Δfi, Δfi+1, ...].
COSMOAccelerators.fill_delta!
— MethodUpdate a single history matrices V = [vgi, vgi+1, ...] .
COSMOAccelerators.initialise_eta!
— MethodAnderson Type 1: Initialise eta = X'f.
COSMOAccelerators.initialise_eta!
— MethodAnderson Type 2: Initialise eta = F'f.
COSMOAccelerators.regularize!
— MethodChoose regularisation parameter based on the frobenius norms of the matrices X, F (Fu, Zhang, Boyd, 2019).
COSMOAccelerators.regularize!
— MethodUse Tikonov - Regularizer for the Least squares matrix M.
COSMOAccelerators.restart!
— MethodPut accelerator into a fresh state, i.e. forget about past history.
COSMOAccelerators.set_prev_vectors!
— MethodStore a copy of the last x, g, f to be able to compute Δx, Δg, Δf at the next step.
COSMOAccelerators.update!
— Methodupdate!(aa, g, x, iter)
- Update history of accelerator
aa
with iterates g = g(xi) - Computes residuals f = x - g
- The iteration
iter
is passed in for logging