runNMDS | R Documentation |
Perform non-metric multi-dimensional scaling (nMDS) on samples, based on the
data in a SingleCellExperiment
object.
calculateNMDS(x, ...)
## S4 method for signature 'ANY'
calculateNMDS(
x,
FUN = vegdist,
nmdsFUN = c("isoMDS", "monoMDS"),
ncomponents = 2,
ntop = 500,
subset_row = NULL,
scale = FALSE,
transposed = FALSE,
keep_dist = FALSE,
...
)
## S4 method for signature 'SummarizedExperiment'
calculateNMDS(
x,
...,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
FUN = vegdist
)
## S4 method for signature 'SingleCellExperiment'
calculateNMDS(
x,
...,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
dimred = NULL,
n_dimred = NULL,
FUN = vegdist
)
runNMDS(x, ..., altexp = NULL, name = "NMDS")
x |
For For |
... |
additional arguments to pass to |
FUN |
a |
nmdsFUN |
a |
ncomponents |
Numeric scalar indicating the number of NMDS dimensions to obtain. |
ntop |
Numeric scalar specifying the number of features with the highest variances to use for dimensionality reduction. |
subset_row |
Vector specifying the subset of features to use for dimensionality reduction. This can be a character vector of row names, an integer vector of row indices or a logical vector. |
scale |
Logical scalar, should the expression values be standardized? |
transposed |
Logical scalar, is x transposed with cells in rows? |
keep_dist |
Logical scalar indicating whether the |
assay.type |
a single |
assay_name |
a single |
exprs_values |
a single |
dimred |
String or integer scalar specifying the existing dimensionality reduction results to use. |
n_dimred |
Integer scalar or vector specifying the dimensions to use if dimred is specified. |
altexp |
String or integer scalar specifying an alternative experiment containing the input data. |
name |
String specifying the name to be used to store the result in the reducedDims of the output. |
Either MASS::isoMDS
or
vegan::monoMDS
are used internally to compute
the NMDS components. If you supply a custom FUN
, make sure that
the arguments of FUN
and nmdsFUN
do not collide.
For calculateNMDS
, a matrix is returned containing the MDS
coordinates for each sample (row) and dimension (column).
Felix Ernst
MASS::isoMDS
,
vegan::monoMDS
for NMDS component calculation.
plotMDS
, to quickly visualize the
results.
# generate some example data
mat <- matrix(1:60, nrow = 6)
df <- DataFrame(n = c(1:6))
tse <- TreeSummarizedExperiment(assays = list(counts = mat),
rowData = df)
#
calculateNMDS(tse)
#
data(esophagus)
esophagus <- runNMDS(esophagus, FUN = vegan::vegdist, name = "BC")
esophagus <- runNMDS(esophagus, FUN = vegan::vegdist, name = "euclidean",
method = "euclidean")
reducedDims(esophagus)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.