Here, we demonstrate how BANKSY can be used with Harmony for integrating multiple spatial omics datasets in the presence of strong batch effects. We use 10x Visium data of the human dorsolateral prefrontal cortex from Maynard et al (2018). The data comprise 12 samples obtained from 3 subjects, with manual annotation of the layers in each sample.
library(Banksy)
library(SummarizedExperiment)
library(SpatialExperiment)
library(Seurat)
library(scran)
library(data.table)
library(harmony)
library(scater)
library(cowplot)
library(ggplot2)
SEED <- 1000
We fetch the data for all 12 DLPFC samples with the spatialLIBD package. This might take awhile.
library(spatialLIBD)
library(ExperimentHub)
ehub <- ExperimentHub::ExperimentHub()
spe <- spatialLIBD::fetch_data(type = "spe", eh = ehub)
After the download is completed, we trim the SpatialExperiment object, retaining only the counts and some metadata such as the sample identifier and pathology annotations. This saves some memory.
#' Remove NA spots
na_id <- which(is.na(spe$layer_guess_reordered_short))
spe <- spe[, -na_id]
#' Trim
imgData(spe) <- NULL
assay(spe, "logcounts") <- NULL
reducedDims(spe) <- NULL
rowData(spe) <- NULL
colData(spe) <- DataFrame(
sample_id = spe$sample_id,
subject_id = factor(spe$sample_id, labels = rep(paste0("Subject", 1:3), each = 4)),
clust_annotation = factor(as.numeric(spe$layer_guess_reordered_short)),
in_tissue = spe$in_tissue,
row.names = colnames(spe)
)
colnames(spe) <- paste0(colnames(spe), "_", spe$sample_id)
invisible(gc())
We analyse the first sample of each subject due to vignette runtime constraints.
spe <- spe[, spe$sample_id %in% c("151507", "151669", "151673")]
sample_names <- unique(spe$sample_id)
Next, stagger the spatial coordinates across the samples so that spots from different samples do not overlap.
#' Stagger spatial coordinates
locs <- spatialCoords(spe)
locs <- cbind(locs, sample_id = factor(spe$sample_id))
locs_dt <- data.table(locs)
colnames(locs_dt) <- c("sdimx", "sdimy", "group")
locs_dt[, sdimx := sdimx - min(sdimx), by = group]
global_max <- max(locs_dt$sdimx) * 1.5
locs_dt[, sdimx := sdimx + group * global_max]
locs <- as.matrix(locs_dt[, 1:2])
rownames(locs) <- colnames(spe)
spatialCoords(spe) <- locs
Find highly variable features and normalize counts. Here we use Seurat, but
other methods may also be used (e.g. scran::getTopHVGs
).
#' Get HVGs
seu <- as.Seurat(spe, data = NULL)
seu <- FindVariableFeatures(seu, nfeatures = 2000)
#' Normalize data
scale_factor <- median(colSums(assay(spe, "counts")))
seu <- NormalizeData(seu, scale.factor = scale_factor, normalization.method = "RC")
#' Add data to SpatialExperiment and subset to HVGs
aname <- "normcounts"
assay(spe, aname) <- GetAssayData(seu)
spe <- spe[VariableFeatures(seu),]
Compute BANKSY neighborhood matrices. We use k_geom=18
corresponding to
first and second-order neighbors in 10x Visium.
compute_agf <- TRUE
k_geom <- 18
spe <- computeBanksy(spe, assay_name = aname, compute_agf = compute_agf, k_geom = k_geom)
Run PCA on the BANKSY matrix:
lambda <- 0.2
npcs <- 20
use_agf <- TRUE
spe <- runBanksyPCA(spe, use_agf = use_agf, lambda = lambda, npcs = npcs, seed = SEED)
We run Harmony on the PCs of the BANKSY matrix:
set.seed(SEED)
harmony_embedding <- HarmonyMatrix(
data_mat = reducedDim(spe, "PCA_M1_lam0.2"),
meta_data = colData(spe),
vars_use = c("sample_id", "subject_id"),
do_pca = FALSE,
max.iter.harmony = 20,
verbose = FALSE
)
reducedDim(spe, "Harmony_BANKSY") <- harmony_embedding
Next, run UMAP on the ‘raw’ and Harmony corrected PCA embeddings:
spe <- runBanksyUMAP(spe, use_agf = TRUE, lambda = lambda, npcs = npcs)
spe <- runBanksyUMAP(spe, dimred = "Harmony_BANKSY")
Visualize the UMAPs annotated by subject ID:
plot_grid(
plotReducedDim(spe, "UMAP_M1_lam0.2",
point_size = 0.1,
point_alpha = 0.5,
color_by = "subject_id") +
theme(legend.position = "none"),
plotReducedDim(spe, "UMAP_Harmony_BANKSY",
point_size = 0.1,
point_alpha = 0.5,
color_by = "subject_id") +
theme(legend.title = element_blank()) +
guides(colour = guide_legend(override.aes = list(size = 5, alpha = 1))),
nrow = 1,
rel_widths = c(1, 1.2)
)
Cluster the Harmony corrected PCA embedding:
spe <- clusterBanksy(spe, dimred = "Harmony_BANKSY", resolution = 0.55, seed = SEED)
spe <- connectClusters(spe, map_to = "clust_annotation")
Generate spatial plots:
cnm <- clusterNames(spe)[2]
spatial_plots <- lapply(sample_names, function(snm) {
x <- spe[, spe$sample_id == snm]
ari <- aricode::ARI(x$clust_annotation, colData(x)[, cnm])
df <- cbind.data.frame(clust=colData(x)[[cnm]], spatialCoords(x))
ggplot(df, aes(x=sdimy, y=sdimx, col=clust)) +
geom_point(size = 0.5) +
scale_color_manual(values = pals::kelly()[-1]) +
theme_classic() +
theme(
legend.position = "none",
axis.text.x=element_blank(),
axis.text.y=element_blank(),
axis.ticks=element_blank(),
axis.title.x=element_blank(),
axis.title.y=element_blank()) +
labs(title = sprintf("Sample %s - ARI: %s", snm, round(ari, 3))) +
coord_equal()
})
plot_grid(plotlist = spatial_plots, ncol = 3, byrow = FALSE)
Vignette runtime:
#> Time difference of 2.795952 mins
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] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] ExperimentHub_2.12.0 AnnotationHub_3.12.0
#> [3] BiocFileCache_2.12.0 dbplyr_2.5.0
#> [5] spatialLIBD_1.15.4 cowplot_1.1.3
#> [7] scater_1.32.0 ggplot2_3.5.1
#> [9] harmony_1.2.0 Rcpp_1.0.12
#> [11] data.table_1.15.4 scran_1.32.0
#> [13] scuttle_1.14.0 Seurat_5.0.3
#> [15] SeuratObject_5.0.1 sp_2.1-4
#> [17] SpatialExperiment_1.14.0 SingleCellExperiment_1.26.0
#> [19] SummarizedExperiment_1.34.0 Biobase_2.64.0
#> [21] GenomicRanges_1.56.0 GenomeInfoDb_1.40.0
#> [23] IRanges_2.38.0 S4Vectors_0.42.0
#> [25] BiocGenerics_0.50.0 MatrixGenerics_1.16.0
#> [27] matrixStats_1.3.0 Banksy_1.0.0
#> [29] BiocStyle_2.32.0
#>
#> loaded via a namespace (and not attached):
#> [1] bitops_1.0-7 spatstat.sparse_3.0-3
#> [3] doParallel_1.0.17 httr_1.4.7
#> [5] RColorBrewer_1.1-3 tools_4.4.0
#> [7] sctransform_0.4.1 DT_0.33
#> [9] utf8_1.2.4 R6_2.5.1
#> [11] lazyeval_0.2.2 uwot_0.2.2
#> [13] withr_3.0.0 gridExtra_2.3
#> [15] progressr_0.14.0 cli_3.6.2
#> [17] spatstat.explore_3.2-7 fastDummies_1.7.3
#> [19] labeling_0.4.3 sass_0.4.9
#> [21] spatstat.data_3.0-4 ggridges_0.5.6
#> [23] pbapply_1.7-2 Rsamtools_2.20.0
#> [25] dbscan_1.1-12 aricode_1.0.3
#> [27] dichromat_2.0-0.1 sessioninfo_1.2.2
#> [29] parallelly_1.37.1 attempt_0.3.1
#> [31] maps_3.4.2 limma_3.60.0
#> [33] pals_1.8 RSQLite_2.3.6
#> [35] BiocIO_1.14.0 generics_0.1.3
#> [37] ica_1.0-3 spatstat.random_3.2-3
#> [39] dplyr_1.1.4 Matrix_1.7-0
#> [41] ggbeeswarm_0.7.2 fansi_1.0.6
#> [43] abind_1.4-5 lifecycle_1.0.4
#> [45] yaml_2.3.8 edgeR_4.2.0
#> [47] SparseArray_1.4.0 Rtsne_0.17
#> [49] paletteer_1.6.0 grid_4.4.0
#> [51] blob_1.2.4 promises_1.3.0
#> [53] dqrng_0.3.2 crayon_1.5.2
#> [55] miniUI_0.1.1.1 lattice_0.22-6
#> [57] beachmat_2.20.0 mapproj_1.2.11
#> [59] KEGGREST_1.44.0 magick_2.8.3
#> [61] pillar_1.9.0 knitr_1.46
#> [63] metapod_1.12.0 rjson_0.2.21
#> [65] future.apply_1.11.2 codetools_0.2-20
#> [67] leiden_0.4.3.1 glue_1.7.0
#> [69] vctrs_0.6.5 png_0.1-8
#> [71] spam_2.10-0 gtable_0.3.5
#> [73] rematch2_2.1.2 cachem_1.0.8
#> [75] xfun_0.43 S4Arrays_1.4.0
#> [77] mime_0.12 survival_3.6-4
#> [79] RcppHungarian_0.3 iterators_1.0.14
#> [81] tinytex_0.50 fields_15.2
#> [83] statmod_1.5.0 bluster_1.14.0
#> [85] fitdistrplus_1.1-11 ROCR_1.0-11
#> [87] nlme_3.1-164 bit64_4.0.5
#> [89] filelock_1.0.3 RcppAnnoy_0.0.22
#> [91] bslib_0.7.0 irlba_2.3.5.1
#> [93] vipor_0.4.7 KernSmooth_2.23-22
#> [95] colorspace_2.1-0 DBI_1.2.2
#> [97] tidyselect_1.2.1 bit_4.0.5
#> [99] compiler_4.4.0 curl_5.2.1
#> [101] BiocNeighbors_1.22.0 DelayedArray_0.30.0
#> [103] plotly_4.10.4 rtracklayer_1.64.0
#> [105] bookdown_0.39 scales_1.3.0
#> [107] lmtest_0.9-40 rappdirs_0.3.3
#> [109] stringr_1.5.1 digest_0.6.35
#> [111] goftest_1.2-3 spatstat.utils_3.0-4
#> [113] rmarkdown_2.26 benchmarkmeData_1.0.4
#> [115] RhpcBLASctl_0.23-42 XVector_0.44.0
#> [117] htmltools_0.5.8.1 pkgconfig_2.0.3
#> [119] sparseMatrixStats_1.16.0 highr_0.10
#> [121] fastmap_1.1.1 rlang_1.1.3
#> [123] htmlwidgets_1.6.4 UCSC.utils_1.0.0
#> [125] shiny_1.8.1.1 DelayedMatrixStats_1.26.0
#> [127] farver_2.1.1 jquerylib_0.1.4
#> [129] zoo_1.8-12 jsonlite_1.8.8
#> [131] BiocParallel_1.38.0 mclust_6.1.1
#> [133] config_0.3.2 RCurl_1.98-1.14
#> [135] BiocSingular_1.20.0 magrittr_2.0.3
#> [137] GenomeInfoDbData_1.2.12 dotCall64_1.1-1
#> [139] patchwork_1.2.0 munsell_0.5.1
#> [141] viridis_0.6.5 reticulate_1.36.1
#> [143] leidenAlg_1.1.3 stringi_1.8.3
#> [145] zlibbioc_1.50.0 MASS_7.3-60.2
#> [147] plyr_1.8.9 parallel_4.4.0
#> [149] listenv_0.9.1 ggrepel_0.9.5
#> [151] deldir_2.0-4 Biostrings_2.72.0
#> [153] sccore_1.0.5 splines_4.4.0
#> [155] tensor_1.5 locfit_1.5-9.9
#> [157] igraph_2.0.3 spatstat.geom_3.2-9
#> [159] RcppHNSW_0.6.0 reshape2_1.4.4
#> [161] ScaledMatrix_1.12.0 XML_3.99-0.16.1
#> [163] BiocVersion_3.19.1 evaluate_0.23
#> [165] golem_0.4.1 BiocManager_1.30.22
#> [167] foreach_1.5.2 httpuv_1.6.15
#> [169] RANN_2.6.1 tidyr_1.3.1
#> [171] purrr_1.0.2 polyclip_1.10-6
#> [173] benchmarkme_1.0.8 future_1.33.2
#> [175] scattermore_1.2 rsvd_1.0.5
#> [177] xtable_1.8-4 restfulr_0.0.15
#> [179] RSpectra_0.16-1 later_1.3.2
#> [181] viridisLite_0.4.2 tibble_3.2.1
#> [183] GenomicAlignments_1.40.0 memoise_2.0.1
#> [185] beeswarm_0.4.0 AnnotationDbi_1.66.0
#> [187] cluster_2.1.6 shinyWidgets_0.8.6
#> [189] globals_0.16.3