Contents

1 Installation

if (!require("BiocManager"))
    install.packages("BiocManager")
BiocManager::install("spicyR")
# load required packages
library(spicyR)
library(lisaClust)
library(ggplot2)

2 Overview

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.

3 Quick start

3.1 Generate toy data

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)

3.2 Create SegmentedCellExperiment object

First we store our data in a SegmentedCells object.


cellExp <- SegmentedCells(cells, cellTypeString = 'cellType')

3.3 Generate LISA curves

We can then 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))

3.4 Perform some clustering

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 region() <- function.


kM <- kmeans(lisaCurves,2)
region(cellExp) <- paste('region',kM$cluster,sep = '_')

3.5 Plot identified regions

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, imageID = c('s1','s2'))

3.6 Alternative hatching plot

We could also create this plot using geom_hatching and scale_region_manual.


df <- region(cellExp, annot = TRUE)

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 = 0.1) +
  theme_minimal() + 
  scale_region_manual(values = 6:7, labels = c('ab','cd'))

p

4 Damond et al. islet data.

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).

4.1 Read in data

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)

4.2 Cluster cell-types

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 = '')

4.3 Generate LISA curves

As before, we can calculate local indicators of spatial association (LISA) functions using the lisa function.


lisaCurves <- lisa(cellExp, Rs = c(10,20,50))

4.4 Perform some clustering

The LISA curves can then be used to cluster the cells. Here we use k-means clustering to cluster the cells into two microenvironments. We can store these cell clusters or cell “regions” in our SegmentedCells object using the region() <- function.


kM <- kmeans(lisaCurves,2)
region(cellExp) <- paste('region',kM$cluster,sep = '_')

4.5 Plot identified regions

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)

5 sessionInfo()

sessionInfo()
## R version 4.2.0 RC (2022-04-19 r82224)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.15-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] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] ggplot2_3.3.5    lisaClust_1.4.0  spicyR_1.8.0     BiocStyle_2.24.0
## 
## loaded via a namespace (and not attached):
##  [1] nlme_3.1-157          spatstat.sparse_2.1-1 RColorBrewer_1.1-3   
##  [4] numDeriv_2016.8-1.1   tools_4.2.0           bslib_0.3.1          
##  [7] utf8_1.2.2            R6_2.5.1              rpart_4.1.16         
## [10] DBI_1.1.2             BiocGenerics_0.42.0   mgcv_1.8-40          
## [13] colorspace_2.0-3      withr_2.5.0           tidyselect_1.1.2     
## [16] curl_4.3.2            compiler_4.2.0        cli_3.3.0            
## [19] labeling_0.4.2        bookdown_0.26         sass_0.4.1           
## [22] scales_1.2.0          spatstat.data_2.2-0   goftest_1.2-3        
## [25] stringr_1.4.0         digest_0.6.29         spatstat.utils_2.3-0 
## [28] fftwtools_0.9-11      minqa_1.2.4           rmarkdown_2.14       
## [31] pkgconfig_2.0.3       htmltools_0.5.2       lme4_1.1-29          
## [34] fastmap_1.1.0         highr_0.9             rlang_1.0.2          
## [37] jquerylib_0.1.4       generics_0.1.2        farver_2.1.0         
## [40] jsonlite_1.8.0        spatstat.random_2.2-0 BiocParallel_1.30.0  
## [43] dplyr_1.0.8           magrittr_2.0.3        scam_1.2-12          
## [46] Matrix_1.4-1          Rcpp_1.0.8.3          munsell_0.5.0        
## [49] S4Vectors_0.34.0      fansi_1.0.3           abind_1.4-5          
## [52] lifecycle_1.0.1       stringi_1.7.6         yaml_2.3.5           
## [55] MASS_7.3-57           grid_4.2.0            parallel_4.2.0       
## [58] crayon_1.5.1          deldir_1.0-6          lattice_0.20-45      
## [61] splines_4.2.0         tensor_1.5            magick_2.7.3         
## [64] knitr_1.38            pillar_1.7.0          boot_1.3-28          
## [67] spatstat.geom_2.4-0   stats4_4.2.0          glue_1.6.2           
## [70] evaluate_0.15         V8_4.1.0              data.table_1.14.2    
## [73] BiocManager_1.30.17   vctrs_0.4.1           nloptr_2.0.0         
## [76] gtable_0.3.0          spatstat.core_2.4-2   purrr_0.3.4          
## [79] polyclip_1.10-0       tidyr_1.2.0           assertthat_0.2.1     
## [82] xfun_0.30             class_7.3-20          tibble_3.1.6         
## [85] pheatmap_1.0.12       lmerTest_3.1-3        IRanges_2.30.0       
## [88] concaveman_1.1.0      ellipsis_0.3.2