These functions perform Non-negative Matrix Factorization on data stored in a
TreeSummarizedExperiment
object.
getNMF(x, ...)
addNMF(x, ...)
# S4 method for class 'SummarizedExperiment'
getNMF(x, k = 2, assay.type = "counts", eval.metric = "evar", ...)
# S4 method for class 'SummarizedExperiment'
addNMF(
x,
k = 2,
assay.type = "counts",
eval.metric = "evar",
name = "NMF",
...
)
a
TreeSummarizedExperiment
object.
optional arguments passed to nmf::NMF
.
numeric vector
. A number of latent vectors/topics.
(Default: 2
)
Character scalar
. Specifies which assay to use for
NMF ordination. (Default: "counts"
)
Character scalar
. Specifies the evaluation metric
that will be used to select the model with the best fit. Must be one of the
following options: "evar"
(explained variance; maximized),
"sparseness.basis"
(degree of sparsity in the basis matrix;
maximized), "sparseness.coef"
(degree of sparsity in the coefficient
matrix; maximized), "rss"
(residual sum of squares; minimized),
"silhouette.coef"
(quality of clustering based on the coefficient
matrix; maximized), "silhouette.basis"
(quality of clustering based
on the basis matrix; maximized), "cophenetic"
(correlation between
cophenetic distances and original distances; maximized), "dispersion"
(spread of data points within clusters; minimized). (Default: "evar"
)
Character scalar
. The name to be used to store the result
in the reducedDims of the output. (Default: "NMF"
)
For getNMF
, the ordination matrix with feature loadings matrix
as attribute "loadings"
.
For addNMF
, a
TreeSummarizedExperiment
object is returned containing the ordination matrix in
reducedDims(x, name)
with the following attributes:
"loadings" which is a matrix containing the feature loadings
"NMF_output" which is the output of function nmf::NMF
"best_fit" which is the result of the best fit if k is a vector of integers
The functions getNMF
and addNMF
internally use nmf::NMF
compute the ordination matrix and
feature loadings.
If k is a vector of integers, NMF output is calculated for all the rank
values contained in k, and the best fit is selected based on
eval.metric
value.
data(GlobalPatterns)
tse <- GlobalPatterns
# Reduce the number of features
tse <- agglomerateByPrevalence(tse, rank = "Phylum")
# Run NMF and add the result to reducedDim(tse, "NMF").
tse <- addNMF(tse, k = 2, name = "NMF")
#> Loading required package: registry
#> Loading required package: rngtools
#> Loading required package: cluster
#> NMF - BioConductor layer [OK] | Shared memory capabilities [NO: bigmemory] | Cores 2/2
#> To enable shared memory capabilities, try: install.extras('
#> NMF
#> ')
#>
#> Attaching package: ‘NMF’
#> The following object is masked from ‘package:S4Vectors’:
#>
#> nrun
# Extract feature loadings
loadings_NMF <- attr(reducedDim(tse, "NMF"), "loadings")
head(loadings_NMF)
#> [,1] [,2]
#> AD3 6.978118e-04 1.568333e-07
#> Acidobacteria 2.589159e-02 1.528508e-05
#> Actinobacteria 7.748740e-02 2.876304e-02
#> Armatimonadetes 3.759120e-04 2.220446e-16
#> BRC1 5.044815e-05 8.456265e-08
#> Bacteroidetes 6.852413e-02 2.170128e-01
# Estimate models with number of topics from 2 to 4. Perform 2 runs.
tse <- addNMF(tse, k = c(2, 3, 4), name = "NMF_4", nrun = 2)
#> NMF - BioConductor layer [OK] | Shared memory capabilities [NO: bigmemory] | Cores 2/2
#> To enable shared memory capabilities, try: install.extras('
#> NMF
#> ')
#>
#> Attaching package: ‘NMF’
#> The following object is masked from ‘package:S4Vectors’:
#>
#> nrun
# Extract feature loadings
loadings_NMF_4 <- attr(reducedDim(tse, "NMF_4"), "loadings")
head(loadings_NMF_4)
#> [,1] [,2] [,3] [,4]
#> AD3 2.220446e-16 1.791067e-07 7.561821e-04 2.220446e-16
#> Acidobacteria 9.403068e-06 2.220446e-16 2.589716e-02 2.504186e-03
#> Actinobacteria 1.516838e-02 2.220446e-16 8.438975e-02 1.275557e-02
#> Armatimonadetes 2.220446e-16 2.220446e-16 2.621654e-04 1.678860e-04
#> BRC1 2.220446e-16 6.012305e-08 4.458859e-05 1.168002e-05
#> Bacteroidetes 1.603669e-02 2.525113e-01 1.893582e-02 1.552504e-02