Skip to contents

scider is an user-friendly R package providing functions to model the global density of cells in a slide of spatial transcriptomics data. All functions in the package are built based on the SpatialExperiment object, allowing integration into various spatial transcriptomics-related packages from Bioconductor. After modelling density, the package allows for serveral downstream analysis, including colocalization analysis, boundary detection analysis and differential density analysis.

Installation

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

BiocManager::install("scider")

The development version of scider can be installed from GitHub:

devtools::install_github("ChenLaboratory/scider")

Load data

In this vignette, we will use a subset of a Xenium Breast Cancer dataset.

data("xenium_bc_spe")

In the data, we have quantification of 541 genes from 10000 cells.

spe
## class: SpatialExperiment 
## dim: 541 10000 
## metadata(0):
## assays(1): counts
## rownames(541): ENSG00000121270 ENSG00000213088 ... BLANK_0444
##   BLANK_0447
## rowData names(3): ID Symbol Type
## colnames(10000): cell_212124 cell_120108 ... cell_252054 cell_568560
## colData names(21): cell_id transcript_counts ... cell_type sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : x_centroid y_centroid
## imgData names(1): sample_id

We also have cell-type annotations of these cells, there are 4 cell types.

table(colData(spe)$cell_type)
## 
##       B cells Breast cancer   Fibroblasts       T cells 
##           643          3550          4234          1573

We can use the function plotSpatial to visualise the cell position and color the cells by cell types.

plotSpatial(spe, group.by = "cell_type", pt.alpha = 0.8)

Grid-based analysis

scider can conduct grid-based density analysis for spatial transcriptomics data.

Density calculation

We can perform density calculation for each cell type using function gridDensity. The calculated density and grid information are saved in the metadata of the SpatialExperimnet object.

spe <- gridDensity(spe)
names(metadata(spe))
## [1] "grid_density" "grid_info"
metadata(spe)$grid_density
## DataFrame with 12319 rows and 10 columns
##          x_grid    y_grid    node_x    node_y        node density_b_cells
##       <numeric> <numeric> <numeric> <numeric> <character>       <numeric>
## 1       280.937   155.873         1         1         1-1      0.00127779
## 2       380.937   155.873         2         1         2-1      0.00170764
## 3       480.937   155.873         3         1         3-1      0.00227219
## 4       580.937   155.873         4         1         4-1      0.00302542
## 5       680.937   155.873         5         1         5-1      0.00405470
## ...         ...       ...       ...       ...         ...             ...
## 12315   9480.94   11067.8        93       127      93-127     0.001167782
## 12316   9580.94   11067.8        94       127      94-127     0.000771056
## 12317   9680.94   11067.8        95       127      95-127     0.000491317
## 12318   9780.94   11067.8        96       127      96-127     0.000303461
## 12319   9880.94   11067.8        97       127      97-127     0.000182427
##       density_breast_cancer density_fibroblasts density_t_cells density_overall
##                   <numeric>           <numeric>       <numeric>       <numeric>
## 1               4.03130e-05          0.00229935     0.000785227      0.00440268
## 2               7.95323e-05          0.00334136     0.001181296      0.00630982
## 3               1.55385e-04          0.00478370     0.001766524      0.00897779
## 4               2.98294e-04          0.00673066     0.002620974      0.01267536
## 5               5.58759e-04          0.00928390     0.003852121      0.01774948
## ...                     ...                 ...             ...             ...
## 12315            0.01224824          0.01392798     2.05759e-05      0.02736458
## 12316            0.00767062          0.00936585     9.86410e-06      0.01781739
## 12317            0.00476467          0.00612204     4.80806e-06      0.01138284
## 12318            0.00294447          0.00390776     2.43691e-06      0.00715813
## 12319            0.00181413          0.00244662     1.30538e-06      0.00444449

We can visualise the overall cell density of the whole tissue using function plotDensity.

We can also visualise the density of individual cell type, e.g., fibroblast cells.

plotDensity(spe, coi = "Fibroblasts")

Find Regions-of-interest (ROIs)

After obtaining grid-based density for each COI, we can then detect regions-of-interest (ROIs) based on density or select by user.

Detected by algorithm

To detect ROIs automatically, we can use the function findROI.

The detected ROIs are saved in the metadata of the SpatialExperiment object.

Here we identify ROIs based on the fibroblasts cell density.

spe <- findROI(spe, coi = "Fibroblasts")
metadata(spe)$roi
## DataFrame with 1833 rows and 6 columns
##      component     members           x           y    xcoord    ycoord
##       <factor> <character> <character> <character> <numeric> <numeric>
## 1            1       38-34          38          34   4030.94   3013.76
## 2            1       39-34          39          34   4130.94   3013.76
## 3            1       40-34          40          34   4230.94   3013.76
## 4            1       41-34          41          34   4330.94   3013.76
## 5            1       42-34          42          34   4430.94   3013.76
## ...        ...         ...         ...         ...       ...       ...
## 1829         8      68-124          68         124   7030.94  10807.99
## 1830         8      69-124          69         124   7130.94  10807.99
## 1831         8      70-124          70         124   7230.94  10807.99
## 1832         8      78-100          78         100   8030.94   8729.52
## 1833         8      69-117          69         117   7080.94  10201.77

We can visualise the ROIs with function plotROI.

plotROI(spe)

Select ROI by user

Alternatively, users can select ROIs based on their own research interest (drawn by hand). This can be done using function selectRegion. This function will open an interactive window with an interactive plot for users to zoom-in/-out and select ROI using either a rectangular or lasso selection tool. Users can also press the Export selected points button to save the ROIs as object in the R environment.

selectRegion(metadata(spe)$grid_density, x.col = "x_grid", y.col = "y_grid")

After closing the interactive window, the selected ROI has been saved as a data.frame object named sel_region in the R environment.

sel_region

We can then use the postSelRegion to save the ROI in the metadata of the SpatialExperiment object.

spe1 <- postSelRegion(spe, sel_region = sel_region)

metadata(spe1)$roi

Similarly, we can plot visualise the user-defined ROI with function plotROI.

plotROI(spe1)

Testing relationship between cell types

After defining ROIs, we can then test the relationship between any two cell types within each ROI or overall but account for ROI variation using a cubic spline or a linear fit.

This can be done with function corrDensity, by setting the celltype1 and celltype2 parameters, the modelling results are saved in the metadata of the SpatialExperiment object.

results <- corDensity(spe)

We can see the correlation between breast cancer and fibroblasts in each ROI.

results$ROI
## DataFrame with 48 rows and 9 columns
##       celltype1     celltype2         ROI     ngrid    cor.coef          t
##     <character>   <character> <character> <numeric>   <numeric>  <numeric>
## 1       B cells Breast cancer           1        31    0.338005   0.755462
## 2       B cells Breast cancer           2       160   -0.913988  -2.483623
## 3       B cells Breast cancer           3       325   -0.482294  -1.722597
## 4       B cells Breast cancer           4       369   -0.425421  -1.203010
## 5       B cells Breast cancer           5        87    0.437659   0.817755
## ...         ...           ...         ...       ...         ...        ...
## 44  Fibroblasts       T cells           4       369  0.06743410  0.2466458
## 45  Fibroblasts       T cells           5        87  0.60889057  2.2109896
## 46  Fibroblasts       T cells           6       173  0.18793427  0.5924011
## 47  Fibroblasts       T cells           7       315 -0.16126966 -0.6808440
## 48  Fibroblasts       T cells           8       373  0.00525121  0.0292494
##            df     p.Pos     p.Neg
##     <numeric> <numeric> <numeric>
## 1     4.42476  0.244103 0.7558973
## 2     1.21560  0.896728 0.1032722
## 3     9.78951  0.941830 0.0581705
## 4     6.54929  0.864675 0.1353253
## 5     2.82248  0.238409 0.7615909
## ...       ...       ...       ...
## 44   13.31707 0.4044718  0.595528
## 45    8.29697 0.0284087  0.971591
## 46    9.58525 0.2836464  0.716354
## 47   17.35982 0.7475245  0.252476
## 48   31.02448 0.4884265  0.511573

Or the correlation between breast cancer and fibroblasts across the whole slide:

results$overall
## DataFrame with 6 rows and 5 columns
##       celltype1     celltype2   cor.coef      p.Pos     p.Neg
##     <character>   <character>  <numeric>  <numeric> <numeric>
## 1       B cells Breast cancer -0.1451480 0.59904768  0.135167
## 2       B cells   Fibroblasts  0.0325429 0.58266567  0.571149
## 3       B cells       T cells  0.4838838 0.00725233  0.996427
## 4 Breast cancer   Fibroblasts  0.0510471 0.31920244  0.876317
## 5 Breast cancer       T cells -0.1847550 0.29867633  0.113621
## 6   Fibroblasts       T cells  0.1756227 0.16256935  0.873096

We can also visualise the fitting using function plotDensCor.

plotDensCor(spe, celltype1 = "Breast cancer", celltype2 = "Fibroblasts")

Or, we can visualise the statistics between each pair of cell types using function plotCorHeatmap in the ROIs:

plotCorHeatmap(results$ROI)

Or the correlation between cell type pairs across the whole slide:

plotCorHeatmap(results$overall)

Cell-based analysis

Based on the grid density, we can ask many biological question about the data. For example, we would like to know if a certain cell type that are located in high density of fibroblast cells are different to the same cell type from a different level of fibroblast region.

cell annotation based on grid density

To address this question, we first need to divide cells into different levels of grid density. This can be done using a contour identification strategy with function getContour.

spe <- getContour(spe, coi = "Fibroblasts", equal.cell = TRUE)

Different level of contour can be visualised with cells using plotContour.

plotContour(spe, coi = "Fibroblasts")

We can then annotate cells by their locations within each contour using function allocateCells.

spe <- allocateCells(spe)
plotSpatial(spe, group.by = "fibroblasts_contour", pt.alpha = 0.5)

We can visualise cell type composition per level.

plotCellCompo(spe, contour = "Fibroblasts")

plotCellCompo(spe, contour = "Fibroblasts", by.roi = TRUE)

## R version 4.4.2 (2024-10-31)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.5 LTS
## 
## Matrix products: default
## BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.20.so;  LAPACK version 3.10.0
## 
## locale:
##  [1] LC_CTYPE=C.UTF-8       LC_NUMERIC=C           LC_TIME=C.UTF-8       
##  [4] LC_COLLATE=C.UTF-8     LC_MONETARY=C.UTF-8    LC_MESSAGES=C.UTF-8   
##  [7] LC_PAPER=C.UTF-8       LC_NAME=C              LC_ADDRESS=C          
## [10] LC_TELEPHONE=C         LC_MEASUREMENT=C.UTF-8 LC_IDENTIFICATION=C   
## 
## time zone: UTC
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] sf_1.0-19                   SpatialExperiment_1.16.0   
##  [3] SingleCellExperiment_1.28.1 SummarizedExperiment_1.36.0
##  [5] Biobase_2.66.0              GenomicRanges_1.58.0       
##  [7] GenomeInfoDb_1.42.0         IRanges_2.40.0             
##  [9] S4Vectors_0.44.0            BiocGenerics_0.52.0        
## [11] MatrixGenerics_1.18.0       matrixStats_1.4.1          
## [13] scider_1.5.3                ggplot2_3.5.1              
## 
## loaded via a namespace (and not attached):
##   [1] DBI_1.2.3               deldir_2.0-4            rlang_1.1.4            
##   [4] magrittr_2.0.3          snakecase_0.11.1        e1071_1.7-16           
##   [7] compiler_4.4.2          spatstat.geom_3.3-4     mgcv_1.9-1             
##  [10] systemfonts_1.1.0       vctrs_0.6.5             stringr_1.5.1          
##  [13] pkgconfig_2.0.3         crayon_1.5.3            fastmap_1.2.0          
##  [16] magick_2.8.5            XVector_0.46.0          lwgeom_0.2-14          
##  [19] labeling_0.4.3          utf8_1.2.4              promises_1.3.0         
##  [22] rmarkdown_2.29          UCSC.utils_1.2.0        ragg_1.3.3             
##  [25] purrr_1.0.2             xfun_0.49               zlibbioc_1.52.0        
##  [28] cachem_1.1.0            jsonlite_1.8.9          goftest_1.2-3          
##  [31] later_1.3.2             DelayedArray_0.32.0     spatstat.utils_3.1-1   
##  [34] R6_2.5.1                RColorBrewer_1.1-3      bslib_0.8.0            
##  [37] stringi_1.8.4           spatstat.data_3.1-4     spatstat.univar_3.1-1  
##  [40] lubridate_1.9.3         jquerylib_0.1.4         Rcpp_1.0.13-1          
##  [43] knitr_1.49              tensor_1.5              splines_4.4.2          
##  [46] igraph_2.1.1            httpuv_1.6.15           Matrix_1.7-1           
##  [49] timechange_0.3.0        tidyselect_1.2.1        abind_1.4-8            
##  [52] yaml_2.3.10             spatstat.random_3.3-2   spatstat.explore_3.3-3 
##  [55] lattice_0.22-6          tibble_3.2.1            shiny_1.9.1            
##  [58] withr_3.0.2             evaluate_1.0.1          desc_1.4.3             
##  [61] units_0.8-5             proxy_0.4-27            polyclip_1.10-7        
##  [64] pillar_1.9.0            KernSmooth_2.23-24      plotly_4.10.4          
##  [67] generics_0.1.3          dbscan_1.2-0            munsell_0.5.1          
##  [70] scales_1.3.0            xtable_1.8-4            class_7.3-22           
##  [73] glue_1.8.0              janitor_2.2.0           pheatmap_1.0.12        
##  [76] lazyeval_0.2.2          tools_4.4.2             hexDensity_1.4.5       
##  [79] hexbin_1.28.5           data.table_1.16.2       fs_1.6.5               
##  [82] grid_4.4.2              tidyr_1.3.1             colorspace_2.1-1       
##  [85] nlme_3.1-166            GenomeInfoDbData_1.2.13 fastmatrix_0.5-7721    
##  [88] cli_3.6.3               spatstat.sparse_3.1-0   textshaping_0.4.0      
##  [91] fansi_1.0.6             S4Arrays_1.6.0          viridisLite_0.4.2      
##  [94] dplyr_1.1.4             gtable_0.3.6            SpatialPack_0.4-1      
##  [97] sass_0.4.9              digest_0.6.37           classInt_0.4-10        
## [100] SparseArray_1.6.0       rjson_0.2.23            htmlwidgets_1.6.4      
## [103] farver_2.1.2            htmltools_0.5.8.1       pkgdown_2.1.1          
## [106] lifecycle_1.0.4         httr_1.4.7              mime_0.12