nmf-class-methods: nmf class methods

subset,nmf-methodR Documentation

nmf class methods

Description

Given an NMF model in the form A = wdh, project projects w onto A to solve for h.

Usage

## S4 method for signature 'nmf'
subset(x, i, ...)

## S4 method for signature 'nmf,ANY,ANY,ANY'
x[i]

## S4 method for signature 'nmf'
head(x, n = getOption("digits"), ...)

## S4 method for signature 'nmf'
show(object)

## S4 method for signature 'nmf'
dimnames(x)

## S4 method for signature 'nmf'
dim(x)

## S4 method for signature 'nmf'
t(x)

## S4 method for signature 'nmf'
sort(x, decreasing = TRUE, ...)

## S4 method for signature 'nmf'
prod(x, ..., na.rm = FALSE)

## S4 method for signature 'nmf'
x$name

## S4 method for signature 'nmf,list'
coerce(from, to)

## S4 method for signature 'nmf'
x[[i]]

## S4 method for signature 'nmf'
predict(object, data, L1 = 0, L2 = 0, mask = NULL, upper_bound = NULL, ...)

Arguments

x

object of class nmf.

i

indices

...

arguments passed to or from other methods

n

number of rows/columns to show

object

fitted model, class nmf, generally the result of calling nmf, with models of equal dimensions as data

decreasing

logical. Should the sort be increasing or decreasing?

na.rm

remove na values

name

name of nmf class slot

from

class which the coerce method should perform coercion from

to

class which the coerce method should perform coercion to

data

dense or sparse matrix of features in rows and samples in columns. Prefer matrix or Matrix::dgCMatrix, respectively

L1

a single LASSO penalty in the range (0, 1]

L2

a single Ridge penalty greater than zero

mask

dense or sparse matrix of values in data to handle as missing. Prefer Matrix::dgCMatrix. Alternatively, specify "zeros" or "NA" to mask either all zeros or NA values.

upper_bound

maximum value permitted in least squares solutions, essentially a bounded-variable least squares problem between 0 and upper_bound

Details

For the alternating least squares matrix factorization update problem A = wh, the updates (or projection) of h is given by the equation:

w^Twh = wA_j

which is in the form ax = b where a = w^Tw x = h and b = wA_j for all columns j in A.

Any L1 penalty is subtracted from b and should generally be scaled to max(b), where b = WA_j for all columns j in A. An easy way to properly scale an L1 penalty is to normalize all columns in w to sum to the same value (e.g. 1). No scaling is applied in this function. Such scaling guarantees that L1 = 1 gives a completely sparse solution.

There are specializations for dense and sparse input matrices, symmetric input matrices, and for rank-1 and rank-2 projections. See documentation for nmf for theoretical details and guidance.

Value

object of class nmf

Author(s)

Zach DeBruine

References

DeBruine, ZJ, Melcher, K, and Triche, TJ. (2021). "High-performance non-negative matrix factorization for large single-cell data." BioRXiv.

See Also

nnls, nmf

Examples

## Not run: 
w <- matrix(runif(1000 * 10), 1000, 10)
h_true <- matrix(runif(10 * 100), 10, 100)
# A is the crossproduct of "w" and "h" with 10% signal dropout
A <- (w %*% h_true) * (r_sparsematrix(1000, 100, 10) > 0)
h <- project(w, A)
cor(as.vector(h_true), as.vector(h))

# alternating projections refine solution (like NMF)
mse_bad <- mse(w, rep(1, ncol(w)), h, A) # mse before alternating updates
h <- project(w, A)
w <- t(project(h, t(A)))
h <- project(w, A)
w <- t(project(h, t(A)))
h <- project(w, A)
w <- t(project(h, t(A)))
mse_better <- mse(w, rep(1, ncol(w)), h, A) # mse after alternating updates
mse_better < mse_bad

## End(Not run)

zdebruine/RcppML documentation built on Sept. 13, 2023, 11:44 p.m.