spicyR 1.10.6
if (!require("BiocManager"))
install.packages("BiocManager")
BiocManager::install("lisaClust")
# load required packages
library(lisaClust)
library(spicyR)
library(ggplot2)
library(SingleCellExperiment)
Clustering local indicators of spatial association (LISA) functions is a
methodology for identifying consistent spatial organisation of multiple
cell-types in an unsupervised way. This can be used to enable the
characterization of interactions between multiple cell-types simultaneously and
can complement traditional pairwise analysis. In our implementation our LISA
curves are a localised summary of an L-function from a Poisson point process
model. Our framework lisaClust
can be used to provide a high-level summary
of cell-type colocalization in high-parameter spatial cytometry data,
facilitating the identification of distinct tissue compartments or
identification of complex cellular microenvironments.
TO illustrate our lisaClust
framework, here we consider a very simple toy
example where two cell-types are completely separated spatially. We simulate
data for two different images.
set.seed(51773)
x <- round(c(runif(200),runif(200)+1,runif(200)+2,runif(200)+3,
runif(200)+3,runif(200)+2,runif(200)+1,runif(200)),4)*100
y <- round(c(runif(200),runif(200)+1,runif(200)+2,runif(200)+3,
runif(200),runif(200)+1,runif(200)+2,runif(200)+3),4)*100
cellType <- factor(paste('c',rep(rep(c(1:2),rep(200,2)),4),sep = ''))
imageID <- rep(c('s1', 's2'),c(800,800))
cells <- data.frame(x, y, cellType, imageID)
ggplot(cells, aes(x,y, colour = cellType)) + geom_point() + facet_wrap(~imageID)
First we store our data in a SegmentedCells
object.
cellExp <- SegmentedCells(cells, cellTypeString = 'cellType')
We can then use a convience function lisaClust
to simultaneously calculate local indicators of spatial association (LISA) functions
using the lisa
function and perform k-means clustering.
cellExp <- lisaClust(cellExp, k = 2)
The hatchingPlot
function can be used to construct a ggplot
object where the
regions are marked by different hatching patterns. This allows us to plot both
regions and cell-types on the same visualization.
hatchingPlot(cellExp, useImages = c('s1','s2'))
While the lisaClust
function is convenient, we have not implemented an exhaustive
suite of clustering methods as it is very easy to do this yourself. There are
just two simple steps.
We can calculate local indicators of spatial association (LISA) functions
using the lisa
function. Here the LISA curves are a
localised summary of an L-function from a Poisson point process model. The radii
that will be calculated over can be set with Rs
.
lisaCurves <- lisa(cellExp, Rs = c(20, 50, 100))
The LISA curves can then be used to cluster the cells. Here we use k-means
clustering, other clustering methods like SOM could be used. We can store these
cell clusters or cell “regions” in our SegmentedCells
object using the
cellAnnotation() <-
function.
kM <- kmeans(lisaCurves,2)
cellAnnotation(cellExp, "region") <- paste('region',kM$cluster,sep = '_')
We could also create this plot using geom_hatching
and scale_region_manual
.
df <- as.data.frame(cellSummary(cellExp))
p <- ggplot(df,aes(x = x,y = y, colour = cellType, region = region)) +
geom_point() +
facet_wrap(~imageID) +
geom_hatching(window = "concave",
line.spacing = 11,
nbp = 50,
line.width = 2,
hatching.colour = "gray20",
window.length = NULL) +
theme_minimal() +
scale_region_manual(values = 6:7, labels = c('ab','cd'))
p
## Faster ploting
The hatchingPlot
can be quite slow for large images and high nbp
or linewidth
.
It is often useful to simply plot the regions without the cell type information.
df <- as.data.frame(cellSummary(cellExp))
df <- df[df$imageID == "s1", ]
p <- ggplot(df,aes(x = x,y = y, colour = region)) +
geom_point() +
theme_classic()
p
The lisaClust
function also works with a SingleCellExperiment
. First lets
create a SingleCellExperiment
object.
sce <- SingleCellExperiment(colData = cellSummary(cellExp))
lisaClust
just needs columns in colData
corresponding to the x and y coordinates of the
cells, a column annotating the cell types of the cells and a column indicating
which image each cell came from.
sce <- lisaClust(sce,
k = 2,
spatialCoords = c("x", "y"),
cellType = "cellType",
imageID = "imageID")
We can then plot the regions using the following.
hatchingPlot(sce)
Here we apply our lisaClust
framework to three images of pancreatic islets
from A Map of Human Type 1 Diabetes Progression by Imaging Mass Cytometry by
Damond et al. (2019).
We will start by reading in the data and storing it as a SegmentedCells
object. Here the data is in a format consistent with that outputted by
CellProfiler.
isletFile <- system.file("extdata","isletCells.txt.gz", package = "spicyR")
cells <- read.table(isletFile, header = TRUE)
cellExp <- SegmentedCells(cells, cellProfiler = TRUE)
This data does not include annotation of the cell-types of each cell. Here we
extract the marker intensities from the SegmentedCells
object using
cellMarks
. We then perform k-means clustering with eight clusters and store
these cell-type clusters in our SegmentedCells
object using cellType() <-
.
markers <- cellMarks(cellExp)
kM <- kmeans(markers,10)
cellType(cellExp) <- paste('cluster', kM$cluster, sep = '')
As before, we can calculate perform k-means clustering on the local indicators
of spatial association (LISA) functions using the lisaClust
function.
cellExp <- lisaClust(cellExp, k = 2, Rs = c(10,20,50))
These regions are stored in cellExp and can be extracted.
cellAnnotation(cellExp, "region") |>
head()
## [1] "region_1" "region_1" "region_1" "region_1" "region_1" "region_1"
We should check to see which cell types appear more frequently in each region than expected by chance.
regionMap(cellExp, type = "bubble")
Finally, we can use hatchingPlot
to construct a ggplot
object where the
regions are marked by different hatching patterns. This allows us to visualize
the two regions and ten cell-types simultaneously.
hatchingPlot(cellExp)
sessionInfo()
## R version 4.2.2 (2022-10-31)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
##
## 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
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] SingleCellExperiment_1.20.0 SummarizedExperiment_1.28.0
## [3] Biobase_2.58.0 GenomicRanges_1.50.2
## [5] GenomeInfoDb_1.34.9 IRanges_2.32.0
## [7] S4Vectors_0.36.1 BiocGenerics_0.44.0
## [9] MatrixGenerics_1.10.0 matrixStats_0.63.0
## [11] ggplot2_3.4.0 spicyR_1.10.6
## [13] lisaClust_1.6.3 BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] SpatialExperiment_1.8.0 minqa_1.2.5
## [3] colorspace_2.1-0 rjson_0.2.21
## [5] deldir_1.0-6 class_7.3-21
## [7] scuttle_1.8.4 XVector_0.38.0
## [9] fftwtools_0.9-11 spatstat.data_3.0-0
## [11] farver_2.1.1 fansi_1.0.4
## [13] codetools_0.2-19 splines_4.2.2
## [15] R.methodsS3_1.8.2 sparseMatrixStats_1.10.0
## [17] cachem_1.0.6 knitr_1.42
## [19] polyclip_1.10-4 jsonlite_1.8.4
## [21] nloptr_2.0.3 R.oo_1.25.0
## [23] pheatmap_1.0.12 ggforce_0.4.1
## [25] HDF5Array_1.26.0 spatstat.sparse_3.0-0
## [27] BiocManager_1.30.19 compiler_4.2.2
## [29] dqrng_0.3.0 Matrix_1.5-3
## [31] fastmap_1.1.0 limma_3.54.1
## [33] cli_3.6.0 tweenr_2.0.2
## [35] htmltools_0.5.4 tools_4.2.2
## [37] lmerTest_3.1-3 gtable_0.3.1
## [39] glue_1.6.2 GenomeInfoDbData_1.2.9
## [41] dplyr_1.1.0 V8_4.2.2
## [43] Rcpp_1.0.10 jquerylib_0.1.4
## [45] vctrs_0.5.2 rhdf5filters_1.10.0
## [47] spatstat.explore_3.0-6 nlme_3.1-162
## [49] DelayedMatrixStats_1.20.0 spatstat.random_3.1-3
## [51] xfun_0.37 beachmat_2.14.0
## [53] lme4_1.1-31 lifecycle_1.0.3
## [55] goftest_1.2-3 scam_1.2-13
## [57] edgeR_3.40.2 zlibbioc_1.44.0
## [59] MASS_7.3-58.2 scales_1.2.1
## [61] spatstat.utils_3.0-1 parallel_4.2.2
## [63] rhdf5_2.42.0 RColorBrewer_1.1-3
## [65] curl_5.0.0 yaml_2.3.7
## [67] sass_0.4.5 highr_0.10
## [69] boot_1.3-28.1 BiocParallel_1.32.5
## [71] rlang_1.0.6 pkgconfig_2.0.3
## [73] bitops_1.0-7 evaluate_0.20
## [75] lattice_0.20-45 purrr_1.0.1
## [77] tensor_1.5 Rhdf5lib_1.20.0
## [79] labeling_0.4.2 tidyselect_1.2.0
## [81] magrittr_2.0.3 bookdown_0.32
## [83] R6_2.5.1 magick_2.7.3
## [85] generics_0.1.3 DelayedArray_0.24.0
## [87] pillar_1.8.1 withr_2.5.0
## [89] mgcv_1.8-41 abind_1.4-5
## [91] RCurl_1.98-1.10 tibble_3.1.8
## [93] DropletUtils_1.18.1 utf8_1.2.3
## [95] spatstat.geom_3.0-6 rmarkdown_2.20
## [97] locfit_1.5-9.7 grid_4.2.2
## [99] data.table_1.14.6 digest_0.6.31
## [101] tidyr_1.3.0 numDeriv_2016.8-1.1
## [103] R.utils_2.12.2 munsell_0.5.0
## [105] concaveman_1.1.0 bslib_0.4.2