runEscape: Enrichment calculation for single-cell workflows

View source: R/runEscape.R

runEscapeR Documentation

Enrichment calculation for single-cell workflows

Description

Run the escape-based gene-set enrichment calculation with Seurat or SingleCellExperiment pipelines

Usage

runEscape(
  input.data,
  gene.sets = NULL,
  method = "ssGSEA",
  groups = 1000,
  min.size = 5,
  normalize = FALSE,
  make.positive = FALSE,
  new.assay.name = "escape",
  BPPARAM = SerialParam(),
  ...
)

Arguments

input.data

The count matrix, Seurat, or Single-Cell Experiment object.

gene.sets

Gene sets can be a list, output from getGeneSets, or the built-in gene sets in the escape package escape.gene.sets.

method

Select the method to calculate enrichment, AUCell, GSVA, ssGSEA or UCell.

groups

The number of cells to separate the enrichment calculation.

min.size

Minimum number of gene necessary to perform the enrichment calculation

normalize

Whether to divide the enrichment score by the number of genes TRUE or report unnormalized FALSE.

make.positive

During normalization shift enrichment values to a positive range TRUE for downstream analysis or not TRUE (default). Will only be applied if normalize = TRUE.

new.assay.name

The new name of the assay to append to the single-cell object containing the enrichment scores.

BPPARAM

A BiocParallel::bpparam() object that for parallelization.

...

pass arguments to AUCell GSVA, ssGSEA or UCell call

Value

Seurat or Single-Cell Experiment object with escape enrichment scores in the assay slot.

Examples

GS <- list(Bcells = c("MS4A1", "CD79B", "CD79A", "IGH1", "IGH2"),
           Tcells = c("CD3E", "CD3D", "CD3G", "CD7","CD8A"))
pbmc_small <- SeuratObject::pbmc_small
pbmc_small <- runEscape(pbmc_small, 
                        gene.sets = GS, 
                        min.size = NULL)


ncborcherding/escape documentation built on Nov. 6, 2024, 1:43 p.m.