Dimensionality reduction and batch effect removal using NewWave

Installation

First of all we need to install NewWave:

if(!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("NewWave")
suppressPackageStartupMessages(
  {library(SingleCellExperiment)
library(splatter)
library(irlba)
library(Rtsne)
library(ggplot2)
library(mclust)
library(NewWave)}
)

Introduction

NewWave is a new package that assumes a Negative Binomial distributions for dimensionality reduction and batch effect removal. In order to reduce the memory consumption it uses a PSOCK cluster combined with the R package SharedObject that allow to share a matrix between different cores without memory duplication. Thanks to that we can massively parallelize the estimation process with huge benefit in terms of time consumption. We can reduce even more the time consumption using some minibatch approaches on the different steps of the optimization.

I am going to show how to use NewWave with example data generated with Splatter.

params <- newSplatParams()
N=500
set.seed(1234)
data <- splatSimulateGroups(params,batchCells=c(N/2,N/2),
                           group.prob = rep(0.1,10),
                           de.prob = 0.2,
                           verbose = FALSE) 

Now we have a dataset with 500 cells and 10000 genes, I will use only the 500 most variable genes. NewWave takes as input raw data, not normalized.

set.seed(12359)
hvg <- rowVars(counts(data))
names(hvg) <- rownames(counts(data))
data <- data[names(sort(hvg,decreasing=TRUE))[1:500],]

As you can see there is a variable called batch in the colData section.

colData(data)
#> DataFrame with 500 rows and 4 columns
#>                Cell       Batch    Group ExpLibSize
#>         <character> <character> <factor>  <numeric>
#> Cell1         Cell1      Batch1  Group5     51275.8
#> Cell2         Cell2      Batch1  Group5     50170.9
#> Cell3         Cell3      Batch1  Group10    60703.1
#> Cell4         Cell4      Batch1  Group3     68564.5
#> Cell5         Cell5      Batch1  Group2     45825.5
#> ...             ...         ...      ...        ...
#> Cell496     Cell496      Batch2  Group8     82344.6
#> Cell497     Cell497      Batch2  Group7     65257.7
#> Cell498     Cell498      Batch2  Group9     69259.2
#> Cell499     Cell499      Batch2  Group6     59976.0
#> Cell500     Cell500      Batch2  Group10    59980.9

IMPORTANT: For batch effecr removal the batch variable must be a factor

data$Batch <- as.factor(data$Batch)

We also have a variable called Group that represent the cell type labels.

We can see the how the cells are distributed between group and batch

pca <- prcomp_irlba(t(counts(data)),n=10)
plot_data <-data.frame(Rtsne(pca$x)$Y)
plot_data$batch <- data$Batch
plot_data$group <- data$Group
ggplot(plot_data, aes(x=X1,y=X2,col=group, shape=batch))+ geom_point()

There is a clear batch effect between the cells.

Let’s try to correct it.

NewWave

I am going to show different implementation and the suggested way to use them with the given hardware.

Some advise:

Standard usage

This is the way to insert the batch variable, in the same manner can be inserted other cell-related variable and if you need some gene related variable those can be inserted in V.

res <- newWave(data,X = "~Batch", K=10, verbose = TRUE)
#> Time of setup
#>    user  system elapsed 
#>   0.009   0.012   1.043 
#> Time of initialization
#>    user  system elapsed 
#>   0.056   0.012   1.034
#> Iteration 1
#> penalized log-likelihood = -1290948.63834991
#> Time of dispersion optimization
#>    user  system elapsed 
#>   1.008   0.028   1.037
#> after optimize dispersion = -1055634.52319307
#> Time of right optimization
#>    user  system elapsed 
#>   0.001   0.000  12.999
#> after right optimization= -1054898.38388464
#> after orthogonalization = -1054898.35300597
#> Time of left optimization
#>    user  system elapsed 
#>   0.001   0.000  11.354
#> after left optimization= -1054515.9551221
#> after orthogonalization = -1054515.95197842
#> Iteration 2
#> penalized log-likelihood = -1054515.95197842
#> Time of dispersion optimization
#>    user  system elapsed 
#>   1.106   0.024   1.318
#> after optimize dispersion = -1054508.01461186
#> Time of right optimization
#>    user  system elapsed 
#>   0.005   0.000  11.866
#> after right optimization= -1054470.05461025
#> after orthogonalization = -1054470.05286745
#> Time of left optimization
#>    user  system elapsed 
#>   0.001   0.000   5.594
#> after left optimization= -1054457.37220533
#> after orthogonalization = -1054457.37218113

In order to make it faster you can increase the number of cores using “children” parameter:

res2 <- newWave(data,X = "~Batch", K=10, verbose = TRUE, children=2)
#> Time of setup
#>    user  system elapsed 
#>   0.020   0.004   0.708 
#> Time of initialization
#>    user  system elapsed 
#>   0.064   0.008   1.209
#> Iteration 1
#> penalized log-likelihood = -1290948.63856983
#> Time of dispersion optimization
#>    user  system elapsed 
#>   1.028   0.032   1.680
#> after optimize dispersion = -1055634.52537363
#> Time of right optimization
#>    user  system elapsed 
#>   0.003   0.000   7.992
#> after right optimization= -1054898.38527306
#> after orthogonalization = -1054898.35439348
#> Time of left optimization
#>    user  system elapsed 
#>   0.003   0.000   6.019
#> after left optimization= -1054515.99132131
#> after orthogonalization = -1054515.98818897
#> Iteration 2
#> penalized log-likelihood = -1054515.98818897
#> Time of dispersion optimization
#>    user  system elapsed 
#>   0.655   0.020   0.675
#> after optimize dispersion = -1054508.04998052
#> Time of right optimization
#>    user  system elapsed 
#>   0.001   0.000   3.548
#> after right optimization= -1054470.11401871
#> after orthogonalization = -1054470.11229962
#> Time of left optimization
#>    user  system elapsed 
#>   0.001   0.000   2.691
#> after left optimization= -1054457.38634829
#> after orthogonalization = -1054457.38632061

Commonwise dispersion and minibatch approaches

If you do not have an high number of cores to run newWave this is the fastest way to run. The optimization process is done by three process itereated until convercence.

Each of these three steps can be accelerated using mini batch, the number of observation is settled with these parameters:

res3 <- newWave(data,X = "~Batch", verbose = TRUE,K=10, children=2,
                n_gene_disp = 100, n_gene_par = 100, n_cell_par = 100)
#> Time of setup
#>    user  system elapsed 
#>   0.009   0.000   0.362 
#> Time of initialization
#>    user  system elapsed 
#>   0.058   0.000   0.491
#> Iteration 1
#> penalized log-likelihood = -1290948.63823182
#> Time of dispersion optimization
#>    user  system elapsed 
#>   0.647   0.004   0.651
#> after optimize dispersion = -1055634.52235217
#> Time of right optimization
#>    user  system elapsed 
#>   0.002   0.000   5.567
#> after right optimization= -1054898.38464066
#> after orthogonalization = -1054898.35377449
#> Time of left optimization
#>    user  system elapsed 
#>   0.002   0.000   8.667
#> after left optimization= -1054515.94268307
#> after orthogonalization = -1054515.93955328
#> Iteration 2
#> penalized log-likelihood = -1054515.93955328
#> Time of dispersion optimization
#>    user  system elapsed 
#>   0.356   0.003   0.364
#> after optimize dispersion = -1054511.551407
#> Time of right optimization
#>    user  system elapsed 
#>   0.000   0.004   2.056
#> after right optimization= -1054503.3307867
#> after orthogonalization = -1054503.3298583
#> Time of left optimization
#>    user  system elapsed 
#>   0.003   0.001   1.108
#> after left optimization= -1054502.99110147
#> after orthogonalization = -1054502.99107411

Genewise dispersion mini-batch

If you have a lot of core disposable or you want to estimate a genewise dispersion parameter this is the fastes configuration:

res3 <- newWave(data,X = "~Batch", verbose = TRUE,K=10, children=2,
                n_gene_par = 100, n_cell_par = 100, commondispersion = FALSE)
#> Time of setup
#>    user  system elapsed 
#>   0.027   0.000   0.982 
#> Time of initialization
#>    user  system elapsed 
#>   0.054   0.004   0.778
#> Iteration 1
#> penalized log-likelihood = -1290948.63834109
#> Time of dispersion optimization
#>    user  system elapsed 
#>   0.626   0.016   0.642
#> after optimize dispersion = -1055634.52177397
#> Time of right optimization
#>    user  system elapsed 
#>   0.001   0.000   3.424
#> after right optimization= -1054898.38242025
#> after orthogonalization = -1054898.35154454
#> Time of left optimization
#>    user  system elapsed 
#>   0.001   0.000   3.668
#> after left optimization= -1054515.98687352
#> after orthogonalization = -1054515.98373779
#> Iteration 2
#> penalized log-likelihood = -1054515.98373779
#> Time of dispersion optimization
#>    user  system elapsed 
#>   0.062   0.012   0.629
#> after optimize dispersion = -1050020.99860375
#> Time of right optimization
#>    user  system elapsed 
#>   0.001   0.000   0.732
#> after right optimization= -1050014.25270341
#> after orthogonalization = -1050014.2518192
#> Time of left optimization
#>    user  system elapsed 
#>   0.001   0.000   0.670
#> after left optimization= -1049975.05924613
#> after orthogonalization = -1049975.05841423
#> Iteration 3
#> penalized log-likelihood = -1049975.05841423
#> Time of dispersion optimization
#>    user  system elapsed 
#>   0.053   0.008   0.241
#> after optimize dispersion = -1049975.07855537
#> Time of right optimization
#>    user  system elapsed 
#>   0.001   0.000   0.612
#> after right optimization= -1049967.28039586
#> after orthogonalization = -1049967.27940269
#> Time of left optimization
#>    user  system elapsed 
#>   0.002   0.000   0.645
#> after left optimization= -1049933.96699264
#> after orthogonalization = -1049933.96608838

NB: do not use n_gene_disp in this case, it will slower the computation.

Now I can use the latent dimension rapresentation for visualization purpose:

latent <- reducedDim(res)

tsne_latent <- data.frame(Rtsne(latent)$Y)
tsne_latent$batch <- data$Batch
tsne_latent$group <- data$Group
ggplot(tsne_latent, aes(x=X1,y=X2,col=group, shape=batch))+ geom_point()

or for clustering:

cluster <- kmeans(latent, 10)

adjustedRandIndex(cluster$cluster, data$Group)
#> [1] 0.8593404

Session Information

sessionInfo()
#> R Under development (unstable) (2023-10-22 r85388)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 22.04.3 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] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] NewWave_1.13.0              mclust_6.0.0               
#>  [3] ggplot2_3.4.4               Rtsne_0.16                 
#>  [5] irlba_2.3.5.1               Matrix_1.6-1.1             
#>  [7] splatter_1.27.0             SingleCellExperiment_1.25.0
#>  [9] SummarizedExperiment_1.33.0 Biobase_2.63.0             
#> [11] GenomicRanges_1.55.0        GenomeInfoDb_1.39.0        
#> [13] IRanges_2.37.0              S4Vectors_0.41.0           
#> [15] BiocGenerics_0.49.0         MatrixGenerics_1.15.0      
#> [17] matrixStats_1.0.0          
#> 
#> loaded via a namespace (and not attached):
#>  [1] gtable_0.3.4            xfun_0.40               bslib_0.5.1            
#>  [4] lattice_0.22-5          vctrs_0.6.4             tools_4.4.0            
#>  [7] bitops_1.0-7            generics_0.1.3          parallel_4.4.0         
#> [10] tibble_3.2.1            fansi_1.0.5             pkgconfig_2.0.3        
#> [13] SharedObject_1.17.0     checkmate_2.2.0         lifecycle_1.0.3        
#> [16] GenomeInfoDbData_1.2.11 farver_2.1.1            compiler_4.4.0         
#> [19] munsell_0.5.0           codetools_0.2-19        htmltools_0.5.6.1      
#> [22] sass_0.4.7              RCurl_1.98-1.12         yaml_2.3.7             
#> [25] pillar_1.9.0            crayon_1.5.2            jquerylib_0.1.4        
#> [28] BiocParallel_1.37.0     DelayedArray_0.29.0     cachem_1.0.8           
#> [31] abind_1.4-5             rsvd_1.0.5              tidyselect_1.2.0       
#> [34] locfit_1.5-9.8          digest_0.6.33           BiocSingular_1.19.0    
#> [37] dplyr_1.1.3             labeling_0.4.3          fastmap_1.1.1          
#> [40] grid_4.4.0              colorspace_2.1-0        cli_3.6.1              
#> [43] SparseArray_1.3.0       magrittr_2.0.3          S4Arrays_1.3.0         
#> [46] utf8_1.2.4              withr_2.5.1             scales_1.2.1           
#> [49] backports_1.4.1         rmarkdown_2.25          XVector_0.43.0         
#> [52] beachmat_2.19.0         ScaledMatrix_1.11.0     evaluate_0.22          
#> [55] knitr_1.44              rlang_1.1.1             Rcpp_1.0.11            
#> [58] glue_1.6.2              jsonlite_1.8.7          R6_2.5.1               
#> [61] zlibbioc_1.49.0