rcur {rsvd}R Documentation

Randomized CUR matrix decomposition.

Description

Randomized CUR matrix decomposition.

Usage

rcur(A, k = NULL, p = 10, q = 0, idx_only = FALSE, rand = TRUE)

Arguments

A

array_like;
numeric (m,n)(m, n) input matrix (or data frame).
If the data contain NANAs na.omit is applied.

k

integer;
target rank of the low-rank approximation, i.e., the number of columns/rows to be selected. It is required that kk is smaller or equal to min(m,n)min(m,n).

p

integer, optional;
oversampling parameter (default p=10p=10).

q

integer, optional;
number of additional power iterations (default q=0q=0).

idx_only

bool, optional;
if (TRUETRUE), only the index set C.idx and R.idx is returned, but not the matrices C and R. This is more memory efficient, when dealing with large-scale data.

rand

bool, optional;
if (TRUETRUE), a probabilistic strategy is used, otherwise a deterministic algorithm is used.

Details

Algorithm for computing the CUR matrix decomposition of a rectangular (m,n)(m, n) matrix AA, with target rank k<<min(m,n)k << min(m,n). The input matrix is factored as

A=CURA = C * U * R

using the rid decomposition. The factor matrix CC is formed using actual columns of AA, also called the partial column skeleton. The factor matrix RR is formed using actual rows of AA, also called the partial row skeleton.

If rand=TRUErand=TRUE a probabilistic strategy is used to compute the decomposition, otherwise a deterministic algorithm is used.

Value

rcur returns a list with class idid containing the following components:

C

array_like;
column subset C=A[,C.idx]C = A[,C.idx]; (m,k)(m, k) dimensional array.

R

Array_like.
row subset R=A[R.idx,]R = A[R.idx, ]; (k,n)(k, n) dimensional array.

U

array_like;
connector matrix; (k,k)(k,k) dimensional array.

C.idx

array_like;
index set of the kk selected columns used to form CC.

R.idx

array_like;
index set of the kk selected rows used to form RR.

C.scores

array_like;
scores of the selected columns.

R.scores

array_like;
scores of the selected rows.

Author(s)

N. Benjamin Erichson, erichson@berkeley.edu

References

See Also

rid

Examples

## Not run: 
# Load test image
data('tiger')

# Compute (column) randomized interpolative decompsition
# Note that the image needs to be transposed for correct plotting
out <- rcur(tiger, k = 150)

# Reconstruct image
tiger.re <- out$C %*% out$U %*% out$R

# Compute relative error
print(norm(tiger-tiger.re, 'F') / norm(tiger, 'F')) 

# Plot approximated image
image(tiger.re, col = gray((0:255)/255))

## End(Not run)

[Package rsvd version 1.0.5 Index]