--- description: Combining and harmonizing samples or datasets from different batches such as experiments or conditions to enable meaningful cross-sample comparisons. subtitle: Bioconductor Toolkit title: Data Integration ---
> **Note** > > Code chunks run R commands unless otherwise specified.
In this tutorial we will look at different ways of integrating multiple single cell RNA-seq datasets. We will explore a few different methods to correct for batch effects across datasets. Seurat uses the data integration method presented in Comprehensive Integration of Single Cell Data, while Scran and Scanpy use a mutual Nearest neighbour method (MNN). Below you can find a list of some methods for single data integration: ----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------- Markdown Language Library Ref ----------------- ----------------- ----------------- ----------------------------------------------------------------------------------------------------------------------------------- CCA R Seurat [Cell](https://www.sciencedirect.com/science/article/pii/S0092867419305598?via%3Dihub) MNN R/Python Scater/Scanpy [Nat. Biotech.](https://www.nature.com/articles/nbt.4091) Conos R conos [Nat. Methods](https://www.nature.com/articles/s41592-019-0466-z?error=cookies_not_supported&code=5680289b-6edb-40ad-9934-415dac4fdb2f) Scanorama Python scanorama [Nat. Biotech.](https://www.nature.com/articles/s41587-019-0113-3) ----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------- ## Data preparation Let's first load necessary libraries and the data saved in the previous lab. ``` {r} # Activate scanorama Python venv reticulate::use_virtualenv("/opt/venv/scanorama") reticulate::py_discover_config() suppressPackageStartupMessages({ library(scater) library(scran) library(patchwork) library(ggplot2) library(batchelor) library(harmony) library(reticulate) }) ``` ``` {r} # download pre-computed data if missing or long compute fetch_data <- TRUE # url for source and intermediate data path_data <- "https://export.uppmax.uu.se/naiss2023-23-3/workshops/workshop-scrnaseq" path_file <- "data/covid/results/bioc_covid_qc_dr.rds" if (!dir.exists(dirname(path_file))) dir.create(dirname(path_file), recursive = TRUE) if (fetch_data && !file.exists(path_file)) download.file(url = file.path(path_data, "covid/results/bioc_covid_qc_dr.rds"), destfile = path_file) sce <- readRDS(path_file) print(reducedDims(sce)) ``` We split the combined object into a list, with each dataset as an element. We perform standard preprocessing (log-normalization), and identify variable features individually for each dataset based on a variance stabilizing transformation (**vst**). ``` {r} sce.list <- lapply(unique(sce$sample), function(x) { x <- sce[, sce$sample == x] }) hvgs_per_dataset <- lapply(sce.list, function(x) { x <- computeSumFactors(x, sizes = c(20, 40, 60, 80)) x <- logNormCounts(x) var.out <- modelGeneVar(x, method = "loess") hvg.out <- var.out[which(var.out$FDR <= 0.05 & var.out$bio >= 0.2), ] hvg.out <- hvg.out[order(hvg.out$bio, decreasing = TRUE), ] return(rownames(hvg.out)) }) names(hvgs_per_dataset) <- unique(sce$sample) # venn::venn(hvgs_per_dataset,opacity = .4,zcolor = scales::hue_pal()(3),cexsn = 1,cexil = 1,lwd=1,col="white",borders = NA) temp <- unique(unlist(hvgs_per_dataset)) overlap <- sapply(hvgs_per_dataset, function(x) { temp %in% x }) ``` ``` {r} #| fig-height: 4 #| fig-width: 8 pheatmap::pheatmap(t(overlap * 1), cluster_rows = F, color = c("grey90", "grey20")) ## MNN ``` The mutual nearest neighbors (MNN) approach within the scran package utilizes a novel approach to adjust for batch effects. The `fastMNN()` function returns a representation of the data with reduced dimensionality, which can be used in a similar fashion to other lower-dimensional representations such as PCA. In particular, this representation can be used for downstream methods such as clustering. The BNPARAM can be used to specify the specific nearest neighbors method to use from the BiocNeighbors package. Here we make use of the [Annoy library](https://github.com/spotify/annoy) via the `BiocNeighbors::AnnoyParam()` argument. We save the reduced-dimension MNN representation into the reducedDims slot of our sce object. ``` {r} mnn_out <- batchelor::fastMNN(sce, subset.row = unique(unlist(hvgs_per_dataset)), batch = factor(sce$sample), k = 20, d = 50) ```
> **Caution** > > `fastMNN()` does not produce a batch-corrected expression matrix.
``` {r} mnn_out <- t(reducedDim(mnn_out, "corrected")) colnames(mnn_out) <- unlist(lapply(sce.list, function(x) { colnames(x) })) mnn_out <- mnn_out[, colnames(sce)] rownames(mnn_out) <- paste0("dim", 1:50) reducedDim(sce, "MNN") <- t(mnn_out) ``` We can observe that a new assay slot is now created under the name `MNN`. ``` {r} reducedDims(sce) ``` Thus, the result from `fastMNN()` should solely be treated as a reduced dimensionality representation, suitable for direct plotting, TSNE/UMAP, clustering, and trajectory analysis that relies on such results. ``` {r} set.seed(42) sce <- runTSNE(sce, dimred = "MNN", n_dimred = 50, perplexity = 30, name = "tSNE_on_MNN") sce <- runUMAP(sce, dimred = "MNN", n_dimred = 50, ncomponents = 2, name = "UMAP_on_MNN") ``` We can now plot the unintegrated and the integrated space reduced dimensions. ``` {r} #| fig-height: 6 #| fig-width: 12 wrap_plots( plotReducedDim(sce, dimred = "PCA", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "PCA"), plotReducedDim(sce, dimred = "tSNE_on_PCA", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "tSNE_on_PCA"), plotReducedDim(sce, dimred = "UMAP_on_PCA", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "UMAP_on_PCA"), plotReducedDim(sce, dimred = "MNN", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "MNN"), plotReducedDim(sce, dimred = "tSNE_on_MNN", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "tSNE_on_MNN"), plotReducedDim(sce, dimred = "UMAP_on_MNN", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "UMAP_on_MNN"), ncol = 3 ) + plot_layout(guides = "collect") ``` Let's plot some marker genes for different cell types onto the embedding. Markers Cell Type -------------------------- ------------------- CD3E T cells CD3E CD4 CD4+ T cells CD3E CD8A CD8+ T cells GNLY, NKG7 NK cells MS4A1 B cells CD14, LYZ, CST3, MS4A7 CD14+ Monocytes FCGR3A, LYZ, CST3, MS4A7 FCGR3A+ Monocytes FCER1A, CST3 DCs ``` {r} #| fig-height: 16 #| fig-width: 13 plotlist <- list() for (i in c("CD3E", "CD4", "CD8A", "NKG7", "GNLY", "MS4A1", "CD14", "LYZ", "MS4A7", "FCGR3A", "CST3", "FCER1A")) { plotlist[[i]] <- plotReducedDim(sce, dimred = "UMAP_on_MNN", colour_by = i, by_exprs_values = "logcounts", point_size = 0.6) + scale_fill_gradientn(colours = colorRampPalette(c("grey90", "orange3", "firebrick", "firebrick", "red", "red"))(10)) + ggtitle(label = i) + theme(plot.title = element_text(size = 20)) } wrap_plots(plotlist = plotlist, ncol = 3) ``` ## Harmony An alternative method for integration is Harmony, for more details on the method, please se their paper [Nat. Methods](https://www.nature.com/articles/s41592-019-0619-0). This method runs the integration on a dimensionality reduction, in most applications the PCA. So first, we will rerun scaling and PCA with the same set of genes that were used for the CCA integration. ``` {r} #| fig-height: 5 #| fig-width: 14 library(harmony) reducedDimNames(sce) sce <- RunHarmony( sce, group.by.vars = "sample", reduction.save = "harmony", reduction = "PCA", dims.use = 1:50 ) # Here we use all PCs computed from Harmony for UMAP calculation sce <- runUMAP(sce, dimred = "harmony", n_dimred = 50, ncomponents = 2, name = "UMAP_on_Harmony") ``` ## Scanorama
> **Important** > > If you are running locally using Docker and you have a Mac with ARM > chip, the Scanorama reticulate module is known to crash. In this case, > you might want to skip this section.
Another integration method is Scanorama (see [Nat. Biotech.](https://www.nature.com/articles/s41587-019-0113-3)). This method is implemented in python, but we can run it through the Reticulate package. We will run it with the same set of variable genes, but first we have to create a list of all the objects per sample. ``` {r} #| fig-height: 5 #| fig-width: 15 hvgs <- unique(unlist(hvgs_per_dataset)) scelist <- list() genelist <- list() for (i in 1:length(sce.list)) { scelist[[i]] <- t(as.matrix(logcounts(sce.list[[i]])[hvgs, ])) genelist[[i]] <- hvgs } lapply(scelist, dim) ``` Now we can run scanorama: ``` {r} #| fig-height: 5 #| fig-width: 15 scanorama <- reticulate::import("scanorama") integrated.data <- scanorama$integrate(datasets_full = scelist, genes_list = genelist) intdimred <- do.call(rbind, integrated.data[[1]]) colnames(intdimred) <- paste0("PC_", 1:100) rownames(intdimred) <- colnames(logcounts(sce)) # Add standard deviations in order to draw Elbow Plots stdevs <- apply(intdimred, MARGIN = 2, FUN = sd) attr(intdimred, "varExplained") <- stdevs reducedDim(sce, "Scanorama_PCA") <- intdimred # Here we use all PCs computed from Scanorama for UMAP calculation sce <- runUMAP(sce, dimred = "Scanorama_PCA", n_dimred = 50, ncomponents = 2, name = "UMAP_on_Scanorama") plotReducedDim(sce, dimred = "UMAP_on_Scanorama", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "UMAP_on_Scanorama") ``` ## Overview all methods Now we will plot UMAPS with all three integration methods side by side. ``` {r} #| fig-height: 8 #| fig-width: 10 p1 <- plotReducedDim(sce, dimred = "UMAP_on_PCA", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "UMAP_on_PCA") p2 <- plotReducedDim(sce, dimred = "UMAP_on_MNN", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "UMAP_on_MNN") p3 <- plotReducedDim(sce, dimred = "UMAP_on_Harmony", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "UMAP_on_Harmony") p4 <- plotReducedDim(sce, dimred = "UMAP_on_Scanorama", colour_by = "sample", point_size = 0.6) + ggplot2::ggtitle(label = "UMAP_on_Scanorama") wrap_plots(p1, p2, p3, p4, nrow = 2) + plot_layout(guides = "collect") ``` Let's save the integrated data for further analysis. ``` {r} saveRDS(sce, "data/covid/results/bioc_covid_qc_dr_int.rds") ``` ## Session info ```{=html}
``` ```{=html} ``` Click here ```{=html} ``` ``` {r} sessionInfo() ``` ```{=html}
```