Contents

1 Introduction

In this vignette, we provide an overview of the basic functionality and usage of the scds package, which interfaces with SingleCellExperiment objects.

2 Installation

Install the scds package using Bioconductor:

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("scds", version = "3.9")

Or from github:

library(devtools)
devtools::install_github('kostkalab/scds')

3 Quick start

scds takes as input a SingleCellExperiment object (see here SingleCellExperiment), where raw counts are stored in a counts assay, i.e. assay(sce,"counts"). An example dataset created by sub-sampling the cell-hashing cell-lines data set (see https://satijalab.org/seurat/hashing_vignette.html) is included with the package and accessible via data("sce").Note that scds is designed to workd with larger datasets, but for the purposes of this vignette, we work with a smaller example dataset. We apply scds to this data and compare/visualize reasults:

3.1 Example data set

Get example data set provided with the package.

library(scds)
library(scater)
library(rsvd)
library(Rtsne)
library(cowplot)
set.seed(30519)
data("sce_chcl")
sce = sce_chcl #- less typing
dim(sce)
## [1] 2000 2000

We see it contains 2,000 genes and 2,000 cells, 216 of which are identified as doublets:

table(sce$hto_classification_global)
## 
##  Doublet Negative  Singlet 
##      216       83     1701

We can visualize cells/doublets after projecting into two dimensions:

logcounts(sce) = log1p(counts(sce))
vrs            = apply(logcounts(sce),1,var)
pc             = rpca(t(logcounts(sce)[order(vrs,decreasing=TRUE)[1:100],]))
ts             = Rtsne(pc$x[,1:10],verb=FALSE)

reducedDim(sce,"tsne") = ts$Y; rm(ts,vrs,pc)
plotReducedDim(sce,"tsne",color_by="hto_classification_global")

3.2 Computational doublet annotation

We now run the scds doublet annotation approaches. Briefly, we identify doublets in two complementary ways: cxds is based on co-expression of gene pairs and works with absence/presence calls only, while bcds uses the full count information and a binary classification approach using artificially generated doublets. cxds_bcds_hybrid combines both approaches, for more details please consult (this manuscript). Each of the three methods returns a doublet score, with higher scores indicating more “doublet-like” barcodes.

#- Annotate doublet using co-expression based doublet scoring:
sce = cxds(sce,retRes = TRUE)
sce = bcds(sce,retRes = TRUE,verb=TRUE)
sce = cxds_bcds_hybrid(sce)
par(mfcol=c(1,3))
boxplot(sce$cxds_score   ~ sce$doublet_true_labels, main="cxds")
boxplot(sce$bcds_score   ~ sce$doublet_true_labels, main="bcds")
boxplot(sce$hybrid_score ~ sce$doublet_true_labels, main="hybrid")

3.3 Visualizing gene pairs

For cxds we can identify and visualize gene pairs driving doublet annoataions, with the expectation that the two genes in a pair might mark different types of cells (see manuscript). In the following we look at the top three pairs, each gene pair is a row in the plot below:

scds =
top3 = metadata(sce)$cxds$topPairs[1:3,]
rs   = rownames(sce)
hb   = rowData(sce)$cxds_hvg_bool
ho   = rowData(sce)$cxds_hvg_ordr[hb]
hgs  = rs[ho]

l1 =  ggdraw() + draw_text("Pair 1", x = 0.5, y = 0.5)
p1 = plotReducedDim(sce,"tsne",color_by=hgs[top3[1,1]])
p2 = plotReducedDim(sce,"tsne",color_by=hgs[top3[1,2]])

l2 =  ggdraw() + draw_text("Pair 2", x = 0.5, y = 0.5)
p3 = plotReducedDim(sce,"tsne",color_by=hgs[top3[2,1]])
p4 = plotReducedDim(sce,"tsne",color_by=hgs[top3[2,2]])

l3 = ggdraw() + draw_text("Pair 3", x = 0.5, y = 0.5)
p5 = plotReducedDim(sce,"tsne",color_by=hgs[top3[3,1]])
p6 = plotReducedDim(sce,"tsne",color_by=hgs[top3[3,2]])

plot_grid(l1,p1,p2,l2,p3,p4,l3,p5,p6,ncol=3, rel_widths = c(1,2,2))

4 Session Info

sessionInfo()
## R version 4.3.1 (2023-06-16)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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] Matrix_1.6-1.1              cowplot_1.1.1              
##  [3] Rtsne_0.16                  rsvd_1.0.5                 
##  [5] scater_1.30.0               ggplot2_3.4.4              
##  [7] scuttle_1.12.0              SingleCellExperiment_1.24.0
##  [9] SummarizedExperiment_1.32.0 Biobase_2.62.0             
## [11] GenomicRanges_1.54.0        GenomeInfoDb_1.38.0        
## [13] IRanges_2.36.0              S4Vectors_0.40.0           
## [15] BiocGenerics_0.48.0         MatrixGenerics_1.14.0      
## [17] matrixStats_1.0.0           scds_1.18.0                
## [19] BiocStyle_2.30.0           
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.0          viridisLite_0.4.2        
##  [3] farver_2.1.1              dplyr_1.1.3              
##  [5] vipor_0.4.5               viridis_0.6.4            
##  [7] bitops_1.0-7              fastmap_1.1.1            
##  [9] RCurl_1.98-1.12           pROC_1.18.4              
## [11] digest_0.6.33             lifecycle_1.0.3          
## [13] magrittr_2.0.3            compiler_4.3.1           
## [15] rlang_1.1.1               sass_0.4.7               
## [17] tools_4.3.1               utf8_1.2.4               
## [19] yaml_2.3.7                data.table_1.14.8        
## [21] knitr_1.44                labeling_0.4.3           
## [23] S4Arrays_1.2.0            xgboost_1.7.5.1          
## [25] DelayedArray_0.28.0       plyr_1.8.9               
## [27] abind_1.4-5               BiocParallel_1.36.0      
## [29] withr_2.5.1               grid_4.3.1               
## [31] fansi_1.0.5               beachmat_2.18.0          
## [33] colorspace_2.1-0          scales_1.2.1             
## [35] cli_3.6.1                 rmarkdown_2.25           
## [37] crayon_1.5.2              generics_0.1.3           
## [39] DelayedMatrixStats_1.24.0 ggbeeswarm_0.7.2         
## [41] cachem_1.0.8              zlibbioc_1.48.0          
## [43] parallel_4.3.1            BiocManager_1.30.22      
## [45] XVector_0.42.0            vctrs_0.6.4              
## [47] jsonlite_1.8.7            bookdown_0.36            
## [49] BiocSingular_1.18.0       BiocNeighbors_1.20.0     
## [51] ggrepel_0.9.4             irlba_2.3.5.1            
## [53] beeswarm_0.4.0            magick_2.8.1             
## [55] jquerylib_0.1.4           glue_1.6.2               
## [57] codetools_0.2-19          gtable_0.3.4             
## [59] ScaledMatrix_1.10.0       munsell_0.5.0            
## [61] tibble_3.2.1              pillar_1.9.0             
## [63] htmltools_0.5.6.1         GenomeInfoDbData_1.2.11  
## [65] R6_2.5.1                  sparseMatrixStats_1.14.0 
## [67] evaluate_0.22             lattice_0.22-5           
## [69] bslib_0.5.1               Rcpp_1.0.11              
## [71] gridExtra_2.3             SparseArray_1.2.0        
## [73] xfun_0.40                 pkgconfig_2.0.3