| Title: | Fast and Accurate Randomized Singular Value Decomposition Algorithms with 'PCAone' |
|---|---|
| Description: | Fast and Accurate Randomized Singular Value Decomposition (RSVD) methods proposed in the 'PCAone' paper by Li (2023) <https://genome.cshlp.org/content/33/9/1599>. |
| Authors: | Zilong Li [aut, cre] |
| Maintainer: | Zilong Li <[email protected]> |
| License: | GPL (>= 3) |
| Version: | 1.5.0 |
| Built: | 2026-05-13 05:55:16 UTC |
| Source: | https://github.com/zilong-li/pcaoner |
dashSVD implements the Dynamic shifted Randomized SVD proposed by Feng et al. 2024
dashSVD(A, k, p = 10, s = 10)dashSVD(A, k, p = 10, s = 10)
A |
array_like; |
k |
integer; |
p |
integer, optional; |
s |
integer, optional; |
Zilong Li [email protected]
Feng, Xu and Yu, Wenjian and Xie, Yuyang and Tang, Jie. "Algorithm 1043: Faster Randomized SVD with Dynamic Shifts" (2024) doi:10.1145/3660629.
mat <- matrix(rnorm(100*20000), 20000, 100) res <- dashSVD(mat, k = 10) str(res)mat <- matrix(rnorm(100*20000), 20000, 100) res <- dashSVD(mat, k = 10) str(res)
The Randomized Singular Value Decomposition (RSVD) computes the near-optimal low-rank approximation of a rectangular matrix using a fast probablistic algorithm.
pcaone( A, k = NULL, p = 10, s = 20, method = "auto", B = 64, shuffle = TRUE, opts = list() ) ## S3 method for class 'matrix' pcaone( A, k = NULL, p = 10, s = 20, method = "auto", B = 64, shuffle = TRUE, opts = list() ) ## S3 method for class 'dgCMatrix' pcaone( A, k = NULL, p = 10, s = 20, method = "auto", B = 64, shuffle = TRUE, opts = list() ) ## S3 method for class 'dgRMatrix' pcaone( A, k = NULL, p = 7, s = 10, method = "winsvd", B = 64, shuffle = TRUE, opts = list() )pcaone( A, k = NULL, p = 10, s = 20, method = "auto", B = 64, shuffle = TRUE, opts = list() ) ## S3 method for class 'matrix' pcaone( A, k = NULL, p = 10, s = 20, method = "auto", B = 64, shuffle = TRUE, opts = list() ) ## S3 method for class 'dgCMatrix' pcaone( A, k = NULL, p = 10, s = 20, method = "auto", B = 64, shuffle = TRUE, opts = list() ) ## S3 method for class 'dgRMatrix' pcaone( A, k = NULL, p = 7, s = 10, method = "winsvd", B = 64, shuffle = TRUE, opts = list() )
A |
array_like; |
k |
integer; |
p |
integer, optional; |
s |
integer, optional; |
method |
string |
B |
integer, optional; |
shuffle |
logical, optional; |
opts |
list, optional; |
The singular value decomposition (SVD) plays an important role in data analysis, and scientific computing.
Given a rectangular matrix , and a target rank ,
the SVD factors the input matrix as
The left singular vectors are the columns of the
real or complex unitary matrix . The right singular vectors are the columns
of the real or complex unitary matrix . The dominant singular values are the
entries of , and non-negative and real numbers.
is an oversampling parameter to improve the approximation.
A value of at least 10 is recommended, and is set by default.
The parameter specifies the number of power (subspace) iterations
to reduce the approximation error. The power scheme is recommended,
especially when the singular values decay slowly. However, computing power iterations increases the
computational costs. Even though most RSVD implementations recommend power iterations by default,
it's always sufficient to run only few power iterations where our window-based power iterations ()
come to play. We recommend using and for pcaone algorithm2. As it is designed for large dataset,
we recommend using when .
If , a deterministic partial or truncated svd
algorithm might be faster.
pcaone returns a list containing the following three components:
array_like;
singular values; vector of length .
array_like;
left singular vectors; or dimensional array.
array_like;
right singular vectors; or dimensional array.
The singular vectors are not unique and only defined up to sign. If a left singular vector has its sign changed, changing the sign of the corresponding right vector gives an equivalent decomposition.
Zilong Li [email protected]
Z. Li, J Meisner, A Albrechtsen. "Fast and accurate out-of-core PCA framework for large scale biobank data" (2023) doi:10.1101/gr.277525.122.
library('pcaone') load(system.file("extdata", "popgen.rda", package="pcaone") ) A <- popgen - rowMeans(popgen) res <- pcaone(A, k = 40) str(res)library('pcaone') load(system.file("extdata", "popgen.rda", package="pcaone") ) A <- popgen - rowMeans(popgen) res <- pcaone(A, k = 40) str(res)
Four populations from East Asia in the 1000 Genome Project
An object of class pcaone.
https://www.internationalgenome.org
winSVD implements the window-based Randomized SVD proposed by Li et al. 2023
winSVD(A, k, p = 10, s = 10, B = 64, perm = 1)winSVD(A, k, p = 10, s = 10, B = 64, perm = 1)
A |
array_like; |
k |
integer; |
p |
integer, optional; |
s |
integer, optional; |
B |
integer, optional; |
perm |
integer, optional; |
winSVD returns a list containing the following three components:
array_like;
singular values; vector of length .
array_like;
left singular vectors; or dimensional array.
array_like;
right singular vectors; or dimensional array.
The singular vectors are not unique and only defined up to sign. If a left singular vector has its sign changed, changing the sign of the corresponding right vector gives an equivalent decomposition.
Zilong Li [email protected]
Z. Li, J Meisner, A Albrechtsen. "Fast and accurate out-of-core PCA framework for large scale biobank data" (2023) doi:10.1101/gr.277525.122.
mat <- matrix(rnorm(100*20000), 20000, 100) res <- winSVD(mat, k = 10) str(res)mat <- matrix(rnorm(100*20000), 20000, 100) res <- winSVD(mat, k = 10) str(res)