miloR 1.4.0
library(miloR)
library(SingleCellExperiment)
library(scater)
library(scran)
library(dplyr)
library(patchwork)
library(MouseThymusAgeing)
library(scuttle)
We have seen how Milo uses graph neighbourhoods to model cell state abundance differences in an experiment, when comparing 2 groups. However, we are often interested in
testing between 2 specific groups in our analysis when our experiment has collected data from \(\gt\) 2 groups. We can focus our analysis to a 2 group comparison and
still make use of all of the data for things like dispersion estimation, by using contrasts. For an in-depth use of contrasts we recommend users refer to the limma
and edgeR
Biostars and Bioconductor community forum threads on the subject. Here I will give an overview of how to use contrasts in the context of a Milo analysis.
We will use the MouseThymusAgeing
data package as there are multiple groups that we can compare.
thy.sce <- MouseSMARTseqData() # this function downloads the full SCE object
## snapshotDate(): 2022-04-19
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
thy.sce <- logNormCounts(thy.sce)
thy.sce
## class: SingleCellExperiment
## dim: 48801 2327
## metadata(0):
## assays(2): counts logcounts
## rownames(48801): ERCC-00002 ERCC-00003 ... ENSMUSG00000064371
## ENSMUSG00000064372
## rowData names(6): Geneid Chr ... Strand Length
## colnames(2327): B13.B002229.1_52.1.32.1_S109
## B13.B002297.1_32.4.52.1_S73 ... P9.B002345.5_52.1.32.1_S93
## P9.B002450.5_4.52.16.1_S261
## colData names(11): CellID ClusterID ... SubType sizeFactor
## reducedDimNames(1): PCA
## mainExpName: NULL
## altExpNames(0):
thy.sce <- runUMAP(thy.sce) # add a UMAP for plotting results later
thy.milo <- Milo(thy.sce) # from the SCE object
reducedDim(thy.milo, "UMAP") <- reducedDim(thy.sce, "UMAP")
plotUMAP(thy.milo, colour_by="SubType") + plotUMAP(thy.milo, colour_by="Age")
These UMAPs shows how the different thymic epithelial cell subtypes and cells from different aged mice are distributed across our single-cell data set. Next we build the KNN graph and define neighbourhoods to quantify cell abundance across our experimental samples.
# we build KNN graph
thy.milo <- buildGraph(thy.milo, k = 10, d = 30)
## Constructing kNN graph with k:10
thy.milo <- makeNhoods(thy.milo, prop = 0.1, k = 10, d=30, refined = TRUE, refinement_scheme="graph") # make nhoods using graph-only as this is faster
## Checking valid object
## Running refined sampling with graph
colData(thy.milo)$Sample <- paste(colData(thy.milo)$SortDay, colData(thy.milo)$Age, sep="_")
thy.milo <- countCells(thy.milo, meta.data = data.frame(colData(thy.milo)), samples="Sample") # make the nhood X sample counts matrix
## Checking meta.data validity
## Counting cells in neighbourhoods
Now we have the pieces in place for DA testing to demonstrate how to use contrasts. We will use these contrasts to explicitly define which groups will be compared to each other.
thy.design <- data.frame(colData(thy.milo))[,c("Sample", "SortDay", "Age")]
thy.design <- distinct(thy.design)
rownames(thy.design) <- thy.design$Sample
## Reorder rownames to match columns of nhoodCounts(milo)
thy.design <- thy.design[colnames(nhoodCounts(thy.milo)), , drop=FALSE]
table(thy.design$Age)
##
## 16wk 1wk 32wk 4wk 52wk
## 5 5 5 5 5
To demonstrate the use of contrasts we will fit the whole model to the whole data set, but we will compare sequential pairs of time points. I’ll start with week 1 vs. week 4 to illustrate the syntax.
rownames(thy.design) <- thy.design$Sample
contrast.1 <- c("Age1wk - Age4wk") # the syntax is <VariableName><ConditionLevel> - <VariableName><ControlLevel>
# we need to use the ~ 0 + Variable expression here so that we have all of the levels of our variable as separate columns in our model matrix
da_results <- testNhoods(thy.milo, design = ~ 0 + Age, design.df = thy.design, model.contrasts = contrast.1, fdr.weighting="graph-overlap")
## Using TMM normalisation
## Warning in max(abs(logR)): no non-missing arguments to max; returning -Inf
## Warning in max(abs(logR)): no non-missing arguments to max; returning -Inf
## Performing spatial FDR correction withgraph-overlap weighting
table(da_results$SpatialFDR < 0.1)
##
## FALSE TRUE
## 175 25
This calculates a Fold-change and corrected P-value for each neighbourhood, which indicates whether there is significant differential abundance between conditions for 25 neighbourhoods.
You will notice that the syntax for the contrasts is quite specific. It starts with the name of the column variable that contains the different group levels; in this case
it is the Age
variable. We then define the comparison levels as level1 - level2
. To understand this syntax we need to consider what we are concretely comparing. In this
case we are asking what is the ratio of the average cell count at week1 compared to the average cell count at week 4, where the averaging is across the replicates. The
reason we express this as a difference rather than a ratio is because we are dealing with the log fold change.
We can also pass multiple comparisons at the same time, for instance if we wished to compare each sequential pair of time points. This will give us a better intuition behind how to use contrasts to compare multiple groups.
contrast.all <- c("Age1wk - Age4wk", "Age4wk - Age16wk", "Age16wk - Age32wk", "Age32wk - Age52wk")
# this is the edgeR code called by `testNhoods`
model <- model.matrix(~ 0 + Age, data=thy.design)
mod.constrast <- makeContrasts(contrasts=contrast.all, levels=model)
mod.constrast
## Contrasts
## Levels Age1wk - Age4wk Age4wk - Age16wk Age16wk - Age32wk Age32wk - Age52wk
## Age16wk 0 -1 1 0
## Age1wk 1 0 0 0
## Age32wk 0 0 -1 1
## Age4wk -1 1 0 0
## Age52wk 0 0 0 -1
This shows the contrast matrix. If we want to test each of these comparisons then we need to pass them sequentially to testNhoods
, then apply an additional
multiple testing correction to the spatial FDR values.
Contrasts are not limited to these simple pair-wise comparisons, we can also group levels together for comparisons. For instance, imagine we want to know what the effect of the cell counts in the week 1 mice is compared to all other time points.
model <- model.matrix(~ 0 + Age, data=thy.design)
ave.contrast <- c("Age1wk - (Age4wk + Age16wk + Age32wk + Age52wk)/4")
mod.constrast <- makeContrasts(contrasts=ave.contrast, levels=model)
mod.constrast
## Contrasts
## Levels Age1wk - (Age4wk + Age16wk + Age32wk + Age52wk)/4
## Age16wk -0.25
## Age1wk 1.00
## Age32wk -0.25
## Age4wk -0.25
## Age52wk -0.25
In this contrasts matrix we can see that we have taken the average effect over the other time points. Now running this using testNhoods
da_results <- testNhoods(thy.milo, design = ~ 0 + Age, design.df = thy.design, model.contrasts = ave.contrast, fdr.weighting="graph-overlap")
## Using TMM normalisation
## Warning in max(abs(logR)): no non-missing arguments to max; returning -Inf
## Warning in max(abs(logR)): no non-missing arguments to max; returning -Inf
## Performing spatial FDR correction withgraph-overlap weighting
table(da_results$SpatialFDR < 0.1)
##
## FALSE TRUE
## 120 80
In this comparison there are 80 DA nhoods - which we can visualise on a superimposed single-cell UMAP.
thy.milo <- buildNhoodGraph(thy.milo)
plotUMAP(thy.milo, colour_by="SubType") + plotNhoodGraphDA(thy.milo, da_results, alpha=0.1) +
plot_layout(guides="auto" )
In these side-by-side UMAPs we can see that there is an enrichment of the Perinatal cTEC and Proliferating TEC populations in the 1 week old compared to the other time points.
For a more extensive description of the uses of contrasts please take a look at the edgeR documentation .
Session Info
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] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] MouseThymusAgeing_1.3.0 patchwork_1.1.1
## [3] dplyr_1.0.8 scran_1.24.0
## [5] scater_1.24.0 ggplot2_3.3.5
## [7] scuttle_1.6.0 SingleCellExperiment_1.18.0
## [9] SummarizedExperiment_1.26.0 Biobase_2.56.0
## [11] GenomicRanges_1.48.0 GenomeInfoDb_1.32.0
## [13] IRanges_2.30.0 S4Vectors_0.34.0
## [15] BiocGenerics_0.42.0 MatrixGenerics_1.8.0
## [17] matrixStats_0.62.0 miloR_1.4.0
## [19] edgeR_3.38.0 limma_3.52.0
## [21] BiocStyle_2.24.0
##
## loaded via a namespace (and not attached):
## [1] AnnotationHub_3.4.0 BiocFileCache_2.4.0
## [3] igraph_1.3.1 splines_4.2.0
## [5] BiocParallel_1.30.0 digest_0.6.29
## [7] htmltools_0.5.2 magick_2.7.3
## [9] viridis_0.6.2 fansi_1.0.3
## [11] magrittr_2.0.3 memoise_2.0.1
## [13] ScaledMatrix_1.4.0 cluster_2.1.3
## [15] Biostrings_2.64.0 graphlayouts_0.8.0
## [17] colorspace_2.0-3 blob_1.2.3
## [19] rappdirs_0.3.3 ggrepel_0.9.1
## [21] xfun_0.30 crayon_1.5.1
## [23] RCurl_1.98-1.6 jsonlite_1.8.0
## [25] glue_1.6.2 polyclip_1.10-0
## [27] gtable_0.3.0 zlibbioc_1.42.0
## [29] XVector_0.36.0 DelayedArray_0.22.0
## [31] BiocSingular_1.12.0 scales_1.2.0
## [33] DBI_1.1.2 Rcpp_1.0.8.3
## [35] viridisLite_0.4.0 xtable_1.8-4
## [37] dqrng_0.3.0 bit_4.0.4
## [39] rsvd_1.0.5 metapod_1.4.0
## [41] httr_1.4.2 FNN_1.1.3
## [43] RColorBrewer_1.1-3 ellipsis_0.3.2
## [45] pkgconfig_2.0.3 farver_2.1.0
## [47] uwot_0.1.11 sass_0.4.1
## [49] dbplyr_2.1.1 locfit_1.5-9.5
## [51] utf8_1.2.2 labeling_0.4.2
## [53] tidyselect_1.1.2 rlang_1.0.2
## [55] later_1.3.0 AnnotationDbi_1.58.0
## [57] munsell_0.5.0 BiocVersion_3.15.2
## [59] tools_4.2.0 cachem_1.0.6
## [61] cli_3.3.0 generics_0.1.2
## [63] RSQLite_2.2.12 ExperimentHub_2.4.0
## [65] evaluate_0.15 stringr_1.4.0
## [67] fastmap_1.1.0 yaml_2.3.5
## [69] knitr_1.38 bit64_4.0.5
## [71] tidygraph_1.2.1 purrr_0.3.4
## [73] KEGGREST_1.36.0 ggraph_2.0.5
## [75] sparseMatrixStats_1.8.0 mime_0.12
## [77] compiler_4.2.0 png_0.1-7
## [79] beeswarm_0.4.0 filelock_1.0.2
## [81] curl_4.3.2 interactiveDisplayBase_1.34.0
## [83] tibble_3.1.6 statmod_1.4.36
## [85] tweenr_1.0.2 bslib_0.3.1
## [87] stringi_1.7.6 highr_0.9
## [89] RSpectra_0.16-1 lattice_0.20-45
## [91] bluster_1.6.0 Matrix_1.4-1
## [93] vctrs_0.4.1 pillar_1.7.0
## [95] lifecycle_1.0.1 BiocManager_1.30.17
## [97] jquerylib_0.1.4 BiocNeighbors_1.14.0
## [99] cowplot_1.1.1 bitops_1.0-7
## [101] irlba_2.3.5 httpuv_1.6.5
## [103] R6_2.5.1 bookdown_0.26
## [105] promises_1.2.0.1 gridExtra_2.3
## [107] vipor_0.4.5 MASS_7.3-57
## [109] gtools_3.9.2 assertthat_0.2.1
## [111] withr_2.5.0 GenomeInfoDbData_1.2.8
## [113] parallel_4.2.0 grid_4.2.0
## [115] beachmat_2.12.0 tidyr_1.2.0
## [117] rmarkdown_2.14 DelayedMatrixStats_1.18.0
## [119] ggforce_0.3.3 shiny_1.7.1
## [121] ggbeeswarm_0.6.0