knitr::opts_chunk$set( collapse = TRUE, comment = "#>" )
The scMerge algorithm allows batch effect removal and normalisation for single cell RNA-Seq data. It comprises of three key components including:
The purpose of this vignette is to illustrate some uses of scMerge
and explain its key components.
We will load the scMerge
package. We designed our package to be consistent with the popular BioConductor's single cell analysis framework, namely the SingleCellExperiment
and scater
package.
suppressPackageStartupMessages({ library(SingleCellExperiment) library(scMerge) library(scater) })
We provided an illustrative mouse embryonic stem cell (mESC) data in our package, as well as a set of pre-computed stably expressed gene (SEG) list to be used as negative control genes.
The full curated, unnormalised mESC data can be found here. The scMerge
package comes with a sub-sampled, two-batches version of this data (named "batch2" and "batch3" to be consistent with the full data) .
library(genefilter) load("~/Downloads/sce_mESC.rda") data("segList_ensemblGeneID", package = "scMerge") set.seed(2019) example_sce = sce_mESC[, sce_mESC$batch %in% c("batch2", "batch3")] example_sce$batch = droplevels(example_sce$batch) batch2Sampled = sample(colnames(example_sce[,example_sce$batch == "batch2"]), 100) batch3Sampled = sample(colnames(example_sce[,example_sce$batch == "batch3"]), 100) countsMat = SingleCellExperiment::counts(example_sce) batchTest = rowFtests(countsMat, fac = example_sce$batch) celltypeTest = rowFtests(countsMat, fac = factor(example_sce$cellTypes)) commonSegGenes = intersect(segList_ensemblGeneID$mouse$mouse_scSEG, rownames(sce_mESC)) keepGenes = unique(c(commonSegGenes, rownames(batchTest)[rank(batchTest$p.value) < 50], rownames(celltypeTest)[rank(celltypeTest$p.value) < 250] )) example_sce = example_sce[keepGenes, c(batch2Sampled, batch3Sampled)] example_sce = example_sce[base::rowSums(counts(example_sce)) != 0, base::colSums(counts(example_sce)) != 0] table(example_sce$batch, example_sce$cellTypes) dim(example_sce) example_sce = runPCA(example_sce, exprs_values = "logcounts") scater::plotPCA(example_sce, colour_by = "cellTypes", shape_by = "batch") save(example_sce, file = "data/example_sce.rda")
## Subsetted mouse ESC data data("example_sce", package = "scMerge")
In this mESC data, we pooled data from 2 different batches from three different cell types. Using a PCA plot, we can see that despite strong separation of cell types, there is also a strong separation due to batch effects. This information is stored in the colData
of example_sce
.
example_sce = runPCA(example_sce, exprs_values = "logcounts") scater::plotPCA( example_sce, colour_by = "cellTypes", shape_by = "batch")
The first major component of scMerge
is to obtain negative controls for our normalisation. In this vignette, we will be using a set of pre-computed SEGs from a single cell mouse data made available through the segList_ensemblGeneID
data in our package. For more information about the selection of negative controls and SEGs, please see Section select SEGs.
## single-cell stably expressed gene list data("segList_ensemblGeneID", package = "scMerge") head(segList_ensemblGeneID$mouse$mouse_scSEG)
The second major component of scMerge
is to compute pseudo-replicates for cells so we can perform normalisation. We offer three major ways of computing this pseudo-replicate information:
scMerge
In unsupervised scMerge
, we will perform a k-means clustering to obtain pseudo-replicates. This requires the users to supply a kmeansK
vector with each element indicating number of clusters in each of the batches. For example, we know "batch2" and "batch3" both contain three cell types. Hence, kmeansK = c(3, 3)
in this case.
t1 = Sys.time()
scMerge_unsupervised <- scMerge( sce_combine = example_sce, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_unsupervised")
t2 = Sys.time()
We now colour construct the PCA plot again on our normalised data. We can observe a much better separation by cell type and less separation by batches.
scMerge_unsupervised = runPCA(scMerge_unsupervised, exprs_values = "scMerge_unsupervised") scater::plotPCA( scMerge_unsupervised, colour_by = "cellTypes", shape_by = "batch")
By default, scMerge
only uses 50% of the cells to perform kmeans clustering. While this is sufficient to perform a satisfactory normalisation in most cases, users can control if they wish all cells be used in the kmeans clustering.
scMerge_unsupervised_all <- scMerge( sce_combine = example_sce, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_unsupervised_all", replicate_prop = 1)
scMerge_unsupervised_all = runPCA(scMerge_unsupervised_all, exprs_values = "scMerge_unsupervised_all") scater::plotPCA( scMerge_unsupervised_all, colour_by = "cellTypes", shape_by = "batch")
scMerge
If all cell type information is available to the user, then it is possible to use this information to create pseudo-replicates. This can be done through the cell_type
argument in the scMerge
function.
scMerge_supervised <- scMerge( sce_combine = example_sce, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_supervised", cell_type = example_sce$cellTypes)
scMerge_supervised = runPCA(scMerge_supervised, exprs_values = "scMerge_supervised") scater::plotPCA( scMerge_supervised, colour_by = "cellTypes", shape_by = "batch")
If the user is only able to access partial cell type information, then it is still possible to use this information to create pseudo-replicates. This can be done through the cell_type
and cell_type_inc
arguments in the scMerge
function. cell_type_inc
should contain a vector of indices indicating which elements in the cell_type
vector should be used to perform semi-supervised scMerge.
scMerge_semisupervised1 <- scMerge( sce_combine = example_sce, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3,3), assay_name = "scMerge_semisupervised1", cell_type = example_sce$cellTypes, cell_type_inc = which(example_sce$cellTypes == "2i"), cell_type_match = FALSE)
scMerge_semisupervised1 = runPCA(scMerge_semisupervised1, exprs_values = "scMerge_semisupervised1") scater::plotPCA( scMerge_semisupervised1, colour_by = "cellTypes", shape_by = "batch")
There is alternative semi-supervised method to create pseudo-replicates for scMerge
. This uses known cell type information to identify mutual nearest clusters and it is achieved via the cell_type
and cell_type_match = TRUE
options in the scMerge
function.
scMerge_semisupervised2 <- scMerge( sce_combine = example_sce, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_semisupervised2", cell_type = example_sce$cellTypes, cell_type_inc = NULL, cell_type_match = TRUE)
scMerge_semisupervised2 = runPCA(scMerge_semisupervised2, exprs_values = "scMerge_semisupervised2") scater::plotPCA( scMerge_semisupervised2, colour_by = "cellTypes", shape_by = "batch")
In simple terms, a negative control is a gene that has expression values relatively constant across these datasets. The concept of using these negative control genes for normalisation was most widely used in the RUV method family (e.g. Gagnon-Bartsch & Speed (2012) and Risso et. al. (2014)) and there exist multiple methods to find these negative controls. In our paper, we recommened the SEGs as negative controls for scRNA-Seq data and SEGs can be found using either a data-adaptive computational method or external knowledge.
scSEGIndex
to calculate the SEG from a data matrix. The output of this function is a data.frame
with a SEG index calculated for each gene. See Lin et. al. (2018) for more details.exprs_mat = SummarizedExperiment::assay(example_sce, 'counts') result = scSEGIndex(exprs_mat = exprs_mat)
## SEG list in ensemblGene ID data("segList_ensemblGeneID", package = "scMerge") ## SEG list in official gene symbols data("segList", package = "scMerge") ## SEG list for human scRNA-Seq data head(segList$human$human_scSEG) ## SEG list for human bulk microarray data head(segList$human$bulkMicroarrayHK) ## SEG list for human bulk RNASeq data head(segList$human$bulkRNAseqHK)
Under most circumstances, scMerge
is fast enough to be used on a personal laptop for a moderately large data. However, we do recognise the difficulties associated with computation when dealing with larger data. To this end, we devised a fast version of scMerge
. The major difference between the two versions lies on the noise estimation component, which utilised singular value decomposition (SVD). In order to speed up scMerge
, we used BiocSingular
package that offers several SVD speed improvements. This computational method is able to speed up scMerge
by obtain a very accurate approximation of the noise structure in the data. This option is achieved via the option BSPARAM = IrlbaParam()
or BSPARAM = RandomParam()
. Additionally, svd_k
is a parameter that controlling the degree of approximations.
We recommend using this option in the case where the number of cells is large in your single cell data. The speed advantage we obtain for large single cell data is much more dramatic than on a smaller dataset like the example mESC data. For example, a single run of normal scMerge
on a human pancreas data (23699 features and 4566 cells) takes about 10 minutes whereas the speed up version takes just under 4 minutes.
t3 = Sys.time()
library(BiocSingular) scMerge_fast <- scMerge( sce_combine = example_sce, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_fast", BSPARAM = IrlbaParam(), svd_k = 20)
t4 = Sys.time()
paste("Normally, scMerge takes ", round(t2 - t1, 2), " seconds") paste("Fast version of scMerge takes ", round(t4 - t3, 2), " seconds") scMerge_fast = runPCA(scMerge_fast, exprs_values = "scMerge_fast") scater::plotPCA( scMerge_fast, colour_by = "cellTypes", shape_by = "batch") + labs(title = "fast scMerge yields similar results to the default version")
scMerge
is implemented with a parallelised computational option via the BiocParallel package. You can enable this option using the BPPARAM
argument with various BiocParallelParam
objects that is suitable for your operating system.
Please note that any parallelisation would incur a small overhead. Hence we recommend you do not use parallelisation for small data.
library(BiocParallel) scMerge_parallel <- scMerge( sce_combine = example_sce, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_parallel", BPPARAM = MulticoreParam(workers = 2) )
scMerge
also supports sparse array input, which could be very helpful in speeding up computations and saving RAM. scMerge
does not perform internal matrix conversion, so you may use the following codes as an example of converting typical matrix
class to sparse matrices before running scMerge
.
library(Matrix) library(DelayedArray) sparse_input = example_sce assay(sparse_input, "counts") = as(counts(sparse_input), "dgeMatrix") assay(sparse_input, "logcounts") = as(logcounts(sparse_input), "dgeMatrix") scMerge_sparse = scMerge( sce_combine = sparse_input, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_sparse")
HDF5Array
)Bioconductor provides an infrastructure for out-of-memory computation through HDF5Array
. In simple terms, we can load an on-disk data into RAM, make computations and write to hard disk. This is particularly helpful when the data is too large for in-RAM computations. You may use the following codes as an example of converting typical matrix
class to HDF5Array
matrices before running scMerge
.
library(HDF5Array) library(DelayedArray) DelayedArray:::set_verbose_block_processing(TRUE) ## To monitor block processing hdf5_input = example_sce assay(hdf5_input, "counts") = as(counts(hdf5_input), "HDF5Array") assay(hdf5_input, "logcounts") = as(logcounts(hdf5_input), "HDF5Array") scMerge_hdf5 = scMerge( sce_combine = sparse_input, ctl = segList_ensemblGeneID$mouse$mouse_scSEG, kmeansK = c(3, 3), assay_name = "scMerge_hdf5")
Please check out our paper for detailed analysis and results on multiple scRNA-Seq data. https://www.biorxiv.org/content/10.1101/393280v2
citation("scMerge")
sessionInfo()
Any scripts or data that you put into this service are public.
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.