Main

After pre-processing, we create a SpaTalk object with st_data and st_meta for spot-based ST data. In practise, we generate simulated spot-based ST from single-cell STARmap data using generate_spot(). For Slide-seq (10um) close to single-cell resolution, we recommend spot_max_cell as 1. For 10X Visum (55um), we recommend spot_max_cell as 30 as stated in the recent review.

library(SpaTalk)
# load starmap data
load(paste0(system.file(package = 'SpaTalk'), "/extdata/starmap_data.rda"))
load(paste0(system.file(package = 'SpaTalk'), "/extdata/starmap_meta.rda"))

# generate spot-based ST from single-cell STARmap data
st_data <- generate_spot(st_data = as.matrix(starmap_data),
                         st_meta = starmap_meta,
                         x_min = 0,
                         x_res = 200,
                         x_max = 3800,
                         y_min = 0,
                         y_res = 500,
                         y_max = 17500)
## Using celltype as value column: use value.var to override.
## Aggregation function missing: defaulting to length
st_meta <- st_data$st_meta
st_data <- st_data$st_data

# create SpaTalk data
obj <- createSpaTalk(st_data = as.matrix(st_data),
                     st_meta = st_meta[, c(1,2,3)],
                     species = "Mouse",
                     if_st_is_sc = F,
                     spot_max_cell = max(st_meta$cell_real))

Use plot_st_pie_generate() to view the cell-type composition

plot_st_pie_generate(st_meta = st_meta, pie_scale = 1.3, xy_ratio = 1.8)

Cell-type decomposition

First, we apply NNLM and spatial mapping to reconstruct the ST data at single-cell resolution using dec_celltype()

# decode the cell-type composition
obj <- dec_celltype(object = obj,
                    sc_data = as.matrix(starmap_data),
                    sc_celltype = starmap_meta$celltype)
## Performing Non-negative regression for each spot 
## Generating single-cell data for each spot 
## ***Done***

Use plot_st_pie() to view predicted cell-type composition

# Scatter pie plot for each spot
plot_st_pie(object = obj, pie_scale = 1.3, xy_ratio = 1.8)

Use plot_st_celltype_percent() to view cell-type percent

# plot cell-type percent across spatial spots
plot_st_celltype_percent(object = obj, celltype = 'Oligo',size = 4)

Use plot_st_gene() to view gene expression

# plot marker gene expression across spatial spots
plot_st_gene(object = obj, gene = 'Plp1',size = 4, if_use_newmeta = F)

Use plot_st_cor_heatmap() to view correlation heatmap

# correlation between marker gene expression and cell type percent across spatial spots
plot_st_cor_heatmap(object = obj,
                    marker_genes = c("Plp1","Vip","Sst","Lamp5","Pcp4","Mfge8","Pvalb"),
                    celltypes = c("Oligo","VIP","SST","eL2_3","eL6","Astro","PVALB"),
                    scale = "none",
                    color_low = 'blue',
                    color_high = 'yellow',
                    color_mid = 'yellow')

ST at single-cell resolution

Use plot_st_celltype() to view cell type in space

# plot cell type in reconstructed ST atlas
plot_st_celltype(object = obj, celltype = 'Oligo', size = 2)

Use plot_st_gene() to view gene expression in space

# plot marker gene expression in reconstructed ST atlas
plot_st_gene(object = obj,gene = 'Plp1', if_use_newmeta = T, size = 2)

Use plot_st_celltype_density() to view cell-type density in space

# plot cell-type density in reconstructed ST atlas
plot_st_celltype_density(object = obj,
                         celltype = 'Oligo',
                         type = 'raster',
                         color_low = 'purple',
                         color_high = 'yellow')

plot_st_celltype_density(object = obj,
                         celltype = 'Oligo',
                         type = 'contour',
                         color_low = 'purple',
                         color_high = 'yellow')

Use plot_st_celltype_all() to view all cell types in space

# plot all cell types in reconstructed ST atlas
plot_st_celltype_all(object = obj, size = 2)

Infer cell-cell communications

Based on the constructed ST atlas at single-cell resolution, we use find_lr_path() and dec_cci() to infer cell-cell communications mediated by ligand-receptor interactions (LRIs) in space using LRIs from CellTalkDB, pathways from KEGG and Reactome, and transcriptional factors (TFs) from AnimalTFDB

# Filter LRIs with downstream targets
obj <- find_lr_path(object = obj, lrpairs = lrpairs, pathways = pathways)
## Checking input data 
## Begin to filter lrpairs and pathways 
## ***Done***
# Infer cell-cell communications from SST to PVALB neurons
obj <- dec_cci(object = obj, 
               celltype_sender = 'SST',
               celltype_receiver = 'PVALB')
## Begin to find LR pairs
# Get LR and downstream pathways
obj_lr_path <- get_lr_path(object = obj,
                           celltype_sender = 'SST',
                           celltype_receiver = 'PVALB',
                           ligand = 'Sst',
                           receptor = 'Sstr2')

obj_lr_path$tf_path
##          src     dest src_tf dest_tf hop  co_ratio    tf
## 909510 Sstr2    Gnai1     NO      NO   1 0.3076923 Smad3
## 268551 Gnai1    Mapk3     NO      NO   2 0.1538462 Smad3
## 300531 Mapk3    Smad3     NO     YES   3 0.1538462 Smad3
## 327207 Smad3 Serpine1    YES      NO   4 0.2307692 Smad3
obj_lr_path$path_pvalue
##   celltype_sender celltype_receiver ligand receptor
## 1             SST             PVALB    Sst    Sstr2
## 2             SST             PVALB    Sst    Sstr2
## 3             SST             PVALB    Sst    Sstr2
##                                receptor_pathways       pvalue gene_count
## 1                         cAMP signaling pathway 1.009418e-05          3
## 2        Neuroactive ligand-receptor interaction 6.039598e-03          2
## 3 Growth hormone synthesis, secretion and action 1.281881e-08          4
##                    gene
## 1       Sst,Sstr2,Gnai1
## 2             Sst,Sstr2
## 3 Sst,Sstr2,Gnai1,Mapk3

Use plot_ccdist() to view distribution of senders and receivers

# Point plot with spatial distribution of celltype_sender and celltype_receiver
plot_ccdist(object = obj, 
            celltype_sender = 'SST', 
            celltype_receiver = 'PVALB',            
            size = 2,
            arrow_length = 0.1)

Use plot_cci_lrpairs() to view all LRIs of senders and receivers

Given the limited LR pairs, we do not show the result here

# Heatmap with LR pairs of celltype_sender and celltype_receiver
# plot_cci_lrpairs(object = obj, celltype_sender = 'SST', celltype_receiver = 'PVALB')

Use plot_lrpair() to view the specific LRI

# Point plot with LR pair from celltype_sender to celltype_receiver
plot_lrpair(object = obj,
            celltype_sender = 'SST',
            ligand = 'Sst',
            celltype_receiver = 'PVALB',
            receptor = 'Sstr2',
            if_plot_density = F,            
            size = 2,
            arrow_length = 0.1)

Use plot_lrpair_vln() to view violin plot with spatial distance of LRI

# Violin plot with spatial distance of LR pair between senders and receivers and between all cell-cell pairs
plot_lrpair_vln(object = obj,
                celltype_sender = 'SST',
                ligand = 'Sst',
                celltype_receiver = 'PVALB',
                receptor = 'Sstr2')
## 
##  Wilcoxon rank sum test with continuity correction
## 
## data:  celltype_dist1 and celltype_dist2
## W = 5685812, p-value = 5.395e-13
## alternative hypothesis: true location shift is greater than 0

Use plot_lr_path() to view network with LR and downstream pathways

# Plot network with LR and downstream pathways
plot_lr_path(object = obj,                
             celltype_sender = 'SST',
             ligand = 'Sst',
             celltype_receiver = 'PVALB',
             receptor = 'Sstr2')

Use plot_path2gene() to view river plot of significantly activated pathways

# River plot of significantly activated pathways and related downstream genes of receptors
plot_path2gene(object = obj,                
               celltype_sender = 'SST',
               ligand = 'Sst',
               celltype_receiver = 'PVALB',
               receptor = 'Sstr2')
## Warning: The `.dots` argument of `group_by()` is deprecated as of dplyr 1.0.0.
## This warning is displayed once every 8 hours.
## Call `lifecycle::last_lifecycle_warnings()` to see where this warning was generated.

Note

To infer all paired cell-cell communications, use dec_cci_all() instead of dec_cci()

# Infer all cell-cell communications
# obj <- dec_cci_all(object = obj)
sessionInfo()
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-w64-mingw32/x64 (64-bit)
## Running under: Windows 10 x64 (build 19044)
## 
## Matrix products: default
## 
## locale:
## [1] LC_COLLATE=Chinese (Simplified)_China.936 
## [2] LC_CTYPE=Chinese (Simplified)_China.936   
## [3] LC_MONETARY=Chinese (Simplified)_China.936
## [4] LC_NUMERIC=C                              
## [5] LC_TIME=Chinese (Simplified)_China.936    
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] SpaTalk_1.0       ggalluvial_0.12.3 ggplot2_3.3.5    
## 
## loaded via a namespace (and not attached):
##   [1] Seurat_4.0.5          Rtsne_0.15            colorspace_2.0-2     
##   [4] ggsignif_0.6.3        deldir_1.0-6          ellipsis_0.3.2       
##   [7] ggridges_0.5.3        spatstat.data_2.1-0   ggpubr_0.4.0         
##  [10] leiden_0.3.9          listenv_0.8.0         farver_2.1.0         
##  [13] ggrepel_0.9.1         fansi_0.5.0           scatterpie_0.1.7     
##  [16] codetools_0.2-18      splines_4.1.1         knitr_1.36           
##  [19] polyclip_1.10-0       jsonlite_1.7.2        broom_0.7.10         
##  [22] ica_1.0-2             cluster_2.1.2         png_0.1-7            
##  [25] uwot_0.1.10           pheatmap_1.0.12       spatstat.sparse_2.0-0
##  [28] ggforce_0.3.3         shiny_1.7.1           sctransform_0.3.2    
##  [31] compiler_4.1.1        httr_1.4.2            backports_1.3.0      
##  [34] assertthat_0.2.1      SeuratObject_4.0.2    Matrix_1.3-4         
##  [37] fastmap_1.1.0         lazyeval_0.2.2        later_1.3.0          
##  [40] tweenr_1.0.2          htmltools_0.5.2       prettyunits_1.1.1    
##  [43] tools_4.1.1           igraph_1.2.7          gtable_0.3.0         
##  [46] glue_1.4.2            RANN_2.6.1            reshape2_1.4.4       
##  [49] dplyr_1.0.7           Rcpp_1.0.7            carData_3.0-4        
##  [52] scattermore_0.7       jquerylib_0.1.4       vctrs_0.3.8          
##  [55] nlme_3.1-152          lmtest_0.9-38         xfun_0.30            
##  [58] stringr_1.4.0         globals_0.14.0        mime_0.12            
##  [61] miniUI_0.1.1.1        lifecycle_1.0.1       irlba_2.3.3          
##  [64] rstatix_0.7.0         goftest_1.2-3         NNLM_0.4.4           
##  [67] future_1.23.0         MASS_7.3-54           zoo_1.8-9            
##  [70] scales_1.1.1          spatstat.core_2.3-0   spatstat.utils_2.2-0 
##  [73] hms_1.1.1             promises_1.2.0.1      parallel_4.1.1       
##  [76] RColorBrewer_1.1-2    yaml_2.2.1            gridExtra_2.3        
##  [79] reticulate_1.22       pbapply_1.5-0         ggfun_0.0.4          
##  [82] sass_0.4.0            rpart_4.1-15          ggExtra_0.9          
##  [85] stringi_1.7.5         highr_0.9             rlang_0.4.12         
##  [88] pkgconfig_2.0.3       matrixStats_0.61.0    evaluate_0.14        
##  [91] lattice_0.20-44       tensor_1.5            ROCR_1.0-11          
##  [94] purrr_0.3.4           labeling_0.4.2        patchwork_1.1.1      
##  [97] htmlwidgets_1.5.4     cowplot_1.1.1         tidyselect_1.1.1     
## [100] ggsci_2.9             parallelly_1.28.1     RcppAnnoy_0.0.19     
## [103] plyr_1.8.6            magrittr_2.0.1        R6_2.5.1             
## [106] generics_0.1.1        DBI_1.1.1             mgcv_1.8-36          
## [109] pillar_1.6.4          withr_2.4.2           fitdistrplus_1.1-6   
## [112] prettydoc_0.4.1       abind_1.4-5           survival_3.2-11      
## [115] tibble_3.1.5          future.apply_1.8.1    car_3.0-12           
## [118] crayon_1.4.2          KernSmooth_2.23-20    utf8_1.2.2           
## [121] spatstat.geom_2.3-0   plotly_4.10.0         rmarkdown_2.13       
## [124] progress_1.2.2        isoband_0.2.5         grid_4.1.1           
## [127] data.table_1.14.2     digest_0.6.28         xtable_1.8-4         
## [130] tidyr_1.1.4           httpuv_1.6.3          munsell_0.5.0        
## [133] viridisLite_0.4.0     bslib_0.3.1