library(demuxSNP)
library(ComplexHeatmap)
library(viridisLite)
library(Seurat)
library(ggpubr)
library(dittoSeq)
library(utils)
library(EnsDb.Hsapiens.v86)
colors <- structure(viridis(n = 3), names = c("-1", "0", "1"))
Multiplexing in scRNAseq involves the sequencing of samples from different patients, treatment types or physiological locations together, resulting in significant cost savings. The cells must then be demultiplexed, or assigned back to their respective groups. A number of experimental and computational methods have been proposed to facilitate this, but a universally robust algorithm remains elusive. Below, we introduce some existing methods, highlight the novel features of our approach and its advantages to the user.
Cells from each group are labelled with a distinct tag (HTO or LMO) which is sequenced to give a counts matrix. Due to non-specific binding, these counts form a bimodal distribution. Such methods are generally computationally efficient. Their classification performance, however, is highly dependent on the tagging quality and many methods do not account for uncertainty in classification (Boggy et al. (2022), Kim et al. (2020) & Stoeckius et al. (2018)).
More recent methods, including demuxmix, assign a probability that a cell is from a particular group, or made up of multiple groups (doublet). This allows users to define a cut-off threshold for the assignment confidence. Accounting for uncertainty is an important feature for these types of algorithms. But, while they give the user greater flexibility in determining which cells to keep, this ultimately results in a trade-off between keeping cells which cannot be confidently called or discarding them - due to issues with tagging quality rather than RNA quality.
The second class of methods exploits natural genetic variation between cells and so can only be used where the groups are genetically distinct. Demuxlet (Kang et al. (2018)) uses genotype information from each group to classify samples. This genotyping incurs additional experimental cost. To address this, Souporcell (Heaton et al. (2020)) and Vireo (Huang, McCarthy, and Stegle (2019)) among other methods were developed to classify cells based on their SNPs in an unsupervised manner. Without prior knowledge of the SNPs associated with each group, these unsupervised methods may confuse groups with lower cell counts for other signals in the data.
Demuxlet remains the standard often used to benchmark other methods but its more widespread adoption has been limited by the requirement of sample genotype information.
With cell hashing, we can confidently demultiplex some but not all cells. Using these high confidence cells, we can learn the SNPs associated with each group. This SNP information can then be used to assign remaining cells (which we could not confidently call using cell hashing) to their most similar group based on their SNP profile.
Novel features:
Uses both cell hashing and SNP data. Current methods are limited to using one or the other.
Selects SNPs based on being located in a gene expressed in a large proportion of cells to reduce noise, computational cost and increase interpretability.
Impact:
Note: the approach used here differs from most SNP methods in that it is supervised. We attain knowledge of which SNPs are associated with which patients from the high confidence cells then use this to train a classifier. It is similar to demuxlet in the sense that the classifier uses group specific SNP information, however our method does not require the expense of genotyping and so may be much more widely applicable.
if (!require("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("demuxSNP")
## or for development version:
devtools::install_github("michaelplynch/demuxSNP")
For full list of arguments and explanation of each function, please refer to relevant documentation.
#Data loading
data(multiplexed_scrnaseq_sce,
commonvariants_1kgenomes_subset,
vartrix_consensus_snps,
package = "demuxSNP")
small_sce<-multiplexed_scrnaseq_sce[,1:100]
ensdb <- EnsDb.Hsapiens.v86::EnsDb.Hsapiens.v86
#Preprocessing
top_genes<-common_genes(small_sce)
vcf_sub<-subset_vcf(commonvariants_1kgenomes_subset, top_genes, ensdb)
small_sce<-high_conf_calls(small_sce)
#Use subsetted vcf with VarTrix in default 'consensus' mode to generate SNPs
#matrix
small_sce<-add_snps(small_sce,vartrix_consensus_snps[,1:100])
small_sce<-reassign(small_sce,k=5)
table(small_sce$knn)
#>
#> Doublet Hashtag1 Hashtag3 Hashtag4 Hashtag5 Hashtag6
#> 22 10 16 7 13 32
top_genes <- common_genes(sce)
: Returns the genes which are expressed
(expression > 0) in the highest proportion of cells. These genes are used below
to subset the .vcf file.
new_vcf <- subset_vcf(vcf = vcf, top_genes = top_genes, ensdb = ensdb)
:
Subsets a supplied .vcf to SNP locations within the genes supplied. The ranges
of the genes are extracted from the EnsDb object.
sce <- high_conf_calls(sce = sce, assay = "HTO")
:
Takes a SingleCellExperiment object with HTO altExp, runs
demuxmix
and returns a factor of assigned labels, a logical vector indicating high
confidence calls and a logical vector indicating which cells to predict (all).
sce <- add_snps(sce = sce, mat = vartrix_consensus_snps, thresh = 0.8)
:
Adds the SNP data from VarTrix (default consensus mode) to the
SingleCellExperiment object as an altExp. Additionally, filters out SNPs with no
reads in less than ‘thresh’ proportion of cells.
sce <- reassign(sce,
k = 10,
d = 10,
train_cells = sce$train,
predict_cells = sce$predict
)
Reassigns cells based on SNP profiles of high confidence cells. Singlet training data is based on high confidence singlet assignment. Doublets are simulated by systematically sampling and combining the SNP profiles of n cells pairs from each grouping combination. Cells to be used for training data are specified by “train_cells” (logical). Cells to be used for prediction are specified by “predict_cells” (logical), this may also include the training data.
We load three data objects: a SingleCellExperiment object containing RNA and HTO counts, a .vcf file of class CollapsedVCF containing SNPs from 1000 Genomes common variants and a matrix containing SNP information for each cell (we will show you how to generate this SNPs matrix using VarTrix outside of R). We have already removed low quality cells (library size<1,000 and percentage of genes mapping to mitochondrial genes>10%).
class(multiplexed_scrnaseq_sce)
#> [1] "SingleCellExperiment"
#> attr(,"package")
#> [1] "SingleCellExperiment"
class(commonvariants_1kgenomes_subset)
#> [1] "CollapsedVCF"
#> attr(,"package")
#> [1] "VariantAnnotation"
class(vartrix_consensus_snps)
#> [1] "matrix" "array"
The HTO or LMO distribution is usually bimodal, with a signal (high counts) and background distribution (low counts) caused by non-specific binding. Ideally, these distributions would be clearly separated with no overlap, but in practice, this is not always the case. In our example data, we see that the signal and noise overlap to varying extents in each group.
We will begin by running Seurat’s HTODemux, a popular HTO demultiplexing algorithm on the data.
logcounts(multiplexed_scrnaseq_sce) <- counts(multiplexed_scrnaseq_sce)
seurat <- as.Seurat(multiplexed_scrnaseq_sce)
seurat <- HTODemux(seurat)
seurat$hash.ID <- factor(as.character(seurat$hash.ID))
multiplexed_scrnaseq_sce$seurat <- seurat$hash.ID
multiplexed_scrnaseq_sce$seurat <- seurat$hash.ID
table(multiplexed_scrnaseq_sce$seurat)
#>
#> Doublet Hashtag1 Hashtag2 Hashtag3 Hashtag4 Hashtag5 Hashtag6 Negative
#> 633 121 29 264 158 177 383 235
Although HTO library size of the Negative group is low, the RNA library size is similar to that of other groups, indicating that they may be misclassified as Negative due to their tagging quality rather than overall RNA quality.
dittoPlot(seurat, "nCount_HTO", group.by = "hash.ID")
dittoPlot(seurat, "nCount_RNA", group.by = "hash.ID")
For the remainder of this vignette we will outline our method of checking whether or not these cells have been called correctly and how to assign them to their appropriate group!
Common variants files, for example from the 1000 Genomes Project, can contain over 7 million SNPs. To reduce computational cost and cell-type effects, we subset our SNPs list to those located within genes expressed across most cells in our data.
We first find the most commonly expressed genes in our RNA data.
top_genes <- common_genes(sce = multiplexed_scrnaseq_sce, n = 100)
top_genes[1:10]
#> [1] "TPT1" "RPL13" "RPL28" "TMSB4X" "RPS27" "EEF1A1" "RPL41" "B2M"
#> [9] "RPLP1" "RPL32"
We have a sample .vcf preloaded, but you can load your .vcf file in using ‘readVcf()’ from VariantAnnotation.
We will subset our .vcf file to SNPs seen in commonly expressed genes from our dataset. Notice that the genome for the vcf and EnsDb object must be compatible!
The returned vcf can be written to file and used with VarTrix later on.
ensdb <- EnsDb.Hsapiens.v86::EnsDb.Hsapiens.v86
genome(commonvariants_1kgenomes_subset)[1] == genome(ensdb)[1]
#> 1
#> TRUE
new_vcf <- subset_vcf(commonvariants_1kgenomes_subset, top_genes = top_genes, ensdb)
commonvariants_1kgenomes_subset
#> class: CollapsedVCF
#> dim: 2609 0
#> rowRanges(vcf):
#> GRanges with 5 metadata columns: paramRangeID, REF, ALT, QUAL, FILTER
#> info(vcf):
#> DataFrame with 1 column: AF
#> info(header(vcf)):
#> Number Type Description
#> AF A Float Estimated allele frequency in the range (0,1)
#> geno(vcf):
#> List of length 0:
new_vcf
#> class: CollapsedVCF
#> dim: 2399 0
#> rowRanges(vcf):
#> GRanges with 5 metadata columns: paramRangeID, REF, ALT, QUAL, FILTER
#> info(vcf):
#> DataFrame with 1 column: AF
#> info(header(vcf)):
#> Number Type Description
#> AF A Float Estimated allele frequency in the range (0,1)
#> geno(vcf):
#> List of length 0:
The subsetted .vcf can be written to disk using ‘writeVcf()’, again from VariantAnnotation package.
Next, we wish to identify cells which we can confidently call to a particular group. There are a number of ways this can be achieved, including probabilistic modelling of the HTO counts, manually setting non-conservative thresholds or using consensus calls. The user may wish to experiment with different approaches. Here we have used demuxmix, a probabilistic model which we have set with a high acceptance threshold to identify high confidence cell calls to use as training data (cells which we can confidently call as a particular singlet group).
multiplexed_scrnaseq_sce <- high_conf_calls(multiplexed_scrnaseq_sce)
table(multiplexed_scrnaseq_sce$train)
#>
#> FALSE TRUE
#> 955 1045
table(multiplexed_scrnaseq_sce$predict)
#>
#> TRUE
#> 2000
table(multiplexed_scrnaseq_sce$labels)
#>
#> Hashtag1 Hashtag2 Hashtag3 Hashtag4 Hashtag5 Hashtag6 multiplet negative
#> 62 15 226 102 348 292 335 12
#> uncertain
#> 608
So, for this particular dataset, we can confidently call 1,045 cells as being from a particular singlet group. 608 cells cannot be called to a group with high confidence.
Variant calling is not done within the package. Instead, we refer the reader to VarTrix, where they can use the subsetted .vcf file along with their .bam, barcodes.tsv and reference genome to call SNPs in each cell.
A sample VarTrix command looks like the following:
./vartrix -v <path_to_input_vcf> -b <path_to_cellranger_bam> -f <path_to_fasta_file> -c <path_to_cell_barcodes_file> -o <path_for_output_matrix>
Using the output matrix from Vartrix and the high confidence classifications from the HTO algorithm, we can reassign cells using k-nearest neighbours.
To keep things tidy, we will add the SNP data to our SingleCellExperiment object as an altExp. We recode the SNP matrix as follows: 0=no read, 1=SNP present, -1=SNP absent. This function also filters out SNPs which are observed at a low frequency in the data, and the frequency threshold can be set manually.
dim(vartrix_consensus_snps)
#> [1] 2542 2000
multiplexed_scrnaseq_sce <- add_snps(multiplexed_scrnaseq_sce, vartrix_consensus_snps, thresh = 0.95)
altExp(multiplexed_scrnaseq_sce, "SNP")
#> class: SingleCellExperiment
#> dim: 85 2000
#> metadata(0):
#> assays(1): counts
#> rownames(85): Snp Snp ... Snp Snp
#> rowData names(0):
#> colnames(2000): AAACCTGAGATCTGCT-1 AAACCTGAGCGTCAAG-1 ...
#> ACTTTCAGTAAGTTCC-1 ACTTTCAGTAGGCATG-1
#> colData names(0):
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
Before we reassign any cells, we will first use the SNPs data to inspect the results from stand-alone algorithms. Splitting the SNP data by Seurat HTODemux classification, we initially see a large number of ‘negative’ cells which appear of good quality (high proportion of reads) which may be assignable to another group. This is consistent with the library size plot we visualised earlier.
hm <- Heatmap(counts(altExp(multiplexed_scrnaseq_sce, "SNP")),
column_split = multiplexed_scrnaseq_sce$seurat,
cluster_rows = FALSE,
show_column_names = FALSE,
cluster_column_slices = FALSE,
column_title_rot = -45,
row_title = "SNPs",
show_row_names = FALSE,
col = colors
)
draw(hm,
column_title = "SNP profiles split by Seurat HTODemux call",
padding = unit(c(2, 15, 2, 2), "mm")
)
We will use our knn method to reassign cells based on their SNP profiles. The training data is the high confidence cells
set.seed(1)
multiplexed_scrnaseq_sce <- reassign(multiplexed_scrnaseq_sce,
k = 10,
d = 10,
train_cells = multiplexed_scrnaseq_sce$train,
predict_cells = multiplexed_scrnaseq_sce$predict
)
table(multiplexed_scrnaseq_sce$knn)
#>
#> Doublet Hashtag1 Hashtag2 Hashtag3 Hashtag4 Hashtag5 Hashtag6
#> 337 119 30 349 174 464 527
hm <- Heatmap(counts(altExp(multiplexed_scrnaseq_sce, "SNP")),
column_split = multiplexed_scrnaseq_sce$knn,
cluster_rows = FALSE,
show_column_names = FALSE,
cluster_column_slices = FALSE,
column_names_rot = 45,
column_title_rot = -45,
row_title = "SNPs",
show_row_names = FALSE,
col = colors
)
draw(hm,
column_title = "SNP profiles split by updated knn classification",
padding = unit(c(2, 15, 2, 2), "mm")
)
Focusing in on the knn Hashtag5 group, we see that a lot of the Negative cells have now been correctly reclassed to this group, as well as a small number of cells from other groups.
hm <- Heatmap(counts(altExp(multiplexed_scrnaseq_sce, "SNP"))[, multiplexed_scrnaseq_sce$knn == "Hashtag5"],
column_split = multiplexed_scrnaseq_sce$seurat[multiplexed_scrnaseq_sce$knn == "Hashtag5"],
cluster_rows = FALSE,
show_column_names = FALSE,
cluster_column_slices = FALSE,
column_names_rot = 45,
column_title_rot = -45,
row_title = "SNPs",
show_row_names = FALSE,
col = colors
)
draw(hm,
column_title = "knn Hashtag5 group split by Seurat HTODemux classification",
padding = unit(c(2, 15, 2, 2), "mm")
)
Next we will run some basic performance checks. We subset our SingleCellExperiment object to only retain cells which we could confidently call as being from a singlet group, then split this into a training and test dataset.
sce_test <- multiplexed_scrnaseq_sce[, multiplexed_scrnaseq_sce$train == TRUE]
sce_test$knn <- NULL
sce_test$labels <- droplevels(sce_test$labels)
sce_test
#> class: SingleCellExperiment
#> dim: 259 1045
#> metadata(0):
#> assays(2): counts logcounts
#> rownames(259): RPL22 CDC42 ... MT-ND5 MT-CYB
#> rowData names(0):
#> colnames(1045): AAACCTGAGCGTCAAG-1 AAACCTGAGGCGTACA-1 ...
#> ACTTTCAGTAAGAGAG-1 ACTTTCAGTAAGTTCC-1
#> colData names(11): orig.ident nCount_RNA ... predict labels
#> reducedDimNames(0):
#> mainExpName: RNA
#> altExpNames(2): HTO SNP
sce_test$train2 <- rep(FALSE, length(sce_test$train))
sce_test$train2[seq_len(500)] <- TRUE
sce_test$test <- sce_test$train2 == FALSE
Comparing the predicted labels in the test dataset with the hidden high confidence labels, we see excellent agreement.
sce_test <- reassign(sce_test, k = 3, train_cells = sce_test$train2, predict_cells = sce_test$test)
table(sce_test$labels[sce_test$test == TRUE], sce_test$knn[sce_test$test == TRUE])
#>
#> Doublet Hashtag1 Hashtag2 Hashtag3 Hashtag4 Hashtag5 Hashtag6
#> Hashtag1 1 35 0 0 0 0 0
#> Hashtag2 1 0 10 0 0 0 0
#> Hashtag3 0 0 0 122 0 0 0
#> Hashtag4 0 0 0 0 46 0 0
#> Hashtag5 1 0 0 0 0 193 1
#> Hashtag6 0 0 0 0 0 0 135
We can also show that that the model can correct misclassified cells when predicted back on the training data. We create a new vector “labels2” which is initially identical to the labels used in training previously.
sce_test$knn <- NULL
sce_test$labels2 <- droplevels(sce_test$labels)
table(sce_test$labels, sce_test$labels2)
#>
#> Hashtag1 Hashtag2 Hashtag3 Hashtag4 Hashtag5 Hashtag6
#> Hashtag1 62 0 0 0 0 0
#> Hashtag2 0 15 0 0 0 0
#> Hashtag3 0 0 226 0 0 0
#> Hashtag4 0 0 0 102 0 0
#> Hashtag5 0 0 0 0 348 0
#> Hashtag6 0 0 0 0 0 292
We then randomly reassign 25 cells from Hashtag5 to Hashtag2 in the training data. Predicting the model back on itself we see that in the new predicted labels, the 25 altered have been correctly reclassified back to Hashtag5
sce_test$labels2[which(sce_test$labels2 == "Hashtag5")[1:25]] <- "Hashtag2"
table(sce_test$labels, sce_test$labels2)
#>
#> Hashtag1 Hashtag2 Hashtag3 Hashtag4 Hashtag5 Hashtag6
#> Hashtag1 62 0 0 0 0 0
#> Hashtag2 0 15 0 0 0 0
#> Hashtag3 0 0 226 0 0 0
#> Hashtag4 0 0 0 102 0 0
#> Hashtag5 0 25 0 0 323 0
#> Hashtag6 0 0 0 0 0 292
sce_test <- reassign(sce_test,
train_cells = sce_test$train,
predict_cells = sce_test$train
)
table(sce_test$labels, sce_test$knn)
#>
#> Doublet Hashtag1 Hashtag2 Hashtag3 Hashtag4 Hashtag5 Hashtag6
#> Hashtag1 0 62 0 0 0 0 0
#> Hashtag2 0 0 15 0 0 0 0
#> Hashtag3 0 0 0 226 0 0 0
#> Hashtag4 0 0 0 0 102 0 0
#> Hashtag5 1 0 0 0 0 346 1
#> Hashtag6 0 0 0 0 0 0 292
For the knn Hashtag6 group, we see mostly good agreement with the original demuxmix labels. However, it appears one cell now being called a Hashtag6 was originally called Hashtag5. Again, we can investigate this by visualising the corresponding SNP profile. In doing so, we can see that this most likely is a true Hashtag6.
hm <- Heatmap(counts(altExp(sce_test, "SNP"))[, sce_test$knn == "Hashtag6"],
column_split = sce_test$labels[sce_test$knn == "Hashtag6"],
cluster_rows = FALSE,
show_column_names = FALSE,
cluster_column_slices = FALSE,
column_names_rot = 45,
column_title_rot = -45,
row_title = "SNPs",
show_row_names = FALSE,
col = colors
)
draw(hm,
column_title = "knn Hashtag6 group split by demuxmix classification",
padding = unit(c(2, 15, 2, 2), "mm")
)
sessionInfo()
#> R version 4.4.0 beta (2024-04-15 r86425)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 22.04.4 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.19-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_US.UTF-8 LC_COLLATE=en_US.UTF-8
#> [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] grid stats4 stats graphics grDevices utils datasets
#> [8] methods base
#>
#> other attached packages:
#> [1] EnsDb.Hsapiens.v86_2.99.0 dittoSeq_1.16.0
#> [3] ggpubr_0.6.0 ggplot2_3.5.1
#> [5] Seurat_5.0.3 SeuratObject_5.0.1
#> [7] sp_2.1-4 viridisLite_0.4.2
#> [9] ComplexHeatmap_2.20.0 demuxSNP_1.2.0
#> [11] ensembldb_2.28.0 AnnotationFilter_1.28.0
#> [13] GenomicFeatures_1.56.0 AnnotationDbi_1.66.0
#> [15] VariantAnnotation_1.50.0 Rsamtools_2.20.0
#> [17] Biostrings_2.72.0 XVector_0.44.0
#> [19] SingleCellExperiment_1.26.0 SummarizedExperiment_1.34.0
#> [21] Biobase_2.64.0 GenomicRanges_1.56.0
#> [23] GenomeInfoDb_1.40.0 IRanges_2.38.0
#> [25] S4Vectors_0.42.0 BiocGenerics_0.50.0
#> [27] MatrixGenerics_1.16.0 matrixStats_1.3.0
#> [29] BiocStyle_2.32.0
#>
#> loaded via a namespace (and not attached):
#> [1] ProtGenerics_1.36.0 spatstat.sparse_3.0-3 bitops_1.0-7
#> [4] httr_1.4.7 RColorBrewer_1.1-3 doParallel_1.0.17
#> [7] ggsci_3.0.3 tools_4.4.0 sctransform_0.4.1
#> [10] backports_1.4.1 utf8_1.2.4 R6_2.5.1
#> [13] lazyeval_0.2.2 uwot_0.2.2 GetoptLong_1.0.5
#> [16] withr_3.0.0 gridExtra_2.3 progressr_0.14.0
#> [19] KernelKnn_1.1.5 cli_3.6.2 Cairo_1.6-2
#> [22] spatstat.explore_3.2-7 fastDummies_1.7.3 labeling_0.4.3
#> [25] sass_0.4.9 spatstat.data_3.0-4 ggridges_0.5.6
#> [28] pbapply_1.7-2 parallelly_1.37.1 BSgenome_1.72.0
#> [31] RSQLite_2.3.6 generics_0.1.3 shape_1.4.6.1
#> [34] BiocIO_1.14.0 ica_1.0-3 spatstat.random_3.2-3
#> [37] car_3.1-2 dplyr_1.1.4 Matrix_1.7-0
#> [40] fansi_1.0.6 abind_1.4-5 lifecycle_1.0.4
#> [43] yaml_2.3.8 carData_3.0-5 SparseArray_1.4.0
#> [46] Rtsne_0.17 blob_1.2.4 promises_1.3.0
#> [49] crayon_1.5.2 miniUI_0.1.1.1 lattice_0.22-6
#> [52] cowplot_1.1.3 KEGGREST_1.44.0 demuxmix_1.6.0
#> [55] magick_2.8.3 pillar_1.9.0 knitr_1.46
#> [58] rjson_0.2.21 future.apply_1.11.2 codetools_0.2-20
#> [61] leiden_0.4.3.1 glue_1.7.0 data.table_1.15.4
#> [64] vctrs_0.6.5 png_0.1-8 spam_2.10-0
#> [67] gtable_0.3.5 cachem_1.0.8 xfun_0.43
#> [70] S4Arrays_1.4.0 mime_0.12 survival_3.6-4
#> [73] pheatmap_1.0.12 iterators_1.0.14 tinytex_0.50
#> [76] fitdistrplus_1.1-11 ROCR_1.0-11 nlme_3.1-164
#> [79] bit64_4.0.5 RcppAnnoy_0.0.22 bslib_0.7.0
#> [82] irlba_2.3.5.1 KernSmooth_2.23-22 colorspace_2.1-0
#> [85] DBI_1.2.2 tidyselect_1.2.1 bit_4.0.5
#> [88] compiler_4.4.0 curl_5.2.1 DelayedArray_0.30.0
#> [91] plotly_4.10.4 bookdown_0.39 rtracklayer_1.64.0
#> [94] scales_1.3.0 lmtest_0.9-40 stringr_1.5.1
#> [97] digest_0.6.35 goftest_1.2-3 spatstat.utils_3.0-4
#> [100] rmarkdown_2.26 htmltools_0.5.8.1 pkgconfig_2.0.3
#> [103] highr_0.10 fastmap_1.1.1 rlang_1.1.3
#> [106] GlobalOptions_0.1.2 htmlwidgets_1.6.4 UCSC.utils_1.0.0
#> [109] shiny_1.8.1.1 farver_2.1.1 jquerylib_0.1.4
#> [112] zoo_1.8-12 jsonlite_1.8.8 BiocParallel_1.38.0
#> [115] RCurl_1.98-1.14 magrittr_2.0.3 GenomeInfoDbData_1.2.12
#> [118] dotCall64_1.1-1 patchwork_1.2.0 munsell_0.5.1
#> [121] Rcpp_1.0.12 reticulate_1.36.1 stringi_1.8.3
#> [124] zlibbioc_1.50.0 MASS_7.3-60.2 plyr_1.8.9
#> [127] parallel_4.4.0 listenv_0.9.1 ggrepel_0.9.5
#> [130] deldir_2.0-4 splines_4.4.0 tensor_1.5
#> [133] circlize_0.4.16 igraph_2.0.3 spatstat.geom_3.2-9
#> [136] ggsignif_0.6.4 RcppHNSW_0.6.0 reshape2_1.4.4
#> [139] XML_3.99-0.16.1 evaluate_0.23 BiocManager_1.30.22
#> [142] foreach_1.5.2 httpuv_1.6.15 RANN_2.6.1
#> [145] tidyr_1.3.1 purrr_1.0.2 polyclip_1.10-6
#> [148] future_1.33.2 clue_0.3-65 scattermore_1.2
#> [151] broom_1.0.5 xtable_1.8-4 restfulr_0.0.15
#> [154] RSpectra_0.16-1 rstatix_0.7.2 later_1.3.2
#> [157] class_7.3-22 tibble_3.2.1 memoise_2.0.1
#> [160] GenomicAlignments_1.40.0 cluster_2.1.6 globals_0.16.3
Boggy, Gregory J, G W McElfresh, Eisa Mahyari, Abigail B Ventura, Scott G Hansen, Louis J Picker, and Benjamin N Bimber. 2022. “BFF and cellhashR: Analysis Tools for Accurate Demultiplexing of Cell Hashing Data.” Bioinformatics 38 (10): 2791–2801. https://doi.org/10.1093/bioinformatics/btac213.
Heaton, Haynes, Arthur M. Talman, Andrew Knights, Maria Imaz, Daniel J. Gaffney, Richard Durbin, Martin Hemberg, and Mara K. N. Lawniczak. 2020. “Souporcell: Robust Clustering of Single-Cell RNA-Seq Data by Genotype Without Reference Genotypes.” Nature Methods 17 (6): 615–20. https://doi.org/10.1038/s41592-020-0820-1.
Huang, Yuanhua, Davis J. McCarthy, and Oliver Stegle. 2019. “Vireo: Bayesian Demultiplexing of Pooled Single-Cell RNA-Seq Data Without Genotype Reference.” Genome Biology 20 (1): 273. https://doi.org/10.1186/s13059-019-1865-2.
Kang, Hyun Min, Meena Subramaniam, Sasha Targ, Michelle Nguyen, Lenka Maliskova, Elizabeth McCarthy, Eunice Wan, et al. 2018. “Multiplexed Droplet Single-Cell RNA-Sequencing Using Natural Genetic Variation.” Nature Biotechnology 36 (1): 89–94. https://doi.org/10.1038/nbt.4042.
Kim, Hani Jieun, Yingxin Lin, Thomas A Geddes, Jean Yee Hwa Yang, and Pengyi Yang. 2020. “CiteFuse Enables Multi-Modal Analysis of CITE-Seq Data.” Bioinformatics 36 (14): 4137–43. https://doi.org/10.1093/bioinformatics/btaa282.
Stoeckius, Marlon, Shiwei Zheng, Brian Houck-Loomis, Stephanie Hao, Bertrand Z. Yeung, William M. Mauck, Peter Smibert, and Rahul Satija. 2018. “Cell Hashing with Barcoded Antibodies Enables Multiplexing and Doublet Detection for Single Cell Genomics.” Genome Biology 19 (1): 224. https://doi.org/10.1186/s13059-018-1603-1.