In this vignette, we provide an overview of the basic functionality and usage of the scds
package, which interfaces with SingleCellExperiment
objects.
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')
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:
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")
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")
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))
sessionInfo()
## R version 4.3.0 RC (2023-04-13 r84269)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.17-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] cowplot_1.1.1 Rtsne_0.16
## [3] rsvd_1.0.5 scater_1.28.0
## [5] ggplot2_3.4.2 scuttle_1.10.0
## [7] SingleCellExperiment_1.22.0 SummarizedExperiment_1.30.0
## [9] Biobase_2.60.0 GenomicRanges_1.52.0
## [11] GenomeInfoDb_1.36.0 IRanges_2.34.0
## [13] S4Vectors_0.38.0 BiocGenerics_0.46.0
## [15] MatrixGenerics_1.12.0 matrixStats_0.63.0
## [17] scds_1.16.0 BiocStyle_2.28.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.0 viridisLite_0.4.1
## [3] dplyr_1.1.2 vipor_0.4.5
## [5] farver_2.1.1 viridis_0.6.2
## [7] bitops_1.0-7 fastmap_1.1.1
## [9] RCurl_1.98-1.12 pROC_1.18.0
## [11] digest_0.6.31 lifecycle_1.0.3
## [13] magrittr_2.0.3 compiler_4.3.0
## [15] rlang_1.1.0 sass_0.4.5
## [17] tools_4.3.0 utf8_1.2.3
## [19] yaml_2.3.7 data.table_1.14.8
## [21] knitr_1.42 labeling_0.4.2
## [23] xgboost_1.7.5.1 DelayedArray_0.26.0
## [25] plyr_1.8.8 BiocParallel_1.34.0
## [27] withr_2.5.0 grid_4.3.0
## [29] fansi_1.0.4 beachmat_2.16.0
## [31] colorspace_2.1-0 scales_1.2.1
## [33] cli_3.6.1 rmarkdown_2.21
## [35] generics_0.1.3 DelayedMatrixStats_1.22.0
## [37] ggbeeswarm_0.7.1 cachem_1.0.7
## [39] zlibbioc_1.46.0 parallel_4.3.0
## [41] BiocManager_1.30.20 XVector_0.40.0
## [43] vctrs_0.6.2 Matrix_1.5-4
## [45] jsonlite_1.8.4 bookdown_0.33
## [47] BiocSingular_1.16.0 BiocNeighbors_1.18.0
## [49] ggrepel_0.9.3 irlba_2.3.5.1
## [51] beeswarm_0.4.0 magick_2.7.4
## [53] jquerylib_0.1.4 glue_1.6.2
## [55] codetools_0.2-19 gtable_0.3.3
## [57] ScaledMatrix_1.8.0 munsell_0.5.0
## [59] tibble_3.2.1 pillar_1.9.0
## [61] htmltools_0.5.5 GenomeInfoDbData_1.2.10
## [63] R6_2.5.1 sparseMatrixStats_1.12.0
## [65] evaluate_0.20 lattice_0.21-8
## [67] highr_0.10 bslib_0.4.2
## [69] Rcpp_1.0.10 gridExtra_2.3
## [71] xfun_0.39 pkgconfig_2.0.3