Title: | Sparse Truncated Singular Value Decomposition (from 'SVDLIBC') |
---|---|
Description: | Wrapper around the 'SVDLIBC' library for (truncated) singular value decomposition of a sparse matrix. Currently, only sparse real matrices in Matrix package format are supported. |
Authors: | Doug Rohde [aut], Michael Berry [aut], Theresa Do [aut], Gavin O'Brien [aut], Vijay Krishna [aut], Sowmini Varadhan [aut], University of Tennessee Research Foundation [cph] (files src/las2.c, src/svdlib.[ch], src/svdutil.[ch]), Stephanie Evert [cre, aut, cph] (copyright holder for files src/main.c, R/*, man/*, tests/*) |
Maintainer: | Stephanie Evert <[email protected]> |
License: | BSD_3_clause + file LICENSE |
Version: | 0.2-2 |
Built: | 2024-11-02 06:18:16 UTC |
Source: | CRAN |
Compute the (usually truncated) singular value decomposition (SVD) of a sparse real matrix. This function is a shallow wrapper around the SVDLIBC implementation of Berry's (1992) single Lanczos algorithm.
sparsesvd(M, rank=0L, tol=1e-15, kappa=1e-6)
sparsesvd(M, rank=0L, tol=1e-15, kappa=1e-6)
M |
a sparse real matrix in Matrix package format.
The preferred format is a |
rank |
an integer specifying the desired number of singular components, i.e. the rank of the truncated SVD.
Specify 0 to return all singular values of magnitude larger than |
tol |
exclude singular values whose magnitude is smaller than |
kappa |
accuracy parameter |
The truncated SVD decomposition
where is the optimal rank
approximation of
.
Note that
may be smaller than the requested number
rank
of singular components.
The returned value is a list with components
d |
a vector containing the first |
u |
a column matrix of the first |
v |
a column matrix of the first |
The SVDLIBC homepage http://tedlab.mit.edu/~dr/SVDLIBC/
seems to be no longer available.
A copy of the source code can be obtained from https://github.com/lucasmaystre/svdlibc.
Berry, Michael~W. (1992). Large scale sparse singular value computations. International Journal of Supercomputer Applications, 6, 13–49.
M <- rbind( c(20, 10, 15, 0, 2), c(10, 5, 8, 1, 0), c( 0, 1, 2, 6, 3), c( 1, 0, 0, 10, 5)) M <- Matrix::Matrix(M, sparse=TRUE) print(M) res <- sparsesvd(M, rank=2L) # compute first 2 singular components print(res, digits=3) M2 <- res$u %*% diag(res$d) %*% t(res$v) # rank-2 approximation print(M2, digits=1) print(as.matrix(M) - M2, digits=2) # approximation error
M <- rbind( c(20, 10, 15, 0, 2), c(10, 5, 8, 1, 0), c( 0, 1, 2, 6, 3), c( 1, 0, 0, 10, 5)) M <- Matrix::Matrix(M, sparse=TRUE) print(M) res <- sparsesvd(M, rank=2L) # compute first 2 singular components print(res, digits=3) M2 <- res$u %*% diag(res$d) %*% t(res$v) # rank-2 approximation print(M2, digits=1) print(as.matrix(M) - M2, digits=2) # approximation error