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.
if (!require("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("scider")
The development version of scider
can be installed from GitHub:
In this vignette, we will use a subset of a Xenium Breast Cancer dataset.
In the data, we have quantification of 541 genes from 10000 cells.
## 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.
##
## 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.
scider
can conduct grid-based density analysis for spatial transcriptomics data.
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.
## [1] "grid_density" "grid_info"
## 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.00161274
## 2 380.937 155.873 2 1 2-1 0.00211721
## 3 480.937 155.873 3 1 3-1 0.00277295
## 4 580.937 155.873 4 1 4-1 0.00364209
## 5 680.937 155.873 5 1 5-1 0.00482578
## ... ... ... ... ... ... ...
## 12315 9480.94 11067.8 93 127 93-127 0.001135422
## 12316 9580.94 11067.8 94 127 94-127 0.000741773
## 12317 9680.94 11067.8 95 127 95-127 0.000466692
## 12318 9780.94 11067.8 96 127 96-127 0.000284044
## 12319 9880.94 11067.8 97 127 97-127 0.000167953
## density_breast_cancer density_fibroblasts density_t_cells density_overall
## <numeric> <numeric> <numeric> <numeric>
## 1 5.08804e-05 0.00290209 0.000991062 0.00555677
## 2 9.86082e-05 0.00414279 0.001464630 0.00782324
## 3 1.89630e-04 0.00583797 0.002155848 0.01095640
## 4 3.59094e-04 0.00810256 0.003155203 0.01525895
## 5 6.65018e-04 0.01104941 0.004584676 0.02112489
## ... ... ... ... ...
## 12315 0.01190884 0.01354203 2.00058e-05 0.02660629
## 12316 0.00737930 0.00901015 9.48948e-06 0.01714072
## 12317 0.00452586 0.00581520 4.56708e-06 0.01081232
## 12318 0.00275608 0.00365773 2.28099e-06 0.00670013
## 12319 0.00167020 0.00225250 1.20181e-06 0.00409186
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.
After obtaining grid-based density for each COI, we can then detect regions-of-interest (ROIs) based on density or select by user.
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.
## 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 66-124 66 124 6830.94 10808
## 1830 8 67-124 67 124 6930.94 10808
## 1831 8 68-124 68 124 7030.94 10808
## 1832 8 69-124 69 124 7130.94 10808
## 1833 8 70-124 70 124 7230.94 10808
We can visualise the ROIs with function plotROI
.
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.
After closing the interactive window, the selected ROI has been saved as a data.frame object named sel_region
in the R environment.
We can then use the postSelRegion
to save the ROI in the metadata of the SpatialExperiment object.
Similarly, we can plot visualise the user-defined ROI with function plotROI
.
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
.
We can see the correlation between each pair of cell types in each fibroblasts 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.913986 -2.483606
## 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.16143840 -0.6815605
## 48 Fibroblasts T cells 8 373 0.00536088 0.0298639
## df p.Pos p.Neg
## <numeric> <numeric> <numeric>
## 1 4.42476 0.244103 0.7558973
## 2 1.21561 0.896728 0.1032725
## 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.35906 0.7477454 0.252255
## 48 31.03191 0.4881834 0.511817
Or the correlation between each pair of cell types across all fibroblasts ROI:
## 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.1451816 0.59921852 0.135175
## 2 B cells Fibroblasts 0.0324573 0.58293635 0.570772
## 3 B cells T cells 0.4838632 0.00725754 0.996425
## 4 Breast cancer Fibroblasts 0.0510354 0.31924365 0.876378
## 5 Breast cancer T cells -0.1847721 0.29873370 0.113617
## 6 Fibroblasts T cells 0.1756014 0.16264804 0.873024
We can also visualise the fitting using function plotDensCor
.
Or, we can visualise the statistics between each pair of cell types using function plotCorHeatmap
in the ROIs:
Or the correlation between cell type pairs across the whole slide:
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.
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
.
Different level of contour can be visualised with cells using plotContour
.
We can then annotate cells by their locations within each contour using function allocateCells
.
We can visualise cell type composition per level.
## R version 4.5.1 Patched (2025-08-23 r88802)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.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] sf_1.0-21 SpatialExperiment_1.19.1
## [3] SingleCellExperiment_1.31.1 SummarizedExperiment_1.39.2
## [5] Biobase_2.69.1 GenomicRanges_1.61.5
## [7] Seqinfo_0.99.2 IRanges_2.43.5
## [9] S4Vectors_0.47.4 BiocGenerics_0.55.1
## [11] generics_0.1.4 MatrixGenerics_1.21.0
## [13] matrixStats_1.5.0 scider_1.7.1
## [15] ggplot2_4.0.0
##
## loaded via a namespace (and not attached):
## [1] DBI_1.2.3 deldir_2.0-4 rlang_1.1.6
## [4] magrittr_2.0.4 snakecase_0.11.1 e1071_1.7-16
## [7] compiler_4.5.1 spatstat.geom_3.6-0 mgcv_1.9-3
## [10] fftwtools_0.9-11 vctrs_0.6.5 stringr_1.5.2
## [13] pkgconfig_2.0.3 crayon_1.5.3 fastmap_1.2.0
## [16] magick_2.9.0 XVector_0.49.1 lwgeom_0.2-14
## [19] labeling_0.4.3 promises_1.3.3 rmarkdown_2.30
## [22] purrr_1.1.0 xfun_0.53 cachem_1.1.0
## [25] jsonlite_2.0.0 goftest_1.2-3 later_1.4.4
## [28] DelayedArray_0.35.3 spatstat.utils_3.2-0 R6_2.6.1
## [31] bslib_0.9.0 stringi_1.8.7 RColorBrewer_1.1-3
## [34] spatstat.data_3.1-8 spatstat.univar_3.1-4 lubridate_1.9.4
## [37] jquerylib_0.1.4 Rcpp_1.1.0 knitr_1.50
## [40] tensor_1.5.1 splines_4.5.1 httpuv_1.6.16
## [43] Matrix_1.7-4 igraph_2.1.4 timechange_0.3.0
## [46] tidyselect_1.2.1 dichromat_2.0-0.1 abind_1.4-8
## [49] yaml_2.3.10 spatstat.random_3.4-2 spatstat.explore_3.5-3
## [52] lattice_0.22-7 tibble_3.3.0 shiny_1.11.1
## [55] withr_3.0.2 S7_0.2.0 evaluate_1.0.5
## [58] units_0.8-7 proxy_0.4-27 polyclip_1.10-7
## [61] pillar_1.11.1 KernSmooth_2.23-26 plotly_4.11.0
## [64] dbscan_1.2.3 scales_1.4.0 xtable_1.8-4
## [67] class_7.3-23 glue_1.8.0 janitor_2.2.1
## [70] pheatmap_1.0.13 lazyeval_0.2.2 tools_4.5.1
## [73] hexDensity_1.4.10 hexbin_1.28.5 data.table_1.17.8
## [76] grid_4.5.1 tidyr_1.3.1 nlme_3.1-168
## [79] fastmatrix_0.6 cli_3.6.5 spatstat.sparse_3.1-0
## [82] S4Arrays_1.9.1 viridisLite_0.4.2 dplyr_1.1.4
## [85] gtable_0.3.6 SpatialPack_0.4-1 sass_0.4.10
## [88] digest_0.6.37 classInt_0.4-11 SparseArray_1.9.1
## [91] rjson_0.2.23 htmlwidgets_1.6.4 farver_2.1.2
## [94] htmltools_0.5.8.1 lifecycle_1.0.4 httr_1.4.7
## [97] mime_0.13