UCell 2.13.3
In single-cell RNA-seq analysis, gene signature (or “module”) scoring constitutes a simple yet powerful approach to evaluate the strength of biological signals, typically associated to a specific cell type or biological process, in a transcriptome.
UCell is an R package for evaluating gene signatures in single-cell datasets. UCell signature scores, based on the Mann-Whitney U statistic, are robust to dataset size and heterogeneity, and their calculation demands less computing time and memory than other available methods, enabling the processing of large datasets in a few minutes even on machines with limited computing power. UCell can be applied to any single-cell data matrix, and includes functions to directly interact with Seurat objects.
To test your installation, load a small sample dataset and run UCell:
library(UCell)
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)##                       Tcell_signature_UCell Myeloid_signature_UCell
## L5_ATTTCTGAGGTCGTGA               0.8991989                       0
## L4_TCACTATTCATCTCTA               0.6900312                       0
## L1_TCCTTCTTCTTTACAC               0.5932354                       0
## L5_AAAGTGAAGGCGCTCT               0.6090343                       0
## E2L3_CCTCAGTAGTGCAGGT             0.8506898                       0
## L5_CCCTCTCGTTCTAAGC               0.6557632                       0For this demo, we will download a single-cell dataset of lung cancer (Zilionis et al. (2019) Immunity) through the scRNA-seq package. This dataset contains >170,000 single cells; for the sake of simplicity, in this demo will we focus on immune cells, according to the annotations by the authors, and downsample to 5000 cells.
library(scRNAseq)
lung <- ZilionisLungData()
immune <- lung$Used & lung$used_in_NSCLC_immune
lung <- lung[,immune]
lung <- lung[,1:5000]
exp.mat <- Matrix::Matrix(counts(lung),sparse = TRUE)
colnames(exp.mat) <- paste0(colnames(exp.mat), seq(1,ncol(exp.mat)))Note: becase UCell scores are based on relative gene ranks, it can be applied both on raw counts or normalized data. As long as the normalization preserves the relative ranks between genes, the results will be equivalent.
Here we define some simple gene sets based on the “Human Cell Landscape” signatures Han et al. (2020) Nature. You may edit existing signatures, or add new one as elements in a list.
signatures <- list(
    Tcell = c("CD3D","CD3E","CD3G","CD2","TRAC"),
    Myeloid = c("CD14","LYZ","CSF1R","FCER1G","SPI1","LCK-"),
    NK = c("KLRD1","NCR1","NKG7","CD3D-","CD3E-"),
    Plasma_cell = c("MZB1","DERL3","CD19-")
)Run ScoreSignatures_UCell and get directly signature scores for all cells
u.scores <- ScoreSignatures_UCell(exp.mat,features=signatures)
head(u.scores)##         Tcell_UCell Myeloid_UCell NK_UCell Plasma_cell_UCell
## bcHTNA1           0     0.5227121        0        0.00000000
## bcHNVA2           0     0.5112892        0        0.00000000
## bcALZN3           0     0.3584502        0        0.07540874
## bcFWBP4           0     0.1546426        0        0.00000000
## bcBJYE5           0     0.4629927        0        0.00000000
## bcGSBJ6           0     0.5452238        0        0.00000000Show the distribution of predicted scores
library(reshape2)
library(ggplot2)
melted <- reshape2::melt(u.scores)
colnames(melted) <- c("Cell","Signature","UCell_score")
p <- ggplot(melted, aes(x=Signature, y=UCell_score)) + 
    geom_violin(aes(fill=Signature), scale = "width") +
    geom_boxplot(width=0.1, outlier.size=0) +
    theme_bw() + theme(axis.text.x=element_blank())
pThe time- and memory-demanding step in UCell is the calculation of gene rankings for each individual cell. If we plan to experiment with signatures, editing them or adding new cell subtypes, it is possible to pre-calculate the gene rankings once and for all and then apply new signatures over these pre-calculated ranks. Run the StoreRankings_UCell function to pre-calculate gene rankings over a dataset:
set.seed(123)
ranks <- StoreRankings_UCell(exp.mat)
ranks[1:5,1:5]## 5 x 5 sparse Matrix of class "dgCMatrix"
##           bcHTNA1 bcHNVA2 bcALZN3 bcFWBP4 bcBJYE5
## 5S_rRNA         .       .       .       .       .
## 5_8S_rRNA       .       .       .       .       .
## 7SK             .       .       .       .       .
## A1BG            .       .       .       .       .
## A1BG-AS1        .       .       .       .       .Then, we can apply our signature set, or any other new signature to the pre-calculated ranks. The calculations will be considerably faster.
set.seed(123)
u.scores.2 <- ScoreSignatures_UCell(features=signatures,
                                    precalc.ranks = ranks)
melted <- reshape2::melt(u.scores.2)
colnames(melted) <- c("Cell","Signature","UCell_score")
p <- ggplot(melted, aes(x=Signature, y=UCell_score)) + 
    geom_violin(aes(fill=Signature), scale = "width") +
    geom_boxplot(width=0.1, outlier.size = 0) + 
    theme_bw() + theme(axis.text.x=element_blank())
pnew.signatures <- list(Mast.cell = c("TPSAB1","TPSB2","CPA3","MS4A2"),
                       Lymphoid = c("LCK"))
u.scores.3 <- ScoreSignatures_UCell(features=new.signatures,
                                    precalc.ranks = ranks)
melted <- reshape2::melt(u.scores.3)
colnames(melted) <- c("Cell","Signature","UCell_score")
p <- ggplot(melted, aes(x=Signature, y=UCell_score)) + 
    geom_violin(aes(fill=Signature), scale = "width") +
    geom_boxplot(width=0.1, outlier.size=0) + 
    theme_bw() + theme(axis.text.x=element_blank())
pIf your machine has multi-core capabilities and enough RAM, running UCell in parallel can speed up considerably your analysis. The example below runs on a single core - you may modify this behavior by setting e.g. workers=4 to parallelize to 4 cores:
BPPARAM <- BiocParallel::MulticoreParam(workers=1)
u.scores <- ScoreSignatures_UCell(exp.mat,features=signatures,
                                  BPPARAM=BPPARAM)SingleCellExperiment and Seurat are popular environments for single-cell analysis. The UCell package implements functions to interact directly with these pipelines, as described in dedicated demos available on the Bioc landing page.
Please report any issues at the UCell GitHub repository.
More demos available on the Bioc landing page and at the UCell demo repository.
If you find UCell useful, you may also check out the scGate package, which relies on UCell scores to automatically purify populations of interest based on gene signatures.
See also SignatuR for easy storing and retrieval of gene signatures.
sessionInfo()## R version 4.5.1 Patched (2025-08-23 r88802)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              LC_COLLATE=C              
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] reshape2_1.4.4              scater_1.37.0              
##  [3] scuttle_1.19.0              patchwork_1.3.2            
##  [5] ggplot2_4.0.0               UCell_2.13.3               
##  [7] Seurat_5.3.0                SeuratObject_5.2.0         
##  [9] sp_2.2-0                    scRNAseq_2.23.0            
## [11] SingleCellExperiment_1.31.1 SummarizedExperiment_1.39.2
## [13] Biobase_2.69.1              GenomicRanges_1.61.5       
## [15] Seqinfo_0.99.2              IRanges_2.43.5             
## [17] S4Vectors_0.47.4            BiocGenerics_0.55.1        
## [19] generics_0.1.4              MatrixGenerics_1.21.0      
## [21] matrixStats_1.5.0           BiocStyle_2.37.1           
## 
## loaded via a namespace (and not attached):
##   [1] ProtGenerics_1.41.0      spatstat.sparse_3.1-0    bitops_1.0-9            
##   [4] httr_1.4.7               RColorBrewer_1.1-3       tools_4.5.1             
##   [7] sctransform_0.4.2        alabaster.base_1.9.5     R6_2.6.1                
##  [10] HDF5Array_1.37.0         lazyeval_0.2.2           uwot_0.2.3              
##  [13] rhdf5filters_1.21.0      withr_3.0.2              gridExtra_2.3           
##  [16] progressr_0.16.0         cli_3.6.5                spatstat.explore_3.5-3  
##  [19] fastDummies_1.7.5        alabaster.se_1.9.0       labeling_0.4.3          
##  [22] sass_0.4.10              S7_0.2.0                 spatstat.data_3.1-8     
##  [25] ggridges_0.5.7           pbapply_1.7-4            Rsamtools_2.25.3        
##  [28] dichromat_2.0-0.1        parallelly_1.45.1        RSQLite_2.4.3           
##  [31] BiocIO_1.19.0            ica_1.0-3                spatstat.random_3.4-2   
##  [34] dplyr_1.1.4              Matrix_1.7-4             ggbeeswarm_0.7.2        
##  [37] abind_1.4-8              lifecycle_1.0.4          yaml_2.3.10             
##  [40] rhdf5_2.53.5             SparseArray_1.9.1        BiocFileCache_2.99.6    
##  [43] Rtsne_0.17               grid_4.5.1               blob_1.2.4              
##  [46] promises_1.3.3           ExperimentHub_2.99.5     crayon_1.5.3            
##  [49] miniUI_0.1.2             lattice_0.22-7           beachmat_2.25.5         
##  [52] cowplot_1.2.0            GenomicFeatures_1.61.6   KEGGREST_1.49.1         
##  [55] magick_2.9.0             pillar_1.11.1            knitr_1.50              
##  [58] rjson_0.2.23             future.apply_1.20.0      codetools_0.2-20        
##  [61] glue_1.8.0               spatstat.univar_3.1-4    data.table_1.17.8       
##  [64] vctrs_0.6.5              png_0.1-8                gypsum_1.5.0            
##  [67] spam_2.11-1              gtable_0.3.6             cachem_1.1.0            
##  [70] xfun_0.53                S4Arrays_1.9.1           mime_0.13               
##  [73] survival_3.8-3           tinytex_0.57             fitdistrplus_1.2-4      
##  [76] ROCR_1.0-11              nlme_3.1-168             bit64_4.6.0-1           
##  [79] alabaster.ranges_1.9.1   filelock_1.0.3           RcppAnnoy_0.0.22        
##  [82] GenomeInfoDb_1.45.12     bslib_0.9.0              irlba_2.3.5.1           
##  [85] vipor_0.4.7              KernSmooth_2.23-26       DBI_1.2.3               
##  [88] ggrastr_1.0.2            tidyselect_1.2.1         bit_4.6.0               
##  [91] compiler_4.5.1           curl_7.0.0               httr2_1.2.1             
##  [94] BiocNeighbors_2.3.1      h5mread_1.1.1            DelayedArray_0.35.3     
##  [97] plotly_4.11.0            bookdown_0.45            rtracklayer_1.69.1      
## [100] scales_1.4.0             lmtest_0.9-40            rappdirs_0.3.3          
## [103] stringr_1.5.2            digest_0.6.37            goftest_1.2-3           
## [106] spatstat.utils_3.2-0     alabaster.matrix_1.9.0   rmarkdown_2.30          
## [109] XVector_0.49.1           htmltools_0.5.8.1        pkgconfig_2.0.3         
## [112] dbplyr_2.5.1             fastmap_1.2.0            ensembldb_2.33.2        
## [115] rlang_1.1.6              htmlwidgets_1.6.4        UCSC.utils_1.5.0        
## [118] shiny_1.11.1             farver_2.1.2             jquerylib_0.1.4         
## [121] zoo_1.8-14               jsonlite_2.0.0           BiocParallel_1.43.4     
## [124] BiocSingular_1.25.0      RCurl_1.98-1.17          magrittr_2.0.4          
## [127] dotCall64_1.2            Rhdf5lib_1.31.0          Rcpp_1.1.0              
## [130] viridis_0.6.5            reticulate_1.43.0        stringi_1.8.7           
## [133] alabaster.schemas_1.9.0  MASS_7.3-65              AnnotationHub_3.99.6    
## [136] plyr_1.8.9               parallel_4.5.1           listenv_0.9.1           
## [139] ggrepel_0.9.6            deldir_2.0-4             Biostrings_2.77.2       
## [142] splines_4.5.1            tensor_1.5.1             igraph_2.1.4            
## [145] spatstat.geom_3.6-0      RcppHNSW_0.6.0           ScaledMatrix_1.17.0     
## [148] BiocVersion_3.22.0       XML_3.99-0.19            evaluate_1.0.5          
## [151] BiocManager_1.30.26      httpuv_1.6.16            RANN_2.6.2              
## [154] tidyr_1.3.1              purrr_1.1.0              polyclip_1.10-7         
## [157] future_1.67.0            scattermore_1.2          alabaster.sce_1.9.0     
## [160] rsvd_1.0.5               xtable_1.8-4             restfulr_0.0.16         
## [163] AnnotationFilter_1.33.0  RSpectra_0.16-2          later_1.4.4             
## [166] viridisLite_0.4.2        tibble_3.3.0             memoise_2.0.1           
## [169] beeswarm_0.4.0           AnnotationDbi_1.71.1     GenomicAlignments_1.45.4
## [172] cluster_2.1.8.1          globals_0.18.0