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",col="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",col=hgs[top3[1,1]])
p2 = plotReducedDim(sce,"tsne",col=hgs[top3[1,2]])

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

l3 = ggdraw() + draw_text("Pair 3", x = 0.5, y = 0.5)
p5 = plotReducedDim(sce,"tsne",col=hgs[top3[3,1]])
p6 = plotReducedDim(sce,"tsne",col=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.0.0 (2020-04-24)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.11-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.11-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        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] parallel  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] cowplot_1.0.0               Rtsne_0.15                 
##  [3] rsvd_1.0.3                  scater_1.16.0              
##  [5] ggplot2_3.3.0               SingleCellExperiment_1.10.0
##  [7] SummarizedExperiment_1.18.0 DelayedArray_0.14.0        
##  [9] matrixStats_0.56.0          Biobase_2.48.0             
## [11] GenomicRanges_1.40.0        GenomeInfoDb_1.24.0        
## [13] IRanges_2.22.0              S4Vectors_0.26.0           
## [15] BiocGenerics_0.34.0         scds_1.4.0                 
## [17] BiocStyle_2.16.0           
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_1.0.4.6              lattice_0.20-41          
##  [3] assertthat_0.2.1          digest_0.6.25            
##  [5] R6_2.4.1                  plyr_1.8.6               
##  [7] evaluate_0.14             pillar_1.4.3             
##  [9] zlibbioc_1.34.0           rlang_0.4.5              
## [11] data.table_1.12.8         irlba_2.3.3              
## [13] magick_2.3                Matrix_1.2-18            
## [15] rmarkdown_2.1             labeling_0.3             
## [17] BiocNeighbors_1.6.0       BiocParallel_1.22.0      
## [19] stringr_1.4.0             RCurl_1.98-1.2           
## [21] munsell_0.5.0             vipor_0.4.5              
## [23] compiler_4.0.0            BiocSingular_1.4.0       
## [25] xfun_0.13                 pkgconfig_2.0.3          
## [27] ggbeeswarm_0.6.0          htmltools_0.4.0          
## [29] tidyselect_1.0.0          gridExtra_2.3            
## [31] tibble_3.0.1              GenomeInfoDbData_1.2.3   
## [33] bookdown_0.18             viridisLite_0.3.0        
## [35] crayon_1.3.4              dplyr_0.8.5              
## [37] withr_2.2.0               bitops_1.0-6             
## [39] grid_4.0.0                gtable_0.3.0             
## [41] lifecycle_0.2.0           magrittr_1.5             
## [43] pROC_1.16.2               scales_1.1.0             
## [45] stringi_1.4.6             farver_2.0.3             
## [47] XVector_0.28.0            viridis_0.5.1            
## [49] DelayedMatrixStats_1.10.0 ellipsis_0.3.0           
## [51] vctrs_0.2.4               xgboost_1.0.0.2          
## [53] tools_4.0.0               glue_1.4.0               
## [55] beeswarm_0.2.3            purrr_0.3.4              
## [57] yaml_2.2.1                colorspace_1.4-1         
## [59] BiocManager_1.30.10       knitr_1.28