SpaTalk: Quick Start Guide
Zaoqu Liu
Maintainerliuzaoqu@163.com
Xin Shao
Original Developerxin_shao@zju.edu.cn
2026-01-23
Source:vignettes/SpaTalk.Rmd
SpaTalk.RmdIntroduction
SpaTalk is a computational framework for inferring spatially resolved cell-cell communications (CCIs) from spatial transcriptomics (ST) data. Published in Nature Communications (2022), SpaTalk integrates graph network modeling and knowledge graph approaches to reconstruct ligand-receptor-target signaling networks between spatially proximal cells.
Citation
Shao, X., Li, C., Yang, H., et al. Knowledge-graph-based cell-cell communication inference for spatially resolved transcriptomic data with SpaTalk. Nature Communications 13, 4429 (2022). https://doi.org/10.1038/s41467-022-32111-8
Installation
# From R-universe (recommended)
install.packages("SpaTalk", repos = "https://zaoqu-liu.r-universe.dev")
# From GitHub
devtools::install_github("Zaoqu-Liu/SpaTalk")Quick Start with STARmap Data
This tutorial uses the built-in STARmap mouse visual cortex data.
Step 1: Load Package and Data
library(SpaTalk)
# Load built-in demo data
load(system.file("extdata", "starmap_data.rda", package = "SpaTalk"))
load(system.file("extdata", "starmap_meta.rda", package = "SpaTalk"))
# Load curated databases
data(lrpairs)
data(pathways)
cat("Expression matrix:", nrow(starmap_data), "genes x", ncol(starmap_data), "cells\n")
#> Expression matrix: 996 genes x 930 cells
cat("Metadata:", nrow(starmap_meta), "cells\n")
#> Metadata: 930 cells
cat("Cell types:", paste(unique(starmap_meta$celltype), collapse = ", "), "\n")
#> Cell types: eL2_3, eL6, Astro, PVALB, Endo, VIP, SST, Smc, eL4, Micro, Oligo, eL5, Reln, HPCStep 2: Create SpaTalk Object
# Prepare spatial metadata
st_meta <- data.frame(
cell = starmap_meta$cell,
x = starmap_meta$x,
y = starmap_meta$y
)
# Create SpaTalk object (single-cell resolution data)
obj <- createSpaTalk(
st_data = starmap_data,
st_meta = st_meta,
species = "Mouse",
if_st_is_sc = TRUE,
spot_max_cell = 1,
celltype = starmap_meta$celltype
)
obj
#> An object of class SpaTalk
#> 996 genes across 930 single-cells (0 lrpair)Step 3: Visualize Spatial Distribution
plot_st_celltype_all(obj, size = 1.2)
Spatial distribution of all cell types
Step 4: Filter LR-Pathway Pairs
# Find LR pairs with downstream pathway targets
obj <- find_lr_path(
object = obj,
lrpairs = lrpairs,
pathways = pathways,
if_doParallel = FALSE
)
#> Checking input data
#> Begin to filter lrpairs and pathways
#> ***Done***
#>
cat("Filtered LR pairs:", nrow(obj@lr_path$lrpairs), "\n")
#> Filtered LR pairs: 6Step 5: Infer Cell-Cell Communications
# Infer CCIs between specific cell types
obj <- dec_cci(
object = obj,
celltype_sender = "eL6",
celltype_receiver = "PVALB",
if_doParallel = FALSE
)
#> Begin to find LR pairs
#>
# View results
if(nrow(obj@lrpair) > 0) {
cat("Found", nrow(obj@lrpair), "significant LR pairs:\n")
print(obj@lrpair[, c("ligand", "receptor", "lr_co_ratio", "lr_co_ratio_pvalue", "score")])
}
#> Found 1 significant LR pairs:
#> ligand receptor lr_co_ratio lr_co_ratio_pvalue score
#> 5 Inhba Acvr1c 0.1666667 0.002 0.8541642Step 6: Visualize Results
Cell-Cell Distance Distribution
plot_ccdist(
object = obj,
celltype_sender = "eL6",
celltype_receiver = "PVALB"
)
Distance distribution between eL6 and PVALB cells
LR Pair Spatial Distribution
if(nrow(obj@lrpair) > 0) {
lr <- obj@lrpair[1, ]
plot_lrpair(
object = obj,
ligand = lr$ligand,
receptor = lr$receptor,
celltype_sender = "eL6",
celltype_receiver = "PVALB",
size = 1.2
)
}
Spatial distribution of ligand-receptor interactions
Next Steps
- See the Algorithm vignette for methodological details
- See the Visualization vignette for all plotting options
- See the Advanced Usage vignette for custom databases and parallel processing
- See the Platforms vignette for platform-specific workflows
Session Info
sessionInfo()
#> R version 4.4.0 (2024-04-24)
#> Platform: aarch64-apple-darwin20
#> Running under: macOS 15.6.1
#>
#> Matrix products: default
#> BLAS: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRblas.0.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.0
#>
#> locale:
#> [1] C
#>
#> time zone: Asia/Shanghai
#> tzcode source: internal
#>
#> attached base packages:
#> [1] parallel stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] SpaTalk_2.0.0 doParallel_1.0.17 iterators_1.0.14 foreach_1.5.2
#> [5] ggalluvial_0.12.5 ggplot2_4.0.1
#>
#> loaded via a namespace (and not attached):
#> [1] RColorBrewer_1.1-3 jsonlite_2.0.0 magrittr_2.0.4
#> [4] spatstat.utils_3.2-1 farver_2.1.2 rmarkdown_2.30
#> [7] fs_1.6.6 ragg_1.5.0 vctrs_0.7.0
#> [10] ROCR_1.0-11 spatstat.explore_3.6-0 rstatix_0.7.3
#> [13] htmltools_0.5.9 progress_1.2.3 broom_1.0.11
#> [16] Formula_1.2-5 sass_0.4.10 sctransform_0.4.3
#> [19] parallelly_1.46.1 KernSmooth_2.23-26 bslib_0.9.0
#> [22] htmlwidgets_1.6.4 desc_1.4.3 ica_1.0-3
#> [25] plyr_1.8.9 plotly_4.11.0 zoo_1.8-15
#> [28] cachem_1.1.0 igraph_2.2.1 mime_0.13
#> [31] lifecycle_1.0.5 pkgconfig_2.0.3 Matrix_1.7-4
#> [34] R6_2.6.1 fastmap_1.2.0 fitdistrplus_1.2-4
#> [37] future_1.69.0 shiny_1.12.1 digest_0.6.39
#> [40] patchwork_1.3.2 Seurat_4.4.0 tensor_1.5.1
#> [43] irlba_2.3.5.1 textshaping_1.0.4 ggpubr_0.6.2
#> [46] labeling_0.4.3 progressr_0.18.0 spatstat.sparse_3.1-0
#> [49] httr_1.4.7 polyclip_1.10-7 abind_1.4-8
#> [52] compiler_4.4.0 withr_3.0.2 backports_1.5.0
#> [55] S7_0.2.1 carData_3.0-5 ggforce_0.5.0
#> [58] ggsignif_0.6.4 MASS_7.3-65 rappdirs_0.3.4
#> [61] ggsci_4.2.0 tools_4.4.0 lmtest_0.9-40
#> [64] otel_0.2.0 scatterpie_0.2.6 httpuv_1.6.16
#> [67] future.apply_1.20.1 goftest_1.2-3 glue_1.8.0
#> [70] nlme_3.1-168 promises_1.5.0 grid_4.4.0
#> [73] Rtsne_0.17 cluster_2.1.8.1 reshape2_1.4.5
#> [76] generics_0.1.4 gtable_0.3.6 spatstat.data_3.1-9
#> [79] tzdb_0.5.0 tidyr_1.3.2 data.table_1.18.0
#> [82] hms_1.1.4 car_3.1-3 sp_2.2-0
#> [85] spatstat.geom_3.6-1 RcppAnnoy_0.0.23 ggrepel_0.9.6
#> [88] RANN_2.6.2 pillar_1.11.1 stringr_1.6.0
#> [91] yulab.utils_0.2.3 ggExtra_0.11.0 later_1.4.5
#> [94] splines_4.4.0 tweenr_2.0.3 dplyr_1.1.4
#> [97] lattice_0.22-7 survival_3.8-3 deldir_2.0-4
#> [100] tidyselect_1.2.1 miniUI_0.1.2 pbapply_1.7-4
#> [103] knitr_1.51 gridExtra_2.3 scattermore_1.2
#> [106] xfun_0.56 matrixStats_1.5.0 pheatmap_1.0.13
#> [109] stringi_1.8.7 ggfun_0.2.0 lazyeval_0.2.2
#> [112] yaml_2.3.12 evaluate_1.0.5 codetools_0.2-20
#> [115] tibble_3.3.1 cli_3.6.5 uwot_0.2.4
#> [118] xtable_1.8-4 reticulate_1.44.1 systemfonts_1.3.1
#> [121] jquerylib_0.1.4 dichromat_2.0-0.1 Rcpp_1.1.1
#> [124] globals_0.18.0 spatstat.random_3.4-3 png_0.1-8
#> [127] spatstat.univar_3.1-6 readr_2.1.6 pkgdown_2.1.3
#> [130] NNLM_0.4.4 prettyunits_1.2.0 listenv_0.10.0
#> [133] viridisLite_0.4.2 scales_1.4.0 ggridges_0.5.7
#> [136] SeuratObject_4.1.4 leiden_0.4.3.1 purrr_1.2.1
#> [139] crayon_1.5.3 rlang_1.1.7 cowplot_1.2.0