AUCell_run: Run AUCell

AUCell_runR Documentation

Run AUCell

Description

Runs AUCell (calculates the ranking + score genesets)

Usage

AUCell_run(
  exprMat,
  geneSets,
  featureType = "genes",
  keepZeroesAsNA = FALSE,
  normAUC = TRUE,
  aucMaxRank = ceiling(0.05 * nrow(exprMat)),
  BPPARAM = NULL,
  ...
)

## S4 method for signature 'dgCMatrix'
AUCell_run(
  exprMat,
  geneSets,
  featureType = "genes",
  keepZeroesAsNA = FALSE,
  normAUC = TRUE,
  aucMaxRank = ceiling(0.05 * nrow(exprMat)),
  BPPARAM = NULL
)

## S4 method for signature 'matrix'
AUCell_run(
  exprMat,
  geneSets,
  featureType = "genes",
  keepZeroesAsNA = FALSE,
  normAUC = TRUE,
  aucMaxRank = ceiling(0.05 * nrow(exprMat)),
  BPPARAM = NULL
)

## S4 method for signature 'SummarizedExperiment'
AUCell_run(
  exprMat,
  geneSets,
  featureType = "genes",
  keepZeroesAsNA = FALSE,
  normAUC = TRUE,
  aucMaxRank = ceiling(0.05 * nrow(exprMat)),
  BPPARAM = NULL,
  assayName = NULL
)

Arguments

exprMat

Expression matrix (genes/regions as rows, cells as columns) The expression matrix can also be provided as one of the R/Bioconductor classes:

  • dgCMatrix-class: Sparse matrix

  • RangedSummarizedExperiment and derived classes (e.g. SingleCellExperiment ): The matrix will be obtained through assay(exprMatrix), -which will extract the first assay (usually the counts)- or the assay name given in 'assayName'

geneSets

List of gene-sets (or signatures) to test in the cells. The gene-sets should be provided as GeneSet, GeneSetCollection or character list (see examples).

featureType

Name for the rows (e.g. "genes"). Only for naming the rankings, not used internally.

keepZeroesAsNA

Convert zeroes to NA instead of locating randomly at the end of the ranking.

normAUC

Whether to normalize the maximum possible AUC to 1 (Default: TRUE).

aucMaxRank

Threshold to calculate the AUC (see 'details' section)

BPPARAM

Set to use multiple cores. Only used if 'splitByBlocks=TRUE'

...

Other arguments

assayName

Name of the assay containing the expression matrix (e.g. in SingleCellExperiment objects)

Details

In a simplified way, the AUC value represents the fraction of genes, within the top X genes in the ranking, that are included in the signature. The parameter 'aucMaxRank' allows to modify the number of genes (maximum ranking) that is used to perform this computation. By default, it is set to 5% of the total number of genes in the rankings. Common values may range from 1 to 20%.

Value

Matrix with the AUC values (gene-sets as rows, cells as columns).

See Also

Includes AUCell_buildRankings and AUCell_calcAUC. Next step in the workflow: AUCell_exploreThresholds.

See the package vignette for examples and more details: vignette("AUCell")

Examples

# This example is run using a fake expression matrix.
# Therefore, the output will be meaningless.

############# Fake expression matrix #############
set.seed(123)
exprMatrix <- matrix(data=sample(c(rep(0, 5000), sample(1:3, 5000, replace=TRUE))),
                     nrow=20, 
                     dimnames=list(paste("Gene", 1:20, sep=""), 
                                   paste("Cell", 1:500, sep="")))
exprMatrix <- as(exprMatrix, "dgCMatrix")

# In this example we use two gene sets: 10 and 5 random genes
# (see other formatting examples at the end)
fewGenes <- sample(rownames(exprMatrix), 10)
otherGenes <- sample(rownames(exprMatrix), 5)

geneSets <- list(geneSet1=fewGenes,
                 geneSet2=otherGenes)
geneSets

# Calculate AUCell score for the genes in the sets
# To be able to run this fake example (which contain only 20 genes),
# we use aucMaxRank=5 (top 25% of the genes in the ranking)
cells_AUC <- AUCell_run(exprMatrix, geneSets, aucMaxRank=5)


## To run in paralell:
# cells_AUC <- AUCell_run(exprMatrix, geneSets, aucMaxRank=5, 
#                   BPPARAM=BiocParallel::MulticoreParam(5))



# Format of the output:
cells_AUC

# To subset & access the AUC slot (as matrix):
cells_AUC[1:2,]
cells_AUC[,3:4]
getAUC(cells_AUC)[,1:5]


# These methods are also available:
dim(cells_AUC)
nrow(cells_AUC)
ncol(cells_AUC)
colnames(cells_AUC)[1:4]
rownames(cells_AUC)

#########################################################
# Alternatives for the input of gene sets:

# a) Character vector (i.e. only one gene-set)
# It will take the default name 'geneSet'
fewGenes
test <- AUCell_run(exprMatrix, fewGenes, aucMaxRank=5)

# b) List
geneSets <- list(geneSet1=fewGenes,
                 geneSet2=otherGenes)
geneSets
test <- AUCell_run(exprMatrix, fewGenes, aucMaxRank=5)

# c) GeneSet object (from GSEABase)
library(GSEABase)
geneSetOne <- GeneSet(fewGenes, setName="geneSetOne")
geneSetOne
test <- AUCell_run(exprMatrix, fewGenes, aucMaxRank=5)


# d) GeneSetCollection object (from GSEABase)
geneSetTwo <- GeneSet(otherGenes, setName="geneSetTwo")
geneSets <- GeneSetCollection(geneSetOne, geneSetTwo)
geneSets
test <- AUCell_run(exprMatrix, fewGenes, aucMaxRank=5)


aertslab/AUCell documentation built on March 12, 2024, 11:40 p.m.