Perform non-metric multi-dimensional scaling (nMDS) on samples, based on the data in a SingleCellExperiment object.

getNMDS(x, ...)

# S4 method for class 'ANY'
getNMDS(
  x,
  FUN = vegdist,
  nmds.fun = nmdsFUN,
  nmdsFUN = c("isoMDS", "monoMDS"),
  ncomponents = 2,
  ntop = 500,
  subset.row = subset_row,
  subset_row = NULL,
  scale = FALSE,
  transposed = FALSE,
  keep.dist = keep_dist,
  keep_dist = FALSE,
  ...
)

# S4 method for class 'SummarizedExperiment'
getNMDS(
  x,
  ...,
  assay.type = assay_name,
  assay_name = exprs_values,
  exprs_values = "counts",
  FUN = vegdist
)

# S4 method for class 'SingleCellExperiment'
getNMDS(
  x,
  ...,
  assay.type = assay_name,
  assay_name = exprs_values,
  exprs_values = "counts",
  dimred = NULL,
  ndimred = n_dimred,
  n_dimred = NULL,
  FUN = vegdist
)

calculateNMDS(x, ...)

addNMDS(x, ..., altexp = NULL, name = "NMDS")

runNMDS(x, ...)

Arguments

x

TreeSummarizedExperiment.

...

additional arguments to pass to FUN and nmds.fun.

FUN

Function or Character scalar. A value with a function name returning a dist object

nmds.fun

Character scalar. A value to choose the scaling implementation, either “isoMDS” for MASS::isoMDS or “monoMDS” for vegan::monoMDS

nmdsFUN

Deprecated. Use nmds.fun instead.

ncomponents

Numeric scalar. Indicates the number of DPCoA dimensions to obtain. (Default: 2)

ntop

Numeric scalar. Specifies the number of features with the highest variances to use for dimensionality reduction. Alternatively NULL, if all features should be used. (Default: NULL)

subset.row

Character Vector. Specifies 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. (Default: NULL)

subset_row

Deprecated. Use subset.row instead.

scale

Logical scalar. Should the expression values be standardized? (Default: FALSE)

transposed

Logical scalar. Specifies if x is transposed with cells in rows. (Default: FALSE)

keep.dist

Logical scalar. Indicates whether the dist object calculated by FUN should be stored as ‘dist’ attribute of the matrix returned/stored by getNMDS/ addNMDS. (Default: FALSE)

keep_dist

Deprecated. Use keep.dist instead.

assay.type

Character scalar. Specifies which assay to use for calculation. (Default: "counts")

assay_name

Deprecated. Use assay.type instead.

exprs_values

Deprecated. Use assay.type instead.

dimred

Character scalar or integer scalar. Specifies the existing dimensionality reduction results to use.

ndimred

integer vector. Specifies the dimensions to use if dimred is specified.

n_dimred

Deprecated. Use ndimred instead.

altexp

Character scalar or integer scalar. Specifies an alternative experiment containing the input data. (Default: NULL)

name

Character scalar. A name for the column of the colData where results will be stored. (Default: "NMDS")

Value

For getNMDS, a matrix is returned containing the MDS coordinates for each sample (row) and dimension (column).

Details

For addNMDS a SingleCellExperiment

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 nmds.fun do not collide.

See also

MASS::isoMDS, vegan::monoMDS for NMDS component calculation.

plotMDS, to quickly visualize the results.

Examples

# generate some example data
mat <- matrix(1:60, nrow = 6)
df <- DataFrame(n = c(1:6))
tse <- TreeSummarizedExperiment(assays = list(counts = mat),
                                rowData = df)
#
getNMDS(tse)
#> initial  value 0.383462 
#> iter   5 value 0.161655
#> iter  10 value 0.113278
#> final  value 0.003270 
#> converged
#>              [,1]         [,2]
#>  [1,] -0.60647893  0.282119175
#>  [2,] -0.38097446 -0.130374816
#>  [3,] -0.17579928 -0.184259003
#>  [4,] -0.02959516 -0.158516980
#>  [5,]  0.09082607 -0.111879947
#>  [6,]  0.14427170 -0.047910688
#>  [7,]  0.20118089 -0.001837781
#>  [8,]  0.22817639  0.061022767
#>  [9,]  0.25717122  0.116809029
#> [10,]  0.27122156  0.174828245
#> attr(,"Stress")
#> [1] 0.003269857

#
data(esophagus)
esophagus <- addNMDS(esophagus, FUN = vegan::vegdist, name = "BC")
#> initial  value 0.000000 
#> final  value 0.000000 
#> converged
esophagus <- addNMDS(esophagus, FUN = vegan::vegdist, name = "euclidean",
                     method = "euclidean")
#> initial  value 0.000000 
#> final  value 0.000000 
#> converged
reducedDims(esophagus)
#> List of length 2
#> names(2): BC euclidean