1 Introduction

Batch effects refer to differences between data sets generated at different times or in different laboratories. These often occur due to uncontrolled variability in experimental factors, e.g., reagent quality, operator skill, atmospheric ozone levels. The presence of batch effects can interfere with downstream analyses if they are not explicitly modelled. For example, differential expression analyses typically use a blocking factor to absorb any batch-to-batch differences.

For single-cell RNA sequencing (scRNA-seq) data analyses, explicit modelling of the batch effect is less relevant. Manny common downstream procedures for exploratory data analysis are not model-based, including clustering and visualization. It is more generally useful to have methods that can remove batch effects to create an corrected expression matrix for further analysis. This follows the same strategy as, e.g., the removeBatchEffect() function in the limma package (Ritchie et al. 2015).

Batch correction methods designed for bulk genomics data usually require knowledge of the other factors of variation. This is usually not known in scRNA-seq experiments where the aim is to explore unknown heterogeneity in cell populations. The batchelor package implements batch correction methods that do not rely on a priori knowledge about the population structure. To demonstrate, we will use a small scRNA-seq data set (Tasic et al. 2016) from the scRNAseq package:

library(scRNAseq)
data(allen)

library(SingleCellExperiment)
sce1 <- as(allen, "SingleCellExperiment")
counts(sce1) <- assay(sce1)

library(scater)
sce1 <- sce1[1:2000,] # reducing the size for demo purposes.
sce1 <- normalize(sce1) # quick and dirty normalization.

We artificially create a batch effect in a separate SingleCellExperiment object:

sce2 <- sce1
logcounts(sce2) <- logcounts(sce2) + rnorm(nrow(sce2), sd=2)

combined <- cbind(sce1, sce2)
combined$batch <- factor(rep(1:2, c(ncol(sce1), ncol(sce2))))

plotPCA(combined, colour_by="batch") # checking there is a batch effect.

2 Mutual nearest neighbors

2.1 Overview

Mutual nearest neighbors (MNNs) are defined as pairs of cells - one from each batch - that are within each other’s set of k nearest neighbors. The idea is that MNN pairs across batches refer to the same cell type, assuming that the batch effect is orthogonal to the biological subspace (Haghverdi et al. 2018). Once MNN pairs are identified, the difference between the paired cells is used to infer the magnitude and direction of the batch effect. It is then straightforward to remove the batch effect and obtain a set of corrected expression values.

2.2 The new, fast method

The fastMNN() function performs a principal components analysis (PCA) to obtain a low-dimensional representation of the input data. MNN identification and correction is performed in this low-dimensional space, which offers some advantages with respect to speed and denoising.
This returns a SingleCellExperiment containing a matrix of corrected PC scores is returned, which can be used directly for downstream analyses such as clustering and visualization.

f.out <- fastMNN(A=sce1, B=sce2)
str(reducedDim(f.out, "corrected"))
##  num [1:758, 1:50] -0.203 -0.206 -0.202 -0.207 -0.202 ...
##  - attr(*, "dimnames")=List of 2
##   ..$ : chr [1:758] "SRR2140028" "SRR2140022" "SRR2140055" "SRR2140083" ...
##   ..$ : NULL

The batch of origin for each row/cell in the output matrix is also returned:

rle(f.out$batch)
## Run Length Encoding
##   lengths: int [1:2] 379 379
##   values : chr [1:2] "A" "B"

Another way to call fastMNN() is to specify the batch manually. This will return a corrected matrix with the cells in the same order as that in the input combined object. In this case, it doesn’t matter as the two batches were concatenated to created combined anyway, but these semantics may be useful when cells from the same batch are not contiguous.

f.out2 <- fastMNN(combined, batch=combined$batch)
str(reducedDim(f.out2, "corrected"))
##  num [1:758, 1:50] -0.203 -0.206 -0.202 -0.207 -0.202 ...
##  - attr(*, "dimnames")=List of 2
##   ..$ : chr [1:758] "SRR2140028" "SRR2140022" "SRR2140055" "SRR2140083" ...
##   ..$ : NULL

As we can see, the batch effect is successfully removed in the PCA plot below.

plotReducedDim(f.out, colour_by="batch", use_dimred="corrected")

We can also obtain per-gene corrected expression values by using the rotation vectors stored in the output. This reverses the original projection used to obtain the initial low-dimensional representation. There are, however, many caveats to using these values for downstream analyses.

cor.exp <- assay(f.out)[1,]
hist(cor.exp, xlab="Corrected expression for gene 1") 

While the default arguments are usually satisfactory, there are many options for running fastMNN(), e.g., to improve speed or to achieve a particular merge order. Refer to the corresponding workflow for more details.

2.3 The old, classic method

The original method described by Haghverdi et al. (2018) is implemented in the mnnCorrect() method. This performs the MNN identification and correction in the gene expression space, and uses a different strategy to overcome high-dimensional noise. mnnCorrect() is called with the same semantics as fastMNN():

classic.out <- mnnCorrect(sce1, sce2)

… but returns the corrected gene expression matrix directly, rather than using a low-dimensional representation1 Again, those readers wanting to use the corrected values for per-gene analyses should consider the caveats mentioned previously.. This is wrapped in a SummarizedExperiment object to store various batch-related metadata.

classic.out
## class: SingleCellExperiment 
## dim: 2000 758 
## metadata(2): merge.order merge.info
## assays(1): corrected
## rownames(2000): 0610007P14Rik 0610009B22Rik ... Awat2 Axdnd1
## rowData names(0):
## colnames(758): SRR2140028 SRR2140022 ... SRR2139341 SRR2139336
## colData names(1): batch
## reducedDimNames(0):
## spikeNames(0):

For scRNA-seq data, fastMNN() tends to be both faster and better at achieving a satisfactory merge. mnnCorrect() is mainly provided here for posterity’s sake, though it is more robust than fastMNN() to certain violations of the orthogonality assumptions.

3 Batch rescaling

rescaleBatches() effectively centers the batches in log-expression space on a per-gene basis. This is conceptually equivalent to running removeBatchEffect() with no covariates other than the batch. However, rescaleBatches() achieves this rescaling by reversing the log-transformation, downscaling the counts so that the average of each batch is equal to the smallest value, and then re-transforming. This preserves sparsity by ensuring that zeroes remain so after correction, and mitigates differences in the variance when dealing with counts of varying size between batches2 Done by downscaling, which increases the shrinkage from the added pseudo-count..

Calling rescaleBatches() returns a corrected matrix of per-gene log-expression values, wrapped in a SummarizedExperiment containin batch-related metadata.

rescale.out <- rescaleBatches(sce1, sce2)
rescale.out
## class: SingleCellExperiment 
## dim: 2000 758 
## metadata(0):
## assays(1): corrected
## rownames(2000): 0610007P14Rik 0610009B22Rik ... Awat2 Axdnd1
## rowData names(0):
## colnames(758): SRR2140028 SRR2140022 ... SRR2139341 SRR2139336
## colData names(1): batch
## reducedDimNames(0):
## spikeNames(0):

While this method is fast and simple, it makes the strong assumption that the population composition of each batch is the same. This is usually not the case for scRNA-seq experiments in real systems that exhibit biological variation. Thus, rescaleBatches() is best suited for merging technical replicates of the same sample, e.g., that have been sequenced separately.

4 Using data subsets

4.1 Selecting genes

The subset.row= argument will only perform the correction on a subset of genes in the data set. This is useful for focusing on highly variable or marker genes during high-dimensional procedures like PCA or neighbor searches, mitigating noise from irrelevant genes. For per-gene methods, this argument provides a convenient alternative to subsetting the input. Functions operating on SingleCellExperiment inputs will also automatically remove spike-in transcripts unless get.spikes=TRUE.

# Just picking a bunch of random genes:
chosen.genes <- sample(nrow(combined), 1000)
f.out3 <- fastMNN(combined, batch=combined$batch,
    subset.row=chosen.genes)
str(reducedDim(f.out2, "corrected"))
##  num [1:758, 1:50] -0.203 -0.206 -0.202 -0.207 -0.202 ...
##  - attr(*, "dimnames")=List of 2
##   ..$ : chr [1:758] "SRR2140028" "SRR2140022" "SRR2140055" "SRR2140083" ...
##   ..$ : NULL

For some functions, it is also possible to set correct.all=TRUE when subset.row= is specified. This will compute corrected values for the unselected genes as well, which is possible once the per-cell statistics are obtained with the gene subset. With this setting, we can guarantee that the output contains all the genes provided in the input.

4.2 Restricted correction

Many functions support the restrict= argument whereby the correction is determined using only a restricted subset of cells in each batch. The effect of the correction is then - for want of a better word - “extrapolated” to all other cells in that batch. This is useful for experimental designs where a control set of cells from the same source population were run on different batches. Any difference in the controls between batches must be artificial in origin, allowing us to estimate and remove the batch effect without making further biological assumptions.

# Pretend the first X cells in each batch are controls.
restrict <- list(1:100, 1:200) 
rescale.out <- rescaleBatches(sce1, sce2, restrict=restrict)

5 Other utilities

5.1 Multi-batch normalization

Differences in sequencing depth between batches are an obvious cause for batch-to-batch differences. These can be removed by multiBatchNorm(), which downscales all batches to match the coverage of the least-sequenced batch. This function returns a list of SingleCellExperiment with log-transformed normalized expression values that can be directly used for further correction.

normed <- multiBatchNorm(A=sce1, B=sce2)
names(normed)
## [1] "A" "B"

Downscaling mitigates differences in variance between batches due to the mean-variance relationship of count data. It is achieved using a median-based estimator to avoid problems with composition biases between batches (Lun, Bach, and Marioni 2016). An example usage of this function is provided in this workflow.

5.2 Multi-batch PCA

Users can perform a PCA across multiple batches using the multiBatchPCA() function. The output of this function is roughly equivalent to cbinding all batches together and performing PCA on the merged matrix. The main difference is that each sample is forced to contribute equally to the identification of the rotation vectors. This allows small batches with unique subpopulations to contribute meaningfully to the definition of the low-dimensional space.

pca.out <- multiBatchPCA(A=sce1, B=sce2)
names(pca.out)
## [1] "A" "B"

This function is used internally in fastMNN() but can be explicitly called by the user to perform hierarchical merges - see here for details.

6 Session information

sessionInfo()
## R version 3.6.0 (2019-04-26)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.9-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.9-bioc/R/lib/libRlapack.so
## 
## 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       
## 
## attached base packages:
## [1] parallel  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] scater_1.12.1               ggplot2_3.1.1              
##  [3] scRNAseq_1.10.0             batchelor_1.0.1            
##  [5] SingleCellExperiment_1.6.0  SummarizedExperiment_1.14.0
##  [7] DelayedArray_0.10.0         BiocParallel_1.18.0        
##  [9] matrixStats_0.54.0          Biobase_2.44.0             
## [11] GenomicRanges_1.36.0        GenomeInfoDb_1.20.0        
## [13] IRanges_2.18.0              S4Vectors_0.22.0           
## [15] BiocGenerics_0.30.0         knitr_1.23                 
## [17] BiocStyle_2.12.0           
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_1.0.1               rsvd_1.0.0              
##  [3] lattice_0.20-38          assertthat_0.2.1        
##  [5] digest_0.6.19            R6_2.4.0                
##  [7] plyr_1.8.4               evaluate_0.13           
##  [9] pillar_1.4.0             zlibbioc_1.30.0         
## [11] rlang_0.3.4              lazyeval_0.2.2          
## [13] irlba_2.3.3              Matrix_1.2-17           
## [15] rmarkdown_1.13           labeling_0.3            
## [17] BiocNeighbors_1.2.0      stringr_1.4.0           
## [19] RCurl_1.95-4.12          munsell_0.5.0           
## [21] compiler_3.6.0           vipor_0.4.5             
## [23] BiocSingular_1.0.0       xfun_0.7                
## [25] pkgconfig_2.0.2          ggbeeswarm_0.6.0        
## [27] htmltools_0.3.6          tidyselect_0.2.5        
## [29] tibble_2.1.1             gridExtra_2.3           
## [31] GenomeInfoDbData_1.2.1   bookdown_0.10           
## [33] viridisLite_0.3.0        crayon_1.3.4            
## [35] dplyr_0.8.1              withr_2.1.2             
## [37] bitops_1.0-6             grid_3.6.0              
## [39] gtable_0.3.0             magrittr_1.5            
## [41] scales_1.0.0             stringi_1.4.3           
## [43] XVector_0.24.0           viridis_0.5.1           
## [45] DelayedMatrixStats_1.6.0 cowplot_0.9.4           
## [47] tools_3.6.0              glue_1.3.1              
## [49] beeswarm_0.2.3           purrr_0.3.2             
## [51] yaml_2.2.0               colorspace_1.4-1        
## [53] BiocManager_1.30.4

References

Haghverdi, L., A. T. L. Lun, M. D. Morgan, and J. C. Marioni. 2018. “Batch effects in single-cell RNA-sequencing data are corrected by matching mutual nearest neighbors.” Nat. Biotechnol. 36 (5):421–27.

Lun, A. T., K. Bach, and J. C. Marioni. 2016. “Pooling across cells to normalize single-cell RNA sequencing data with many zero counts.” Genome Biol. 17 (April):75.

Ritchie, M. E., B. Phipson, D. Wu, Y. Hu, C. W. Law, W. Shi, and G. K. Smyth. 2015. “limma powers differential expression analyses for RNA-sequencing and microarray studies.” Nucleic Acids Res. 43 (7):e47.

Tasic, B., V. Menon, T. N. Nguyen, T. K. Kim, T. Jarsky, Z. Yao, B. Levi, et al. 2016. “Adult mouse cortical cell taxonomy revealed by single cell transcriptomics.” Nat. Neurosci. 19 (2):335–46.