UCell 2.10.1
SingleCellExperiment is Bioconductor’s data structure of choice for storing single-cell experiment data. The function ScoreSignatures_UCell() allows performing signature scoring with UCell directly on sce objects. UCell scores are returned in a altExp object: altExp(sce, 'UCell')
For 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)))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-")
)library(UCell)
library(SingleCellExperiment)
sce <- SingleCellExperiment(list(counts=exp.mat))
sce <- ScoreSignatures_UCell(sce, features=signatures, 
                                 assay = 'counts', name=NULL)
altExp(sce, 'UCell')## class: SummarizedExperiment 
## dim: 4 5000 
## metadata(0):
## assays(1): UCell
## rownames(4): Tcell Myeloid NK Plasma_cell
## rowData names(0):
## colnames(5000): bcHTNA1 bcHNVA2 ... bcGVZB4999 bcHGKL5000
## colData names(0):Dimensionality reduction and visualization
library(scater)
library(patchwork)
#PCA
sce <- logNormCounts(sce)
sce <- runPCA(sce, scale=TRUE, ncomponents=10)
#UMAP
set.seed(1234)
sce <- runUMAP(sce, dimred="PCA")Visualize UCell scores on low-dimensional representation (UMAP)
pll <- lapply(names(signatures), function(x) {
    plotUMAP(sce, colour_by = x, by_exprs_values = "UCell",
             point_size=0.2) + theme(aspect.ratio = 1)
})
wrap_plots(pll)Single-cell data are sparse. It can be useful to ‘impute’ scores by neighboring cells and partially correct this sparsity. The function SmoothKNN performs smoothing of single-cell scores by weighted average of the k-nearest neighbors in a given dimensionality reduction. It can be applied directly on SingleCellExperiment objects to smooth UCell scores:
sce <- SmoothKNN(sce, signature.names = names(signatures), reduction="PCA")a <- plotUMAP(sce, colour_by="Myeloid", by_exprs_values="UCell",
         point_size=0.2) + ggtitle("UCell") + theme(aspect.ratio = 1)
b <- plotUMAP(sce, colour_by="Myeloid_kNN", by_exprs_values="UCell_kNN",
         point_size=0.2) + ggtitle("Smoothed UCell") + theme(aspect.ratio = 1)
a | bPlease 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.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
## 
## Matrix products: default
## BLAS:   /media/volume/teran2_disk/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.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] scater_1.34.0               scuttle_1.16.0             
##  [3] patchwork_1.3.0             ggplot2_3.5.1              
##  [5] Seurat_5.1.0                SeuratObject_5.0.2         
##  [7] sp_2.1-4                    UCell_2.10.1               
##  [9] scRNAseq_2.19.1             SingleCellExperiment_1.28.0
## [11] SummarizedExperiment_1.36.0 Biobase_2.66.0             
## [13] GenomicRanges_1.58.0        GenomeInfoDb_1.42.0        
## [15] IRanges_2.40.0              S4Vectors_0.44.0           
## [17] BiocGenerics_0.52.0         MatrixGenerics_1.18.0      
## [19] matrixStats_1.4.1           BiocStyle_2.34.0           
## 
## loaded via a namespace (and not attached):
##   [1] ProtGenerics_1.38.0      spatstat.sparse_3.1-0    bitops_1.0-9            
##   [4] httr_1.4.7               RColorBrewer_1.1-3       tools_4.4.1             
##   [7] sctransform_0.4.1        alabaster.base_1.6.0     utf8_1.2.4              
##  [10] R6_2.5.1                 HDF5Array_1.34.0         lazyeval_0.2.2          
##  [13] uwot_0.2.2               rhdf5filters_1.18.0      withr_3.0.2             
##  [16] gridExtra_2.3            progressr_0.15.0         cli_3.6.3               
##  [19] spatstat.explore_3.3-3   fastDummies_1.7.4        alabaster.se_1.6.0      
##  [22] labeling_0.4.3           sass_0.4.9               spatstat.data_3.1-2     
##  [25] ggridges_0.5.6           pbapply_1.7-2            Rsamtools_2.22.0        
##  [28] parallelly_1.38.0        RSQLite_2.3.7            generics_0.1.3          
##  [31] BiocIO_1.16.0            ica_1.0-3                spatstat.random_3.3-2   
##  [34] dplyr_1.1.4              Matrix_1.7-1             ggbeeswarm_0.7.2        
##  [37] fansi_1.0.6              abind_1.4-8              lifecycle_1.0.4         
##  [40] yaml_2.3.10              rhdf5_2.50.0             SparseArray_1.6.0       
##  [43] BiocFileCache_2.14.0     Rtsne_0.17               grid_4.4.1              
##  [46] blob_1.2.4               promises_1.3.0           ExperimentHub_2.14.0    
##  [49] crayon_1.5.3             miniUI_0.1.1.1           lattice_0.22-6          
##  [52] beachmat_2.22.0          cowplot_1.1.3            GenomicFeatures_1.58.0  
##  [55] KEGGREST_1.46.0          magick_2.8.5             pillar_1.9.0            
##  [58] knitr_1.48               rjson_0.2.23             future.apply_1.11.3     
##  [61] codetools_0.2-20         leiden_0.4.3.1           glue_1.8.0              
##  [64] spatstat.univar_3.0-1    data.table_1.16.2        vctrs_0.6.5             
##  [67] png_0.1-8                gypsum_1.2.0             spam_2.11-0             
##  [70] gtable_0.3.6             cachem_1.1.0             xfun_0.49               
##  [73] S4Arrays_1.6.0           mime_0.12                survival_3.7-0          
##  [76] tinytex_0.53             fitdistrplus_1.2-1       ROCR_1.0-11             
##  [79] nlme_3.1-166             bit64_4.5.2              alabaster.ranges_1.6.0  
##  [82] filelock_1.0.3           RcppAnnoy_0.0.22         bslib_0.8.0             
##  [85] irlba_2.3.5.1            vipor_0.4.7              KernSmooth_2.23-24      
##  [88] colorspace_2.1-1         DBI_1.2.3                tidyselect_1.2.1        
##  [91] bit_4.5.0                compiler_4.4.1           curl_5.2.3              
##  [94] httr2_1.0.5              BiocNeighbors_2.0.0      DelayedArray_0.32.0     
##  [97] plotly_4.10.4            bookdown_0.41            rtracklayer_1.66.0      
## [100] scales_1.3.0             lmtest_0.9-40            rappdirs_0.3.3          
## [103] stringr_1.5.1            digest_0.6.37            goftest_1.2-3           
## [106] spatstat.utils_3.1-0     alabaster.matrix_1.6.0   rmarkdown_2.28          
## [109] XVector_0.46.0           htmltools_0.5.8.1        pkgconfig_2.0.3         
## [112] highr_0.11               dbplyr_2.5.0             fastmap_1.2.0           
## [115] ensembldb_2.30.0         rlang_1.1.4              htmlwidgets_1.6.4       
## [118] UCSC.utils_1.2.0         shiny_1.9.1              farver_2.1.2            
## [121] jquerylib_0.1.4          zoo_1.8-12               jsonlite_1.8.9          
## [124] BiocParallel_1.40.0      BiocSingular_1.22.0      RCurl_1.98-1.16         
## [127] magrittr_2.0.3           GenomeInfoDbData_1.2.13  dotCall64_1.2           
## [130] Rhdf5lib_1.28.0          munsell_0.5.1            Rcpp_1.0.13             
## [133] viridis_0.6.5            reticulate_1.39.0        stringi_1.8.4           
## [136] alabaster.schemas_1.6.0  zlibbioc_1.52.0          MASS_7.3-61             
## [139] AnnotationHub_3.14.0     plyr_1.8.9               parallel_4.4.1          
## [142] listenv_0.9.1            ggrepel_0.9.6            deldir_2.0-4            
## [145] Biostrings_2.74.0        splines_4.4.1            tensor_1.5              
## [148] igraph_2.1.1             spatstat.geom_3.3-3      RcppHNSW_0.6.0          
## [151] ScaledMatrix_1.14.0      reshape2_1.4.4           BiocVersion_3.20.0      
## [154] XML_3.99-0.17            evaluate_1.0.1           BiocManager_1.30.25     
## [157] httpuv_1.6.15            RANN_2.6.2               tidyr_1.3.1             
## [160] purrr_1.0.2              polyclip_1.10-7          future_1.34.0           
## [163] scattermore_1.2          alabaster.sce_1.6.0      rsvd_1.0.5              
## [166] xtable_1.8-4             restfulr_0.0.15          AnnotationFilter_1.30.0 
## [169] RSpectra_0.16-2          later_1.3.2              viridisLite_0.4.2       
## [172] tibble_3.2.1             beeswarm_0.4.0           memoise_2.0.1           
## [175] AnnotationDbi_1.68.0     GenomicAlignments_1.42.0 cluster_2.1.6           
## [178] globals_0.16.3