Perform multi-dimensional scaling (MDS) also know as Principal Coordinate Analysis (PCoA). These functions are wrappers for scater::calculateMDS.

getMDS(x, ...)

addMDS(x, ...)

# S4 method for class 'SingleCellExperiment'
addMDS(x, name = "MDS", ...)

# S4 method for class 'SingleCellExperiment'
getMDS(x, assay.type = "counts", ...)

# S4 method for class 'TreeSummarizedExperiment'
getMDS(x, assay.type = "counts", ...)

Arguments

x

a SummarizedExperiment object.

...

additional arguments.

  • FUN: Function. A function that is applied to calculate dissimilarity. (Default: getDissimilarity)

    \item \code{subset.result}: \code{Logical result}. Specifies whether to
    subset \code{x} to match the result if some samples were removed during
    calculation. (Default: \code{TRUE})
    
    \item \code{ntop}: \code{Numeric scalar}. Forwarded to
    \code{scater::calculateMDS}. Controls how many of the most variable
    features are kept before MDS is run. (Default: \code{500}). Set
    \code{ntop = Inf} (or supply \code{subset.row}) when you need the full
    assay for operations such as rarefaction where absolute library sizes
    matter.

name

Character scalar. A name for the reducedDim() where results will be stored. (Default: "MDS")

assay.type

Character scalar. Specifies the name of assay used in calculation. (Default: "counts")

Value

getMDS returns a MDS results. addMDS returns a x with MDS results added to its reducedDim(x, name).

Details

These functions are wrappers for scater::calculateMDS and scater::runMDS. While getMDS returns the results, addMDS adds them to reducedDim(x). The difference is that these functions apply microbiome-specific options such as getDissimilarity function by default.

By default scater::calculateMDS keeps only the ntop = 500 most variable features. When rarefaction is enabled (niter > 0) the sampling depth is evaluated on this filtered matrix, so samples whose depth falls below the requested sample are dropped. Pass ntop = Inf (or specify subset.row) if you need to retain all taxa so that the minimum depth of assay(x, assay.type) remains reachable.

See scater::calculateMDS for details.

Examples

library(mia)
library(scater)
#> Loading required package: scuttle
#> Loading required package: ggplot2
library(patchwork)

data(GlobalPatterns)
tse <- GlobalPatterns

# Calculate PCoA with Bray-Curtis dissimilarity
tse <- transformAssay(tse, method = "relabundance")
tse <- addMDS(tse, assay.type = "relabundance", method = "bray")

# Calculate PCoA with Unifrac and rarefaction. (Note: increase iterations)
tse <- addMDS(tse, method = "unifrac", name = "unifrac")
#> Warning: Pruning tree...

# Calculate PCoA with Unifrac and rarefaction. (Note: increase iterations)
tse <- addMDS(tse, method = "unifrac", name = "unifrac_rare", niter = 2L)
#> Warning: Pruning tree...
#> Warning: Pruning tree...

# Disable feature filtering when rarefying at the minimum depth
min_depth <- min(colSums(assay(tse, "counts")))
tse <- addMDS(
    tse,
    assay.type = "counts",
    method = "bray",
    niter = 2L,
    sample = min_depth,
    ntop = Inf
)
reducedDim(tse, "MDS")[1:3, 1:2]
#>            [,1]        [,2]
#> CL3 -0.10001156  0.02224679
#> CC1 -0.12328430  0.02173183
#> SV1 -0.08561559 -0.01413502

# Visualize results
p1 <- plotReducedDim(tse, "unifrac", colour_by = "SampleType") +
    labs(title = "Not rarefied")
p2 <- plotReducedDim(tse, "unifrac_rare", colour_by = "SampleType") +
    labs(title = "Rarefied")
p1 + p2