svds | R Documentation |
Given an m
by n
matrix A
,
function svds()
can find its largest k
singular values and the corresponding singular vectors.
It is also called the Truncated SVD or Partial SVD
since it only calculates a subset of the whole singular triplets.
Currently svds()
supports matrices of the following classes:
matrix | The most commonly used matrix type, defined in the base package. |
dgeMatrix | General matrix, equivalent to matrix ,
defined in the Matrix package. |
dgCMatrix | Column oriented sparse matrix, defined in the Matrix package. |
dgRMatrix | Row oriented sparse matrix, defined in the Matrix package. |
dsyMatrix | Symmetrix matrix, defined in the Matrix package. |
dsCMatrix | Symmetric column oriented sparse matrix, defined in the Matrix package. |
dsRMatrix | Symmetric row oriented sparse matrix, defined in the Matrix package. |
function | Implicitly specify the matrix through two
functions that calculate
f(x)=Ax and
g(x)=A'x . See section
Function Interface for details.
|
Note that when A
is symmetric and positive semi-definite,
SVD reduces to eigen decomposition, so you may consider using
eigs()
instead. When A
is symmetric but
not necessarily positive semi-definite, the left
and right singular vectors are the same as the left and right
eigenvectors, but the singular values and eigenvalues will
not be the same. In particular, if \lambda
is a negative
eigenvalue of A
, then |\lambda|
will be the
corresponding singular value.
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class 'matrix'
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class 'dgeMatrix'
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class 'dgCMatrix'
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class 'dgRMatrix'
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class 'dsyMatrix'
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class 'dsCMatrix'
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class 'dsRMatrix'
svds(A, k, nu = k, nv = k, opts = list(), ...)
## S3 method for class ''function''
svds(A, k, nu = k, nv = k, opts = list(), ..., Atrans, dim, args = NULL)
A |
The matrix whose truncated SVD is to be computed. |
k |
Number of singular values requested. |
nu |
Number of left singular vectors to be computed. This must
be between 0 and |
nv |
Number of right singular vectors to be computed. This must
be between 0 and |
opts |
Control parameters related to the computing algorithm. See Details below. |
... |
Arguments for specialized S3 function calls, for example
|
Atrans |
Only used when |
dim |
Only used when |
args |
Only used when |
The opts
argument is a list that can supply any of the
following parameters:
ncv
Number of Lanzcos basis vectors to use. More vectors
will result in faster convergence, but with greater
memory use. ncv
must be satisfy
k < ncv \le p
where
p = min(m, n)
.
Default is min(p, max(2*k+1, 20))
.
tol
Precision parameter. Default is 1e-10.
maxitr
Maximum number of iterations. Default is 1000.
center
Either a logical value (TRUE
/FALSE
), or a numeric
vector of length n
. If a vector c
is supplied, then
SVD is computed on the matrix A - 1c'
,
in an implicit way without actually forming this matrix.
center = TRUE
has the same effect as
center = colMeans(A)
. Default is FALSE
.
scale
Either a logical value (TRUE
/FALSE
), or a numeric
vector of length n
. If a vector s
is supplied, then
SVD is computed on the matrix (A - 1c')S
,
where c
is the centering vector and S = diag(1/s)
.
If scale = TRUE
, then the vector s
is computed as
the column norm of A - 1c'
.
Default is FALSE
.
A list with the following components:
d |
A vector of the computed singular values. |
u |
An |
v |
An |
nconv |
Number of converged singular values. |
niter |
Number of iterations used. |
nops |
Number of matrix-vector multiplications used. |
The matrix A
can be specified through two functions with
the following definitions
A <- function(x, args) { ## should return A %*% x } Atrans <- function(x, args) { ## should return t(A) %*% x }
They receive a vector x
as an argument and returns a vector
of the proper dimension. These two functions should have the effect of
calculating Ax
and A'x
respectively, and extra
arguments can be passed in through the
args
parameter. In svds()
, user should also provide
the dimension of the implicit matrix through the argument dim
.
The function interface does not support the center
and scale
parameters
in opts
.
Yixuan Qiu <https://statr.me>
eigen()
, svd()
,
eigs()
.
m = 100
n = 20
k = 5
set.seed(111)
A = matrix(rnorm(m * n), m)
svds(A, k)
svds(t(A), k, nu = 0, nv = 3)
## Sparse matrices
library(Matrix)
A[sample(m * n, m * n / 2)] = 0
Asp1 = as(A, "dgCMatrix")
Asp2 = as(A, "dgRMatrix")
svds(Asp1, k)
svds(Asp2, k, nu = 0, nv = 0)
## Function interface
Af = function(x, args)
{
as.numeric(args %*% x)
}
Atf = function(x, args)
{
as.numeric(crossprod(args, x))
}
svds(Af, k, Atrans = Atf, dim = c(m, n), args = Asp1)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.