scCB2 1.4.0
Droplet-based single-cell RNA sequencing (scRNA-seq) is a powerful and widely-used approach for profiling genome-wide gene expression
in individual cells. Current commercial droplet-based technologies such as 10X Genomics utilize gel beads, each containing
oligonucleotide indexes made up of bead-specific barcodes combined with unique molecular identifiers (UMIs) and oligo-dT tags
to prime polyadenylated RNA. Single cells of interest are combined with reagents in one channel of a microfluidic chip,
and gel beads in another, to form gel-beads in emulsion, or GEMs. Oligonucleotide indexes bind polyadenylated RNA within each GEM
reaction vesicle before gel beads are dissolved releasing the bound oligos into solution for reverse transcription.
By design, each resulting cDNA molecule contains a UMI and a GEM-specific barcode that, ideally,
tags mRNA from an individual cell, but this is often not the case in practice. To distinguish true cells from background barcodes in
droplet-based single cell RNA-seq experiments, we introduce CB2 and scCB2
, its corresponding R package.
CB2 extends the EmptyDrops approach by introducing a clustering step that groups similar barcodes and
then conducts a statistical test to identify groups with expression distributions that vary from the background.
While advantages are expected in many settings, users will benefit from
noting that CB2 does not test for doublets or multiplets and, consequently, some of the high count identifications
may consist of two or more cells. Methods for identifying multiplets may prove useful after applying CB2.
It is also important to note that any method for distinguishing cells from background barcodes is technically correct in
identifying low-quality cells given that damaged cells exhibit expression profiles that differ from the background.
Specifically, mitochondrial gene expression is often high in damaged cells.
Such cells are typically not of interest in downstream analysis and should therefore be removed.
The GetCellMat function in scCB2
may be used toward this end.
Install from Bioconductor:
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("scCB2")
QuickCB2
is an all-in-one function to apply CB2 on 10x Cell Ranger raw data and get a matrix of real cells identified by CB2 under default settings. By specifying AsSeurat = TRUE
, a Seurat object is returned so that users can directly apply the Seurat pipeline for downstream analyses.
Usage:
library(scCB2)
# If raw data has three separate files within one directory
# and you want to control FDR at the default 1%:
RealCell <- QuickCB2(dir = "/path/to/raw/data/directory")
# If raw data is in HDF5 format and
# you'd like a Seurat object under default FDR threshold:
RealCell_S <- QuickCB2(h5file = "/path/to/raw/data/HDF5",
AsSeurat = TRUE)
An example illustrating how it works and what the final output looks like can be found at the end of Detailed Steps.
The computational speed is related to the size and structure of input datasets. As a reference, scCB2 running on a medium-size dataset (around 30,000 genes and 8,000 cells after cell calling) using 6 cores takes less than 10 minutes.
For users who would like to save the CB2 output cell matrix to 10x format (e.g. “matrix.mtx”, “barcodes.tsv” and “genes.tsv”), there are existing packages to help. For example in package DropletUtils
:
DropletUtils::write10xCounts(path = "/path/to/save/data",
x = RealCell)
Currently, the most widely-used droplet-based protocol is 10x Chromium. Our package provides functions to directly read 10x Cell Ranger output files and generate a feature-by-barcode count matrix that may be read into R. Public 10x datasets can be found here.
Our package contains a small subset of 10x data, mbrainSub
, corresponding to the first 50,000 barcodes of 1k Brain Cells from an E18 Mouse.
We first generate 10x output files of mbrainSub
, then read it using our built-in functions.
library(scCB2)
library(SummarizedExperiment)
data(mbrainSub)
data.dir <- file.path(tempdir(),"CB2_example")
DropletUtils::write10xCounts(data.dir,
mbrainSub,
version = "3")
list.files(data.dir)
## [1] "barcodes.tsv.gz" "features.tsv.gz" "matrix.mtx.gz"
For Cell Ranger version <3, the raw data from 10x Cell Ranger output contains “barcodes.tsv”, “genes.tsv” and “matrix.mtx”. For Cell Ranger version >=3, the output files are “barcodes.tsv.gz”, “features.tsv.gz” and “matrix.mtx.gz”. We now read these files back into R and compare with original data matrix.
mbrainSub_2 <- Read10xRaw(data.dir)
identical(mbrainSub, mbrainSub_2)
## [1] TRUE
If raw data is not from the 10x Chromium pipeline, a user may manually create the feature-by-barcode count matrix with rows representing genes and columns representing barcodes. Gene and barcode IDs should be unique. The format of the count matrix can be either a sparse matrix or standard matrix.
The main function CB2FindCell
takes a raw count matrix as input and returns real cells, test statistics, and p-values. Now we apply CB2FindCell
on mbrainSub
, controlling FDR at 0.01 level (Default), assuming all barcodes with total count less than or equal to 100 are background empty droplets (Default), using 2 cores parallel computation (Default: number of total cores in the machine minus 2). For detailed information, see ?CB2FindCell
.
CBOut <- CB2FindCell(mbrainSub, Ncores = 2)
## Time difference of 1.162474 mins
str(assay(CBOut)) # cell matrix
## Formal class 'dgCMatrix' [package "Matrix"] with 6 slots
## ..@ i : int [1:273471] 62 469 538 663 758 903 1251 1266 1396 1532 ...
## ..@ p : int [1:173] 0 219 2503 4608 6288 6512 6796 7194 7335 7789 ...
## ..@ Dim : int [1:2] 27998 172
## ..@ Dimnames:List of 2
## .. ..$ : chr [1:27998] "Xkr4" "Gm1992" "Gm37381" "Rp1" ...
## .. ..$ : chr [1:172] "AAACCTGGTGCTCTTC-1" "AAACGGGAGCCACGTC-1" "AAACGGGAGCGAGAAA-1" "AAACGGGCAGTTTACG-1" ...
## ..@ x : num [1:273471] 1 1 1 1 1 1 1 1 1 1 ...
## ..@ factors : list()
str(metadata(CBOut)) # test statistics, p-values, etc
## List of 4
## $ ClusterStat:'data.frame': 5 obs. of 4 variables:
## ..$ corr : num [1:5] 0.6715 0.6773 0.0509 0.38 0.5547
## ..$ count: num [1:5] 260 315 271 113 4679
## ..$ pval : num [1:5] 0.273726 0.081918 0.000999 0.005994 0.000999
## ..$ padj : num [1:5] 0.27373 0.1024 0.0025 0.00999 0.0025
## $ Cluster :List of 5
## ..$ : chr [1:85] "AACACGTAGGCAAAGA-1" "ACAGCCGCAAACGCGA-1" "AACACGTGTATAGGTA-1" "AACTCAGTCTCGGACG-1" ...
## ..$ : chr [1:41] "AACTTTCAGGCCCGTT-1" "AACTCTTTCGCCGTGA-1" "AACTTTCGTACCGTTA-1" "AACCATGGTTGCCTCT-1" ...
## ..$ : chr [1:3] "AAGGCAGGTCTTGCGG-1" "ACATCAGCACATTAGC-1" "AAGACCTAGACTGGGT-1"
## ..$ : chr [1:6] "AATCCAGAGATCTGAA-1" "AACCATGTCTCGCATC-1" "AACACGTTCTCTAGGA-1" "AACTCCCGTGCTTCTC-1" ...
## ..$ : chr [1:23] "ACCAGTAAGTGCGATG-1" "AAGTCTGAGGACAGCT-1" "ACATACGCAAGCCGTC-1" "AAGGCAGGTCTAGAGG-1" ...
## $ BarcodeStat:'data.frame': 203 obs. of 4 variables:
## ..$ logLH: num [1:203] -762 -770 -886 -454 -9992 ...
## ..$ count: num [1:203] 337 304 312 207 5717 ...
## ..$ pval : num [1:203] 0.9428 0.4289 0.0011 1 0.0001 ...
## ..$ padj : num [1:203] 1 0.66828 0.002069 1 0.000203 ...
## $ background : Named num [1:11121] 8 3 2 0 2 2 0 1 3 1 ...
## ..- attr(*, "names")= chr [1:11121] "Mrpl15" "Lypla1" "Tcea1" "Rgs20" ...
If readers are not interested in the output testing information, GetCellMat
can extract the real cell matrix directly from CB2FindCell
output. It also provides a filtering option to remove broken cells based on the proportion of mitochondrial gene expressions. Now we apply GetCellMat
on CBOut
, filtering out cells whose mitochondrial proportions are greater than 0.25 (Default: 1, No filtering).
RealCell <- GetCellMat(CBOut, MTfilter = 0.25)
str(RealCell)
## Formal class 'dgCMatrix' [package "Matrix"] with 6 slots
## ..@ i : int [1:272453] 62 469 538 663 758 903 1251 1266 1396 1532 ...
## ..@ p : int [1:170] 0 219 2503 4608 6288 6512 6796 7194 7335 7789 ...
## ..@ Dim : int [1:2] 27998 169
## ..@ Dimnames:List of 2
## .. ..$ : chr [1:27998] "Xkr4" "Gm1992" "Gm37381" "Rp1" ...
## .. ..$ : chr [1:169] "AAACCTGGTGCTCTTC-1" "AAACGGGAGCCACGTC-1" "AAACGGGAGCGAGAAA-1" "AAACGGGCAGTTTACG-1" ...
## ..@ x : num [1:272453] 1 1 1 1 1 1 1 1 1 1 ...
## ..@ factors : list()
After CB2
pre-processing, the real cell matrix is still in matrix format, so it can be directly used in downstream statistical analyses. For example, if we want to use the Seurat pipeline, we can easily create a Seurat object using
SeuratObj <- Seurat::CreateSeuratObject(counts = RealCell,
project = "mbrain_example")
SeuratObj
## An object of class Seurat
## 27998 features across 169 samples within 1 assay
## Active assay: RNA (27998 features, 0 variable features)
Under default parameters, we can directly use the all-in-one function QuickCB2
to get the real cell matrix from 10x raw data.
RealCell_Quick <- QuickCB2(dir = data.dir, Ncores = 2)
## Time difference of 1.395012 mins
str(RealCell_Quick)
## Formal class 'dgCMatrix' [package "Matrix"] with 6 slots
## ..@ i : int [1:274054] 62 469 538 663 758 903 1251 1266 1396 1532 ...
## ..@ p : int [1:175] 0 219 2503 4608 6288 6512 6796 7194 7335 7789 ...
## ..@ Dim : int [1:2] 27998 174
## ..@ Dimnames:List of 2
## .. ..$ : chr [1:27998] "Xkr4" "Gm1992" "Gm37381" "Rp1" ...
## .. ..$ : chr [1:174] "AAACCTGGTGCTCTTC-1" "AAACGGGAGCCACGTC-1" "AAACGGGAGCGAGAAA-1" "AAACGGGCAGTTTACG-1" ...
## ..@ x : num [1:274054] 1 1 1 1 1 1 1 1 1 1 ...
## ..@ factors : list()
Now it’s ready for downstream analysis such as normalization and clustering. Example Seurat tutorial: https://satijalab.org/seurat/v3.1/pbmc3k_tutorial.html
sessionInfo()
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
##
## 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
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] SummarizedExperiment_1.24.0 Biobase_2.54.0
## [3] GenomicRanges_1.46.0 GenomeInfoDb_1.30.0
## [5] IRanges_2.28.0 S4Vectors_0.32.0
## [7] BiocGenerics_0.40.0 MatrixGenerics_1.6.0
## [9] matrixStats_0.61.0 scCB2_1.4.0
## [11] knitr_1.36 BiocStyle_2.22.0
##
## loaded via a namespace (and not attached):
## [1] plyr_1.8.6 igraph_1.2.7
## [3] lazyeval_0.2.2 splines_4.1.1
## [5] BiocParallel_1.28.0 listenv_0.8.0
## [7] scattermore_0.7 ggplot2_3.3.5
## [9] digest_0.6.28 foreach_1.5.1
## [11] htmltools_0.5.2 fansi_0.5.0
## [13] magrittr_2.0.1 tensor_1.5
## [15] cluster_2.1.2 doParallel_1.0.16
## [17] ROCR_1.0-11 limma_3.50.0
## [19] globals_0.14.0 R.utils_2.11.0
## [21] spatstat.sparse_2.0-0 colorspace_2.0-2
## [23] ggrepel_0.9.1 xfun_0.27
## [25] dplyr_1.0.7 crayon_1.4.1
## [27] RCurl_1.98-1.5 jsonlite_1.7.2
## [29] spatstat.data_2.1-0 survival_3.2-13
## [31] zoo_1.8-9 iterators_1.0.13
## [33] glue_1.4.2 polyclip_1.10-0
## [35] gtable_0.3.0 zlibbioc_1.40.0
## [37] XVector_0.34.0 leiden_0.3.9
## [39] DelayedArray_0.20.0 DropletUtils_1.14.0
## [41] Rhdf5lib_1.16.0 future.apply_1.8.1
## [43] SingleCellExperiment_1.16.0 HDF5Array_1.22.0
## [45] abind_1.4-5 scales_1.1.1
## [47] DBI_1.1.1 edgeR_3.36.0
## [49] miniUI_0.1.1.1 Rcpp_1.0.7
## [51] viridisLite_0.4.0 xtable_1.8-4
## [53] reticulate_1.22 spatstat.core_2.3-0
## [55] dqrng_0.3.0 htmlwidgets_1.5.4
## [57] httr_1.4.2 RColorBrewer_1.1-2
## [59] ellipsis_0.3.2 Seurat_4.0.5
## [61] ica_1.0-2 pkgconfig_2.0.3
## [63] R.methodsS3_1.8.1 scuttle_1.4.0
## [65] uwot_0.1.10 deldir_1.0-6
## [67] sass_0.4.0 locfit_1.5-9.4
## [69] utf8_1.2.2 tidyselect_1.1.1
## [71] rlang_0.4.12 reshape2_1.4.4
## [73] later_1.3.0 munsell_0.5.0
## [75] tools_4.1.1 generics_0.1.1
## [77] ggridges_0.5.3 evaluate_0.14
## [79] stringr_1.4.0 fastmap_1.1.0
## [81] goftest_1.2-3 yaml_2.2.1
## [83] fitdistrplus_1.1-6 purrr_0.3.4
## [85] RANN_2.6.1 nlme_3.1-153
## [87] pbapply_1.5-0 future_1.22.1
## [89] sparseMatrixStats_1.6.0 mime_0.12
## [91] R.oo_1.24.0 compiler_4.1.1
## [93] plotly_4.10.0 png_0.1-7
## [95] spatstat.utils_2.2-0 tibble_3.1.5
## [97] bslib_0.3.1 stringi_1.7.5
## [99] lattice_0.20-45 Matrix_1.3-4
## [101] vctrs_0.3.8 pillar_1.6.4
## [103] lifecycle_1.0.1 rhdf5filters_1.6.0
## [105] BiocManager_1.30.16 spatstat.geom_2.3-0
## [107] lmtest_0.9-38 jquerylib_0.1.4
## [109] RcppAnnoy_0.0.19 data.table_1.14.2
## [111] cowplot_1.1.1 bitops_1.0-7
## [113] irlba_2.3.3 httpuv_1.6.3
## [115] patchwork_1.1.1 R6_2.5.1
## [117] bookdown_0.24 promises_1.2.0.1
## [119] KernSmooth_2.23-20 gridExtra_2.3
## [121] parallelly_1.28.1 codetools_0.2-18
## [123] MASS_7.3-54 assertthat_0.2.1
## [125] rhdf5_2.38.0 SeuratObject_4.0.2
## [127] sctransform_0.3.2 GenomeInfoDbData_1.2.7
## [129] mgcv_1.8-38 parallel_4.1.1
## [131] rpart_4.1-15 grid_4.1.1
## [133] beachmat_2.10.0 tidyr_1.1.4
## [135] rmarkdown_2.11 DelayedMatrixStats_1.16.0
## [137] Rtsne_0.15 shiny_1.7.1
Ni, Z., Chen, S., Brown, J., & Kendziorski, C. (2020). CB2 improves power of cell detection in droplet-based single-cell RNA sequencing data. Genome Biology, 21(1), 137. https://doi.org/10.1186/s13059-020-02054-8