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.5.1 Patched (2025-08-23 r88802)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 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] cowplot_1.2.0 Rtsne_0.17
## [3] rsvd_1.0.5 scater_1.37.0
## [5] ggplot2_4.0.0 scuttle_1.19.0
## [7] SingleCellExperiment_1.31.1 SummarizedExperiment_1.39.2
## [9] Biobase_2.69.1 GenomicRanges_1.61.5
## [11] Seqinfo_0.99.2 IRanges_2.43.5
## [13] S4Vectors_0.47.4 BiocGenerics_0.55.1
## [15] generics_0.1.4 MatrixGenerics_1.21.0
## [17] matrixStats_1.5.0 scds_1.25.0
## [19] BiocStyle_2.37.1
##
## loaded via a namespace (and not attached):
## [1] beeswarm_0.4.0 gtable_0.3.6 xfun_0.53
## [4] bslib_0.9.0 ggrepel_0.9.6 lattice_0.22-7
## [7] vctrs_0.6.5 tools_4.5.1 parallel_4.5.1
## [10] tibble_3.3.0 pkgconfig_2.0.3 BiocNeighbors_2.3.1
## [13] Matrix_1.7-4 data.table_1.17.8 RColorBrewer_1.1-3
## [16] S7_0.2.0 lifecycle_1.0.4 compiler_4.5.1
## [19] farver_2.1.2 tinytex_0.57 codetools_0.2-20
## [22] vipor_0.4.7 htmltools_0.5.8.1 sass_0.4.10
## [25] yaml_2.3.10 pillar_1.11.1 crayon_1.5.3
## [28] jquerylib_0.1.4 BiocParallel_1.43.4 DelayedArray_0.35.3
## [31] cachem_1.1.0 magick_2.9.0 viridis_0.6.5
## [34] abind_1.4-8 tidyselect_1.2.1 digest_0.6.37
## [37] BiocSingular_1.25.0 dplyr_1.1.4 bookdown_0.45
## [40] labeling_0.4.3 fastmap_1.2.0 grid_4.5.1
## [43] cli_3.6.5 SparseArray_1.9.1 magrittr_2.0.4
## [46] S4Arrays_1.9.1 dichromat_2.0-0.1 withr_3.0.2
## [49] scales_1.4.0 xgboost_1.7.11.1 ggbeeswarm_0.7.2
## [52] rmarkdown_2.30 XVector_0.49.1 gridExtra_2.3
## [55] ScaledMatrix_1.17.0 beachmat_2.25.5 evaluate_1.0.5
## [58] knitr_1.50 viridisLite_0.4.2 irlba_2.3.5.1
## [61] rlang_1.1.6 Rcpp_1.1.0 glue_1.8.0
## [64] BiocManager_1.30.26 pROC_1.19.0.1 jsonlite_2.0.0
## [67] R6_2.6.1