API Reference
Complete API documentation for SketchySVD.jl.
SketchySVD.SRFTMatrixSketchySVD.SketchySketchySVD._compute_boundSketchySVD.gaussian_rngSketchySVD.init_sketchySketchySVD.optimal_rhoSketchySVD.rademacher_rngSketchySVD.rsvdSketchySVD.rsvd_adaptiveSketchySVD.rsvd_standardSketchySVD.rsvd_transposeSketchySVD.sparse_gaussian_rngSketchySVD.sparse_rngSketchySVD.srft_rngSketchySVD.tail_energySketchySVD.tail_energy_squaredSketchySVD.truncated_approx_error_boundSketchySVD.truncated_approx_error_boundSketchySVD.uniform_rng
Public API
These are the main functions and types you'll use:
SketchySVD.Sketchy — Type
SketchySketchy is a randomized linear dimension reduction map that is used to reduce the dimension of the data matrix in the data-streaming setting. It is based on the sketching technique that is used to approximate the Singular Value Decomposition of the data matrix. The Sketchy algorithm is used to compute the SVD of the data matrix in an incremental fashion. For more details see [TYUC2019].
Fields
V::AbstractArray{T}: Left singular vector matrix.Σ::AbstractArray{T}: Singular value matrix.W::AbstractArray{T}: Right singular vector matrix.Ξ::AbstractArray{T}: Test matrix for range (k x m).Ω::AbstractArray{T}: Test matrix for corange (k x n).Φ::AbstractArray{T}: Test matrix for core (s x m).Ψ::AbstractArray{T}: Test matrix for core (s x n).Θ::AbstractArray{T}: Gaussian test matrix for error (q x m).X::AbstractArray{T}: Corange sketch (k x n).Y::AbstractArray{T}: Range sketch (m x k).Z::AbstractArray{T}: Core sketch (s x s).E::AbstractArray{T}: Error sketch (q x n).ct::Integer: Counter for the number of data streams.dims::Dict{<:Symbol,<:Integer}: Dictionary to store dimensions of the data matrix and sketching parameters.α::Integer: Constant depending on data type (1 for real, 0 for complex).β::Integer: Constant depending on data type (1 for real, 2 for complex).verbose::Bool: Flag for verbosity.ErrorEstimate::Bool: Flag to estimate the error.SpectralDecay::Bool: Flag to compute the scree plot.
References
- [TYUC2019] J. A. Tropp, A. Yurtsever, M. Udell, and V. Cevher, “Streaming Low-Rank Matrix Approximation with an Application to Scientific Simulation,” SIAM J. Sci. Comput., vol. 41, no. 4, pp. A2430-A2463, Jan. 2019, doi: 10.1137/18M1201068.
- [TYUC2019SUP] J. A. Tropp, A. Yurtsever, M. Udell, and V. Cevher, “Suplementary Materials: Streaming Low-Rank Matrix Approximation with an Application to Scientific Simulation,” SIAM J. Sci. Comput., vol. 41, no. 4, pp. A2430-A2463, Jan. 2019, doi: 10.1137/18M1201068.
SketchySVD.gaussian_rng — Method
gaussian_rng(m, n)Standard Gaussian random matrix generator (default).
SketchySVD.init_sketchy — Method
init_sketchy(; m, n, r, k=0, s=0, q=0, T=0, ReduxMap=:sparse, dtype="real",
verbose=false, ErrorEstimate=false, SpectralDecay=false, ζ=min(k,8))Initialize the Sketchy algorithm with the given parameters.
Arguments
m: Row size of the data matrix.n: Column size of the data matrix.r: Selected rank to truncate the SVD.k: Range dimension for sketching (default: 0, will be set based on r).s: Core dimension for sketching (default: 0, will be set based on k).q: Error dimension for sketching (default: 0, will be set based on m).T: Storage budget for sketching (default: 0, will be ignored if k and s are provided).ReduxMap: Type of reduction map for sketching (default: :sparse).dtype: Data type for the data (default: "real").verbose: Verbosity flag (default: false).ErrorEstimate: Flag to estimate the error (default: false).SpectralDecay: Flag to compute the scree plot (default: false).ζ: Sparsity parameter for sparse reduction map (default: min(k,8)).
Returns
An instance of the Sketchy struct initialized with the specified parameters and random matrices.
SketchySVD.rademacher_rng — Method
rademacher_rng(m, n)Rademacher random matrix (entries are ±1 with equal probability). More efficient than Gaussian for some applications.
SketchySVD.rsvd — Method
rsvd(A::AbstractArray, k::Int; p::Int=5, q::Int=0,
rng::Union{DimRedux,Function}=randn,
transpose_trick::Bool=true)Compute a rank-k randomized SVD approximation of matrix A using power iterations and orthonormalization for improved accuracy.
Arguments
A::AbstractArray: Input matrix (m x n)k::Int: Target rank for the approximation
Keyword Arguments
p::Int=5: Oversampling parameter (samples k+p random vectors)q::Int=0: Number of power iterations (0 means no power iterations)rng::Union{DimRedux,Function}=Gauss: Random matrix generator function with signature (nrows, ncols)transpose_trick::Bool=true: Use transpose trick when m >> n for efficiency
Returns
V: Left singular vectors (m x k)Σ: Singular values (vector of length k)W: Right singular vectors (n x k)
Examples
# Basic usage with default Gaussian random matrix
V, Σ, W = rsvd(A, 50)
# With power iterations for better accuracy
V, Σ, W = rsvd(A, 50, q=2)
# With custom random matrix (e.g., sparse random)
using SparseArrays
sparse_rng(m, n) = sprandn(m, n, 0.1) # 10% density
V, Σ, W = rsvd(A, 50, rng=sparse_rng)
# With Rademacher (±1) random matrix (faster generation)
V, Σ, W = rsvd(A, 50, rng=rademacher_rng)
# With SRFT (Subsampled Randomized Fourier Transform) - very efficient
V, Σ, W = rsvd(A, 50, rng=srft_rng, q=1)SketchySVD.rsvd_adaptive — Method
rsvd_adaptive(A, k; tol=1e-10, max_iter=10, p=5, q=0, rng=randn)Adaptive randomized SVD that automatically determines the rank based on singular value decay.
Arguments
A: Input matrixk: Initial target rank estimate
Keyword Arguments
tol::Float64=1e-10: Tolerance for singular value cutoffmax_iter::Int=10: Maximum number of adaptive iterationsp::Int=5: Oversampling parameterq::Int=0: Number of power iterationsrng::Function=randn: Random matrix generator
Returns
U, S, V: SVD factors where rank is automatically determined
SketchySVD.sparse_gaussian_rng — Method
sparse_gaussian_rng(density)Create a sparse Gaussian random matrix generator with specified density.
SketchySVD.sparse_rng — Method
sparse_rng(m, n; zeta=8, field="real")Sparse random matrix generator based on the Sparse dimension reduction structure. Creates a sparse matrix with approximately zeta non-zero entries per column, where each entry is randomly ±1 (real) or a complex sign.
This is more efficient than dense Gaussian matrices and is particularly useful for large-scale randomized SVD computations.
Arguments
m::Int: Number of rows (output dimension)n::Int: Number of columns (input dimension)
Keyword Arguments
zeta::Int=8: Number of non-zero entries per column (default: 8, clamped to min(m,8))field::String="real": "real" for ±1 entries, "complex" for complex unit entries
Returns
Sparse: sparse dimension reduction object
Notes
- Returns a dense matrix for compatibility with existing rsvd code
- For very large problems, consider using the
Sparseobject directly withLeftApplyorRightApplyfor memory efficiency - The sparsity pattern is random with
zetanon-zeros per column
SketchySVD.srft_rng — Method
srft_rng(m, n)Subsampled Randomized Fourier Transform (SRFT) matrix generator. More efficient than Gaussian random matrices, especially for large problems.
The SRFT is an implicit structured random matrix of the form: SRFT = √(m/k) · R · F · D
where:
- D is a diagonal matrix with random ±1 entries (size m × m)
- F is the DFT matrix (applied via FFT)
- R is a random row sampling matrix (selects n rows from m)
The resulting matrix has size (n, m).
This implementation doesn't form the matrix explicitly but returns a lazy operator that applies the SRFT efficiently.
SketchySVD.uniform_rng — Method
uniform_rng(m, n)Uniform random matrix generator on [-1, 1].
Internal Implementation
These are internal functions used by the package. Most users won't need these directly.
SketchySVD.SRFTMatrix — Type
SRFTMatrixLazy representation of a Subsampled Randomized Fourier Transform matrix. The matrix has size (mout, nin): takes nin-dimensional vectors and produces mout-dimensional outputs.
SketchySVD._compute_bound — Method
Internal function to compute the bound given singular values and parameters.
SketchySVD.optimal_rho — Method
optimal_rho(σ::AbstractVector, k::Int; field::Symbol=:real)Find the optimal ρ that minimizes the term in the error bound.
Returns (ρopt, minvalue).
SketchySVD.rsvd_standard — Method
rsvd_standard(A, k, p, q, rng)Standard randomized SVD implementation (for m ≤ n or when transpose trick is disabled).
SketchySVD.rsvd_transpose — Method
rsvd_transpose(A, k, p, q, rng)Transpose trick for tall-thin matrices (m >> n). More efficient when the matrix has many more rows than columns.
SketchySVD.tail_energy — Method
tail_energy(σ::AbstractVector, r::Int)Compute the tail energy τ²ᵣ₊₁(A) = Σⱼ₌ᵣ₊₁ σⱼ², where σ contains the singular values of A.
Returns τᵣ₊₁(A) (the square root of the sum).
SketchySVD.tail_energy_squared — Method
tail_energy_squared(σ::AbstractVector, r::Int)Compute τ²ᵣ₊₁(A) = Σⱼ₌ᵣ₊₁ σⱼ².
SketchySVD.truncated_approx_error_bound — Method
truncated_approx_error_bound(A::AbstractMatrix, r::Int, k::Int, s::Int;
field::Symbol=:real)Compute the error bound from Corollary 5.5 for the rank-r truncated approximation.
Arguments
A: The input matrixr: Target rank for the truncated approximationk: Range sketch size parameter (number of columns in the random test matrix Ω)s: Core sketch size parameter (total sketch dimension, s > k)field: Either:real(default) or:complex
Returns
The expected spectral norm error bound: E‖A - [[Â]]ᵣ‖₂
Mathematical formulation
The bound is: τᵣ₊₁(A) + 2[((s-α)/(s-k-α)) · min{ρ < k-α} ((k+ρ-α)/(k-ρ-α)) · τ²ᵨ₊₁(A)]^(1/2)
where α = 1 for real matrices, α = 0 for complex matrices.
SketchySVD.truncated_approx_error_bound — Method
truncated_approx_error_bound(σ::AbstractVector, r::Int, k::Int, s::Int;
field::Symbol=:real)Compute the error bound directly from singular values.