KernelDowndaters

A KernelDowndater is an abstract type used for computing kernel vectors of a subselection of the rows of the transpose of the matrix $V$.

Denote $V_{S*}$ to be the matrix formed by subselecting the indices specified by $S = \{s_1,\ldots,s_m\}$, and similarly, $x_{S}$ to be the subselection of elements of $x$. This is typically done by forming QR-factorizations of $V_{S*}$, and pulling a trailing column, $q$, from the Q-factor, such that ${V_{S*}}^T q = 0$. The reason we are interested in forming kernel vectors is that then, we can update a set of weights without changing the moments.

\[{V_{S*}}^T q = 0,\quad {V_{S*}}^T w_{S} = \eta \implies {V_{S*}}^T (w_{S} + \alpha q) = \eta,\]

where $\alpha$ is a constant chosen to zero out one of the weights while keeping all others nonnegative. We then update the $S$-indices of $w$ as $w_{S} = w_{S} + \alpha q$.

However, in Carathéodory pruning, at each step, $S$ is typically only changed by removing (or sometimes adding) a few elements at a time. Thus, it can be wasteful to fully recompute the QR decomposition at each step.

CaratheodoryPruning.jl comes with several built-in Downdater options. They can be easily used by calling caratheodory_pruning(V, w_in, kernel=:KERNEL), replacing :KERNEL with the appropriate symbol.

To implement your own KernelDowndater type, create a struct that inherits from KernelDowndater, and implement the necessary methods.

struct MyKernelDowndater <: KernelDowndater
    V::AbstractMatrix
    # Other necessary components
end

function get_inds(kd::MyKernelDowndater)
    error("Still to be implemented")
end
function get_kernel_vectors(kd::MyKernelDowndater)
    error("Still to be implemented")
end
function downdate!(kd::MyKernelDowndater, idx::Int)
    error("Still to be implemented")
end

You can then use that downdater, along with say the prune_weights_first! pruning method, as follows.

kd = MyKernelDowndater(V, additional_args)
caratheodory_pruning(V, w_in, kd, prune_weights_first!)

Below are the available KernelDowndater options implemented in CaratheodoryPruning.jl.

FullQRDowndater

Access with the kernel symbols :FullQRDowndater or :FullQR.

CaratheodoryPruning.FullQRDowndaterType

FullQRDowndater <: KernelDowndater

A mutable struct to hold the Q-factor of the QR decomposition of the matrix V[inds,:] for generating vectors in the kernel of its transpose. Downdates by forming a new QR factor, which takes O(M N²) flops where V is an M x N matrix.

Form with FullQRDowndater(V[; k=1]) where k is the (maximum) number of kernel vectors returned each time get_kernel_vectors is called.

source

GivensDowndater

Access with the kernel symbols :GivensDowndater or :Givens.

CaratheodoryPruning.GivensDowndaterType

GivensDowndater <: KernelDowndater

A mutable struct to hold the Q-factor of the QR decomposition of the matrix V[inds,:] for generating vectors in the kernel of its transpose. Downdates by applying Givens rotations to the old, full, Q-factor, which takes O(M²) flops where V is an M x N matrix.

Form with GivensDowndater(V[; k=1]) where k is the (maximum) number of kernel vectors returned each time get_kernel_vectors is called.

source

CholeskyDowndater

Access with the kernel symbols :CholeskyDowndater or :Cholesky.

CaratheodoryPruning.CholeskyDowndaterType

CholeskyDowndater <: KernelDowndater

A mutable struct to hold the Q-factor of the QR decomposition of the matrix V[inds,:] for generating vectors in the kernel of its transpose. Downdates by reorthogonalizing the old Q-factor, with a row removed, by multiplication by the inverse transpose of its Cholesky factor, which takes O(N³ + MN) flops where V is an M x N matrix.

Form with CholeskyDowndater(V[; k=1, pct_full_qr=10.0, SM_tol=1e-6, full_Q=false)]). k is the (maximum) number of kernel vectors returned each time get_kernel_vectors is called. pct_full_qr is the percentage (between 0 and 100), of times, logarithmically spaced, that a full QR reset will be done to prevent accumulation of error. SM_tol is a tolerance on the denominator of the Sherman Morrison formula to prevent error from division close to zero. full_Q determines whether or not the full Q matrix is updated or just its Cholesky factor; if set to true, will take O(N³ + MN²) flops instead of O(N³ + MN).

From testing, seems to have minimal error accumulation if pct_full_qr ≥ 10.0.

source

FullQRUpDowndater

Access with the kernel symbols :FullQRUpDowndater or :FullQRUpDown.

CaratheodoryPruning.FullQRUpDowndaterType

FullQRUpDowndater <: KernelDowndater

A mutable struct to hold the QR decomposition of the matrix V[inds,:] for generating vectors in the kernel of its transpose. Only acts on N+k indices at a time. When downdate is called, it removes that index, and adds one of the remaining index, calling a new full QR factorization to complete the down and update. Takes O((N+k)³) flops.

Form with FullQRUpDowndater(V[; ind_order=randperm(size(V,1)), k=1]). ind_order is the order in which the indices are added. k is the (maximum) number of kernel vectors returned each time get_kernel_vectors is called.

source

GivensUpDowndater

Access with the kernel symbols :GivensUpDowndater or :GivensUpDown.

CaratheodoryPruning.GivensUpDowndaterType

GivensUpDowndater <: KernelDowndater

A mutable struct to hold the QR decomposition of the matrix V[inds,:] for generating vectors in the kernel of its transpose. Only acts on N+k indices at a time. When downdate is called, it removes that index, and adds one of the remaining index, using Givens rotations to complete the down and update. Takes O((N+k)²) flops.

Form with GivensUpDowndater(V[; ind_order=randperm(size(V,1)), k=1, pct_full_qr=2.0]). ind_order is the order in which the indices are added. k is the (maximum) number of kernel vectors returned each time get_kernel_vectors is called. pct_full_qr is the percent of times, linearly spaced, that full QR factorizations are performed to prevent accumulation error in Q.

source