library(miloR)
library(SingleCellExperiment)
library(scater)
library(scran)
library(dplyr)
library(patchwork)
library(MouseGastrulationData)

1 Load data

For this vignette we will use the mouse gastrulation single-cell data from Pijuan-Sala et al. 2019. The dataset can be downloaded as a SingleCellExperiment object from the MouseGastrulationData package on Bioconductor. To make computations faster, here we will download just a subset of samples, 4 samples at stage E7 and 4 samples at stage E7.5.

This dataset has already been pre-processed and contains a pca.corrected dimensionality reduction, which was built after batch correction using fastMNN.

select_samples <- c(2,  3,  6, 4, #15,
                    # 19, 
                    10, 14#, 20 #30
                    #31, 32
                    )
embryo_data = EmbryoAtlasData(samples = select_samples)
embryo_data
## class: SingleCellExperiment 
## dim: 29452 7558 
## metadata(0):
## assays(1): counts
## rownames(29452): ENSMUSG00000051951 ENSMUSG00000089699 ...
##   ENSMUSG00000096730 ENSMUSG00000095742
## rowData names(2): ENSEMBL SYMBOL
## colnames(7558): cell_361 cell_362 ... cell_29013 cell_29014
## colData names(17): cell barcode ... colour sizeFactor
## reducedDimNames(2): pca.corrected umap
## mainExpName: NULL
## altExpNames(0):

2 Visualize the data

We recompute the UMAP embedding for this subset of cells to visualize the data.

embryo_data <- embryo_data[,apply(reducedDim(embryo_data, "pca.corrected"), 1, function(x) !all(is.na(x)))]
embryo_data <- runUMAP(embryo_data, dimred = "pca.corrected", name = 'umap')

plotReducedDim(embryo_data, colour_by="stage", dimred = "umap") 

We will test for significant differences in abundance of cells between these stages of development, and the associated gene signatures.

3 Differential abundance testing

3.1 Create a Milo object

For differential abundance analysis on graph neighbourhoods we first construct a Milo object. This extends the SingleCellExperiment class to store information about neighbourhoods on the KNN graph.

embryo_milo <- Milo(embryo_data)
embryo_milo
## class: Milo 
## dim: 29452 6875 
## metadata(0):
## assays(1): counts
## rownames(29452): ENSMUSG00000051951 ENSMUSG00000089699 ...
##   ENSMUSG00000096730 ENSMUSG00000095742
## rowData names(2): ENSEMBL SYMBOL
## colnames(6875): cell_361 cell_362 ... cell_29013 cell_29014
## colData names(17): cell barcode ... colour sizeFactor
## reducedDimNames(2): pca.corrected umap
## mainExpName: NULL
## altExpNames(0):
## nhoods dimensions(2): 1 1
## nhoodCounts dimensions(2): 1 1
## nhoodDistances dimension(1): 0
## graph names(0):
## nhoodIndex names(1): 0
## nhoodExpression dimension(2): 1 1
## nhoodReducedDim names(0):
## nhoodGraph names(0):
## nhoodAdjacency dimension(2): 1 1

3.2 Construct KNN graph

We need to add the KNN graph to the Milo object. This is stored in the graph slot, in igraph format. The miloR package includes functionality to build and store the graph from the PCA dimensions stored in the reducedDim slot. In this case, we specify that we want to build the graph from the MNN corrected PCA dimensions.

For graph building you need to define a few parameters:

  • d: the number of reduced dimensions to use for KNN refinement. We recommend using the same \(d\) used for KNN graph building, or to select PCs by inspecting the scree plot.
  • k: this affects the power of DA testing, since we need to have enough cells from each sample represented in a neighbourhood to estimate the variance between replicates. On the other side, increasing \(k\) too much might lead to over-smoothing. We suggest to start by using the same value for \(k\) used for KNN graph building for clustering and UMAP visualization. We will later use some heuristics to evaluate whether the value of \(k\) should be increased.
embryo_milo <- buildGraph(embryo_milo, k = 30, d = 30, reduced.dim = "pca.corrected")

Alternatively, one can add a precomputed KNN graph (for example constructed with Seurat or scanpy) to the graph slot using the adjacency matrix, through the helper function buildFromAdjacency.

3.3 Defining representative neighbourhoods on the KNN graph

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.

As well as \(d\) and \(k\), for sampling we need to define a few additional parameters:

  • prop: the proportion of cells to randomly sample to start with. We suggest using prop=0.1 for datasets of less than 30k cells. For bigger datasets using prop=0.05 should be sufficient (and makes computation faster).
  • refined: indicates whether you want to use the sampling refinement algorith, 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, but the results of DA testing will be suboptimal.
embryo_milo <- makeNhoods(embryo_milo, prop = 0.1, k = 30, d=30, refined = TRUE, reduced_dims = "pca.corrected")

Once we have defined neighbourhoods, we plot the distribution of neighbourhood sizes (i.e. how many cells form each neighbourhood) to evaluate whether the value of \(k\) used for graph building was appropriate. We can check this out using the plotNhoodSizeHist function.

As a rule of thumb we want to have an average neighbourhood size over 5 x N_samples. If the mean is lower, or if the distribution is

plotNhoodSizeHist(embryo_milo)

3.4 Counting cells in neighbourhoods

Milo leverages the variation in cell numbers between replicates for the same experimental condition to test for differential abundance. Therefore 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.

embryo_milo <- countCells(embryo_milo, meta.data = as.data.frame(colData(embryo_milo)), sample="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.

head(nhoodCounts(embryo_milo))
## 6 x 6 sparse Matrix of class "dgCMatrix"
##    2  3  6 4 10 14
## 1  1  2 14 1 38 17
## 2  8 10 42 4  4  1
## 3  4  9 70 3  3  4
## 4  6  3 50 6  . 16
## 5 13  5 24 7  1  2
## 6  .  .  5 . 89 19

3.5 Defining experimental design

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.

We first need to think about our experimental design. The design matrix should match each sample to the experimental condition of interest for DA testing. In this case, we want to detect DA between embryonic stages, stored in the stage column of the dataset colData. We also include the sequencing.batch column in the design matrix. This represents a known technical covariate that we want to account for in DA testing.

embryo_design <- data.frame(colData(embryo_milo))[,c("sample", "stage", "sequencing.batch")]

## Convert batch info from integer to factor
embryo_design$sequencing.batch <- as.factor(embryo_design$sequencing.batch) 
embryo_design <- distinct(embryo_design)
rownames(embryo_design) <- embryo_design$sample

embryo_design
##    sample stage sequencing.batch
## 2       2  E7.5                1
## 3       3  E7.5                1
## 6       6  E7.5                1
## 4       4  E7.5                1
## 10     10  E7.0                1
## 14     14  E7.0                2

3.6 Computing neighbourhood connectivity

Milo uses an adaptation of the Spatial FDR correction introduced by cydar, where we correct p-values accounting for the amount of 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. This is done by the calcNhoodDistance function (N.B. this step is the most time consuming of the analysis workflow and might take a couple of minutes for large datasets).

embryo_milo <- calcNhoodDistance(embryo_milo, d=30, reduced.dim = "pca.corrected")

3.7 Testing

Now we can do the DA test, explicitly defining our experimental design. In this case, we want to test for differences between experimental stages, while accounting for the variability between technical batches (You can find more info on how to use formulas to define a testing design in R here)

da_results <- testNhoods(embryo_milo, design = ~ sequencing.batch + stage, design.df = embryo_design, reduced.dim="pca.corrected")
head(da_results)
##       logFC   logCPM         F       PValue          FDR Nhood   SpatialFDR
## 1 -2.721173 11.46736 10.457176 1.250647e-03 2.448075e-03     1 2.497085e-03
## 2  2.565522 11.35377  5.190424 2.286396e-02 3.414747e-02     2 3.521091e-02
## 3  2.917284 11.41216  6.401897 1.151030e-02 1.843178e-02     3 1.895083e-02
## 4  6.963673 11.55476 13.147819 2.983308e-04 8.367814e-04     4 8.077549e-04
## 5  4.443449 11.34353 10.020401 1.581997e-03 3.032161e-03     5 3.100921e-03
## 6 -5.816565 11.93378 36.826365 1.665432e-09 5.078089e-08     6 4.149905e-08

This calculates a Fold-change and corrected P-value for each neighbourhood, which indicates wheather there is significant differential abundance between developmental stages. The main statistics we consider here are:

  • logFC: indicates the log-Fold change in cell numbers between samples from E7.5 and samples from E7.0
  • PValue: reports P-values before FDR correction
  • SpatialFDR: reports P-values corrected for multiple testing accounting for overlap between neighbourhoods
da_results %>%
  arrange(SpatialFDR) %>%
  head() 
##         logFC   logCPM        F       PValue          FDR Nhood   SpatialFDR
## 239 -6.959240 11.77269 45.97094 1.774069e-11 8.160718e-09   239 6.385240e-09
## 248 -7.640327 12.09614 43.54613 5.886980e-11 1.354005e-08   248 1.076701e-08
## 38  -5.945704 11.88412 41.36196 1.739368e-10 2.000273e-08    38 1.575844e-08
## 383 -6.418828 11.83716 41.40096 1.706000e-10 2.000273e-08   383 1.575844e-08
## 41  -5.736706 11.77298 39.63585 4.103068e-10 3.133913e-08    41 2.509064e-08
## 78  -6.190863 11.67516 39.33374 4.768998e-10 3.133913e-08    78 2.509064e-08

4 Inspecting DA testing results

We can start inspecting the results of our DA analysis from a couple of standard diagnostic plots. We first inspect the distribution of uncorrected P values, to verify that the test was balanced.

ggplot(da_results, aes(PValue)) + geom_histogram(bins=50)

Then we visualize the test results with a volcano plot (remember that each point here represents a neighbourhood, not a cell).

ggplot(da_results, aes(logFC, -log10(SpatialFDR))) + 
  geom_point() +
  geom_hline(yintercept = 1) ## Mark significance threshold (10% FDR)

Looks like we have detected several neighbourhoods were there is a significant difference in cell abundances between developmental stages.

To visualize DA results relating them to the embedding of single cells, we can build an abstracted graph of neighbourhoods that we can superimpose on the single-cell embedding. Here each node represents a neighbourhood, while edges indicate how many cells two neighbourhoods have in common. Here the layout of nodes is determined by the position of the index cell in the UMAP embedding of all single-cells. The neighbourhoods displaying singificant DA are colored by their log-Fold Change.

embryo_milo <- buildNhoodGraph(embryo_milo)

## Plot single-cell UMAP
umap_pl <- plotReducedDim(embryo_milo, dimred = "umap", colour_by="stage", text_by = "celltype", 
                          text_size = 3, point_size=0.5) +
  guides(fill="none")

## Plot neighbourhood graph
nh_graph_pl <- plotNhoodGraphDA(embryo_milo, da_results, layout="umap",alpha=0.1) 
  
umap_pl + nh_graph_pl +
  plot_layout(guides="collect")
## Warning: ggrepel: 1 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps

We might also be interested in visualizing wheather DA is particularly evident in certain cell types. To do this, we assign a cell type label to each neighbourhood by finding the most abundant cell type within cells in each neighbourhood. We can label neighbourhoods in the results data.frame using the function annotateNhoods. This also saves the fraction of cells harbouring the label.

da_results <- annotateNhoods(embryo_milo, da_results, coldata_col = "celltype")
head(da_results)
##       logFC   logCPM         F       PValue          FDR Nhood   SpatialFDR
## 1 -2.721173 11.46736 10.457176 1.250647e-03 2.448075e-03     1 2.497085e-03
## 2  2.565522 11.35377  5.190424 2.286396e-02 3.414747e-02     2 3.521091e-02
## 3  2.917284 11.41216  6.401897 1.151030e-02 1.843178e-02     3 1.895083e-02
## 4  6.963673 11.55476 13.147819 2.983308e-04 8.367814e-04     4 8.077549e-04
## 5  4.443449 11.34353 10.020401 1.581997e-03 3.032161e-03     5 3.100921e-03
## 6 -5.816565 11.93378 36.826365 1.665432e-09 5.078089e-08     6 4.149905e-08
##                         celltype celltype_fraction
## 1                   ExE endoderm         1.0000000
## 2                 Mixed mesoderm         0.7246377
## 3           Rostral neurectoderm         0.4946237
## 4 Haematoendothelial progenitors         0.6049383
## 5                   ExE endoderm         1.0000000
## 6                       Epiblast         0.7522124

While neighbourhoods tend to be homogeneous, we can define a threshold for celltype_fraction to exclude neighbourhoods that are a mix of cell types.

ggplot(da_results, aes(celltype_fraction)) + geom_histogram(bins=50)

da_results$celltype <- ifelse(da_results$celltype_fraction < 0.7, "Mixed", da_results$celltype)

Now we can visualize the distribution of DA Fold Changes in different cell types

plotDAbeeswarm(da_results, group.by = "celltype")

This is already quite informative: we can see that certain early development cell types, such as epiblast and primitive streak, are enriched in the earliest time stage, while others are enriched later in development, such as ectoderm cells. Interestingly, we also see plenty of DA neighbourhood with a mixed label. This could indicate that transitional states show changes in abundance in time.

5 Finding markers of DA populations

Once you have found your neighbourhoods showindg significant DA between conditions, you might want to find gene signatures specific to the cells in those neighbourhoods. The function findNhoodGroupMarkers runs a one-VS-all differential gene expression test to identify marker genes for a group of neighbourhoods of interest. Before running this function you will need to define your neighbourhood groups depending on your biological question, that need to be stored as a NhoodGroup column in the da_results data.frame.

5.0.1 Custom grouping

In a case where all the DA neighbourhoods seem to belong to the same region of the graph, you might just want to test the significant DA neighbourhoods with the same logFC against all the rest (N.B. for illustration purposes, here I am testing on a randomly selected set of 10 genes).

## Add log normalized count to Milo object
embryo_milo <- logNormCounts(embryo_milo)

da_results$NhoodGroup <- as.numeric(da_results$SpatialFDR < 0.1 & da_results$logFC < 0)
da_nhood_markers <- findNhoodGroupMarkers(embryo_milo, da_results, subset.row = rownames(embryo_milo)[1:10])
## Warning: Zero sample variances detected, have been offset away from zero

## Warning: Zero sample variances detected, have been offset away from zero
head(da_nhood_markers)
##                          logFC_1  adj.P.Val_1       logFC_0  adj.P.Val_0
## ENSMUSG00000025902 -0.0973140537 0.0004260754  0.0973140537 0.0004260754
## ENSMUSG00000025903  0.0626711968 0.0088091897 -0.0626711968 0.0088091897
## ENSMUSG00000033845  0.0639820071 0.0145189433 -0.0639820071 0.0145189433
## ENSMUSG00000033813  0.0408076160 0.1258138161 -0.0408076160 0.1258138161
## ENSMUSG00000104217 -0.0008970000 0.4137276081  0.0008970000 0.4137276081
## ENSMUSG00000102343 -0.0001962896 1.0000000000  0.0001962896 1.0000000000
##                                GeneID
## ENSMUSG00000025902 ENSMUSG00000025902
## ENSMUSG00000025903 ENSMUSG00000025903
## ENSMUSG00000033845 ENSMUSG00000033845
## ENSMUSG00000033813 ENSMUSG00000033813
## ENSMUSG00000104217 ENSMUSG00000104217
## ENSMUSG00000102343 ENSMUSG00000102343

For this analysis we recommend aggregating the neighbourhood expression profiles by experimental samples (the same used for DA testing), by setting aggregate.samples=TRUE. This way single-cells will not be considered as “replicates” during DGE testing, and dispersion will be estimated between true biological replicates. Like so:

da_nhood_markers <- findNhoodGroupMarkers(embryo_milo, da_results, subset.row = rownames(embryo_milo)[1:10], 
                                          aggregate.samples = TRUE, sample_col = "sample")
## Warning: Zero sample variances detected, have been offset away from zero

## Warning: Zero sample variances detected, have been offset away from zero
head(da_nhood_markers)
##                         logFC_1 adj.P.Val_1      logFC_0 adj.P.Val_0
## ENSMUSG00000033845  0.173794674   0.9419092 -0.173794674   0.9419092
## ENSMUSG00000033813  0.143303604   0.9419092 -0.143303604   0.9419092
## ENSMUSG00000025900 -0.001704791   0.9419092  0.001704791   0.9419092
## ENSMUSG00000025902  0.238330971   0.9419092 -0.238330971   0.9419092
## ENSMUSG00000025903  0.124027237   0.9419092 -0.124027237   0.9419092
## ENSMUSG00000104217 -0.000357727   0.9419092  0.000357727   0.9419092
##                                GeneID
## ENSMUSG00000033845 ENSMUSG00000033845
## ENSMUSG00000033813 ENSMUSG00000033813
## ENSMUSG00000025900 ENSMUSG00000025900
## ENSMUSG00000025902 ENSMUSG00000025902
## ENSMUSG00000025903 ENSMUSG00000025903
## ENSMUSG00000104217 ENSMUSG00000104217

(Notice the difference in p values)

5.1 Automatic grouping of neighbourhoods

In many cases, such as this example, DA neighbourhoods are found in different areas of the KNN graph, and grouping together all significant DA populations might not be ideal, as they might include cells of very different celltypes. For this kind of scenario, we have implemented a neighbourhood function that uses community detection to partition neighbourhoods into groups on the basis of (1) the number of shared cells between 2 neighbourhoods; (2) the direction of fold-change for DA neighbourhoods; (3) the difference in fold change.

## Run buildNhoodGraph to store nhood adjacency matrix
embryo_milo <- buildNhoodGraph(embryo_milo)

## Find groups
da_results <- groupNhoods(embryo_milo, da_results, max.lfc.delta = 10)
head(da_results)
##       logFC   logCPM         F       PValue          FDR Nhood   SpatialFDR
## 1 -2.721173 11.46736 10.457176 1.250647e-03 2.448075e-03     1 2.497085e-03
## 2  2.565522 11.35377  5.190424 2.286396e-02 3.414747e-02     2 3.521091e-02
## 3  2.917284 11.41216  6.401897 1.151030e-02 1.843178e-02     3 1.895083e-02
## 4  6.963673 11.55476 13.147819 2.983308e-04 8.367814e-04     4 8.077549e-04
## 5  4.443449 11.34353 10.020401 1.581997e-03 3.032161e-03     5 3.100921e-03
## 6 -5.816565 11.93378 36.826365 1.665432e-09 5.078089e-08     6 4.149905e-08
##         celltype celltype_fraction NhoodGroup
## 1   ExE endoderm         1.0000000          1
## 2 Mixed mesoderm         0.7246377          2
## 3          Mixed         0.4946237          3
## 4          Mixed         0.6049383          4
## 5   ExE endoderm         1.0000000          1
## 6       Epiblast         0.7522124          5

Let’s have a look at the detected groups

plotNhoodGroups(embryo_milo, da_results, layout="umap") 

plotDAbeeswarm(da_results, "NhoodGroup")

We can see that some groups are made up of single neighbourhoods.

We can easily check how changing the grouping parameters changes the groups we obtain, starting with the LFC delta by plotting with different values of max.lfc.delta (not executed here).

# code not run - uncomment to run.
# plotDAbeeswarm(groupNhoods(embryo_milo, da_results, max.lfc.delta = 1) , group.by = "NhoodGroup") + ggtitle("max LFC delta=1")
# plotDAbeeswarm(groupNhoods(embryo_milo, da_results, max.lfc.delta = 2)   , group.by = "NhoodGroup") + ggtitle("max LFC delta=2")
# plotDAbeeswarm(groupNhoods(embryo_milo, da_results, max.lfc.delta = 3)   , group.by = "NhoodGroup") + ggtitle("max LFC delta=3")

…and we can do the same for the minimum overlap between neighbourhoods… (code not executed).

# code not run - uncomment to run.
# plotDAbeeswarm(groupNhoods(embryo_milo, da_results, max.lfc.delta = 5, overlap=1) , group.by = "NhoodGroup") + ggtitle("overlap=5")
# plotDAbeeswarm(groupNhoods(embryo_milo, da_results, max.lfc.delta = 5, overlap=5)   , group.by = "NhoodGroup") + ggtitle("overlap=10")
# plotDAbeeswarm(groupNhoods(embryo_milo, da_results, max.lfc.delta = 5, overlap=10)   , group.by = "NhoodGroup") + ggtitle("overlap=20")

In these examples we settle for overlap=5 and max.lfc.delta=5, as we need at least 2 neighbourhoods assigned to each group.

set.seed(42)
da_results <- groupNhoods(embryo_milo, da_results, max.lfc.delta = 5, overlap=5)
plotNhoodGroups(embryo_milo, da_results, layout="umap")

plotDAbeeswarm(da_results, group.by = "NhoodGroup")

5.2 Finding gene signatures for neighbourhoods

Once we have grouped neighbourhoods using groupNhoods we are now all set to identifying gene signatures between neighbourhood groups.

Let’s restrict the testing to highly variable genes in this case

## Exclude zero counts genes
keep.rows <- rowSums(logcounts(embryo_milo)) != 0
embryo_milo <- embryo_milo[keep.rows, ]

## Find HVGs
set.seed(101)
dec <- modelGeneVar(embryo_milo)
hvgs <- getTopHVGs(dec, n=2000)
head(hvgs)
## [1] "ENSMUSG00000032083" "ENSMUSG00000095180" "ENSMUSG00000061808"
## [4] "ENSMUSG00000002985" "ENSMUSG00000024990" "ENSMUSG00000024391"

We run findNhoodGroupMarkers to test for one-vs-all differential gene expression for each neighbourhood group

set.seed(42)
nhood_markers <- findNhoodGroupMarkers(embryo_milo, da_results, subset.row = hvgs, 
                                       aggregate.samples = TRUE, sample_col = "sample")

head(nhood_markers)
##                     logFC_1  adj.P.Val_1   logFC_2  adj.P.Val_2   logFC_3
## ENSMUSG00000043501 2.948745 2.220034e-37 0.9693472 1.723093e-19 0.4103104
## ENSMUSG00000035041 1.264302 6.126952e-33 1.0287601 1.237509e-17 0.7620470
## ENSMUSG00000031883 2.761650 5.276126e-31 2.1489992 1.189918e-16 0.6042184
## ENSMUSG00000054422 1.276137 3.329111e-30 0.2096194 1.158245e-15 0.5412801
## ENSMUSG00000028435 1.828046 1.174155e-29 0.2420920 5.706957e-15 0.5503038
## ENSMUSG00000005681 1.871771 2.275887e-27 0.3619340 2.493124e-14 0.3695746
##                     adj.P.Val_3   logFC_4  adj.P.Val_4    logFC_5  adj.P.Val_5
## ENSMUSG00000043501 1.460601e-11 3.1439261 2.439132e-44  0.5387108 0.0006730315
## ENSMUSG00000035041 2.779088e-04 1.6199004 2.560820e-34 -0.9460564 0.0016831686
## ENSMUSG00000031883 3.360353e-03 2.1633663 4.577294e-29 -1.1081364 0.0016831686
## ENSMUSG00000054422 4.489286e-03 1.6483094 4.577294e-29 -1.0007066 0.0023191663
## ENSMUSG00000028435 1.796440e-02 0.8458409 2.786815e-21 -0.7865266 0.0028215534
## ENSMUSG00000005681 1.796440e-02 1.8037589 3.216965e-19  0.8796745 0.0035658194
##                      logFC_6  adj.P.Val_6   logFC_7  adj.P.Val_7   logFC_8
## ENSMUSG00000043501 2.6310877 1.141566e-50 0.4810104 1.375964e-06 0.7679110
## ENSMUSG00000035041 1.5816460 9.906948e-38 0.3660816 1.375964e-06 0.2772358
## ENSMUSG00000031883 1.8730565 1.449204e-37 0.6355904 4.453254e-05 0.4205827
## ENSMUSG00000054422 0.9415344 2.132318e-36 1.8820759 7.414995e-05 1.0350622
## ENSMUSG00000028435 1.4970619 2.611012e-35 0.4747850 4.957546e-04 0.7366430
## ENSMUSG00000005681 2.4500122 2.611012e-35 0.3750758 5.491657e-04 0.4446317
##                     adj.P.Val_8   logFC_9  adj.P.Val_9 logFC_10 adj.P.Val_10
## ENSMUSG00000043501 0.0005813141 0.8340250 5.777537e-32 3.326014 8.259188e-55
## ENSMUSG00000035041 0.0005813141 1.4668628 7.693912e-25 3.939964 3.610358e-51
## ENSMUSG00000031883 0.0005813141 0.6190077 1.338025e-23 3.689324 1.104441e-45
## ENSMUSG00000054422 0.0005813141 0.5922399 1.338025e-23 4.431534 2.084097e-37
## ENSMUSG00000028435 0.0005813141 1.4579824 4.412035e-22 3.293152 3.789401e-37
## ENSMUSG00000005681 0.0005813141 0.7231761 5.924794e-22 3.173665 5.521996e-35
##                                GeneID
## ENSMUSG00000043501 ENSMUSG00000043501
## ENSMUSG00000035041 ENSMUSG00000035041
## ENSMUSG00000031883 ENSMUSG00000031883
## ENSMUSG00000054422 ENSMUSG00000054422
## ENSMUSG00000028435 ENSMUSG00000028435
## ENSMUSG00000005681 ENSMUSG00000005681

Let’s check out the markers for group 5 for example

gr5_markers <- nhood_markers[c("logFC_5", "adj.P.Val_5")] 
colnames(gr5_markers) <- c("logFC", "adj.P.Val")

head(gr5_markers[order(gr5_markers$adj.P.Val), ])
##                         logFC    adj.P.Val
## ENSMUSG00000043501  0.5387108 0.0006730315
## ENSMUSG00000035041 -0.9460564 0.0016831686
## ENSMUSG00000031883 -1.1081364 0.0016831686
## ENSMUSG00000054422 -1.0007066 0.0023191663
## ENSMUSG00000028435 -0.7865266 0.0028215534
## ENSMUSG00000005681  0.8796745 0.0035658194

If you already know you are interested only in the markers for group 2, you might want to test just 8-VS-all using the subset.groups parameter:

nhood_markers <- findNhoodGroupMarkers(embryo_milo, da_results, subset.row = hvgs, 
                                       aggregate.samples = TRUE, sample_col = "sample",
                                       subset.groups = c("5")
                                       )
## Warning in all(structure(c(0.538710819617361, -0.94605639493378,
## -1.10813638270555, : coercing argument of type 'double' to logical
head(nhood_markers)
##                       logFC_5  adj.P.Val_5             GeneID
## ENSMUSG00000019796  0.5387108 0.0006730315 ENSMUSG00000019796
## ENSMUSG00000005575 -0.9460564 0.0016831686 ENSMUSG00000005575
## ENSMUSG00000034349 -1.1081364 0.0016831686 ENSMUSG00000034349
## ENSMUSG00000027496 -1.0007066 0.0023191663 ENSMUSG00000027496
## ENSMUSG00000029177 -0.7865266 0.0028215534 ENSMUSG00000029177
## ENSMUSG00000031586  0.8796745 0.0035658194 ENSMUSG00000031586

You might also want to compare a subset of neighbourhoods between each other. You can specify the neighbourhoods to use for testing by setting the parameter subset.nhoods.

For example, you might want to compare just one pair of neighbourhood groups against each other:

nhood_markers <- findNhoodGroupMarkers(embryo_milo, da_results, subset.row = hvgs,
                                       subset.nhoods = da_results$NhoodGroup %in% c('5','6'),
                                       aggregate.samples = TRUE, sample_col = "sample")
## Warning: Zero sample variances detected, have been offset away from zero

## Warning: Zero sample variances detected, have been offset away from zero
head(nhood_markers)
##                      logFC_5  adj.P.Val_5   logFC_6  adj.P.Val_6
## ENSMUSG00000095180 -5.094314 4.614495e-12  5.094314 4.614495e-12
## ENSMUSG00000028832  2.892024 4.164581e-10 -2.892024 4.164581e-10
## ENSMUSG00000024406  3.985672 4.164581e-10 -3.985672 4.164581e-10
## ENSMUSG00000028773 -3.745348 4.486631e-10  3.745348 4.486631e-10
## ENSMUSG00000027239  2.708278 1.165526e-09 -2.708278 1.165526e-09
## ENSMUSG00000047139  2.834740 1.165526e-09 -2.834740 1.165526e-09
##                                GeneID
## ENSMUSG00000095180 ENSMUSG00000095180
## ENSMUSG00000028832 ENSMUSG00000028832
## ENSMUSG00000024406 ENSMUSG00000024406
## ENSMUSG00000028773 ENSMUSG00000028773
## ENSMUSG00000027239 ENSMUSG00000027239
## ENSMUSG00000047139 ENSMUSG00000047139

or you might use subset.nhoods to exclude singleton neighbourhoods, or to subset to the neighbourhoods that show significant DA.

5.3 Visualize detected markers

Lets select marker genes for group 10 at FDR 1% and log-fold-Change > 1.

ggplot(nhood_markers, aes(logFC_5, -log10(adj.P.Val_5 ))) + 
  geom_point(alpha=0.5, size=0.5) +
  geom_hline(yintercept = 3)

markers <- rownames(nhood_markers)[nhood_markers$adj.P.Val_5 < 0.01 & nhood_markers$logFC_5 > 0]

We can visualize the expression in neighbourhoods using plotNhoodExpressionGroups.

set.seed(42)
plotNhoodExpressionGroups(embryo_milo, da_results, features=intersect(rownames(embryo_milo), markers[1:10]),
                          subset.nhoods = da_results$NhoodGroup %in% c('6','5'), 
                          scale=TRUE,
                          grid.space = "fixed")
## Warning in plotNhoodExpressionGroups(embryo_milo, da_results, features =
## intersect(rownames(embryo_milo), : Nothing in nhoodExpression(x): computing for
## requested features...

5.4 DGE testing within a group

In some cases you might want to test for differential expression between cells in different conditions within the same neighbourhood group. You can do that using testDiffExp:

dge_6 <- testDiffExp(embryo_milo, da_results, design = ~ stage, meta.data = data.frame(colData(embryo_milo)),
                     subset.row = rownames(embryo_milo)[1:5], subset.nhoods=da_results$NhoodGroup=="6")

dge_6
## $`6`
##                          logFC     AveExpr         t     P.Value  adj.P.Val
## ENSMUSG00000025902 -0.01121445 0.008496462 -2.803917 0.005111052 0.02555526
## ENSMUSG00000033845  0.03969262 2.452590660  1.282031 0.200022200 0.50005550
## ENSMUSG00000051951  0.00000000 0.000000000  0.000000 1.000000000 1.00000000
## ENSMUSG00000102343  0.00000000 0.000000000  0.000000 1.000000000 1.00000000
## ENSMUSG00000025900  0.00000000 0.000000000  0.000000 1.000000000 1.00000000
##                            B Nhood.Group
## ENSMUSG00000025902  -7.31175           6
## ENSMUSG00000033845 -10.41348           6
## ENSMUSG00000051951 -11.23537           6
## ENSMUSG00000102343 -11.23537           6
## ENSMUSG00000025900 -11.23537           6

Session Info

sessionInfo()
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-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] MouseGastrulationData_1.7.0 SpatialExperiment_1.4.0    
##  [3] patchwork_1.1.1             dplyr_1.0.7                
##  [5] scran_1.22.0                scater_1.22.0              
##  [7] ggplot2_3.3.5               scuttle_1.4.0              
##  [9] SingleCellExperiment_1.16.0 SummarizedExperiment_1.24.0
## [11] Biobase_2.54.0              GenomicRanges_1.46.0       
## [13] GenomeInfoDb_1.30.0         IRanges_2.28.0             
## [15] S4Vectors_0.32.0            BiocGenerics_0.40.0        
## [17] MatrixGenerics_1.6.0        matrixStats_0.61.0         
## [19] miloR_1.2.0                 edgeR_3.36.0               
## [21] limma_3.50.0                BiocStyle_2.22.0           
## 
## loaded via a namespace (and not attached):
##   [1] AnnotationHub_3.2.0           BiocFileCache_2.2.0          
##   [3] igraph_1.2.7                  splines_4.1.1                
##   [5] BiocParallel_1.28.0           digest_0.6.28                
##   [7] BumpyMatrix_1.2.0             htmltools_0.5.2              
##   [9] viridis_0.6.2                 magick_2.7.3                 
##  [11] fansi_0.5.0                   memoise_2.0.0                
##  [13] magrittr_2.0.1                ScaledMatrix_1.2.0           
##  [15] cluster_2.1.2                 Biostrings_2.62.0            
##  [17] graphlayouts_0.7.1            R.utils_2.11.0               
##  [19] colorspace_2.0-2              rappdirs_0.3.3               
##  [21] blob_1.2.2                    ggrepel_0.9.1                
##  [23] xfun_0.27                     crayon_1.4.1                 
##  [25] RCurl_1.98-1.5                jsonlite_1.7.2               
##  [27] glue_1.4.2                    polyclip_1.10-0              
##  [29] gtable_0.3.0                  zlibbioc_1.40.0              
##  [31] XVector_0.34.0                DelayedArray_0.20.0          
##  [33] BiocSingular_1.10.0           DropletUtils_1.14.0          
##  [35] Rhdf5lib_1.16.0               HDF5Array_1.22.0             
##  [37] scales_1.1.1                  DBI_1.1.1                    
##  [39] Rcpp_1.0.7                    xtable_1.8-4                 
##  [41] viridisLite_0.4.0             dqrng_0.3.0                  
##  [43] bit_4.0.4                     rsvd_1.0.5                   
##  [45] httr_1.4.2                    metapod_1.2.0                
##  [47] FNN_1.1.3                     RColorBrewer_1.1-2           
##  [49] ellipsis_0.3.2                pkgconfig_2.0.3              
##  [51] R.methodsS3_1.8.1             farver_2.1.0                 
##  [53] dbplyr_2.1.1                  sass_0.4.0                   
##  [55] uwot_0.1.10                   locfit_1.5-9.4               
##  [57] utf8_1.2.2                    AnnotationDbi_1.56.0         
##  [59] later_1.3.0                   tidyselect_1.1.1             
##  [61] labeling_0.4.2                rlang_0.4.12                 
##  [63] BiocVersion_3.14.0            cachem_1.0.6                 
##  [65] munsell_0.5.0                 tools_4.1.1                  
##  [67] generics_0.1.1                RSQLite_2.2.8                
##  [69] ExperimentHub_2.2.0           evaluate_0.14                
##  [71] stringr_1.4.0                 fastmap_1.1.0                
##  [73] yaml_2.2.1                    knitr_1.36                   
##  [75] bit64_4.0.5                   tidygraph_1.2.0              
##  [77] purrr_0.3.4                   KEGGREST_1.34.0              
##  [79] ggraph_2.0.5                  sparseMatrixStats_1.6.0      
##  [81] mime_0.12                     R.oo_1.24.0                  
##  [83] compiler_4.1.1                png_0.1-7                    
##  [85] interactiveDisplayBase_1.32.0 filelock_1.0.2               
##  [87] curl_4.3.2                    beeswarm_0.4.0               
##  [89] tibble_3.1.5                  statmod_1.4.36               
##  [91] tweenr_1.0.2                  bslib_0.3.1                  
##  [93] stringi_1.7.5                 highr_0.9                    
##  [95] RSpectra_0.16-0               lattice_0.20-45              
##  [97] bluster_1.4.0                 Matrix_1.3-4                 
##  [99] vctrs_0.3.8                   pillar_1.6.4                 
## [101] lifecycle_1.0.1               rhdf5filters_1.6.0           
## [103] BiocManager_1.30.16           jquerylib_0.1.4              
## [105] RcppAnnoy_0.0.19              BiocNeighbors_1.12.0         
## [107] cowplot_1.1.1                 bitops_1.0-7                 
## [109] irlba_2.3.3                   httpuv_1.6.3                 
## [111] R6_2.5.1                      promises_1.2.0.1             
## [113] bookdown_0.24                 gridExtra_2.3                
## [115] vipor_0.4.5                   codetools_0.2-18             
## [117] MASS_7.3-54                   gtools_3.9.2                 
## [119] assertthat_0.2.1              rhdf5_2.38.0                 
## [121] rjson_0.2.20                  withr_2.4.2                  
## [123] GenomeInfoDbData_1.2.7        parallel_4.1.1               
## [125] grid_4.1.1                    beachmat_2.10.0              
## [127] tidyr_1.1.4                   rmarkdown_2.11               
## [129] DelayedMatrixStats_1.16.0     ggforce_0.3.3                
## [131] shiny_1.7.1                   ggbeeswarm_0.6.0