library(miloR)
library(SingleCellExperiment)
library(scater)
library(scran)
library(dplyr)
library(patchwork)
library(MouseThymusAgeing)
library(scuttle)

1 Introduction

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.

2 Load data

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
## see ?MouseThymusAgeing and browseVignettes('MouseThymusAgeing') for documentation
## loading from cache
## field not found in version - adding
## 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):

3 Define cell neighbourhoods

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 = 11, d = 40)
## Constructing kNN graph with k:11
thy.milo <- makeNhoods(thy.milo, prop = 0.2, k = 11, d=40, 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
plotNhoodSizeHist(thy.milo)

4 Differential abundance testing with contrasts

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", norm.method="TMM")
## Using TMM normalisation
## Running with model contrasts
## Performing spatial FDR correction with graph-overlap weighting
table(da_results$SpatialFDR < 0.1)
## 
## FALSE  TRUE 
##   331    29

This calculates a Fold-change and corrected P-value for each neighbourhood, which indicates whether there is significant differential abundance between conditions for 29 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")
# contrast.all <- c("Age1wk - Age4wk", "Age16wk - Age32wk")

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

contrast1.res <- testNhoods(thy.milo, design=~0+ Age, design.df=thy.design, fdr.weighting="graph-overlap", model.contrasts = contrast.all)
## Using TMM normalisation
## Running with model contrasts
## Performing spatial FDR correction with graph-overlap weighting
head(contrast1.res)
##   logFC.Age1wk...Age4wk logFC.Age4wk...Age16wk logFC.Age16wk...Age32wk
## 1            -0.9567742             1.50884755             -0.07591883
## 2             1.6776047             2.15244207              0.00000000
## 3             0.5244291             0.63724176              1.08850822
## 4            -4.8234046             2.67860956              0.68113243
## 5            -0.9156193            -0.89601759             -0.57589780
## 6            -2.6505102            -0.02183302             -0.40124094
##   logFC.Age32wk...Age52wk   logCPM         F       PValue          FDR Nhood
## 1              -0.5079649 13.61263 0.5694729 6.848051e-01 0.7539138508     1
## 2               0.0000000 13.42815 4.1180552 2.464861e-03 0.0089631292     2
## 3              -1.0285693 13.66500 0.9535196 4.318199e-01 0.5342102467     3
## 4               1.4636626 13.64878 6.4145787 3.781243e-05 0.0007562485     4
## 5              -0.1800067 13.95991 1.4954560 2.006340e-01 0.2936107645     5
## 6              -0.8446143 13.63399 2.2335148 6.287931e-02 0.1143260108     6
##     SpatialFDR
## 1 0.7362177461
## 2 0.0094986344
## 3 0.5169025712
## 4 0.0007888425
## 5 0.2696852700
## 6 0.1295117334

This matrix of contrasts will perform a quasi-likelihood F-test over all 5 contrasts, hence a single p-value and spatial FDR are returned. Log fold changes are returned for each contrast of the Age variable, which gives 1 log-fold change column for each - this is the default behaviour of glmQLFTest in the edgeR package which is what Milo uses for hypothesis testing. In general, and to avoid confusion, we recommend testing each pair of contrasts separately if these are the comparisons of interest, as shown below.

# compare weeks 4 and 16, with week 4 as the reference.
cont.4vs16.res <- testNhoods(thy.milo, design=~0+ Age, design.df=thy.design, fdr.weighting="graph-overlap", model.contrasts = c("Age4wk - Age16wk"))
## Using TMM normalisation
## Running with model contrasts
## Performing spatial FDR correction with graph-overlap weighting
head(cont.4vs16.res)
##         logFC   logCPM            F     PValue       FDR Nhood SpatialFDR
## 1  1.50884755 13.61263 1.5524215933 0.21281884 0.5393339     1  0.4542207
## 2  2.15244207 13.42815 2.7304188735 0.09849706 0.4098461     2  0.2917372
## 3  0.63724176 13.66500 0.3068607839 0.57963016 0.8730831     3  0.8444921
## 4  2.67860956 13.64878 5.3741636320 0.02046501 0.2389419     4  0.1492485
## 5 -0.89601759 13.95991 0.5700749070 0.45025385 0.7682056     5  0.7305036
## 6 -0.02183302 13.63399 0.0000849153 0.99264789 1.0000000     6  1.0000000

Now we have a single logFC which compares nhood abundance between week 4 and week 16 - as we can see the LFC estimates should be the same, but the SpatialFDR will be different.

par(mfrow=c(1, 2))
plot(contrast1.res$logFC.Age4wk...Age16wk, cont.4vs16.res$logFC,
     xlab="4wk vs. 16wk LFC\nsingle contrast", ylab="4wk vs. 16wk LFC\nmultiple contrast")

plot(contrast1.res$SpatialFDR, cont.4vs16.res$SpatialFDR,
     xlab="Spatial FDR\nsingle contrast", ylab="Spatial FDR\nmultiple contrast")

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
## Running with model contrasts
## Performing spatial FDR correction with graph-overlap weighting
table(da_results$SpatialFDR < 0.1)
## 
## FALSE  TRUE 
##   253   107
head(da_results)
##          logFC   logCPM            F       PValue         FDR Nhood  SpatialFDR
## 1  0.009910809 13.61263 9.092468e-04 9.759453e-01 0.999702249     1 0.999702249
## 2  3.291936253 13.42815 1.553354e+01 8.182174e-05 0.001132916     2 0.001066515
## 3  1.289472212 13.66500 2.566105e+00 1.092201e-01 0.251808769     3 0.293599509
## 4 -2.107965537 13.64878 3.063789e-07 9.995584e-01 0.999702249     4 0.999702249
## 5 -1.920583073 13.95991 3.169485e+00 7.506825e-02 0.194421375     5 0.223518569
## 6 -3.078659063 13.63399 6.778060e+00 9.247452e-03 0.044987604     6 0.051234971

The results table In this comparison there are 107 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" )
## Adding nhood effect sizes to neighbourhood graph attributes

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.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.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] MouseThymusAgeing_1.13.0    patchwork_1.3.0            
##  [3] dplyr_1.1.4                 scran_1.34.0               
##  [5] scater_1.34.0               ggplot2_3.5.1              
##  [7] scuttle_1.16.0              SingleCellExperiment_1.28.0
##  [9] SummarizedExperiment_1.36.0 Biobase_2.66.0             
## [11] GenomicRanges_1.58.0        GenomeInfoDb_1.42.0        
## [13] IRanges_2.40.0              S4Vectors_0.44.0           
## [15] BiocGenerics_0.52.0         MatrixGenerics_1.18.0      
## [17] matrixStats_1.4.1           miloR_2.2.0                
## [19] edgeR_4.4.0                 limma_3.62.0               
## [21] BiocStyle_2.34.0           
## 
## loaded via a namespace (and not attached):
##   [1] RColorBrewer_1.1-3      jsonlite_1.8.9          magrittr_2.0.3         
##   [4] magick_2.8.5            ggbeeswarm_0.7.2        farver_2.1.2           
##   [7] rmarkdown_2.28          zlibbioc_1.52.0         vctrs_0.6.5            
##  [10] memoise_2.0.1           tinytex_0.53            htmltools_0.5.8.1      
##  [13] S4Arrays_1.6.0          AnnotationHub_3.14.0    curl_5.2.3             
##  [16] BiocNeighbors_2.0.0     SparseArray_1.6.0       sass_0.4.9             
##  [19] pracma_2.4.4            bslib_0.8.0             cachem_1.1.0           
##  [22] igraph_2.1.1            mime_0.12               lifecycle_1.0.4        
##  [25] pkgconfig_2.0.3         rsvd_1.0.5              Matrix_1.7-1           
##  [28] R6_2.5.1                fastmap_1.2.0           GenomeInfoDbData_1.2.13
##  [31] digest_0.6.37           numDeriv_2016.8-1.1     colorspace_2.1-1       
##  [34] AnnotationDbi_1.68.0    dqrng_0.4.1             irlba_2.3.5.1          
##  [37] ExperimentHub_2.14.0    RSQLite_2.3.7           beachmat_2.22.0        
##  [40] labeling_0.4.3          filelock_1.0.3          fansi_1.0.6            
##  [43] httr_1.4.7              polyclip_1.10-7         abind_1.4-8            
##  [46] compiler_4.4.1          bit64_4.5.2             withr_3.0.2            
##  [49] BiocParallel_1.40.0     viridis_0.6.5           DBI_1.2.3              
##  [52] highr_0.11              ggforce_0.4.2           MASS_7.3-61            
##  [55] rappdirs_0.3.3          DelayedArray_0.32.0     bluster_1.16.0         
##  [58] gtools_3.9.5            tools_4.4.1             vipor_0.4.7            
##  [61] beeswarm_0.4.0          glue_1.8.0              grid_4.4.1             
##  [64] cluster_2.1.6           generics_0.1.3          gtable_0.3.6           
##  [67] tidyr_1.3.1             BiocSingular_1.22.0     tidygraph_1.3.1        
##  [70] ScaledMatrix_1.14.0     metapod_1.14.0          utf8_1.2.4             
##  [73] XVector_0.46.0          ggrepel_0.9.6           BiocVersion_3.20.0     
##  [76] pillar_1.9.0            stringr_1.5.1           splines_4.4.1          
##  [79] tweenr_2.0.3            BiocFileCache_2.14.0    lattice_0.22-6         
##  [82] FNN_1.1.4.1             bit_4.5.0               tidyselect_1.2.1       
##  [85] locfit_1.5-9.10         Biostrings_2.74.0       knitr_1.48             
##  [88] gridExtra_2.3           bookdown_0.41           xfun_0.48              
##  [91] graphlayouts_1.2.0      statmod_1.5.0           stringi_1.8.4          
##  [94] UCSC.utils_1.2.0        yaml_2.3.10             evaluate_1.0.1         
##  [97] codetools_0.2-20        ggraph_2.2.1            tibble_3.2.1           
## [100] BiocManager_1.30.25     cli_3.6.3               uwot_0.2.2             
## [103] munsell_0.5.1           jquerylib_0.1.4         Rcpp_1.0.13            
## [106] dbplyr_2.5.0            png_0.1-8               parallel_4.4.1         
## [109] blob_1.2.4              viridisLite_0.4.2       scales_1.3.0           
## [112] purrr_1.0.2             crayon_1.5.3            rlang_1.1.4            
## [115] cowplot_1.1.3           KEGGREST_1.46.0