As the number of cells that are routinely interrogated in single-cell experiments increases rapidly and can easily reach hundreds of thousands, the computational resource requirements also grow larger. Several approaches have been proposed to either subsample or aggregate cells in order to reduce the size of the data and enable the application of standard analysis procedures. One such approach is geometric sketching - subsampling in a density-aware manner in such a way that densely populated regions of the gene expression space are subsampled more aggressively, while a larger fraction of cells are retained in sparsely populated regions. In addition to reducing the size of the data set, this often also increases the relative representation of rare cell types in the subsampled data set.
Several tools have been developed for applying sketching to
single-cell (or other) data sets, but not all of them are easily
applicable from R. The sketchR
package implements an
R/Bioconductor interface to some of the most popular python packages for
geometric sketching, allowing them to be directly incorporated into
Bioconductor-based single-cell analysis workflows. The interaction with
python is managed using the basilisk
package, which
automatically takes care of generating and activating a suitable conda
environment with the required packages.
This vignette showcases the main functionalities of the
sketchR
package, and illustrates how it can be used to
generate a subsample of a dataset using the geometric
sketching/subsampling algorithms and implementations proposed by Hie et al. (2019) and Song et al. (2022), as well as create a set of
diagnostic plots.
sketchR
can be installed from Bioconductor using the
following code:
if (!require("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("sketchR")
We start by loading the required packages and preparing an example data set.
suppressPackageStartupMessages({
library(sketchR)
library(TENxPBMCData)
library(scuttle)
library(scran)
library(scater)
library(SingleR)
library(celldex)
library(cowplot)
library(SummarizedExperiment)
library(SingleCellExperiment)
library(beachmat.hdf5)
})
We will use the PBMC3k data set from the TENxPBMCData Bioconductor package for illustration. The chunk below prepares the data set by calculating log-transformed normalized counts, finding highly variable genes, performing dimensionality reduction and predicting cell type labels using the SingleR package.
## Load data
pbmc3k <- TENxPBMCData::TENxPBMCData(dataset = "pbmc3k")
#> see ?TENxPBMCData and browseVignettes('TENxPBMCData') for documentation
#> loading from cache
## Set row and column names
colnames(pbmc3k) <- paste0("Cell", seq_len(ncol(pbmc3k)))
rownames(pbmc3k) <- scuttle::uniquifyFeatureNames(
ID = SummarizedExperiment::rowData(pbmc3k)$ENSEMBL_ID,
names = SummarizedExperiment::rowData(pbmc3k)$Symbol_TENx
)
## Normalize and log-transform counts
pbmc3k <- scuttle::logNormCounts(pbmc3k)
## Find highly variable genes
dec <- scran::modelGeneVar(pbmc3k)
top.hvgs <- scran::getTopHVGs(dec, n = 2000)
## Perform dimensionality reduction
set.seed(100)
pbmc3k <- scater::runPCA(pbmc3k, subset_row = top.hvgs)
pbmc3k <- scater::runTSNE(pbmc3k, dimred = "PCA")
## Predict cell type labels
ref_monaco <- celldex::MonacoImmuneData()
pred_monaco_main <- SingleR::SingleR(test = pbmc3k, ref = ref_monaco,
labels = ref_monaco$label.main)
pbmc3k$labels_main <- pred_monaco_main$labels
dim(pbmc3k)
#> [1] 32738 2700
The geosketch()
function performs geometric sketching by
calling the geosketch
python package. The
output is a vector of indices that can be used to subset the full
dataset. The provided seed will be propagated to the python code to
achieve reproducibility.
idx800gs <- geosketch(SingleCellExperiment::reducedDim(pbmc3k, "PCA"),
N = 800, seed = 123)
head(idx800gs)
#> [1] 3 6 9 11 15 16
length(idx800gs)
#> [1] 800
Similarly, the scsampler()
function calls the
scSampler
python package
to perform subsampling.
idx800scs <- scsampler(SingleCellExperiment::reducedDim(pbmc3k, "PCA"),
N = 800, seed = 123)
head(idx800scs)
#> [1] 1079 644 1494 1278 5 391
length(idx800scs)
#> [1] 800
To illustrate the result of the subsampling, we plot the tSNE representation of the original data as well as the two subsets (using the tSNE coordinates derived from the full dataset).
cowplot::plot_grid(
scater::plotTSNE(pbmc3k, colour_by = "labels_main"),
scater::plotTSNE(pbmc3k[, idx800gs], colour_by = "labels_main"),
scater::plotTSNE(pbmc3k[, idx800scs], colour_by = "labels_main")
)
We can also illustrate the relative abundance of each cell type in the full data and in the subsets, respectively.
compareCompositionPlot(SummarizedExperiment::colData(pbmc3k),
idx = list(geosketch = idx800gs,
scSampler = idx800scs),
column = "labels_main")
sketchR
provides a convenient function to plot the
Hausdorff distance between the full dataset and the subsample, for a
range of sketch sizes (to make this plot reproducible, we use
set.seed
before the call).
set.seed(123)
hausdorffDistPlot(mat = SingleCellExperiment::reducedDim(pbmc3k, "PCA"),
Nvec = c(400, 800, 2000),
Nrep = 3, methods = c("geosketch", "scsampler", "uniform"))
sessionInfo()
#> R Under development (unstable) (2024-10-21 r87258)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.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] beachmat.hdf5_1.3.3 cowplot_1.1.3
#> [3] celldex_1.15.0 SingleR_2.7.6
#> [5] scater_1.33.4 ggplot2_3.5.1
#> [7] scran_1.33.2 scuttle_1.15.4
#> [9] TENxPBMCData_1.23.0 HDF5Array_1.33.8
#> [11] rhdf5_2.49.0 DelayedArray_0.31.14
#> [13] SparseArray_1.5.45 S4Arrays_1.5.11
#> [15] abind_1.4-8 Matrix_1.7-1
#> [17] SingleCellExperiment_1.27.2 SummarizedExperiment_1.35.5
#> [19] Biobase_2.65.1 GenomicRanges_1.57.2
#> [21] GenomeInfoDb_1.41.2 IRanges_2.39.2
#> [23] S4Vectors_0.43.2 BiocGenerics_0.51.3
#> [25] MatrixGenerics_1.17.1 matrixStats_1.4.1
#> [27] sketchR_1.1.3 BiocStyle_2.33.1
#>
#> loaded via a namespace (and not attached):
#> [1] jsonlite_1.8.9 magrittr_2.0.3
#> [3] ggbeeswarm_0.7.2 gypsum_1.1.6
#> [5] farver_2.1.2 rmarkdown_2.28
#> [7] zlibbioc_1.51.2 vctrs_0.6.5
#> [9] memoise_2.0.1 DelayedMatrixStats_1.27.3
#> [11] htmltools_0.5.8.1 AnnotationHub_3.13.3
#> [13] curl_5.2.3 BiocNeighbors_1.99.3
#> [15] Rhdf5lib_1.27.0 sass_0.4.9
#> [17] alabaster.base_1.5.10 bslib_0.8.0
#> [19] basilisk_1.17.2 httr2_1.0.5
#> [21] cachem_1.1.0 igraph_2.1.1
#> [23] mime_0.12 lifecycle_1.0.4
#> [25] pkgconfig_2.0.3 rsvd_1.0.5
#> [27] R6_2.5.1 fastmap_1.2.0
#> [29] GenomeInfoDbData_1.2.13 digest_0.6.37
#> [31] colorspace_2.1-1 AnnotationDbi_1.67.0
#> [33] dqrng_0.4.1 irlba_2.3.5.1
#> [35] ExperimentHub_2.13.1 RSQLite_2.3.7
#> [37] beachmat_2.21.9 labeling_0.4.3
#> [39] filelock_1.0.3 fansi_1.0.6
#> [41] httr_1.4.7 compiler_4.5.0
#> [43] bit64_4.5.2 withr_3.0.1
#> [45] BiocParallel_1.39.0 viridis_0.6.5
#> [47] DBI_1.2.3 highr_0.11
#> [49] alabaster.ranges_1.5.2 alabaster.schemas_1.5.0
#> [51] rappdirs_0.3.3 bluster_1.15.1
#> [53] tools_4.5.0 vipor_0.4.7
#> [55] beeswarm_0.4.0 glue_1.8.0
#> [57] rhdf5filters_1.17.0 grid_4.5.0
#> [59] Rtsne_0.17 cluster_2.1.6
#> [61] generics_0.1.3 gtable_0.3.6
#> [63] BiocSingular_1.21.4 ScaledMatrix_1.13.0
#> [65] metapod_1.13.0 utf8_1.2.4
#> [67] XVector_0.45.0 ggrepel_0.9.6
#> [69] BiocVersion_3.20.0 pillar_1.9.0
#> [71] limma_3.61.12 dplyr_1.1.4
#> [73] BiocFileCache_2.13.2 lattice_0.22-6
#> [75] bit_4.5.0 tidyselect_1.2.1
#> [77] locfit_1.5-9.10 Biostrings_2.73.2
#> [79] knitr_1.48 gridExtra_2.3
#> [81] edgeR_4.3.21 xfun_0.48
#> [83] statmod_1.5.0 UCSC.utils_1.1.0
#> [85] yaml_2.3.10 evaluate_1.0.1
#> [87] codetools_0.2-20 tibble_3.2.1
#> [89] alabaster.matrix_1.5.10 BiocManager_1.30.25
#> [91] cli_3.6.3 reticulate_1.39.0
#> [93] munsell_0.5.1 jquerylib_0.1.4
#> [95] Rcpp_1.0.13 dir.expiry_1.13.1
#> [97] dbplyr_2.5.0 png_0.1-8
#> [99] parallel_4.5.0 blob_1.2.4
#> [101] basilisk.utils_1.17.3 sparseMatrixStats_1.17.2
#> [103] alabaster.se_1.5.4 viridisLite_0.4.2
#> [105] scales_1.3.0 purrr_1.0.2
#> [107] crayon_1.5.3 rlang_1.1.4
#> [109] KEGGREST_1.45.1