runUMAP: Perform UMAP on cell-level data

View source: R/runUMAP.R

calculateUMAPR Documentation

Perform UMAP on cell-level data

Description

Perform uniform manifold approximation and projection (UMAP) for the cells, based on the data in a SingleCellExperiment object.

Usage

calculateUMAP(x, ...)

## S4 method for signature 'ANY'
calculateUMAP(
  x,
  ncomponents = 2,
  ntop = 500,
  subset_row = NULL,
  scale = FALSE,
  transposed = FALSE,
  pca = if (transposed) NULL else 50,
  n_neighbors = 15,
  n_threads = NULL,
  ...,
  external_neighbors = FALSE,
  BNPARAM = KmknnParam(),
  BPPARAM = SerialParam(),
  use_densvis = FALSE,
  dens_frac = 0.3,
  dens_lambda = 0.1
)

## S4 method for signature 'SummarizedExperiment'
calculateUMAP(x, ..., exprs_values = "logcounts", assay.type = exprs_values)

## S4 method for signature 'SingleCellExperiment'
calculateUMAP(
  x,
  ...,
  pca = if (!is.null(dimred)) NULL else 50,
  exprs_values = "logcounts",
  dimred = NULL,
  n_dimred = NULL,
  assay.type = exprs_values
)

runUMAP(x, ..., altexp = NULL, name = "UMAP")

Arguments

x

For calculateUMAP, a numeric matrix of log-expression values where rows are features and columns are cells. Alternatively, a SummarizedExperiment or SingleCellExperiment containing such a matrix.

For runTSNE, a SingleCellExperiment object containing such a matrix.

...

For the calculateUMAP generic, additional arguments to pass to specific methods. For the ANY method, additional arguments to pass to umap. For the SummarizedExperiment and SingleCellExperiment methods, additional arguments to pass to the ANY method.

For runUMAP, additional arguments to pass to calculateUMAP.

ncomponents

Numeric scalar indicating the number of UMAP 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?

pca

Integer scalar specifying how many PCs should be used as input into the UMAP algorithm. By default, no PCA is performed if the input is a dimensionality reduction result.

n_neighbors

Integer scalar, number of nearest neighbors to identify when constructing the initial graph.

n_threads

Integer scalar specifying the number of threads to use in umap. If NULL and BPPARAM is a MulticoreParam, it is set to the number of workers in BPPARAM; otherwise, the umap defaults are used.

external_neighbors

Logical scalar indicating whether a nearest neighbors search should be computed externally with findKNN.

BNPARAM

A BiocNeighborParam object specifying the neighbor search algorithm to use when external_neighbors=TRUE.

BPPARAM

A BiocParallelParam object specifying whether the PCA should be parallelized.

use_densvis

Logical scalar indicating whether densne should be used to perform density-preserving t-SNE.

dens_frac, dens_lambda

See densne

exprs_values

Alias to assay.type.

assay.type

Integer scalar or string indicating which assay of x contains the expression values.

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.

Details

The function umap is used internally to compute the UMAP. Note that the algorithm is not deterministic, so different runs of the function will produce differing results. Users are advised to test multiple random seeds, and then use set.seed to set a random seed for replicable results.

If external_neighbors=TRUE, the nearest neighbor search is conducted using a different algorithm to that in the umap function. This can be parallelized or approximate to achieve greater speed for large data sets. The neighbor search results are then used directly to create the UMAP embedding.

Value

For calculateUMAP, a matrix is returned containing the UMAP coordinates for each cell (row) and dimension (column).

For runUMAP, a modified x is returned that contains the UMAP coordinates in reducedDim(x, name).

Feature selection

This section is relevant if x is a numeric matrix of (log-)expression values with features in rows and cells in columns; or if x is a SingleCellExperiment and dimred=NULL. In the latter, the expression values are obtained from the assay specified by assay.type.

The subset_row argument specifies the features to use for dimensionality reduction. The aim is to allow users to specify highly variable features to improve the signal/noise ratio, or to specify genes in a pathway of interest to focus on particular aspects of heterogeneity.

If subset_row=NULL, the ntop features with the largest variances are used instead. We literally compute the variances from the expression values without considering any mean-variance trend, so often a more considered choice of genes is possible, e.g., with scran functions. Note that the value of ntop is ignored if subset_row is specified.

If scale=TRUE, the expression values for each feature are standardized so that their variance is unity. This will also remove features with standard deviations below 1e-8.

Using reduced dimensions

If x is a SingleCellExperiment, the method can be applied on existing dimensionality reduction results in x by setting the dimred argument. This is typically used to run slower non-linear algorithms (t-SNE, UMAP) on the results of fast linear decompositions (PCA). We might also use this with existing reduced dimensions computed from a priori knowledge (e.g., gene set scores), where further dimensionality reduction could be applied to compress the data.

The matrix of existing reduced dimensions is taken from reducedDim(x, dimred). By default, all dimensions are used to compute the second set of reduced dimensions. If n_dimred is also specified, only the first n_dimred columns are used. Alternatively, n_dimred can be an integer vector specifying the column indices of the dimensions to use.

When dimred is specified, no additional feature selection or standardization is performed. This means that any settings of ntop, subset_row and scale are ignored.

If x is a numeric matrix, setting transposed=TRUE will treat the rows as cells and the columns as the variables/diemnsions. This allows users to manually pass in dimensionality reduction results without needing to wrap them in a SingleCellExperiment. As such, no feature selection or standardization is performed, i.e., ntop, subset_row and scale are ignored.

Using alternative Experiments

This section is relevant if x is a SingleCellExperiment and altexp is not NULL. In such cases, the method is run on data from an alternative SummarizedExperiment nested within x. This is useful for performing dimensionality reduction on other features stored in altExp(x, altexp), e.g., antibody tags.

Setting altexp with assay.type will use the specified assay from the alternative SummarizedExperiment. If the alternative is a SingleCellExperiment, setting dimred will use the specified dimensionality reduction results from the alternative. This option will also interact as expected with n_dimred.

Note that the output is still stored in the reducedDims of the output SingleCellExperiment. It is advisable to use a different name to distinguish this output from the results generated from the main experiment's assay values.

Author(s)

Aaron Lun

References

McInnes L, Healy J, Melville J (2018). UMAP: uniform manifold approximation and projection for dimension reduction. arXiv.

See Also

umap, for the underlying calculations.

plotUMAP, to quickly visualize the results.

Examples

example_sce <- mockSCE() 
example_sce <- logNormCounts(example_sce)

example_sce <- runUMAP(example_sce)
reducedDimNames(example_sce)
head(reducedDim(example_sce))

LTLA/scater documentation built on July 21, 2024, 5:43 p.m.