Beyond Sequence-based Spatially-Resolved Data

Starting from Version 1.2.0, escheR package supports additional two data structures as input, including SpatialExperiment and data.frame from base R. In addition, escheR supports in-situ visualization of image-based spatially resolved data, which will be the focus of future development.

Visualized Dimensionality Reduced Embedding with SingleCellExperiment

SpatialExperiment inherits SingleCellExperiment

Following the same syntax, one can also visualize dimensionality reduced embeddings of a SpatialExperiment object by providing the argument dimred with a non-null value. Hence, the first 2 columns of the corresponding reducedDim(spe) assay will be used as the x-y coordinate of the plot, replacing spatialCoords(spe).

library(escheR)
library(STexampleData)
library(scater)
library(scran)

spe <- Visium_humanDLPFC() |> 
  logNormCounts()
spe <- spe[, spe$in_tissue == 1]
spe <- spe[, !is.na(spe$ground_truth)]
top.gene <- getTopHVGs(spe, n=500)

set.seed(100) # See below.
spe <- runPCA(spe, subset_row = top.gene) 

make_escheR(
  spe,
  dimred = "PCA"
) |> 
  add_fill(var = "ground_truth") +
  theme_minimal()

Hex Binning

spe$counts_MOBP <- counts(spe)[which(rowData(spe)$gene_name=="MOBP"),]
spe$ground_truth <- factor(spe$ground_truth)

# Point Binning version
make_escheR(
  spe,
  dimred = "PCA"
) |> 
  add_ground_bin(
    var = "ground_truth"
    ) |> 
  add_fill_bin(
    var = "counts_MOBP"
    ) +
  # Customize aesthetics
   scale_fill_gradient(low = "white", high = "black", name = "MOBP Count")+
  scale_color_discrete(name = "Spatial Domains") +
  theme_minimal()

Note 1: The strategy of binning to avoid overplotting is previously proposed in schex. While we provide an implementation in escheR, we would caution our users that the binning strategy could lead to intermixing of cluster memberships. In our implementation, the majority membership of the data points belonging to a bin is selected as the label of the bin. Users should use the binning strategy under their own discretion, and interpret the visualization carefully.

Note 2: add_fill_bin() shoudl be applied after add_ground_bin() for the better visualization outcome.

Image-based SpatialExperiment Object

To demonstrate the principle that escheR can be used to visualize image-based spatially-resolved data pending optimization, we include two image-based spatially resolved transcriptomics data generated via seqFish platform and Slide-seq V2 platform respectively. The two datasets have been previously curated in the STexampleData package

seqFISH

library(STexampleData)
library(escheR)
spe_seqFISH <- seqFISH_mouseEmbryo()

make_escheR(spe_seqFISH) |>
  add_fill(var = "embryo")

NOTE: trimming down the colData(spe) before piping into make-escheR could reduce the computation time to make the plots, specifically when colData(spe) contains extremely large number of irrelavent features/columns.

SlideSeqV2

library(STexampleData)
library(escheR)
spe_slideseq <- SlideSeqV2_mouseHPC()

make_escheR(spe_slideseq) |>
  add_fill(var = "celltype")

Beyond Bioconductor Eco-system

We aim to provide accessibility to all users regardless of their programming background and preferred single-cell analysis pipelines. Nevertheless , with limited resource, our sustaining efforts will prioritize towards the maintenance of the established functionality and the optimization for image-based spatially resolved data. We regret we are not be able to provide seamless interface to other R pipelines such as Seurat and Giotto in foreseeable future.

Instead, we provide a generic function that works with a data.frame object as input. For example, relevant features in Suerat can be easily exported as a data.frame object manually or via tidyseurat[https://github.com/stemangiola/tidyseurat]. The exported data frame can be pipe into escheR.

library(escheR)
library(Seurat)
pbmc_small <- SeuratObject::pbmc_small
pbmc_2pc <- pbmc_small@reductions$pca@cell.embeddings[,1:2]
pbmc_meta <- pbmc_small@meta.data

#> Call generic function for make_escheR.data.frame
make_escheR(
  object = pbmc_meta,
  .x = pbmc_2pc[,1],
  .y = pbmc_2pc[,2]) |> 
  add_fill(var = "groups")

Session information

utils::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_GB              LC_COLLATE=C              
#>  [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] BumpyMatrix_1.11.0          scran_1.31.3               
#>  [3] scater_1.31.2               scuttle_1.13.1             
#>  [5] ggpubr_0.6.0                STexampleData_1.11.3       
#>  [7] SpatialExperiment_1.13.2    SingleCellExperiment_1.25.1
#>  [9] SummarizedExperiment_1.33.3 Biobase_2.63.1             
#> [11] GenomicRanges_1.55.4        GenomeInfoDb_1.39.14       
#> [13] IRanges_2.37.1              S4Vectors_0.41.6           
#> [15] MatrixGenerics_1.15.1       matrixStats_1.3.0          
#> [17] ExperimentHub_2.11.3        AnnotationHub_3.11.5       
#> [19] BiocFileCache_2.11.2        dbplyr_2.5.0               
#> [21] BiocGenerics_0.49.1         escheR_1.3.2               
#> [23] ggplot2_3.5.0               BiocStyle_2.31.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] RColorBrewer_1.1-3        jsonlite_1.8.8           
#>   [3] magrittr_2.0.3            ggbeeswarm_0.7.2         
#>   [5] magick_2.8.3              farver_2.1.1             
#>   [7] rmarkdown_2.26            zlibbioc_1.49.3          
#>   [9] vctrs_0.6.5               memoise_2.0.1            
#>  [11] DelayedMatrixStats_1.25.3 rstatix_0.7.2            
#>  [13] tinytex_0.50              htmltools_0.5.8.1        
#>  [15] S4Arrays_1.3.7            curl_5.2.1               
#>  [17] BiocNeighbors_1.21.2      broom_1.0.5              
#>  [19] SparseArray_1.3.5         sass_0.4.9               
#>  [21] bslib_0.7.0               cachem_1.0.8             
#>  [23] igraph_2.0.3              mime_0.12                
#>  [25] lifecycle_1.0.4           pkgconfig_2.0.3          
#>  [27] rsvd_1.0.5                Matrix_1.7-0             
#>  [29] R6_2.5.1                  fastmap_1.1.1            
#>  [31] GenomeInfoDbData_1.2.12   digest_0.6.35            
#>  [33] colorspace_2.1-0          AnnotationDbi_1.65.2     
#>  [35] dqrng_0.3.2               irlba_2.3.5.1            
#>  [37] RSQLite_2.3.6             beachmat_2.19.4          
#>  [39] filelock_1.0.3            labeling_0.4.3           
#>  [41] fansi_1.0.6               httr_1.4.7               
#>  [43] abind_1.4-5               compiler_4.4.0           
#>  [45] bit64_4.0.5               withr_3.0.0              
#>  [47] backports_1.4.1           BiocParallel_1.37.1      
#>  [49] carData_3.0-5             viridis_0.6.5            
#>  [51] DBI_1.2.2                 hexbin_1.28.3            
#>  [53] highr_0.10                ggsignif_0.6.4           
#>  [55] rappdirs_0.3.3            DelayedArray_0.29.9      
#>  [57] rjson_0.2.21              bluster_1.13.0           
#>  [59] tools_4.4.0               vipor_0.4.7              
#>  [61] beeswarm_0.4.0            glue_1.7.0               
#>  [63] grid_4.4.0                cluster_2.1.6            
#>  [65] generics_0.1.3            gtable_0.3.4             
#>  [67] tidyr_1.3.1               metapod_1.11.1           
#>  [69] BiocSingular_1.19.0       ScaledMatrix_1.11.1      
#>  [71] car_3.1-2                 utf8_1.2.4               
#>  [73] XVector_0.43.1            ggrepel_0.9.5            
#>  [75] BiocVersion_3.19.1        pillar_1.9.0             
#>  [77] limma_3.59.8              dplyr_1.1.4              
#>  [79] lattice_0.22-6            bit_4.0.5                
#>  [81] tidyselect_1.2.1          locfit_1.5-9.9           
#>  [83] Biostrings_2.71.5         knitr_1.46               
#>  [85] gridExtra_2.3             bookdown_0.39            
#>  [87] edgeR_4.1.23              xfun_0.43                
#>  [89] statmod_1.5.0             UCSC.utils_0.99.7        
#>  [91] yaml_2.3.8                evaluate_0.23            
#>  [93] codetools_0.2-20          tibble_3.2.1             
#>  [95] BiocManager_1.30.22       cli_3.6.2                
#>  [97] munsell_0.5.1             jquerylib_0.1.4          
#>  [99] Rcpp_1.0.12               png_0.1-8                
#> [101] parallel_4.4.0            blob_1.2.4               
#> [103] sparseMatrixStats_1.15.1  viridisLite_0.4.2        
#> [105] scales_1.3.0              purrr_1.0.2              
#> [107] crayon_1.5.2              rlang_1.1.3              
#> [109] cowplot_1.1.3             KEGGREST_1.43.0