Droplet-based microfluidic devices have become widely used to perform single-cell RNA sequencing (scRNA-seq). However, ambient RNA present in the cell suspension can be aberrantly counted along with a cell’s native mRNA and result in cross-contamination of transcripts between different cell populations. DecontX is a Bayesian method to estimate and remove contamination in individual cells. DecontX assumes the observed expression of a cell is a mixture of counts from two multinomial distributions: (1) a distribution of native transcript counts from the cell’s actual population and (2) a distribution of contaminating transcript counts from all other cell populations captured in the assay. Overall, computational decontamination of single cell counts can aid in downstream clustering and visualization.
DecontX Package can be installed from Bioconductor:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("decontX")
Then the package can be loaded in R using the following command:
To see the latest updates and releases or to post a bug, see our GitHub page at https://github.com/campbio/decontX.
DecontX can take either a SingleCellExperiment
object or a counts matrix as input. decontX
will attempt to
convert any input matrix to class dgCMatrix
from package Matrix
before starting the analysis.
To import datasets directly into an SCE object, the singleCellTK package has several importing functions for different preprocessing tools including CellRanger, STARsolo, BUStools, Optimus, DropEST, SEQC, and Alevin/Salmon. For example, the following code can be used as a template to read in the filtered and raw matrices for multiple samples processed with CellRanger:
library(singleCellTK)
sce <- importCellRanger(sampleDirs = c("path/to/sample1/", "path/to/sample2/"))
Within each sample directory, there should be subfolders called
"outs/filtered_feature_bc_matrix/"
or
"outs/raw_feature_bc_matrix/"
with files called
matrix.mtx.gz
, features.tsv.gz
and
barcodes.tsv.gz
. If these files are in different
subdirectories, the importCellRangerV3Sample
function can
be used to import data from a different directory instead.
Optionally, the “raw” or “droplet” matrix can also be easily imported
by setting the dataType
argument to “raw”:
sce.raw <- importCellRanger(sampleDirs = c("path/to/sample1/", "path/to/sample2/"), dataType = "raw")
The raw matrix can be passed to the background
parameter
in decontX
as described below. If using Seurat, go to the
Working with Seurat section for details on how to
convert between SCE and Seurat objects.
We will utilize the 10X PBMC 4K dataset as an example in this vignette. This data can be easily retrieved from the package TENxPBMCData. Make sure the the column names are set before running decontX.
A SingleCellExperiment (SCE) object or a sparse matrix containing the
counts for filtered cells can be passed to decontX via the
x
parameter. The matrix to use in an SCE object can be
specified with the assayName
parameter, which is set to
"counts"
by default. There are two major ways to run
decontX: with and without the raw/droplet matrix containing empty
droplets. Here is an example of running decontX without supplying the
background:
In this scenario, decontX
will estimate the
contamination distribution for each cell cluster based on the profiles
of the other cell clusters in the filtered dataset. The estimated
contamination results can be found in the
colData(sce)$decontX_contamination
and the decontaminated
counts can be accessed with decontXcounts(sce)
.
decontX
will perform heuristic clustering to quickly define
major cell clusters. However if you have your own cell cluster labels,
they can be specified with the z
parameter. These results
will be used throughout the rest of the vignette.
The raw/droplet matrix can be used to empirically estimate the
distribution of ambient RNA, which is especially useful when cells that
contributed to the ambient RNA are not accurately represented in the
filtered count matrix containing the cells. For example, cells that were
removed via flow cytometry or that were more sensitive to lysis during
dissociation may have contributed to the ambient RNA but were not
measured in the filtered/cell matrix. The raw/droplet matrix can be
input as an SCE object or a sparse matrix using the
background
parameter:
Only empty droplets in the background matrix should be used to
estimate the ambient RNA. If any cell ids (i.e. colnames
)
in the raw/droplet matrix supplied to the background
parameter are also found in the filtered counts matrix (x
),
decontX will automatically remove them from the raw matrix. However, if
the cell ids are not available for the input matrices, decontX will
treat the entire background
input as empty droplets. All of
the outputs are the same as when running decontX without setting the
background
parameter.
Note: If the input object is just a matrix and not an SCE object, make sure to save the output into a variable with a different name (e.g.
result <- decontX(mat)
). The result object will be a list with contamination inresult$contamination
and the decontaminated counts inresult$decontXcounts
.
DecontX creates a UMAP which we can use to plot the cluster labels automatically identified in the analysis. Note that the clustering approach used here is designed to find “broad” cell types rather than individual cell subpopulations within a cell type.
##
## Attaching package: 'celda'
## The following object is masked from 'package:S4Vectors':
##
## params
## The following objects are masked from 'package:decontX':
##
## decontX, decontXcounts, decontXcounts<-, plotDecontXContamination,
## plotDecontXMarkerExpression, plotDecontXMarkerPercentage,
## retrieveFeatureIndex, simulateContamination
The percentage of contamination in each cell can be plotting on the UMAP to visualize what what clusters may have higher levels of ambient RNA.
Known marker genes can also be plotted on the UMAP to identify the cell types for each cluster. We will use CD3D and CD3E for T-cells, LYZ, S100A8, and S100A9 for monocytes, CD79A, CD79B, and MS4A1 for B-cells, GNLY for NK-cells, and PPBP for megakaryocytes.
library(scater)
sce <- logNormCounts(sce)
plotDimReduceFeature(as.matrix(logcounts(sce)),
dim1 = umap[, 1],
dim2 = umap[, 2],
features = c("CD3D", "CD3E", "GNLY",
"LYZ", "S100A8", "S100A9",
"CD79A", "CD79B", "MS4A1"),
exactMatch = TRUE)
## Warning in asMethod(object): sparse->dense coercion: allocating vector of size
## 1.1 GiB
The percetage of cells within a cluster that have detectable
expression of marker genes can be displayed in a barplot. Markers for
cell types need to be supplied in a named list. First, the detection of
marker genes in the original counts
assay is shown:
markers <- list(Tcell_Markers = c("CD3E", "CD3D"),
Bcell_Markers = c("CD79A", "CD79B", "MS4A1"),
Monocyte_Markers = c("S100A8", "S100A9", "LYZ"),
NKcell_Markers = "GNLY")
cellTypeMappings <- list(Tcells = 2, Bcells = 5, Monocytes = 1, NKcells = 6)
plotDecontXMarkerPercentage(sce,
markers = markers,
groupClusters = cellTypeMappings,
assayName = "counts")
We can then look to see how much decontX removed aberrant expression
of marker genes in each cell type by changing the assayName
to decontXcounts
:
plotDecontXMarkerPercentage(sce,
markers = markers,
groupClusters = cellTypeMappings,
assayName = "decontXcounts")
Percentages of marker genes detected in other cell types were reduced
or completely removed. For example, the percentage of cells that
expressed Monocyte marker genes was greatly reduced in T-cells, B-cells,
and NK-cells. The original counts and decontamined counts can be plotted
side-by-side by listing multiple assays in the assayName
parameter. This option is only available if the data is stored in
SingleCellExperiment
object.
plotDecontXMarkerPercentage(sce,
markers = markers,
groupClusters = cellTypeMappings,
assayName = c("counts", "decontXcounts"))
Some helpful hints when using
plotDecontXMarkerPercentage
:
groupCluster
parameter, which also needs to be a named
list. If groupCluster
is used, cell clusters not included
in the list will be excluded in the barplot. For example, if we wanted
to group T-cells and NK-cells together, we could set
cellTypeMappings <- list(NK_Tcells = c(2,6), Bcells = 5, Monocytes = 1)
threshold
parameter.SingleCellExperiment
, then you
will need to supply the original counts matrix or the decontaminated
counts matrix as the first argument to generate the barplots.Another useful way to assess the amount of decontamination is to view
the expression of marker genes before and after decontX
across cell types. Here we view the monocyte markers in each cell type.
The violin plot shows that the markers have been removed from T-cells,
B-cells, and NK-cells, but are largely unaffected in monocytes.
plotDecontXMarkerExpression(sce,
markers = markers[["Monocyte_Markers"]],
groupClusters = cellTypeMappings,
ncol = 3)
Some helpful hints when using
plotDecontXMarkerExpression
:
groupClusters
works the same way as in
plotDecontXMarkerPercentage
.groupClusters
). Therefore, you may want
to keep the number of markers small in each plot and call the function
multiple times for different sets of marker genes.plotDots = TRUE
and/or log transform the points on the fly
by setting log1p = TRUE
.SingleCellExperiment
. Therefore you could also examine
normalized expression of the original and decontaminated counts. For
example:The ability of DecontX to accurately identify contamination is
dependent on the cell cluster labels. DecontX assumes that contamination
for a cell cluster comes from combination of counts from all other
clusters. The default clustering approach used by DecontX tends to
select fewer clusters that represent broader cell types. For example,
all T-cells tend to be clustered together rather than splitting naive
and cytotoxic T-cells into separate clusters. Custom cell type labels
can be suppled via the z
parameter if some cells are not
being clustered appropriately by the default method.
There are ways to force decontX
to estimate more or less
contamination across a dataset by manipulating the priors. The
delta
parameter is a numeric vector of length two. It is
the concentration parameter for the Dirichlet distribution which serves
as the prior for the proportions of native and contamination counts in
each cell. The first element is the prior for the proportion of native
counts while the second element is the prior for the proportion of
contamination counts. These essentially act as pseudocounts for the
native and contamination in each cell. If
estimateDelta = TRUE
, delta
is only used to
produce a random sample of proportions for an initial value of
contamination in each cell. Then delta
is updated in each
iteration. If estimateDelta = FALSE
, then
delta
is fixed with these values for the entire inference
procedure. Fixing delta
and setting a high number in the
second element will force decontX
to be more aggressive and
estimate higher levels of contamination in each cell at the expense of
potentially removing native expression. For example, in the previous
PBMC example, we can see what the estimated delta
was by
looking in the estimates:
## [1] 9.287164 1.038217
Setting a higher value in the second element of delta and
estimateDelta = FALSE
will force decontX
to
estimate higher levels of contamination per cell:
If you are using the Seurat package for downstream analysis, the following code can be used to read in a matrix and convert between Seurat and SCE objects:
# Read counts from CellRanger output
library(Seurat)
counts <- Read10X("sample/outs/filtered_feature_bc_matrix/")
# Create a SingleCellExperiment object and run decontX
sce <- SingleCellExperiment(list(counts = counts))
sce <- decontX(sce)
# Create a Seurat object from a SCE with decontX results
seuratObject <- CreateSeuratObject(round(decontXcounts(sce)))
Optionally, the “raw” matrix can be also be imported and used as the background:
counts.raw <- Read10X("sample/outs/raw_feature_bc_matrix/")
sce.raw <- SingleCellExperiment(list(counts = counts.raw))
sce <- decontX(sce, background = sce.raw)
Note that the decontaminated matrix of decontX consists of floating point numbers and must be rounded to integers before adding it to a Seurat object. If you already have a Seurat object containing the counts matrix and would like to run decontX, you can retrieve the count matrix, create a SCE object, and run decontX, and then add it back to the Seurat object:
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## 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
##
## time zone: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] scater_1.33.4 ggplot2_3.5.1
## [3] scuttle_1.15.5 celda_1.21.0
## [5] TENxPBMCData_1.23.0 HDF5Array_1.33.8
## [7] rhdf5_2.49.0 DelayedArray_0.31.14
## [9] SparseArray_1.5.45 S4Arrays_1.5.11
## [11] abind_1.4-8 Matrix_1.7-1
## [13] SingleCellExperiment_1.27.2 dplyr_1.1.4
## [15] Seurat_5.1.0 SeuratObject_5.0.2
## [17] sp_2.1-4 SingleCellMultiModal_1.17.4
## [19] MultiAssayExperiment_1.31.5 SummarizedExperiment_1.35.5
## [21] Biobase_2.67.0 GenomicRanges_1.57.2
## [23] GenomeInfoDb_1.41.2 IRanges_2.39.2
## [25] S4Vectors_0.43.2 BiocGenerics_0.53.0
## [27] MatrixGenerics_1.17.1 matrixStats_1.4.1
## [29] decontX_1.5.0 BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] spatstat.sparse_3.1-0 doParallel_1.0.17 httr_1.4.7
## [4] RColorBrewer_1.1-3 tools_4.4.1 sctransform_0.4.1
## [7] utf8_1.2.4 R6_2.5.1 lazyeval_0.2.2
## [10] uwot_0.2.2 rhdf5filters_1.17.0 withr_3.0.2
## [13] gridExtra_2.3 progressr_0.15.0 cli_3.6.3
## [16] spatstat.explore_3.3-3 fastDummies_1.7.4 enrichR_3.2
## [19] labeling_0.4.3 sass_0.4.9 spatstat.data_3.1-2
## [22] ggridges_0.5.6 pbapply_1.7-2 QuickJSR_1.4.0
## [25] StanHeaders_2.32.10 dbscan_1.2-0 WriteXLS_6.7.0
## [28] parallelly_1.38.0 RSQLite_2.3.7 generics_0.1.3
## [31] combinat_0.0-8 ica_1.0-3 spatstat.random_3.3-2
## [34] inline_0.3.19 loo_2.8.0 ggbeeswarm_0.7.2
## [37] fansi_1.0.6 lifecycle_1.0.4 yaml_2.3.10
## [40] BiocFileCache_2.15.0 Rtsne_0.17 grid_4.4.1
## [43] blob_1.2.4 promises_1.3.0 ExperimentHub_2.13.1
## [46] crayon_1.5.3 miniUI_0.1.1.1 lattice_0.22-6
## [49] beachmat_2.23.0 cowplot_1.1.3 KEGGREST_1.45.1
## [52] magick_2.8.5 sys_3.4.3 maketools_1.3.1
## [55] pillar_1.9.0 knitr_1.48 rjson_0.2.23
## [58] future.apply_1.11.3 codetools_0.2-20 leiden_0.4.3.1
## [61] glue_1.8.0 V8_6.0.0 spatstat.univar_3.0-1
## [64] data.table_1.16.2 vctrs_0.6.5 png_0.1-8
## [67] spam_2.11-0 gtable_0.3.6 cachem_1.1.0
## [70] xfun_0.48 mime_0.12 RcppEigen_0.3.4.0.2
## [73] survival_3.7-0 iterators_1.0.14 fitdistrplus_1.2-1
## [76] ROCR_1.0-11 nlme_3.1-166 bit64_4.5.2
## [79] filelock_1.0.3 RcppAnnoy_0.0.22 rstan_2.32.6
## [82] bslib_0.8.0 irlba_2.3.5.1 vipor_0.4.7
## [85] KernSmooth_2.23-24 colorspace_2.1-1 DBI_1.2.3
## [88] tidyselect_1.2.1 bit_4.5.0 compiler_4.4.1
## [91] curl_5.2.3 BiocNeighbors_2.1.0 plotly_4.10.4
## [94] scales_1.3.0 lmtest_0.9-40 rappdirs_0.3.3
## [97] stringr_1.5.1 SpatialExperiment_1.15.1 digest_0.6.37
## [100] goftest_1.2-3 spatstat.utils_3.1-0 rmarkdown_2.28
## [103] XVector_0.45.0 htmltools_0.5.8.1 pkgconfig_2.0.3
## [106] highr_0.11 dbplyr_2.5.0 fastmap_1.2.0
## [109] rlang_1.1.4 htmlwidgets_1.6.4 UCSC.utils_1.1.0
## [112] shiny_1.9.1 farver_2.1.2 jquerylib_0.1.4
## [115] zoo_1.8-12 jsonlite_1.8.9 BiocParallel_1.39.0
## [118] BiocSingular_1.23.0 magrittr_2.0.3 GenomeInfoDbData_1.2.13
## [121] dotCall64_1.2 patchwork_1.3.0 Rhdf5lib_1.27.0
## [124] munsell_0.5.1 Rcpp_1.0.13 viridis_0.6.5
## [127] reticulate_1.39.0 stringi_1.8.4 MCMCprecision_0.4.0
## [130] zlibbioc_1.51.2 MASS_7.3-61 AnnotationHub_3.15.0
## [133] plyr_1.8.9 pkgbuild_1.4.5 parallel_4.4.1
## [136] listenv_0.9.1 ggrepel_0.9.6 deldir_2.0-4
## [139] Biostrings_2.75.0 splines_4.4.1 tensor_1.5
## [142] igraph_2.1.1 spatstat.geom_3.3-3 RcppHNSW_0.6.0
## [145] buildtools_1.0.0 ScaledMatrix_1.13.0 reshape2_1.4.4
## [148] rstantools_2.4.0 BiocVersion_3.21.1 evaluate_1.0.1
## [151] RcppParallel_5.1.9 BiocManager_1.30.25 foreach_1.5.2
## [154] httpuv_1.6.15 RANN_2.6.2 tidyr_1.3.1
## [157] purrr_1.0.2 polyclip_1.10-7 future_1.34.0
## [160] scattermore_1.2 BiocBaseUtils_1.9.0 rsvd_1.0.5
## [163] xtable_1.8-4 RSpectra_0.16-2 later_1.3.2
## [166] viridisLite_0.4.2 tibble_3.2.1 beeswarm_0.4.0
## [169] memoise_2.0.1 AnnotationDbi_1.69.0 cluster_2.1.6
## [172] globals_0.16.3