--- title: "Differential abundance testing with Milo" author: - Emma Dann - Mike Morgan output: BiocStyle::html_document: toc_float: true BiocStyle::pdf_document: default package: miloR vignette: | %\VignetteIndexEntry{Differential abundance testing with Milo} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} --- ```{r, include = FALSE} knitr::opts_chunk$set( collapse = FALSE ) ``` ```{r setup, message=FALSE, warning=FALSE} library(miloR) library(SingleCellExperiment) library(scater) library(scran) library(dplyr) library(patchwork) ``` # Introduction Milo is a tool for analysis of complex single cell datasets generated from replicated multi-condition experiments, which detects changes in composition between conditions. While differential abundance (DA) is commonly quantified in discrete cell clusters, Milo uses partially overlapping neighbourhoods of cells on a KNN graph. Starting from a graph that faithfully recapitulates the biology of the cell population, Milo analysis consists of 3 steps: 1. Sampling of representative neighbourhoods 2. Testing for differential abundance of conditions in all neighbourhoods 3. Accounting for multiple hypothesis testing using a weighted FDR procedure that accounts for the overlap of neighbourhoods In this vignette we will elaborate on how these steps are implemented in the `miloR` package. # Load data For this demo we will use a synthetic dataset simulating a developmental trajectory, generated using [dyntoy](https://github.com/dynverse/dyntoy). ```{r} data("sim_trajectory", package = "miloR") ## Extract SingleCellExperiment object traj_sce <- sim_trajectory[['SCE']] ## Extract sample metadata to use for testing traj_meta <- sim_trajectory[["meta"]] ## Add metadata to colData slot colData(traj_sce) <- DataFrame(traj_meta) colnames(traj_sce) <- colData(traj_sce)$cell_id redim <- reducedDim(traj_sce, "PCA") dimnames(redim) <- list(colnames(traj_sce), paste0("PC", c(1:50))) reducedDim(traj_sce, "PCA") <- redim ``` # Pre-processing For DA analysis we need to construct an undirected KNN graph of single-cells. Standard single-cell analysis pipelines usually do this from distances in PCA. We normalize and calculate principal components using `scater`. I also run UMAP for visualization purposes. ```{r} logcounts(traj_sce) <- log(counts(traj_sce) + 1) traj_sce <- runPCA(traj_sce, ncomponents=30) traj_sce <- runUMAP(traj_sce) plotUMAP(traj_sce) ``` # Create a Milo object For differential abundance analysis on graph neighbourhoods we first construct a `Milo` object. This extends the [`SingleCellExperiment`](https://bioconductor.org/packages/release/bioc/html/SingleCellExperiment.html) class to store information about neighbourhoods on the KNN graph. ## From SingleCellExperiment object The `Milo` constructor takes as input a `SingleCellExperiment` object. ```{r} traj_milo <- Milo(traj_sce) reducedDim(traj_milo, "UMAP") <- reducedDim(traj_sce, "UMAP") traj_milo ``` ## From AnnData object (.h5ad) We can use the [`zellkonverter`](https://theislab.github.io/zellkonverter/articles/zellkonverter.html) package to make a `SingleCellExperiment` object from an `AnnData` object stored as `h5ad` file. ```{r, eval=FALSE} library(zellkonverter) # Obtaining an example H5AD file. example_h5ad <- system.file("extdata", "krumsiek11.h5ad", package = "zellkonverter") example_h5ad_sce <- readH5AD(example_h5ad) example_h5ad_milo <- Milo(example_h5ad_sce) ``` ## From Seurat object The `Seurat` package includes a converter to `SingleCellExperiment`. ```{r, eval=FALSE} library(Seurat) data("pbmc_small") pbmc_small_sce <- as.SingleCellExperiment(pbmc_small) pbmc_small_milo <- Milo(pbmc_small_sce) ``` # Construct KNN graph We need to add the KNN graph to the Milo object. This is stored in the `graph` slot, in [`igraph`](https://igraph.org/r/) format. The `miloR` package includes functionality to build and store the graph from the PCA dimensions stored in the `reducedDim` slot. ```{r} traj_milo <- buildGraph(traj_milo, k = 10, d = 30) ``` **In progress:** we are perfecting the functionality to add a precomputed KNN graph (for example constructed with Seurat or scanpy) to the `graph` slot using the adjacency matrix. # 1. Defining representative neighbourhoods We define the neighbourhood of a cell, the index, as the group of cells connected by an edge in the KNN graph to the index cell. For efficiency, we don't test for DA in the neighbourhood of every cell, but we sample as indices a subset of representative cells, using a KNN sampling algorithm used by [Gut et al. 2015](https://www.nature.com/articles/nmeth.3545). For sampling you need to define a few parameters: - `prop`: the proportion of cells to randomly sample to start with (usually 0.1 - 0.2 is sufficient) - `k`: the k to use for KNN refinement (we recommend using the same k used for KNN graph building) - `d`: the number of reduced dimensions to use for KNN refinement (we recommend using the same d used for KNN graph building) - `refined` indicated whether you want to use the sampling refinement algorithm, or just pick cells at random. The default and recommended way to go is to use refinement. The only situation in which you might consider using random instead, is if you have batch corrected your data with a graph based correction algorithm, such as [BBKNN](https://github.com/Teichlab/bbknn), but the results of DA testing will be suboptimal. ```{r} traj_milo <- makeNhoods(traj_milo, prop = 0.1, k = 10, d=30, refined = TRUE) ``` Once we have defined neighbourhoods, it's good to take a look at how big the neighbourhoods are (i.e. how many cells form each neighbourhood). This affects the power of DA testing. We can check this out using the `plotNhoodSizeHist` function. Empirically, we found it's best to have a distribution peaking between 50 and 100. Otherwise you might consider rerunning `makeNhoods` increasing `k` and/or `prop` (here the distribution looks ludicrous because it's a small dataset). ```{r} plotNhoodSizeHist(traj_milo) ``` # Counting cells in neighbourhoods Now we have to count how many cells from each sample are in each neighbourhood. We need to use the cell metadata and specify which column contains the sample information. ```{r} traj_milo <- countCells(traj_milo, meta.data = data.frame(colData(traj_milo)), samples="Sample") ``` This adds to the `Milo` object a `n \times m` matrix, where n is the number of neighbourhoods and $m$ is the number of experimental samples. Values indicate the number of cells from each sample counted in a neighbourhood. This count matrix will be used for DA testing. ```{r} head(nhoodCounts(traj_milo)) ``` # Differential abundance testing Now we are all set to test for differential abundance in neighbourhoods. We implement this hypothesis testing in a generalized linear model (GLM) framework, specifically using the Negative Binomial GLM implementation in [`edgeR`](https://bioconductor.org/packages/release/bioc/html/edgeR.html). We first need to think about our experimental design. The design matrix should match samples to a condition of interest. In this case the `Condition` is the covariate we are going to test for. ```{r} traj_design <- data.frame(colData(traj_milo))[,c("Sample", "Condition")] traj_design <- distinct(traj_design) rownames(traj_design) <- traj_design$Sample ## Reorder rownames to match columns of nhoodCounts(milo) traj_design <- traj_design[colnames(nhoodCounts(traj_milo)), , drop=FALSE] traj_design ``` Milo uses an adaptation of the Spatial FDR correction introduced by [cydar](https://bioconductor.org/packages/release/bioc/html/cydar.html), which accounts for the overlap between neighbourhoods. Specifically, each hypothesis test P-value is weighted by the reciprocal of the kth nearest neighbour distance. To use this statistic we first need to store the distances between nearest neighbors in the Milo object. ```{r} traj_milo <- calcNhoodDistance(traj_milo, d=30) ``` Now we can do the test, explicitly defining our experimental design. ```{r} rownames(traj_design) <- traj_design$Sample da_results <- testNhoods(traj_milo, design = ~ Condition, design.df = traj_design) ``` This calculates a Fold-change and corrected P-value for each neighbourhood, which indicates whether there is significant differential abundance between conditions. ```{r} da_results %>% arrange(- SpatialFDR) %>% head() ``` # Visualize neighbourhoods displaying DA To visualize DA results, we build an abstracted graph of neighbourhoods that we can superimpose on the single-cell embedding. ```{r, fig.width=10, fig.height=6} traj_milo <- buildNhoodGraph(traj_milo) plotUMAP(traj_milo) + plotNhoodGraphDA(traj_milo, da_results, alpha=0.05) + plot_layout(guides="collect") ```
**Session Info** ```{r} sessionInfo() ```