1 Introduction

SingleR is an automatic annotation method for single-cell RNA sequencing (scRNAseq) data (Aran et al. 2019). Given a reference dataset of samples (single-cell or bulk) with known labels, it labels new cells from a test dataset based on similarity to the reference. Thus, the burden of manually interpreting clusters and defining marker genes only has to be done once, for the reference dataset, and this biological knowledge can be propagated to new datasets in an automated manner.

To keep things brief, this vignette only provides a brief summary of the basic capabilities of SingleR. However, the package also provides more advanced functionality that includes the use of multiple references simultaneously, manipulating the cell ontology and improving performance on big datasets. Readers are referred to the book for more details.

2 Using built-in references

The easiest way to use SingleR is to annotate cells against built-in references. In particular, the celldex package provides access to several reference datasets (mostly derived from bulk RNA-seq or microarray data) through dedicated retrieval functions. Here, we will use the Human Primary Cell Atlas (Mabbott et al. 2013), represented as a SummarizedExperiment object containing a matrix of log-expression values with sample-level labels.

library(celldex)
hpca.se <- HumanPrimaryCellAtlasData()
hpca.se
## class: SummarizedExperiment 
## dim: 19363 713 
## metadata(0):
## assays(1): logcounts
## rownames(19363): A1BG A1BG-AS1 ... ZZEF1 ZZZ3
## rowData names(0):
## colnames(713): GSM112490 GSM112491 ... GSM92233 GSM92234
## colData names(3): label.main label.fine label.ont

Our test dataset consists of some human embryonic stem cells (La Manno et al. 2016) from the scRNAseq package. For the sake of speed, we will only label the first 100 cells from this dataset.

library(scRNAseq)
hESCs <- LaMannoBrainData('human-es')
hESCs <- hESCs[,1:100]

We use our hpca.se reference to annotate each cell in hESCs via the SingleR() function. This identifies marker genes from the reference and uses them to compute assignment scores (based on the Spearman correlation across markers) for each cell in the test dataset against each label in the reference. The label with the highest score is the assigned to the test cell, possibly with further fine-tuning to resolve closely related labels.

library(SingleR)
pred.hesc <- SingleR(test = hESCs, ref = hpca.se, assay.type.test=1,
    labels = hpca.se$label.main)

Each row of the output DataFrame contains prediction results for a single cell. Labels are shown before fine-tuning (first.labels), after fine-tuning (labels) and after pruning (pruned.labels), along with the associated scores.

pred.hesc
## DataFrame with 100 rows and 5 columns
##                                          scores         first.labels
##                                        <matrix>          <character>
## 1772122_301_C02  0.347652:0.109547:0.123901:... Neuroepithelial_cell
## 1772122_180_E05  0.361187:0.134934:0.148672:... Neuroepithelial_cell
## 1772122_300_H02  0.446411:0.190084:0.222594:... Neuroepithelial_cell
## 1772122_180_B09  0.373512:0.143537:0.164743:... Neuroepithelial_cell
## 1772122_180_G04  0.357341:0.126511:0.141987:... Neuroepithelial_cell
## ...                                         ...                  ...
## 1772122_299_E07 0.371989:0.169379:0.1986877:... Neuroepithelial_cell
## 1772122_180_D02 0.353314:0.115864:0.1374981:... Neuroepithelial_cell
## 1772122_300_D09 0.348789:0.136732:0.1303042:... Neuroepithelial_cell
## 1772122_298_F09 0.332361:0.141439:0.1437860:... Neuroepithelial_cell
## 1772122_302_A11 0.324928:0.101609:0.0949826:... Neuroepithelial_cell
##                       tuning.scores               labels        pruned.labels
##                         <DataFrame>          <character>          <character>
## 1772122_301_C02 0.1824402:0.0991116 Neuroepithelial_cell Neuroepithelial_cell
## 1772122_180_E05 0.1375484:0.0647134              Neurons              Neurons
## 1772122_300_H02 0.2757982:0.1369690 Neuroepithelial_cell Neuroepithelial_cell
## 1772122_180_B09 0.0851623:0.0819878 Neuroepithelial_cell Neuroepithelial_cell
## 1772122_180_G04 0.1988415:0.1016622 Neuroepithelial_cell Neuroepithelial_cell
## ...                             ...                  ...                  ...
## 1772122_299_E07 0.1760025:0.0922504 Neuroepithelial_cell Neuroepithelial_cell
## 1772122_180_D02 0.1967609:0.1124805 Neuroepithelial_cell Neuroepithelial_cell
## 1772122_300_D09 0.0816424:0.0221368 Neuroepithelial_cell Neuroepithelial_cell
## 1772122_298_F09 0.1872499:0.0671893 Neuroepithelial_cell Neuroepithelial_cell
## 1772122_302_A11 0.1560800:0.1051322            Astrocyte            Astrocyte
# Summarizing the distribution:
table(pred.hesc$labels)
## 
##            Astrocyte Neuroepithelial_cell              Neurons 
##                   14                   81                    5

At this point, it is worth noting that SingleR is workflow/package agnostic. The above example uses SummarizedExperiment objects, but the same functions will accept any (log-)normalized expression matrix.

3 Using single-cell references

Here, we will use two human pancreas datasets from the scRNAseq package. The aim is to use one pre-labelled dataset to annotate the other unlabelled dataset. First, we set up the Muraro et al. (2016) dataset to be our reference.

library(scRNAseq)
sceM <- MuraroPancreasData()

# One should normally do cell-based quality control at this point, but for
# brevity's sake, we will just remove the unlabelled libraries here.
sceM <- sceM[,!is.na(sceM$label)]

# SingleR() expects reference datasets to be normalized and log-transformed.
library(scuttle)
sceM <- logNormCounts(sceM)

We then set up our test dataset from Grun et al. (2016). To speed up this demonstration, we will subset to the first 100 cells.

sceG <- GrunPancreasData()
sceG <- sceG[,colSums(counts(sceG)) > 0] # Remove libraries with no counts.
sceG <- logNormCounts(sceG) 
sceG <- sceG[,1:100]

We then run SingleR() as described previously but with a marker detection mode that considers the variance of expression across cells. Here, we will use the Wilcoxon ranked sum test to identify the top markers for each pairwise comparison between labels. This is slower but more appropriate for single-cell data compared to the default marker detection algorithm (which may fail for low-coverage data where the median is frequently zero).

pred.grun <- SingleR(test=sceG, ref=sceM, labels=sceM$label, de.method="wilcox")
table(pred.grun$labels)
## 
## acinar   beta  delta   duct 
##     53      4      2     41

4 Annotation diagnostics

plotScoreHeatmap() displays the scores for all cells across all reference labels, which allows users to inspect the confidence of the predicted labels across the dataset. Ideally, each cell (i.e., column of the heatmap) should have one score that is obviously larger than the rest, indicating that it is unambiguously assigned to a single label. A spread of similar scores for a given cell indicates that the assignment is uncertain, though this may be acceptable if the uncertainty is distributed across similar cell types that cannot be easily resolved.

plotScoreHeatmap(pred.grun)

Another diagnostic is based on the per-cell “deltas”, i.e., the difference between the score for the assigned label and the median across all labels for each cell. Low deltas indicate that the assignment is uncertain, which is especially relevant if the cell’s true label does not exist in the reference. We can inspect these deltas across cells for each label using the plotDeltaDistribution() function.

plotDeltaDistribution(pred.grun, ncol = 3)

The pruneScores() function will remove potentially poor-quality or ambiguous assignments based on the deltas. The minimum threshold on the deltas is defined using an outlier-based approach that accounts for differences in the scale of the correlations in various contexts - see ?pruneScores for more details. SingleR() will also report the pruned scores automatically in the pruned.labels field where low-quality assignments are replaced with NA.

summary(is.na(pred.grun$pruned.labels))
##    Mode   FALSE    TRUE 
## logical      96       4

Finally, a simple yet effective diagnostic is to examine the expression of the marker genes for each label in the test dataset. We extract the identity of the markers from the metadata of the SingleR() results and use them in the plotHeatmap() function from scater, as shown below for beta cell markers. If a cell in the test dataset is confidently assigned to a particular label, we would expect it to have strong expression of that label’s markers. At the very least, it should exhibit upregulation of those markers relative to cells assigned to other labels.

all.markers <- metadata(pred.grun)$de.genes
sceG$labels <- pred.grun$labels

# Beta cell-related markers
library(scater)
plotHeatmap(sceG, order_columns_by="labels",
    features=unique(unlist(all.markers$beta))) 

5 FAQs

How can I use this with my Seurat, SingleCellExperiment, or cell_data_set object?

SingleR is workflow agnostic - all it needs is normalized counts. An example showing how to map its results back to common single-cell data objects is available in the README.

Where can I find reference sets appropriate for my data?

celldex contains many built-in references that were previously in SingleR but have been migrated into a separate package for more general use by other Bioconductor packages. scRNAseq contains many single-cell datasets, many of which contain the authors’ manual annotations. ArrayExpress and GEOquery can be used to download any of the bulk or single-cell datasets in ArrayExpress or GEO, respectively.

Where can I get more help?

It is likely that your questions is already answered by the function documentation (e.g., ?SingleR). Further explanations on the reasoning behind certain functions can be found in the book. If this is not sufficient, we recommend posting an issue on the Bioconductor support site or the GitHub repository for this package. Be sure to include your session information and a minimal reproducible example.

6 Session information

sessionInfo()
## R version 4.2.0 RC (2022-04-19 r82224)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
## 
## 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       
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] scater_1.24.0               ggplot2_3.3.5              
##  [3] scuttle_1.6.0               SingleR_1.10.0             
##  [5] scRNAseq_2.9.2              SingleCellExperiment_1.18.0
##  [7] celldex_1.5.0               SummarizedExperiment_1.26.0
##  [9] Biobase_2.56.0              GenomicRanges_1.48.0       
## [11] GenomeInfoDb_1.32.0         IRanges_2.30.0             
## [13] S4Vectors_0.34.0            BiocGenerics_0.42.0        
## [15] MatrixGenerics_1.8.0        matrixStats_0.62.0         
## [17] BiocStyle_2.24.0           
## 
## loaded via a namespace (and not attached):
##   [1] AnnotationHub_3.4.0           BiocFileCache_2.4.0          
##   [3] igraph_1.3.1                  lazyeval_0.2.2               
##   [5] BiocParallel_1.30.0           digest_0.6.29                
##   [7] ensembldb_2.20.0              htmltools_0.5.2              
##   [9] magick_2.7.3                  viridis_0.6.2                
##  [11] fansi_1.0.3                   magrittr_2.0.3               
##  [13] memoise_2.0.1                 ScaledMatrix_1.4.0           
##  [15] cluster_2.1.3                 limma_3.52.0                 
##  [17] Biostrings_2.64.0             prettyunits_1.1.1            
##  [19] colorspace_2.0-3              ggrepel_0.9.1                
##  [21] blob_1.2.3                    rappdirs_0.3.3               
##  [23] xfun_0.30                     dplyr_1.0.8                  
##  [25] crayon_1.5.1                  RCurl_1.98-1.6               
##  [27] jsonlite_1.8.0                glue_1.6.2                   
##  [29] gtable_0.3.0                  zlibbioc_1.42.0              
##  [31] XVector_0.36.0                DelayedArray_0.22.0          
##  [33] BiocSingular_1.12.0           scales_1.2.0                 
##  [35] pheatmap_1.0.12               DBI_1.1.2                    
##  [37] edgeR_3.38.0                  Rcpp_1.0.8.3                 
##  [39] viridisLite_0.4.0             xtable_1.8-4                 
##  [41] progress_1.2.2                dqrng_0.3.0                  
##  [43] bit_4.0.4                     rsvd_1.0.5                   
##  [45] metapod_1.4.0                 httr_1.4.2                   
##  [47] RColorBrewer_1.1-3            ellipsis_0.3.2               
##  [49] farver_2.1.0                  pkgconfig_2.0.3              
##  [51] XML_3.99-0.9                  sass_0.4.1                   
##  [53] dbplyr_2.1.1                  locfit_1.5-9.5               
##  [55] utf8_1.2.2                    labeling_0.4.2               
##  [57] tidyselect_1.1.2              rlang_1.0.2                  
##  [59] later_1.3.0                   AnnotationDbi_1.58.0         
##  [61] munsell_0.5.0                 BiocVersion_3.15.2           
##  [63] tools_4.2.0                   cachem_1.0.6                 
##  [65] cli_3.3.0                     generics_0.1.2               
##  [67] RSQLite_2.2.12                ExperimentHub_2.4.0          
##  [69] evaluate_0.15                 stringr_1.4.0                
##  [71] fastmap_1.1.0                 yaml_2.3.5                   
##  [73] knitr_1.38                    bit64_4.0.5                  
##  [75] purrr_0.3.4                   KEGGREST_1.36.0              
##  [77] AnnotationFilter_1.20.0       sparseMatrixStats_1.8.0      
##  [79] mime_0.12                     scran_1.24.0                 
##  [81] xml2_1.3.3                    biomaRt_2.52.0               
##  [83] compiler_4.2.0                beeswarm_0.4.0               
##  [85] filelock_1.0.2                curl_4.3.2                   
##  [87] png_0.1-7                     interactiveDisplayBase_1.34.0
##  [89] tibble_3.1.6                  statmod_1.4.36               
##  [91] bslib_0.3.1                   stringi_1.7.6                
##  [93] highr_0.9                     GenomicFeatures_1.48.0       
##  [95] lattice_0.20-45               bluster_1.6.0                
##  [97] ProtGenerics_1.28.0           Matrix_1.4-1                 
##  [99] vctrs_0.4.1                   pillar_1.7.0                 
## [101] lifecycle_1.0.1               BiocManager_1.30.17          
## [103] jquerylib_0.1.4               BiocNeighbors_1.14.0         
## [105] bitops_1.0-7                  irlba_2.3.5                  
## [107] httpuv_1.6.5                  rtracklayer_1.56.0           
## [109] R6_2.5.1                      BiocIO_1.6.0                 
## [111] bookdown_0.26                 promises_1.2.0.1             
## [113] gridExtra_2.3                 vipor_0.4.5                  
## [115] assertthat_0.2.1              rjson_0.2.21                 
## [117] withr_2.5.0                   GenomicAlignments_1.32.0     
## [119] Rsamtools_2.12.0              GenomeInfoDbData_1.2.8       
## [121] parallel_4.2.0                hms_1.1.1                    
## [123] grid_4.2.0                    beachmat_2.12.0              
## [125] rmarkdown_2.14                DelayedMatrixStats_1.18.0    
## [127] shiny_1.7.1                   ggbeeswarm_0.6.0             
## [129] restfulr_0.0.13

References

Aran, D., A. P. Looney, L. Liu, E. Wu, V. Fong, A. Hsu, S. Chak, et al. 2019. “Reference-based analysis of lung single-cell sequencing reveals a transitional profibrotic macrophage.” Nat. Immunol. 20 (2): 163–72.

Grun, D., M. J. Muraro, J. C. Boisset, K. Wiebrands, A. Lyubimova, G. Dharmadhikari, M. van den Born, et al. 2016. “De Novo Prediction of Stem Cell Identity using Single-Cell Transcriptome Data.” Cell Stem Cell 19 (2): 266–77.

La Manno, G., D. Gyllborg, S. Codeluppi, K. Nishimura, C. Salto, A. Zeisel, L. E. Borm, et al. 2016. “Molecular Diversity of Midbrain Development in Mouse, Human, and Stem Cells.” Cell 167 (2): 566–80.

Mabbott, Neil A., J. K. Baillie, Helen Brown, Tom C. Freeman, and David A. Hume. 2013. “An expression atlas of human primary cells: Inference of gene function from coexpression networks.” BMC Genomics 14. https://doi.org/10.1186/1471-2164-14-632.

Muraro, M. J., G. Dharmadhikari, D. Grun, N. Groen, T. Dielen, E. Jansen, L. van Gurp, et al. 2016. “A Single-Cell Transcriptome Atlas of the Human Pancreas.” Cell Syst 3 (4): 385–94.