ScoreSignatures_UCell: Calculate module enrichment scores from single-cell data

View source: R/ScoreSignatures_UCell.R

ScoreSignatures_UCellR Documentation

Calculate module enrichment scores from single-cell data

Description

Given a gene vs. cell matrix, calculates module/signature enrichment scores on single-cell level using Mann-Whitney U statistic. UCell scores are normalized U statistics (between 0 and 1), and they are mathematically related to the Area under the ROC curve (see Mason and Graham) These scores only depend on the gene expression ranks of individual cell, and therefore they are robust across datasets regardless of dataset composition.

Usage

ScoreSignatures_UCell(
  matrix = NULL,
  features,
  precalc.ranks = NULL,
  maxRank = 1500,
  w_neg = 1,
  name = "_UCell",
  assay = "counts",
  chunk.size = 100,
  BPPARAM = NULL,
  ncores = 1,
  ties.method = "average",
  force.gc = FALSE
)

Arguments

matrix

Input matrix, either stored in a SingleCellExperiment object or as a raw matrix. dgCMatrix format supported.

features

A list of signatures, for example: list(Tcell_signature = c("CD2","CD3E","CD3D"), Myeloid_signature = c("SPI1","FCER1G","CSF1R")) You can also specify positive and negative gene sets by adding a + or - sign to genes in the signature; see an example below

precalc.ranks

If you have pre-calculated ranks using StoreRankings_UCell, you can specify the pre-calculated ranks instead of the gene vs. cell matrix.

maxRank

Maximum number of genes to rank per cell; above this rank, a given gene is considered as not expressed. Note: this parameter is ignored if precalc.ranks are specified

w_neg

Weight on negative genes in signature. e.g. w_neg=1 weighs equally up- and down-regulated genes, w_neg=0.5 gives 50% less importance to negative genes

name

Name suffix appended to signature names

assay

The sce object assay where the data is to be found

chunk.size

Number of cells to be processed simultaneously (lower size requires slightly more computation but reduces memory demands)

BPPARAM

A BiocParallel::bpparam() object that tells UCell how to parallelize. If provided, it overrides the ncores parameter.

ncores

Number of processors to parallelize computation. If BPPARAM = NULL, the function uses BiocParallel::MulticoreParam(workers=ncores)

ties.method

How ranking ties should be resolved - passed on to data.table::frank

force.gc

Explicitly call garbage collector to reduce memory footprint

Value

Returns input SingleCellExperiment object with UCell scores added to altExp

Examples

library(UCell)
# Using sparse matrix
data(sample.matrix)
gene.sets <- list( Tcell_signature = c("CD2","CD3E","CD3D"),
                 Myeloid_signature = c("SPI1","FCER1G","CSF1R"))
scores <- ScoreSignatures_UCell(sample.matrix, features=gene.sets)
head(scores)

# Using sce object
library(SingleCellExperiment)
data(sample.matrix)
my.sce <- SingleCellExperiment(list(counts=sample.matrix))
gene.sets <- list( Tcell_signature = c("CD2","CD3E","CD3D"),
                 Myeloid_signature = c("SPI1","FCER1G","CSF1R"))
my.sce <- ScoreSignatures_UCell(my.sce, features=gene.sets)
altExp(my.sce, 'UCell')


carmonalab/UCell documentation built on Nov. 4, 2024, 5:32 p.m.