Introduction
Eigenvalue decomposition is a commonly used technique in numerous statistical problems. For example, principal component analysis (PCA) basically conducts eigenvalue decomposition on the sample covariance of a data matrix: the eigenvalues are the component variances, and eigenvectors are the variable loadings.
In R, the standard way to compute eigenvalues is the
eigen() function. However, when the matrix becomes large,
eigen() can be very time-consuming: the complexity to
calculate all eigenvalues of a \(n \times
n\) matrix is \(O(n^3)\).
While in real applications, we usually only need to compute a few
eigenvalues or eigenvectors, for example to visualize high dimensional
data using PCA, we may only use the first two or three components to
draw a scatterplot. Unfortunately in eigen(), there is no
option to limit the number of eigenvalues to be computed. This means
that we always need to do the full eigen decomposition, which can cause
a huge waste in computation.
The same thing happens in Singular Value Decomposition (SVD). It is often the case that only a Partial SVD or Truncated SVD is needed, and moreover the matrix is usually stored in sparse format. Base R does not provide functions suitable for these special needs.
And this is why the RSpectra package was developed. RSpectra provides an R interface to the Spectra library, which is used to solve large scale eigenvalue problems. The core part of Spectra and RSpectra was written in efficient C++ code, and they can handle large scale matrices in either dense or sparse formats.
Eigenvalue Problem
Basic Usage
The RSpectra package provides functions
eigs() and eigs_sym() to calculate eigenvalues
of general and symmetric matrices respectively. If the matrix is known
to be symmetric, eigs_sym() is preferred since it
guarantees that the eigenvalues are real.
To obtain eigenvalues of a square matrix A, simply call
the eigs() or eigs_sym() function, tell it how
many eigenvalues are requested (argument k), and which ones
are going to be selected (argument which). By default,
which = "LM" means to pick the eigenvalues with the largest
magnitude (modulus for complex numbers and absolute value for real
numbers).
Below we first generate some random matrices and display some of the eigenvalues and eigenvectors:
set.seed(123)
n = 100 # matrix size
k = 5 # number of eigenvalues to calculate
# Some random data
M = matrix(rnorm(n^2), n)
# Make it symmetric
A = crossprod(M)
# Show its largest 5 eigenvalues and associated eigenvectors
head(eigen(A)$values, 5)## [1] 391.3649 367.1143 353.5425 322.6301 315.4341
head(eigen(A)$vectors[, 1:5])## [,1] [,2] [,3] [,4] [,5]
## [1,] 0.12981428 -0.011305031 -0.0001022712 0.088485329 0.05035607
## [2,] -0.06690188 0.001012109 0.0103448033 -0.009411417 -0.07992121
## [3,] -0.05333064 -0.057256548 0.0998913793 -0.171712454 0.01334678
## [4,] 0.13423510 -0.165303749 0.0076458165 -0.006789840 0.10246642
## [5,] -0.13023611 -0.019154947 0.0425667905 -0.160273074 -0.06246803
## [6,] -0.06632982 -0.053457170 0.0323196492 -0.074827181 -0.16365506
Now we use RSpectra to directly obtain the largest 5 eigenvalues:
library(RSpectra)
res = eigs_sym(A, k, which = "LM") # "LM" is the default
res$values## [1] 391.3649 367.1143 353.5425 322.6301 315.4341
head(res$vectors)## [,1] [,2] [,3] [,4] [,5]
## [1,] -0.12981428 0.011305031 0.0001022712 0.088485329 0.05035607
## [2,] 0.06690188 -0.001012109 -0.0103448033 -0.009411417 -0.07992121
## [3,] 0.05333064 0.057256548 -0.0998913793 -0.171712454 0.01334678
## [4,] -0.13423510 0.165303749 -0.0076458165 -0.006789840 0.10246642
## [5,] 0.13023611 0.019154947 -0.0425667905 -0.160273074 -0.06246803
## [6,] 0.06632982 0.053457170 -0.0323196492 -0.074827181 -0.16365506
If only eigenvalues are requested, the retvec option can
be used:
eigs_sym(A, k, opts = list(retvec = FALSE))## $values
## [1] 391.3649 367.1143 353.5425 322.6301 315.4341
##
## $vectors
## NULL
##
## $nconv
## [1] 5
##
## $niter
## [1] 4
##
## $nops
## [1] 59
Matrix in Sparse Format
For really large data, the matrix is usually stored in sparse format.
RSpectra supports the dgCMatrix and
dgRMatrix matrix types defined in the
Matrix package.
library(Matrix)
Msp = as(M, "dgCMatrix")
Asp = as(A, "dgRMatrix")
eigs(Msp, k, which = "LR", opts = list(retvec = FALSE))$values # largest real part## [1] 9.685399-1.964917i 9.685399+1.964917i 9.696766+0.000000i 8.270271+1.751652i
## [5] 8.270271-1.751652i
eigs_sym(Asp, k, opts = list(retvec = FALSE))$values## [1] 391.3649 367.1143 353.5425 322.6301 315.4341
An even more general way to specify the matrix A is to
define a function that calculates A %*% x for any vector
x. This representation is called the Function
Interface, which can also be regarded as a sparse format, since
users do not need to store the matrix A, but only the
operator x -> Ax. For example, to represent a diagonal
matrix \(diag(1, 2, \ldots, 10)\) and
calculate its eigenvalues, we can define the function f and
call the eigs_sym() function:
# Implicitly define the matrix by a function that calculates A %*% x
# Below represents a diagonal matrix whose elements are stored in the `args` parameter
f = function(x, args)
{
# diag(args) %*% x == x * args
return(x * args)
}
eigs_sym(f, k = 3, n = 10, args = 1:10)## $values
## [1] 10 9 8
##
## $vectors
## [,1] [,2] [,3]
## [1,] 2.307346e-16 1.087717e-16 -1.186185e-17
## [2,] 9.737491e-17 8.384949e-17 5.104424e-16
## [3,] -6.923647e-17 -7.659888e-17 5.290907e-17
## [4,] -6.640738e-18 -7.632783e-17 -7.285839e-17
## [5,] -8.239937e-18 8.673617e-18 2.775558e-17
## [6,] 1.368195e-16 2.836273e-16 -3.677614e-16
## [7,] 3.652440e-16 5.412879e-16 2.547875e-16
## [8,] -8.036987e-17 9.469422e-16 1.000000e+00
## [9,] 2.618602e-17 -1.000000e+00 9.074772e-16
## [10,] 1.000000e+00 -2.263780e-16 3.543410e-16
##
## $nconv
## [1] 3
##
## $niter
## [1] 1
##
## $nops
## [1] 10
n gives the dimension of the matrix, and
args contains user data that will be passed to
f.
Smallest Eigenvalues
Sometimes you may need to calculate the smallest (in magnitude)
k eigenvalues of a matrix. A direct way to do so is to use
eigs(..., which = "SM"). However, the algorithm of
RSpectra is good at finding large eigenvalues rather
than small ones, so which = "SM" tends to require much more
iterations or even may not converge.
The recommended way to retrieve the smallest eigenvalues is to utilize the spectral transformation: If \(A\) has eigenvalues \(\lambda_i\), then \((A-\sigma I)^{-1}\) has eigenvalues \(1/(\lambda_i-\sigma)\). Therefore, we can set \(\sigma = 0\) and calculate the largest eigenvalues of \(A^{-1}\), whose reciprocals are exactly the smallest eigenvalues of \(A\).
eigs() implements the spectral transformation via the
parameter sigma, so the following code returns the smallest
eigenvalues of A. Note that eigs() always
returns eigenvalues in the original scale (i.e., \(\lambda_i\) instead of \(1/(\lambda_i-\sigma)\)).
eigs_sym(A, k, which = "LM", sigma = 0)$values # recommended way## [1] 0.39278178 0.21130356 0.11411209 0.02151974 0.00473068
eigs_sym(A, k, which = "SM")$values # not recommended## [1] 0.39278178 0.21130356 0.11411209 0.02151974 0.00473068
More generally, the option which = "LM", sigma = s
selects eigenvalues that are closest to s.
Truncated SVD
Truncated SVD (or Partial SVD) is frequently used in text mining and image compression, which computes the leading singular values and singular vectors of a rectangular matrix.
RSpectra has the svds() function to
compute Truncated SVD:
set.seed(123)
m = 100
n = 20
k = 5
A = matrix(rnorm(m * n), m)
str(svds(A, k, nu = k, nv = k))## List of 5
## $ d : num [1:5] 14.1 13.8 13 11.8 11.3
## $ u : num [1:100, 1:5] 0.1616 0.0871 0.0536 0.0537 0.1271 ...
## $ v : num [1:20, 1:5] -0.186 0.116 0.27 -0.273 0.1 ...
## $ niter: num 1
## $ nops : num 45
Similar to eigs(), svds() supports sparse
matrices:
Asp = as(A, "dgCMatrix")
svds(Asp, k, nu = 0, nv = 0)## $d
## [1] 14.13079 13.76937 12.95764 11.82205 11.30081
##
## $u
## NULL
##
## $v
## NULL
##
## $niter
## [1] 1
##
## $nops
## [1] 40
svds() also has the function interface: users provide
functions A and Atrans that calcualtes \(Ax\) and \(A'x\) respectively, and
svds() will compute the Truncated SVD.
Af = function(x, args) as.numeric(args %*% x)
Atransf = function(x, args) as.numeric(crossprod(args, x))
str(svds(Af, k, Atrans = Atransf, dim = c(m, n), args = Asp))## List of 5
## $ d : num [1:5] 14.1 13.8 13 11.8 11.3
## $ u : num [1:100, 1:5] 0.1616 0.0871 0.0536 0.0537 0.1271 ...
## $ v : num [1:20, 1:5] -0.186 0.116 0.27 -0.273 0.1 ...
## $ niter: num 1
## $ nops : num 45
Performance
RSpectra is written in efficient C++ code. This page gives some
benchmarking results of svds() with other similar functions
available in R and extension packages.
### Introduction
**RSpectra** is an R interface to the
[Spectra library](https://spectralib.org/).
It is typically used to compute a few eigenvalues/vectors of an `n` by `n`
matrix, e.g., the `k` largest eigen values, which
is usually more efficient than `eigen()` if `k << n`.
Currently this package provides the function `eigs()` for eigenvalue/eigenvector
problems, and `svds()` for truncated SVD. Different matrix types in R,
including sparse matrices, are supported. Below is a list of implemented ones:
- `matrix` (defined in base R)
- `dgeMatrix` (defined in **Matrix** package, for general matrices)
- `dgCMatrix` (defined in **Matrix** package, for column oriented sparse matrices)
- `dgRMatrix` (defined in **Matrix** package, for row oriented sparse matrices)
- `dsyMatrix` (defined in **Matrix** package, for symmetric matrices)
- `dsCMatrix` (defined in **Matrix** package, for symmetric column oriented sparse matrices)
- `dsRMatrix` (defined in **Matrix** package, for symmetric row oriented sparse matrices)
- `function` (implicitly specify the matrix by providing a function that calculates matrix product `A %*% x`)
### Examples
We first generate some matrices:
```r
library(Matrix)
n = 20
k = 5
set.seed(111)
A1 = matrix(rnorm(n^2), n) ## class "matrix"
A2 = Matrix(A1) ## class "dgeMatrix"
```
General matrices have complex eigenvalues:
```r
eigs(A1, k)
eigs(A2, k, opts = list(retvec = FALSE)) ## eigenvalues only
```
**RSpectra** also works on sparse matrices:
```r
A1[sample(n^2, n^2 / 2)] = 0
A3 = as(A1, "dgCMatrix")
A4 = as(A1, "dgRMatrix")
eigs(A3, k)
eigs(A4, k)
```
Function interface is also supported:
```r
f = function(x, args)
{
as.numeric(args %*% x)
}
eigs(f, k, n = n, args = A3)
```
Symmetric matrices have real eigenvalues.
```r
A5 = crossprod(A1)
eigs_sym(A5, k)
```
To find the smallest (in absolute value) `k` eigenvalues of `A5`,
we have two approaches:
```r
eigs_sym(A5, k, which = "SM")
eigs_sym(A5, k, sigma = 0)
```
The results should be the same, but the latter method is far more
stable on large matrices.
For SVD problems, you can specify the number of singular values
(`k`), number of left singular vectors (`nu`) and number of right
singular vectors(`nv`).
```r
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)
```
Similar to `eigs()`, `svds()` supports sparse matrices:
```r
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)
```
and function interface
```r
f = function(x, args)
{
as.numeric(args %*% x)
}
g = function(x, args)
{
as.numeric(crossprod(args, x))
}
svds(f, k, Atrans = g, dim = c(m, n), args = Asp1)
```
RSpectra/man/ 0000755 0001762 0000144 00000000000 14231257703 012554 5 ustar ligges users RSpectra/man/svds.Rd 0000644 0001762 0000144 00000017231 14231257703 014026 0 ustar ligges users % Generated by roxygen2: do not edit by hand
% Please edit documentation in R/30_svds.R
\name{svds}
\alias{svds}
\alias{svds.matrix}
\alias{svds.dgeMatrix}
\alias{svds.dgCMatrix}
\alias{svds.dgRMatrix}
\alias{svds.dsyMatrix}
\alias{svds.dsCMatrix}
\alias{svds.dsRMatrix}
\alias{svds.function}
\title{Find the Largest k Singular Values/Vectors of a Matrix}
\usage{
svds(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{matrix}(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{dgeMatrix}(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{dgCMatrix}(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{dgRMatrix}(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{dsyMatrix}(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{dsCMatrix}(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{dsRMatrix}(A, k, nu = k, nv = k, opts = list(), ...)
\method{svds}{`function`}(A, k, nu = k, nv = k, opts = list(), ..., Atrans, dim, args = NULL)
}
\arguments{
\item{A}{The matrix whose truncated SVD is to be computed.}
\item{k}{Number of singular values requested.}
\item{nu}{Number of left singular vectors to be computed. This must
be between 0 and \code{k}.}
\item{nv}{Number of right singular vectors to be computed. This must
be between 0 and \code{k}.}
\item{opts}{Control parameters related to the computing
algorithm. See \strong{Details} below.}
\item{\dots}{Arguments for specialized S3 function calls, for example
\code{Atrans}, \code{dim} and \code{args}.}
\item{Atrans}{Only used when \code{A} is a function. \code{A} is a function
that calculates the matrix multiplication \eqn{Ax}{A * x}, and
\code{Atrans} is a function that calculates the transpose
multiplication \eqn{A'x}{A' * x}.}
\item{dim}{Only used when \code{A} is a function, to specify the
dimension of the implicit matrix. A vector of length two.}
\item{args}{Only used when \code{A} is a function. This argument
will be passed to the \code{A} and \code{Atrans} functions.}
}
\value{
A list with the following components:
\item{d}{A vector of the computed singular values.}
\item{u}{An \code{m} by \code{nu} matrix whose columns contain
the left singular vectors. If \code{nu == 0}, \code{NULL}
will be returned.}
\item{v}{An \code{n} by \code{nv} matrix whose columns contain
the right singular vectors. If \code{nv == 0}, \code{NULL}
will be returned.}
\item{nconv}{Number of converged singular values.}
\item{niter}{Number of iterations used.}
\item{nops}{Number of matrix-vector multiplications used.}
}
\description{
Given an \eqn{m} by \eqn{n} matrix \eqn{A},
function \code{svds()} can find its largest \eqn{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 \code{svds()} supports matrices of the following classes:
\tabular{ll}{
\code{matrix} \tab The most commonly used matrix type,
defined in the \strong{base} package.\cr
\code{dgeMatrix} \tab General matrix, equivalent to \code{matrix},
defined in the \strong{Matrix} package.\cr
\code{dgCMatrix} \tab Column oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{dgRMatrix} \tab Row oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{dsyMatrix} \tab Symmetrix matrix, defined in the \strong{Matrix}
package.\cr
\code{dsCMatrix} \tab Symmetric column oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{dsRMatrix} \tab Symmetric row oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{function} \tab Implicitly specify the matrix through two
functions that calculate
\eqn{f(x)=Ax}{f(x) = A * x} and
\eqn{g(x)=A'x}{g(x) = A' * x}. See section
\strong{Function Interface} for details.
}
Note that when \eqn{A} is symmetric and positive semi-definite,
SVD reduces to eigen decomposition, so you may consider using
\code{\link{eigs}()} instead. When \eqn{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 \eqn{\lambda} is a negative
eigenvalue of \eqn{A}, then \eqn{|\lambda|} will be the
corresponding singular value.
}
\details{
The \code{opts} argument is a list that can supply any of the
following parameters:
\describe{
\item{\code{ncv}}{Number of Lanzcos basis vectors to use. More vectors
will result in faster convergence, but with greater
memory use. \code{ncv} must be satisfy
\eqn{k < ncv \le p}{k < ncv <= p} where
\code{p = min(m, n)}.
Default is \code{min(p, max(2*k+1, 20))}.}
\item{\code{tol}}{Precision parameter. Default is 1e-10.}
\item{\code{maxitr}}{Maximum number of iterations. Default is 1000.}
\item{\code{center}}{Either a logical value (\code{TRUE}/\code{FALSE}), or a numeric
vector of length \eqn{n}. If a vector \eqn{c} is supplied, then
SVD is computed on the matrix \eqn{A - 1c'}{A - 1 * c'},
in an implicit way without actually forming this matrix.
\code{center = TRUE} has the same effect as
\code{center = colMeans(A)}. Default is \code{FALSE}.}
\item{\code{scale}}{Either a logical value (\code{TRUE}/\code{FALSE}), or a numeric
vector of length \eqn{n}. If a vector \eqn{s} is supplied, then
SVD is computed on the matrix \eqn{(A - 1c')S}{(A - 1 * c')S},
where \eqn{c} is the centering vector and \eqn{S = diag(1/s)}.
If \code{scale = TRUE}, then the vector \eqn{s} is computed as
the column norm of \eqn{A - 1c'}{A - 1 * c'}.
Default is \code{FALSE}.}
}
}
\section{Function Interface}{
The matrix \eqn{A} can be specified through two functions with
the following definitions
\preformatted{A <- function(x, args)
{
## should return A \%*\% x
}
Atrans <- function(x, args)
{
## should return t(A) \%*\% x
}}
They receive a vector \code{x} as an argument and returns a vector
of the proper dimension. These two functions should have the effect of
calculating \eqn{Ax}{A * x} and \eqn{A'x}{A' * x} respectively, and extra
arguments can be passed in through the
\code{args} parameter. In \code{svds()}, user should also provide
the dimension of the implicit matrix through the argument \code{dim}.
The function interface does not support the \code{center} and \code{scale} parameters
in \code{opts}.
}
\examples{
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)
}
\seealso{
\code{\link[base]{eigen}()}, \code{\link[base]{svd}()},
\code{\link[RSpectra]{eigs}()}.
}
\author{
Yixuan Qiu <\url{https://statr.me}>
}
\keyword{array}
RSpectra/man/eigs.Rd 0000644 0001762 0000144 00000023514 14231265526 014001 0 ustar ligges users % Generated by roxygen2: do not edit by hand
% Please edit documentation in R/00_eigs.R
\name{eigs}
\alias{eigs}
\alias{eigs.matrix}
\alias{eigs.dgeMatrix}
\alias{eigs.dsyMatrix}
\alias{eigs.dgCMatrix}
\alias{eigs.dsCMatrix}
\alias{eigs.dgRMatrix}
\alias{eigs.dsRMatrix}
\alias{eigs.function}
\alias{eigs_sym}
\alias{eigs_sym.function}
\title{Find a Specified Number of Eigenvalues/vectors of a Square Matrix}
\usage{
eigs(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{matrix}(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{dgeMatrix}(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{dsyMatrix}(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{dgCMatrix}(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{dsCMatrix}(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{dgRMatrix}(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{dsRMatrix}(A, k, which = "LM", sigma = NULL, opts = list(), ...)
\method{eigs}{`function`}(
A,
k,
which = "LM",
sigma = NULL,
opts = list(),
...,
n = NULL,
args = NULL
)
eigs_sym(A, k, which = "LM", sigma = NULL, opts = list(),
lower = TRUE, ...)
\method{eigs_sym}{`function`}(
A,
k,
which = "LM",
sigma = NULL,
opts = list(),
lower = TRUE,
...,
n = NULL,
args = NULL
)
}
\arguments{
\item{A}{The matrix whose eigenvalues/vectors are to be computed.
It can also be a function which receives a vector \eqn{x}
and calculates \eqn{Ax}{A * x}.
See section \strong{Function Interface} for details.}
\item{k}{Number of eigenvalues requested.}
\item{which}{Selection criterion. See \strong{Details} below.}
\item{sigma}{Shift parameter. See section \strong{Shift-And-Invert Mode}.}
\item{opts}{Control parameters related to the computing
algorithm. See \strong{Details} below.}
\item{\dots}{Arguments for specialized S3 function calls, for example
\code{lower}, \code{n} and \code{args}.}
\item{n}{Only used when \code{A} is a function, to specify the
dimension of the implicit matrix. See section
\strong{Function Interface} for details.}
\item{args}{Only used when \code{A} is a function. This argument
will be passed to the \code{A} function when it is called.
See section \strong{Function Interface} for details.}
\item{lower}{For symmetric matrices, should the lower triangle
or upper triangle be used.}
}
\value{
A list of converged eigenvalues and eigenvectors.
\item{values}{Computed eigenvalues.}
\item{vectors}{Computed eigenvectors. \code{vectors[, j]} corresponds to \code{values[j]}.}
\item{nconv}{Number of converged eigenvalues.}
\item{niter}{Number of iterations used in the computation.}
\item{nops}{Number of matrix operations used in the computation.}
}
\description{
Given an \eqn{n} by \eqn{n} matrix \eqn{A},
function \code{eigs()} can calculate a specified
number of eigenvalues and eigenvectors of \eqn{A}.
Users can specify the selection criterion by argument
\code{which}, e.g., choosing the \eqn{k} largest or smallest
eigenvalues and the corresponding eigenvectors.
Currently \code{eigs()} supports matrices of the following classes:
\tabular{ll}{
\code{matrix} \tab The most commonly used matrix type,
defined in the \strong{base} package.\cr
\code{dgeMatrix} \tab General matrix, equivalent to \code{matrix},
defined in the \strong{Matrix} package.\cr
\code{dgCMatrix} \tab Column oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{dgRMatrix} \tab Row oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{dsyMatrix} \tab Symmetric matrix, defined in the \strong{Matrix}
package.\cr
\code{dsCMatrix} \tab Symmetric column oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{dsRMatrix} \tab Symmetric row oriented sparse matrix, defined in
the \strong{Matrix} package.\cr
\code{function} \tab Implicitly specify the matrix through a
function that has the effect of calculating
\eqn{f(x)=Ax}{f(x) = A * x}. See section
\strong{Function Interface} for details.
}
\code{eigs_sym()} assumes the matrix is symmetric,
and only the lower triangle (or upper triangle, which is
controlled by the argument \code{lower}) is used for
computation, which guarantees that the eigenvalues and eigenvectors are
real, and in general results in faster and more stable computation.
One exception is when \code{A} is a function, in which case the user is
responsible for the symmetry of the operator.
\code{eigs_sym()} supports "matrix", "dgeMatrix", "dgCMatrix", "dgRMatrix"
and "function" typed matrices.
}
\details{
The \code{which} argument is a character string
that specifies the type of eigenvalues to be computed.
Possible values are:
\tabular{ll}{
"LM" \tab The \eqn{k} eigenvalues with largest magnitude. Here the
magnitude means the Euclidean norm of complex numbers.\cr
"SM" \tab The \eqn{k} eigenvalues with smallest magnitude.\cr
"LR" \tab The \eqn{k} eigenvalues with largest real part.\cr
"SR" \tab The \eqn{k} eigenvalues with smallest real part.\cr
"LI" \tab The \eqn{k} eigenvalues with largest imaginary part.\cr
"SI" \tab The \eqn{k} eigenvalues with smallest imaginary part.\cr
"LA" \tab The \eqn{k} largest (algebraic) eigenvalues, considering any
negative sign.\cr
"SA" \tab The \eqn{k} smallest (algebraic) eigenvalues, considering any
negative sign.\cr
"BE" \tab Compute \eqn{k} eigenvalues, half from each end of the
spectrum. When \eqn{k} is odd, compute more from the high
and then from the low end.
}
\code{eigs()} with matrix types "matrix", "dgeMatrix", "dgCMatrix"
and "dgRMatrix" can use "LM", "SM", "LR", "SR", "LI" and "SI".
\code{eigs_sym()} with all supported matrix types,
and \code{eigs()} with symmetric matrix types
("dsyMatrix", "dsCMatrix", and "dsRMatrix") can use "LM", "SM", "LA", "SA" and "BE".
The \code{opts} argument is a list that can supply any of the
following parameters:
\describe{
\item{\code{ncv}}{Number of Lanzcos basis vectors to use. More vectors
will result in faster convergence, but with greater
memory use. For general matrix, \code{ncv} must satisfy
\eqn{k+2\le ncv \le n}{k+2 <= ncv <= n}, and
for symmetric matrix, the constraint is
\eqn{k < ncv \le n}{k < ncv <= n}.
Default is \code{min(n, max(2*k+1, 20))}.}
\item{\code{tol}}{Precision parameter. Default is 1e-10.}
\item{\code{maxitr}}{Maximum number of iterations. Default is 1000.}
\item{\code{retvec}}{Whether to compute eigenvectors. If FALSE,
only calculate and return eigenvalues.}
\item{\code{initvec}}{Initial vector of length \eqn{n} supplied to the
Arnoldi/Lanczos iteration. It may speed up the convergence
if \code{initvec} is close to an eigenvector of \eqn{A}.}
}
}
\section{Shift-And-Invert Mode}{
The \code{sigma} argument is used in the shift-and-invert mode.
When \code{sigma} is not \code{NULL}, the selection criteria specified
by argument \code{which} will apply to
\deqn{\frac{1}{\lambda-\sigma}}{1/(\lambda-\sigma)}
where \eqn{\lambda}'s are the eigenvalues of \eqn{A}. This mode is useful
when user wants to find eigenvalues closest to a given number.
For example, if \eqn{\sigma=0}, then \code{which = "LM"} will select the
largest values of \eqn{1/|\lambda|}, which turns out to select
eigenvalues of \eqn{A} that have the smallest magnitude. The result of
using \code{which = "LM", sigma = 0} will be the same as
\code{which = "SM"}, but the former one is preferable
in that \code{eigs()} is good at finding large
eigenvalues rather than small ones. More explanation of the
shift-and-invert mode can be found in the SciPy document,
\url{https://docs.scipy.org/doc/scipy/tutorial/arpack.html}.
}
\section{Function Interface}{
The matrix \eqn{A} can be specified through a function with
the definition
\preformatted{function(x, args)
{
## should return A \%*\% x
}}
which receives a vector \code{x} as an argument and returns a vector
of the same length. The function should have the effect of calculating
\eqn{Ax}{A * x}, and extra arguments can be passed in through the
\code{args} parameter. In \code{eigs()}, user should also provide
the dimension of the implicit matrix through the argument \code{n}.
}
\examples{
library(Matrix)
n = 20
k = 5
## general matrices have complex eigenvalues
set.seed(111)
A1 = matrix(rnorm(n^2), n) ## class "matrix"
A2 = Matrix(A1) ## class "dgeMatrix"
eigs(A1, k)
eigs(A2, k, opts = list(retvec = FALSE)) ## eigenvalues only
## Sparse matrices
A1[sample(n^2, n^2 / 2)] = 0
A3 = as(A1, "dgCMatrix")
A4 = as(A1, "dgRMatrix")
eigs(A3, k)
eigs(A4, k)
## Function interface
f = function(x, args)
{
as.numeric(args \%*\% x)
}
eigs(f, k, n = n, args = A3)
## Symmetric matrices have real eigenvalues
A5 = crossprod(A1)
eigs_sym(A5, k)
## Find the smallest (in absolute value) k eigenvalues of A5
eigs_sym(A5, k, which = "SM")
## Another way to do this: use the sigma argument
eigs_sym(A5, k, sigma = 0)
## The results should be the same,
## but the latter method is far more stable on large matrices
}
\seealso{
\code{\link[base]{eigen}()}, \code{\link[base]{svd}()},
\code{\link[RSpectra]{svds}()}
}
\author{
Yixuan Qiu \url{https://statr.me}
Jiali Mei \email{vermouthmjl@gmail.com}
}
\keyword{array}
RSpectra/DESCRIPTION 0000644 0001762 0000144 00000003730 14231324752 013511 0 ustar ligges users Package: RSpectra
Type: Package
Title: Solvers for Large-Scale Eigenvalue and SVD Problems
Version: 0.16-1
Date: 2022-04-24
Authors@R: c(
person("Yixuan", "Qiu", , "yixuan.qiu@cos.name", c("aut", "cre")),
person("Jiali", "Mei", , "vermouthmjl@gmail.com", "aut",
comment = "Function interface of matrix operation"),
person("Gael", "Guennebaud", , "gael.guennebaud@inria.fr", "ctb",
comment = "Eigenvalue solvers from the 'Eigen' library"),
person("Jitse", "Niesen", , "jitse@maths.leeds.ac.uk", "ctb",
comment = "Eigenvalue solvers from the 'Eigen' library")
)
Description: R interface to the 'Spectra' library