spicyR 1.10.7
# load required packages
library(spicyR)
library(ggplot2)
library(SingleCellExperiment)
if (!require("BiocManager")) {
install.packages("BiocManager")
}
BiocManager::install("spicyR")
This guide will provide a step-by-step guide on how mixed effects models can be applied to multiple segmented and labelled images to identify how the localisation of different cell types can change across different conditions. Here, the subject is modelled as a random effect, and the different conditions are modelled as a fixed effect.
Here, we use a subset of the Damond et al 2019 imaging mass cytometry dataset. We will compare the spatial distributions of cells in the pancreatic islets of individuals with early onset diabetes and healthy controls.
diabetesData_SCE
is a SingleCellExperiment
object containing single-cell data of 160 images
from 8 subjects, with 20 images per subjects.
cellSummary()
returns a DataFrame
object providing the location (x
and y
)
and cell type (cellType
) of each cell and the image it belongs to (imageID
).
imagePheno()
returns a tibble
object providing the corresponding subject
(subject
) and condition (condition
) for each image.
data("diabetesData_SCE")
diabetesData_SCE
#> class: SingleCellExperiment
#> dim: 0 253777
#> metadata(0):
#> assays(0):
#> rownames: NULL
#> rowData names(0):
#> colnames: NULL
#> colData names(11): imageID cellID ... group stage
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
In this data set, cell types include immune cell types (B cells, naive T cells, T Helper cells, T cytotoxic cells, neutrophils, macrophages) and pancreatic islet cells (alpha, beta, gamma, delta).
To investigate changes in colocalisation between two different cell types, we
measure the level of colocalisation between two cell types by modelling with the
Lcross()
function in the spatstat.explore
package. Specifically, the mean difference
between the obtained function and the theoretical function is used as a measure
for the level of colocalisation. Differences of this statistic between two
conditions is modelled using a weighted mixed effects model, with condition as
the fixed effect and subject as the random effect.
spicyTestBootstrap
## Testing for change in colocalisation for a specific pair of cells
Firstly, we can see whether one cell type tends to be around another cell type
in one condition compared to the other. This can be done using the spicy()
function, where we include condition
, and subject
. In this example, we want
to see whether or not Delta cells (to
) tend to be found around Beta cells (from
)
in onset diabetes images compared to non-diabetic images.
spicyTestPair <- spicy(
diabetesData_SCE,
condition = "stage",
subject = "case",
from = "beta",
to = "delta"
)
#> Warning: The dim() method for DataFrameList objects is deprecated. Please use
#> dims() on these objects instead.
#> Warning: The nrow() method for DataFrameList objects is deprecated. Please use
#> nrows() on these objects instead.
#> Warning: The ncol() method for CompressedSplitDataFrameList objects is
#> deprecated. Please use ncols() on these objects instead.
#> Warning: The dim() method for DataFrameList objects is deprecated. Please use
#> dims() on these objects instead.
#> Warning: The nrow() method for DataFrameList objects is deprecated. Please use
#> nrows() on these objects instead.
#> Warning: The ncol() method for CompressedSplitDataFrameList objects is
#> deprecated. Please use ncols() on these objects instead.
topPairs(spicyTestPair)
#> intercept coefficient p.value adj.pvalue from to
#> beta__delta 179.729 -58.24478 0.000109702 0.000109702 beta delta
We obtain a spicy
object which details the results of the mixed effects
modelling performed. As the coefficient
in spicyTest
is positive, we find
that Th cells cells are more likely to be found around beta cells in the onset
diabetes group compared to the non-diabetic control.
Here, we can perform what we did above for all pairwise combinations of cell
types by excluding the from
and to
parameters from spicy()
.
spicyTest <- spicy(
diabetesData_SCE,
condition = "stage",
subject = "case"
)
#> Warning: The dim() method for DataFrameList objects is deprecated. Please use
#> dims() on these objects instead.
#> Warning: The nrow() method for DataFrameList objects is deprecated. Please use
#> nrows() on these objects instead.
#> Warning: The ncol() method for CompressedSplitDataFrameList objects is
#> deprecated. Please use ncols() on these objects instead.
#> Warning: The dim() method for DataFrameList objects is deprecated. Please use
#> dims() on these objects instead.
#> Warning: The nrow() method for DataFrameList objects is deprecated. Please use
#> nrows() on these objects instead.
#> Warning: The ncol() method for CompressedSplitDataFrameList objects is
#> deprecated. Please use ncols() on these objects instead.
spicyTest
#> conditionOnset conditionLong-duration
#> 0 15
topPairs(spicyTest)
#> intercept coefficient p.value adj.pvalue
#> beta__delta 1.815458e+02 -48.735693 0.0005033247 0.07169649
#> delta__beta 1.817943e+02 -48.166076 0.0005601288 0.07169649
#> B__unknown 7.212593e-16 11.770938 0.0052338392 0.42051606
#> delta__delta 2.089550e+02 -52.061196 0.0125129422 0.42051606
#> unknown__macrophage 1.007337e+01 -15.826919 0.0207410908 0.42051606
#> unknown__B 2.212257e-15 12.142848 0.0225855418 0.42051606
#> macrophage__unknown 1.004424e+01 -14.471666 0.0244668075 0.42051606
#> B__Th 6.285916e-15 26.847934 0.0245039854 0.42051606
#> otherimmune__naiveTc -9.292508e+00 33.584755 0.0255812944 0.42051606
#> ductal__ductal 1.481580e+01 -8.632569 0.0266935703 0.42051606
#> from to
#> beta__delta beta delta
#> delta__beta delta beta
#> B__unknown B unknown
#> delta__delta delta delta
#> unknown__macrophage unknown macrophage
#> unknown__B unknown B
#> macrophage__unknown macrophage unknown
#> B__Th B Th
#> otherimmune__naiveTc otherimmune naiveTc
#> ductal__ductal ductal ductal
Again, we obtain a spicy
object which outlines the result of the mixed effects
models performed for each pairwise combination if cell types.
We can represent this as a heatmap using the spatialMEMMultiPlot()
function by
providing it the spicy
object obtained.
signifPlot(
spicyTest,
breaks = c(-3, 3, 1),
marksToPlot = c(
"alpha", "beta", "gamma", "delta",
"B", "naiveTc", "Th", "Tc", "neutrophil", "macrophage"
)
)
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.1 spicyR_1.10.7
#> [13] BiocStyle_2.26.0
#>
#> loaded via a namespace (and not attached):
#> [1] nlme_3.1-162 spatstat.sparse_3.0-0
#> [3] bitops_1.0-7 RColorBrewer_1.1-3
#> [5] numDeriv_2016.8-1.1 tools_4.2.2
#> [7] bslib_0.4.2 utf8_1.2.3
#> [9] R6_2.5.1 HDF5Array_1.26.0
#> [11] mgcv_1.8-41 colorspace_2.1-0
#> [13] rhdf5filters_1.10.0 withr_2.5.0
#> [15] tidyselect_1.2.0 compiler_4.2.2
#> [17] cli_3.6.0 spatstat.explore_3.0-6
#> [19] DelayedArray_0.24.0 labeling_0.4.2
#> [21] bookdown_0.32 sass_0.4.5
#> [23] scales_1.2.1 spatstat.data_3.0-0
#> [25] goftest_1.2-3 digest_0.6.31
#> [27] SpatialExperiment_1.8.0 spatstat.utils_3.0-1
#> [29] minqa_1.2.5 rmarkdown_2.20
#> [31] R.utils_2.12.2 XVector_0.38.0
#> [33] pkgconfig_2.0.3 htmltools_0.5.4
#> [35] lme4_1.1-31 sparseMatrixStats_1.10.0
#> [37] highr_0.10 fastmap_1.1.0
#> [39] limma_3.54.1 rlang_1.0.6
#> [41] DelayedMatrixStats_1.20.0 jquerylib_0.1.4
#> [43] generics_0.1.3 farver_2.1.1
#> [45] jsonlite_1.8.4 spatstat.random_3.1-3
#> [47] BiocParallel_1.32.5 dplyr_1.1.0
#> [49] R.oo_1.25.0 RCurl_1.98-1.10
#> [51] magrittr_2.0.3 GenomeInfoDbData_1.2.9
#> [53] scuttle_1.8.4 scam_1.2-13
#> [55] Matrix_1.5-3 Rcpp_1.0.10
#> [57] munsell_0.5.0 Rhdf5lib_1.20.0
#> [59] fansi_1.0.4 abind_1.4-5
#> [61] lifecycle_1.0.3 R.methodsS3_1.8.2
#> [63] yaml_2.3.7 edgeR_3.40.2
#> [65] MASS_7.3-58.2 zlibbioc_1.44.0
#> [67] rhdf5_2.42.0 grid_4.2.2
#> [69] parallel_4.2.2 dqrng_0.3.0
#> [71] deldir_1.0-6 lattice_0.20-45
#> [73] beachmat_2.14.0 splines_4.2.2
#> [75] tensor_1.5 locfit_1.5-9.7
#> [77] magick_2.7.3 knitr_1.42
#> [79] pillar_1.8.1 spatstat.geom_3.0-6
#> [81] boot_1.3-28.1 rjson_0.2.21
#> [83] codetools_0.2-19 glue_1.6.2
#> [85] evaluate_0.20 data.table_1.14.6
#> [87] BiocManager_1.30.19 vctrs_0.5.2
#> [89] nloptr_2.0.3 tweenr_2.0.2
#> [91] purrr_1.0.1 tidyr_1.3.0
#> [93] gtable_0.3.1 polyclip_1.10-4
#> [95] cachem_1.0.6 xfun_0.37
#> [97] ggforce_0.4.1 DropletUtils_1.18.1
#> [99] tibble_3.1.8 pheatmap_1.0.12
#> [101] lmerTest_3.1-3 concaveman_1.1.0