Contents

1 Introduction

The goal of escheR is to create an unified multi-dimensional spatial visualizations for spatially-resolved transcriptomics data following Gestalt principles.

Our preprint describing the innovative visualization is available from bioRxiv.

1.1 Installation

You can install the latest release version of escheR from Bioconductor. with using the following code will install version of the nnSVG package from Bioconductor. Additional details are shown on the Bioconductor page.

if (!require("BiocManager", quietly = TRUE)) {
    install.packages("BiocManager")
}
BiocManager::install("escheR")

The latest development version can also be installed from the devel version of Bioconductor or from GitHub following

if (!require("devtools")) install.packages("devtools")
devtools::install_github("boyiguo1/escheR")

2 Input data format

In the examples below, we assume the input data are provided as a SpatialExperiment Bioconductor object. For people whose data are stored as a Seurat object, we advise to convert to a SpatialExperiment object before applying the workflow below.

3 Tutorial

3.1 Load Packages

To run the demonstration, there are two necessary packages to load, escheR and spatialLIBD. spatialLIBD contains a pre-processed 10x Visium dataset.

To note, escheR will automatically load ggplot2 package. Hence, explicitly loading ggplot2 is not required.

library(escheR)
#> Warning: replacing previous import 'utils::findMatches' by
#> 'S4Vectors::findMatches' when loading 'AnnotationDbi'
library(STexampleData)
library(spatialLIBD)

3.2 Preparing example data

In this step, we will find one 10xVisium sample from STexampleData package, indexed by brain number of “151673”. For more information, please see the vignettes of STexampleData.

spe <- Visium_humanDLPFC()

# Subset in-tissue spots
spe <- spe[, spe$in_tissue == 1]

Here is a summary of the SpatialExperiment object called spe.

spe
#> class: SpatialExperiment 
#> dim: 33538 3639 
#> metadata(0):
#> assays(1): counts
#> rownames(33538): ENSG00000243485 ENSG00000237613 ... ENSG00000277475
#>   ENSG00000268674
#> rowData names(3): gene_id gene_name feature_type
#> colnames(3639): AAACAAGTATCTCCCA-1 AAACAATCTACTAGCA-1 ...
#>   TTGTTTGTATTACACG-1 TTGTTTGTGTAAATTC-1
#> colData names(7): barcode_id sample_id ... ground_truth cell_count
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
#> imgData names(4): sample_id image_id data scaleFactor

3.3 Setting up for escheR plot

Similar to ggplot2::ggplot(), we first use the function make_escheR() to create an empty plot. The input of make_escheR() is a SpatialExperiment object. The output of the function is a ggplot object with no layer in it.

p <- make_escheR(spe)

3.4 Adding layers

Unlike ggplot2, we use piping |> instead of + to apply layers the figure. Mainly, we have three functions add_fill, add_ground, add_symbol. The inputs of these add_* functions include the plots created using make_scheR() and the variable name for the layer. Currently, the variable name should be in the the column data of the spe object, i.e. colData(spe).

Here we first apply the add_fill to add the spots color-coded by the total number of cells all spots(sum_umi).

(p1 <- p |>
    add_fill(var = "cell_count"))


p1 +
    scale_fill_viridis_c()

It is okay to use any combination of the add_* functions. For example, we want to show the spatial domains of the samples as the ground of the figure and use symbols to denote if each spot is within the outline of the tissue slice. In this example, all plotted spots are in the outlines of the tissue slice and hence marked with dots.

(p2 <- p |>
    add_ground(var = "ground_truth")) # round layer

p2 |>
    add_symbol(var = "ground_truth", size = 0.2) # Symbol layer
#> Warning: The shape palette can deal with a maximum of 6 discrete values because
#> more than 6 becomes difficult to discriminate; you have 7. Consider
#> specifying shapes manually if you must have them.
#> Warning: Removed 541 rows containing missing values (`geom_point()`).

It is okay to change the ordering of these add_* functions. However, we advise to always have the add_fill as the first step to achieve the best visual effect due to the laying mechanism.

3.5 Adjusting aesthetics

To change the aesthetics of each layer, one can simply use the scale_* from ggplot2 to optimize the final visual presentation. For example, to optimize add_fill, one can use scale_fill_*; to optimize add_ground, one can use scale_color_*; to optimize add_sumbol, one use scale_shape_*. Here, we demonstrate how to change the color for the ground layer ( add_ground) using scale_color_manual.

# Currated color pallette from spatialLIBD
spatialLIBD::libd_layer_colors
#>        Layer1        Layer2        Layer3        Layer4        Layer5 
#>     "#F0027F"     "#377EB8"     "#4DAF4A"     "#984EA3"     "#FFD700" 
#>        Layer6            WM            NA           WM2 
#>     "#FF7F00"     "#1A1A1A" "transparent"     "#666666"

p2 +
    scale_color_manual(
        name = "", # No legend name
        values = spatialLIBD::libd_layer_colors
    ) +
    labs(title = "Example Title")

4 Session information

sessionInfo()
#> R version 4.3.0 RC (2023-04-13 r84269)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.2 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.17-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] spatialLIBD_1.11.14         STexampleData_1.7.0        
#>  [3] SpatialExperiment_1.10.0    SingleCellExperiment_1.22.0
#>  [5] SummarizedExperiment_1.30.0 Biobase_2.60.0             
#>  [7] GenomicRanges_1.52.0        GenomeInfoDb_1.36.0        
#>  [9] IRanges_2.34.0              S4Vectors_0.38.0           
#> [11] MatrixGenerics_1.12.0       matrixStats_0.63.0         
#> [13] ExperimentHub_2.8.0         AnnotationHub_3.8.0        
#> [15] BiocFileCache_2.8.0         dbplyr_2.3.2               
#> [17] BiocGenerics_0.46.0         escheR_1.0.0               
#> [19] ggplot2_3.4.2               BiocStyle_2.28.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] later_1.3.0                   BiocIO_1.10.0                
#>   [3] bitops_1.0-7                  filelock_1.0.2               
#>   [5] fields_14.1                   tibble_3.2.1                 
#>   [7] R.oo_1.25.0                   XML_3.99-0.14                
#>   [9] lifecycle_1.0.3               edgeR_3.42.0                 
#>  [11] doParallel_1.0.17             lattice_0.21-8               
#>  [13] magrittr_2.0.3                limma_3.56.0                 
#>  [15] plotly_4.10.1                 sass_0.4.5                   
#>  [17] rmarkdown_2.21                jquerylib_0.1.4              
#>  [19] yaml_2.3.7                    httpuv_1.6.9                 
#>  [21] spam_2.9-1                    sessioninfo_1.2.2            
#>  [23] cowplot_1.1.1                 DBI_1.1.3                    
#>  [25] RColorBrewer_1.1-3            golem_0.4.0                  
#>  [27] maps_3.4.1                    zlibbioc_1.46.0              
#>  [29] purrr_1.0.1                   R.utils_2.12.2               
#>  [31] RCurl_1.98-1.12               rappdirs_0.3.3               
#>  [33] GenomeInfoDbData_1.2.10       ggrepel_0.9.3                
#>  [35] irlba_2.3.5.1                 dqrng_0.3.0                  
#>  [37] DelayedMatrixStats_1.22.0     codetools_0.2-19             
#>  [39] DropletUtils_1.20.0           DelayedArray_0.26.0          
#>  [41] DT_0.27                       scuttle_1.10.0               
#>  [43] tidyselect_1.2.0              farver_2.1.1                 
#>  [45] ScaledMatrix_1.8.0            viridis_0.6.2                
#>  [47] shinyWidgets_0.7.6            GenomicAlignments_1.36.0     
#>  [49] jsonlite_1.8.4                BiocNeighbors_1.18.0         
#>  [51] ellipsis_0.3.2                scater_1.28.0                
#>  [53] iterators_1.0.14              foreach_1.5.2                
#>  [55] tools_4.3.0                   Rcpp_1.0.10                  
#>  [57] glue_1.6.2                    gridExtra_2.3                
#>  [59] xfun_0.39                     dplyr_1.1.2                  
#>  [61] HDF5Array_1.28.0              withr_2.5.0                  
#>  [63] BiocManager_1.30.20           fastmap_1.1.1                
#>  [65] rhdf5filters_1.12.0           fansi_1.0.4                  
#>  [67] digest_0.6.31                 rsvd_1.0.5                   
#>  [69] R6_2.5.1                      mime_0.12                    
#>  [71] colorspace_2.1-0              RSQLite_2.3.1                
#>  [73] R.methodsS3_1.8.2             config_0.3.1                 
#>  [75] utf8_1.2.3                    tidyr_1.3.0                  
#>  [77] generics_0.1.3                data.table_1.14.8            
#>  [79] rtracklayer_1.60.0            httr_1.4.5                   
#>  [81] htmlwidgets_1.6.2             pkgconfig_2.0.3              
#>  [83] gtable_0.3.3                  blob_1.2.4                   
#>  [85] XVector_0.40.0                htmltools_0.5.5              
#>  [87] dotCall64_1.0-2               bookdown_0.33                
#>  [89] scales_1.2.1                  png_0.1-8                    
#>  [91] attempt_0.3.1                 knitr_1.42                   
#>  [93] rjson_0.2.21                  curl_5.0.0                   
#>  [95] cachem_1.0.7                  rhdf5_2.44.0                 
#>  [97] BiocVersion_3.17.1            vipor_0.4.5                  
#>  [99] parallel_4.3.0                AnnotationDbi_1.62.0         
#> [101] restfulr_0.0.15               pillar_1.9.0                 
#> [103] grid_4.3.0                    vctrs_0.6.2                  
#> [105] promises_1.2.0.1              BiocSingular_1.16.0          
#> [107] beachmat_2.16.0               xtable_1.8-4                 
#> [109] beeswarm_0.4.0                paletteer_1.5.0              
#> [111] evaluate_0.20                 magick_2.7.4                 
#> [113] cli_3.6.1                     locfit_1.5-9.7               
#> [115] compiler_4.3.0                Rsamtools_2.16.0             
#> [117] rlang_1.1.0                   crayon_1.5.2                 
#> [119] labeling_0.4.2                rematch2_2.1.2               
#> [121] ggbeeswarm_0.7.1              viridisLite_0.4.1            
#> [123] BiocParallel_1.34.0           munsell_0.5.0                
#> [125] Biostrings_2.68.0             lazyeval_0.2.2               
#> [127] Matrix_1.5-4                  benchmarkme_1.0.8            
#> [129] sparseMatrixStats_1.12.0      bit64_4.0.5                  
#> [131] Rhdf5lib_1.22.0               statmod_1.5.0                
#> [133] KEGGREST_1.40.0               shiny_1.7.4                  
#> [135] highr_0.10                    interactiveDisplayBase_1.38.0
#> [137] memoise_2.0.1                 bslib_0.4.2                  
#> [139] benchmarkmeData_1.0.4         bit_4.0.5