Approximators
Abstract Types
RLinearAlgebra.Approximator
— TypeApproximator
An abstract supertype for structures that store user-controlled parameters corresponding to techniques that form low-rank approximations of the matrix A
.
RLinearAlgebra.ApproximatorRecipe
— TypeApproximatorRecipe
An abstract supertype for structures that store user-controlled parameters, linear system dependent parameters and preallocated memory corresponding to techniques that form low-rank approximations of the matrix A
.
RLinearAlgebra.ApproximatorAdjoint
— TypeApproximatorAdjoint{S<:ApproximatorRecipe} <: ApproximatorRecipe
A structure for the adjoint of an ApproximatorRecipe
.
Fields
Parent::ApproximatorRecipe
, the approximator that we compute the adjoint of.
RLinearAlgebra.RangeApproximator
— TypeRangeApproximator
An abstract type for the structures that contain the user-controlled parameters corresponding to the Approximator methods that produce an orthogonal approximation to the range of a matrix A. This includes methods like the RandomizedSVD and randomized rangefinder.
RLinearAlgebra.RangeApproximatorRecipe
— TypeRangeApproximatorRecipe
An abstract type for the structures that contain the user-controlled parameters, linear system information, and preallocated memory for methods corresponding to the Approximator methods that produce an orthogonal approximation to the range of a matrix A. This includes methods like the RandomizedSVD and randomized rangefinder.
Approximator Structures
RLinearAlgebra.RandSVD
— TypeRandSVD
A struct that implements the Randomized SVD. The Randomized SVD technique compresses a matrix from the right to compute a rank $k$ estimate to the truncated SVD of a matrix $A$. See Algorithm 5.1 in [9] for additional details.
Mathematical Description
Suppose we have a matrix $A \in \mathbb{R}^{m \times n}$ for which we wish to form a low-rank approximation with the form of an SVD. Specifically, we wish to find an orthogonal matrix, $U$, a diagonal matrix, $\Sigma$, and an orthogonal matrix, $V$ such that $U\Sigma V^\top \approx A$. This can be done effectively using the Randomized SVD presented in Algorithm 5.1 of [9]. This algorithm progresses by first having the user select a $k > 2$ to be the compression_dim
of the compressor, which will also correspond to the desired rank of the approximation. With this $k$, RandomizedSVD generates a compression matrix $S\in\mathbb{R}^{n \times k}$ and computes $Q = \text{qr}(AS)$ as in the RangeFinder. With the $Q$, the RandomizedSVD concludes by computing $W,\Sigma,V = \text{SVD}(Q^\top A)$ and setting $U = QW$. With high probability the approximation generated by the RandomizedSVD will satisfy $\mathbb{E} \|A - U \Sigma V^\top \|_F \leq \sqrt{k+1} (\sum_{i=k+1}^{\min{(m,n)}}\sigma_{i})^{1/2}$ , where $\sigma_{k+1}$ is the $k+1^\text{th}$ singular value of (see Theorem 10.5 of [9]). This bound is often conservative as long as the singular values of $A$ decay quickly.
When the singular values decay slowly, we can improve the quality of the approximation using the power iteration, which applies $A$ and $A^\top$, $q$ times and take the qr factorization of $(AA^\top)^q AS$. Using these power iterations increases the relative gap between the singular values leading to better RandomizedSVD performance.
Performing power iterations in floating points can destroy all information related to the smallest singular values of $A$ (see Remark 4.3 in [9]). We can preserve this information by orthogonalizing inbetween the products of $AS$ with $A$ or $A^\top$ in the power iteration. These steps are known as the orthogonalized power iteration (see Algorithm 4.4 of [9]). Orthogonalized power iterations progress according to the following steps:
- $\tilde{A}_1 = AS$
- $Q_1,R_1 = \textbf{qr}(\tilde{A}_1)$
- $\tilde{A}_2 = A^\top Q_1$
- $Q_2,R_2 = \textbf{qr}(\tilde{A}_2)$
- $\tilde{A}_1 = A Q_2$
- $Q_1, R_1 = \textbf{qr}(\tilde{A}_1)$
- Repeat Steps 3 through 6 for the desired number of power iterations set $Q = Q_1$.
Fields
compressor::Compressor
, the technique for compressing the matrix from the right.power_its::Int64
, the number of power iterations that should be performed.orthogonalize::Bool
, a boolean indicating whether thepower_its
should be performed with orthogonalization.block_size::Int64
, the size of the tile when performing matrix multiplication. By default,block_size = 0
, this will be set to be the number of columns in the original matrix, $A$.
Constructor
RandSVD(;
compressor::Compressor = SparseSign(cardinality = Right()),
orthogonalize::Bool = false,
power_its::Int64 = 0,
block_size::Int64 = 0,
)
Keywords
compressor::Compressor
, the technique for compressing the matrix from the right. By default this isSparseSign
with a cardinalityRight()
.power_its::Int64
, the number of power iterations that should be performed. By default this is zero.orthogonalize::Bool
, a boolean indicating whether thepower_its
should be performed with orthogonalization. By default is false.block_size::Int64
, the size of the tile when performing matrix multiplication. By default,block_size = 0
, this will be set to be the number of columns in the original matrix, $A$. By default this is zero.
Returns
- A
RandSVD
object.
Throws
ArgumentError
ifpower_its
orblock_size
are negative.
RLinearAlgebra.RandSVDRecipe
— TypeRandSVDRecipe
A struct that contains the preallocated memory and completed compressor to form a RandSVD approximation to the matrix $A$.
Fields
n_rows::Int64
, the number of rows in the approximation.n_cols::Int64
, the number of columns in the approximation.compressor::CompressorRecipe
, the compressor to be applied from the right to $A$.power_its::Int64
, the number of power iterations that should be performed.orthogonalize::Bool
, a boolean indicating whether thepower_its
should be performed with orthogonalization.U::AbstractArray
, the orthogonal matrix that approximates the topcompressor_dim
left singular vectors of $A$.S::AbstractVector
, a vector containing the topcompressor_dim
singular values of $A$.V::AbstractArray
, the orthogonal matrix that approximates the topcompressor_dim
right singular vectors of $A$.buffer::AbstractArray
, the storage for matrix multiplication with this low-rank approximation. Will have thecompression_dim
number of rows andblock_size
number of columns.
RLinearAlgebra.RangeFinder
— TypeRangeFinder
A struct that implements the Randomized Range Finder technique which uses compression from the right to form an low-dimensional orthogonal matrix $Q$ that approximates the range of $A$. See [9] for additional details.
Mathematical Description
Suppose we have a matrix $A \in \mathbb{R}^{m \times n}$ of which we wish to form a low rank approximation that approximately captures the range of $A$. Specifically, we wish to find an Orthogonal matrix $Q$ such that $QQ^\top A \approx A$.
A simple way to find such a matrix is to choose a $k$ representing the number of vectors we wish to have in the subspace. Then we can generate a compression matrix $S\in\mathbb{R}^{n \times k}$ and compute $Q = \text{qr}(AS)$. With high probability we will have $\|A - QQ^\top A\|_2 \leq \sqrt{k+1} (\sum_{i=k+1}^{\min{(m,n)}}\sigma_{i})^{1/2}$, where $\sigma_{k+1}$ is the $k+1^\text{th}$ singular value of A (see Theorem 10.5 of [9]). This bound is often conservative as long as the singular values of $A$ decay quickly.
When the singular values decay slowly, we can improve the quality of the approximation using the power iteration, which applies $A$ and $A^\top$, $q$ times and take the qr factorization of $(AA^\top)^q AS$. Using these power iterations increases the relative gap between the singular values leading to better Rangefinder performance.
Performing power iterations in floating points can destroy all information related to the smallest singular values of $A$ (see Remark 4.3 in [9]). We can preserve this information by orthogonalizing inbetween the products of $AS$ with $A$ or $A^\top$ in the power iteration. These steps are known as the orthogonalized power iteration (see Algorithm 4.4 of [9]). Orthogonalized power iterations progress according to the following steps:
- $\tilde{A}_1 = AS$
- $Q_1,R_1 = \textbf{qr}(\tilde{A}_1)$
- $\tilde{A}_2 = A^\top Q_1$
- $Q_2,R_2 = \textbf{qr}(\tilde{A}_2)$
- $\tilde{A}_1 = A Q_2$
- $Q_1, R_1 = \textbf{qr}(\tilde{A}_1)$
- Repeat Steps 3 through 6 for the desired number of power iterations set $Q = Q_1$.
Fields
compressor::Compressor
, the technique that will compress the matrix from the right.power_its::Int64
, the number of power iterations that should be performed.orthogonalize::Bool
, a boolean indicating whether thepower_its
should be performed with orthogonalization.
Constructor
RangeFinder(;
compressor = SparseSign(),
orthogonalize = false,
power_its = 0
)
Keywords
compressor::Compressor
, the technique that will compress the matrix from the right.power_its::Int64
, the number of power iterations that should be performed. Default is zero.orthogonalize::Bool
, a boolean indicating whether thepower_its
should be performed with orthogonalization. Default is false.
Returns
- A
RangeFinder
object.
Throws
ArgumentError
ifpower_its
is negative.
RLinearAlgebra.RangeFinderRecipe
— TypeRangeFinderRecipe
A struct that contains the preallocated memory and completed compressor to form a RangeFinder approximation to the matrix $A$.
Fields
compressor::CompressorRecipe
, the compressor to be applied from the right to $A$.power_its::Int64
, the number of power iterations that should be performed.orthogonalize::Bool
, a boolean indicating whether thepower_its
should be performed with orthogonalization.range::AbstractMatrix
, the orthogonal matrix that approximates the range of $A$.
Exported Functions
RLinearAlgebra.complete_approximator
— Functioncomplete_approximator(approximator::Approximator, A::AbstractMatrix)
A function that generates an ApproximatorRecipe
given arguments.
Arguments
approximator::Approximator
, a data structure containing the user-defined parameters associated with a particular low-rank approximation.A::AbstractMatrix
, a target matrix for approximation.
Outputs
- An
ApproximatorRecipe
object.
RLinearAlgebra.rapproximate
— Functionrapproximate(approximator::Approximator, A::AbstractMatrix)
A function that computes a low-rank approximation of the matrix A
using the information in the provided Approximator
data structure.
Arguments
approximator::Approximator
, a data structure containing the user-defined parameters associated with a particular low-rank approximation.A::AbstractMatrix
, a target matrix for approximation.
Outputs
- An
ApproximatorRecipe
object.
RLinearAlgebra.rapproximate!
— Functionrapproximate!(approximator::ApproximatorRecipe, A::AbstractMatrix)
A function that computes a low-rank approximation of the matrix A
using the information in the provided Approximator
data structure.
Arguments
approximator::ApproximatorRecipe
, a fully initialized realization for a low rank approximation method for a particular matrix.A::AbstractMatrix
, a target matrix for approximation.
Outputs
- An
ApproximatorRecipe
object.
Internal Functions
RLinearAlgebra.rand_power_it
— Functionrand_power_it(approx::RangeFinderRecipe, A::AbstractMatrix)
Function that performs the randomized rangefinder procedure presented in Algortihm 4.3 of [9].
Arguments
approx::RangeFinderRecipe
, aRangeFinderRecipe
structure that contains the compressor
recipe.
A::AbstractMatrix
, the matrix being approximated.
Returns
Q::AbstractMatrix
, an economicalQ
approximating the range of A.
RLinearAlgebra.rand_ortho_it
— Functionrand_ortho_it(approx::RangeFinderRecipe, A::AbstractMatrix)
Function that performs the randomized rangefinder procedure presented in Algortihm 4.4 of [9].
Arguments
approx::RangeFinderRecipe
, aRangeFinderRecipe
structure that contains the compressor
recipe.
A::AbstractMatrix
, the matrix being approximated.
Returns
Q::AbstractMatrix
, an economicalQ
approximating the range of A.