escape.matrix: Calculate Single-Cell Gene-Set Enrichment Scores

View source: R/runEscape.R

escape.matrixR Documentation

Calculate Single-Cell Gene-Set Enrichment Scores

Description

'escape.matrix()' computes per-cell enrichment for arbitrary gene-set collections using one of four scoring back-ends and returns a dense numeric matrix (cells × gene-sets). The expression matrix is processed in user-defined *chunks* ('groups') so that memory use remains predictable; each chunk is dispatched in parallel via a BiocParallel 'BPPARAM' backend. Heavy engines (GSVA, UCell, AUCell) are loaded lazily, keeping them in the package’s Suggests field.

Usage

escape.matrix(
  input.data,
  gene.sets = NULL,
  method = "ssGSEA",
  groups = 1000,
  min.size = 5,
  normalize = FALSE,
  make.positive = FALSE,
  min.expr.cells = 0,
  min.filter.by = NULL,
  BPPARAM = NULL,
  ...
)

Arguments

input.data

A raw‐counts matrix ('genes × cells'), a Seurat object, or a SingleCellExperiment. Gene identifiers must match those in 'gene.sets'.

gene.sets

A named list of character vectors, the result of [getGeneSets()], or the built-in data object [escape.gene.sets]. List names become column names in the result.

method

Scoring algorithm (case-insensitive). One of '"GSVA"', '"ssGSEA"', '"UCell"', or '"AUCell"'. Default **'"ssGSEA"'**.

groups

Integer >= 1. Number of cells per processing chunk. Larger values reduce overhead but increase memory usage. Default **1000**.

min.size

Minimum number of genes from a set that must be detected in the expression matrix for that set to be scored. Default **5**. Use 'NULL' to disable filtering.

normalize

Logical. If 'TRUE', the score matrix is passed to [performNormalization()] (drop-out scaling and optional log transform). Default **FALSE**.

make.positive

Logical. If 'TRUE' *and* 'normalize = TRUE', shifts every gene-set column so its global minimum is zero, facilitating downstream log-ratio analyses. Default **FALSE**.

min.expr.cells

Numeric. Gene-expression filter threshold (see details above). Default **0** (no gene filtering).

min.filter.by

Character or 'NULL'. Column name in 'meta.data' (Seurat) or 'colData' (SCE) defining groups within which the 'min.expr.cells' rule is applied. Default **'NULL'**.

BPPARAM

A BiocParallel parameter object describing the parallel backend.

...

Extra arguments passed verbatim to the chosen back-end scoring function ('gsva()', 'ScoreSignatures_UCell()', or 'AUCell_calcAUC()').

Value

A numeric matrix with one row per cell and one column per gene set, ordered as in 'gene.sets'.

Supported methods

'"GSVA"'

Gene-set variation analysis (Poisson kernel).

'"ssGSEA"'

Single-sample GSEA.

'"UCell"'

Rank-based UCell scoring.

'"AUCell"'

Area-under-the-curve ranking score.

Author(s)

Nick Borcherding, Jared Andrews

See Also

[runEscape()] to attach scores to a single-cell object; [getGeneSets()] for MSigDB retrieval; [performNormalization()] for the optional normalization workflow.

Examples

gs <- list(Bcells = c("MS4A1", "CD79B", "CD79A", "IGH1", "IGH2"),
           Tcells = c("CD3E", "CD3D", "CD3G", "CD7","CD8A"))

pbmc <- SeuratObject::pbmc_small
es   <- escape.matrix(pbmc, 
                      gene.sets = gs,
                      method = "ssGSEA", 
                      groups = 500, 
                      min.size = 3)


ncborcherding/escape documentation built on June 12, 2025, 1 p.m.