library(miloR)
library(SingleCellExperiment)
library(scater)
library(scran)
library(dplyr)
library(patchwork)
library(MouseGastrulationData)
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):
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.
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.
## 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
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.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
.
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 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, 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
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 × 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.
## 6 x 6 sparse Matrix of class "dgCMatrix"
## 2 3 6 4 10 14
## 1 1 10 48 3 25 18
## 2 . . 5 . 70 17
## 3 3 13 29 9 1 1
## 4 6 9 33 . 1 1
## 5 . . 17 1 60 13
## 6 . . 20 1 58 21
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
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).
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 -0.7237126 11.90797 0.6365293 4.251116e-01 4.884959e-01 1 4.914081e-01
## 2 -5.7598806 11.86201 34.7361033 4.769964e-09 6.461133e-08 2 5.382106e-08
## 3 4.1534123 11.31723 8.8186169 3.034725e-03 5.404471e-03 3 5.479943e-03
## 4 3.6061367 10.89608 7.0073522 8.212665e-03 1.349655e-02 4 1.374992e-02
## 5 -3.9437343 11.73461 22.5565767 2.259575e-06 8.782871e-06 5 7.825692e-06
## 6 -3.7417990 11.91113 18.3847843 1.934110e-05 6.219763e-05 6 5.661026e-05
This calculates a Fold-change and corrected P-value for each neighbourhood, which indicates whether 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.0PValue
: reports P-values before FDR correctionSpatialFDR
: reports P-values corrected for multiple
testing accounting for overlap between neighbourhoods## logFC logCPM F PValue FDR Nhood SpatialFDR
## 367 -7.454803 11.54279 54.30281 2.993680e-13 1.338175e-10 367 1.048427e-10
## 331 -6.663261 11.50129 49.77912 2.750012e-12 6.146276e-10 331 4.914137e-10
## 12 -7.159055 11.54040 48.89025 4.257272e-12 6.343335e-10 12 5.184586e-10
## 172 -8.491592 11.36422 47.65750 7.810205e-12 8.727904e-10 172 7.098436e-10
## 125 -6.388135 11.62419 46.95398 1.104635e-11 8.915607e-10 125 7.157460e-10
## 239 -8.856344 11.82424 46.72405 1.237230e-11 8.915607e-10 239 7.157460e-10
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.
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 significant 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: 8 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps
We might also be interested in visualizing whether 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.
## logFC logCPM F PValue FDR Nhood SpatialFDR
## 1 -0.7237126 11.90797 0.6365293 4.251116e-01 4.884959e-01 1 4.914081e-01
## 2 -5.7598806 11.86201 34.7361033 4.769964e-09 6.461133e-08 2 5.382106e-08
## 3 4.1534123 11.31723 8.8186169 3.034725e-03 5.404471e-03 3 5.479943e-03
## 4 3.6061367 10.89608 7.0073522 8.212665e-03 1.349655e-02 4 1.374992e-02
## 5 -3.9437343 11.73461 22.5565767 2.259575e-06 8.782871e-06 5 7.825692e-06
## 6 -3.7417990 11.91113 18.3847843 1.934110e-05 6.219763e-05 6 5.661026e-05
## celltype celltype_fraction
## 1 ExE endoderm 1.0000000
## 2 Primitive Streak 0.7826087
## 3 ExE endoderm 1.0000000
## 4 Rostral neurectoderm 0.6400000
## 5 ExE ectoderm 1.0000000
## 6 Primitive Streak 1.0000000
While neighbourhoods tend to be homogeneous, we can define a
threshold for celltype_fraction
to exclude neighbourhoods
that are a mix of cell types.
Now we can visualize the distribution of DA Fold Changes in different cell types
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.
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.
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
## GeneID logFC_0 adj.P.Val_0 logFC_1 adj.P.Val_1
## 1 ENSMUSG00000025900 -0.0001177587 1.000000e+00 0.0001177587 1.000000e+00
## 2 ENSMUSG00000025902 0.1167223707 3.821486e-06 -0.1167223707 3.821486e-06
## 3 ENSMUSG00000025903 -0.0767830322 5.702953e-04 0.0767830322 5.702953e-04
## 4 ENSMUSG00000033813 -0.0398163352 1.297670e-01 0.0398163352 1.297670e-01
## 5 ENSMUSG00000033845 -0.0557989158 4.014127e-02 0.0557989158 4.014127e-02
## 6 ENSMUSG00000051951 0.0005291873 5.177551e-01 -0.0005291873 5.177551e-01
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
## GeneID logFC_0 adj.P.Val_0 logFC_1 adj.P.Val_1
## 1 ENSMUSG00000025900 0.0008840373 1 -0.0008840373 1
## 2 ENSMUSG00000025902 -0.0768508384 1 0.0768508384 1
## 3 ENSMUSG00000025903 -0.1378004915 1 0.1378004915 1
## 4 ENSMUSG00000033813 -0.1082667762 1 0.1082667762 1
## 5 ENSMUSG00000033845 -0.0723378847 1 0.0723378847 1
## 6 ENSMUSG00000051951 0.0001224230 1 -0.0001224230 1
(Notice the difference in p values)
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 -0.7237126 11.90797 0.6365293 4.251116e-01 4.884959e-01 1 4.914081e-01
## 2 -5.7598806 11.86201 34.7361033 4.769964e-09 6.461133e-08 2 5.382106e-08
## 3 4.1534123 11.31723 8.8186169 3.034725e-03 5.404471e-03 3 5.479943e-03
## 4 3.6061367 10.89608 7.0073522 8.212665e-03 1.349655e-02 4 1.374992e-02
## 5 -3.9437343 11.73461 22.5565767 2.259575e-06 8.782871e-06 5 7.825692e-06
## 6 -3.7417990 11.91113 18.3847843 1.934110e-05 6.219763e-05 6 5.661026e-05
## celltype celltype_fraction NhoodGroup
## 1 ExE endoderm 1.0000000 1
## 2 Primitive Streak 0.7826087 2
## 3 ExE endoderm 1.0000000 1
## 4 Mixed 0.6400000 3
## 5 ExE ectoderm 1.0000000 4
## 6 Primitive Streak 1.0000000 3
Let’s have a look at the detected groups
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.
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)
# this vignette randomly fails to identify HVGs for some reason
if(!length(hvgs)){
set.seed(42)
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)
## GeneID logFC_1 adj.P.Val_1 logFC_2 adj.P.Val_2
## 1 ENSMUSG00000000031 3.114645697 9.503699e-08 -1.37296901 0.1337551
## 2 ENSMUSG00000000078 0.002827669 9.907720e-01 -0.44805855 0.1125936
## 3 ENSMUSG00000000088 0.749520703 3.367626e-03 -0.15473683 0.6375747
## 4 ENSMUSG00000000125 0.028639091 8.049136e-01 -0.13564108 0.3355058
## 5 ENSMUSG00000000149 0.249528009 6.234720e-03 -0.17099410 0.1740285
## 6 ENSMUSG00000000184 -1.019165296 2.919158e-02 -0.06536514 0.9074434
## logFC_3 adj.P.Val_3 logFC_4 adj.P.Val_4 logFC_5 adj.P.Val_5
## 1 -1.34889993 0.21103043 1.2076313 0.17151372 -1.57106689 0.6517954
## 2 -0.04870682 0.87929475 -0.1492375 0.58604515 -0.22423915 0.6616877
## 3 -0.31015113 0.44129745 -0.1928773 0.57307724 -0.33827949 0.6572689
## 4 0.31217904 0.04053752 -0.1530002 0.26441829 0.09458576 0.6995623
## 5 -0.10471826 0.46328943 -0.1267350 0.28724341 0.04878042 0.8216152
## 6 0.92296507 0.19806599 -1.0160922 0.06705351 0.82300344 0.6572689
## logFC_6 adj.P.Val_6 logFC_7 adj.P.Val_7 logFC_8 adj.P.Val_8
## 1 -0.54441683 0.6414584395 1.2069374 0.1720793 0.2563076 8.546921e-01
## 2 0.03138011 0.9403742330 0.2551751 0.3947142 1.0809088 4.377328e-05
## 3 -0.32250015 0.4360183655 0.4160891 0.2342569 0.6036689 1.442214e-01
## 4 0.14667718 0.3899163019 -0.1302541 0.3895379 -0.2005313 2.447164e-01
## 5 -0.07856908 0.6091230945 0.1200438 0.3560700 0.2678906 5.842356e-02
## 6 1.79248682 0.0001653437 -0.4390981 0.5145155 -0.9371240 1.967578e-01
## logFC_9 adj.P.Val_9
## 1 -1.25665952 0.2989986
## 2 -0.04711928 0.9185011
## 3 -0.31533884 0.5492755
## 4 -0.06432338 0.7821122
## 5 -0.14043397 0.4244043
## 6 -0.59534929 0.5082156
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
## 1941 0.8134631 0.04170456
## 1607 0.5352858 0.12799373
## 1351 1.3770657 0.14523135
## 584 0.3325138 0.24397031
## 309 0.3027518 0.49260948
## 565 0.3194266 0.49260948
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")
)
head(nhood_markers)
## logFC_5 adj.P.Val_5 GeneID
## ENSMUSG00000079042 0.8134631 0.04170456 ENSMUSG00000079042
## ENSMUSG00000050288 0.5352858 0.12799373 ENSMUSG00000050288
## ENSMUSG00000038192 1.3770657 0.14523135 ENSMUSG00000038192
## ENSMUSG00000024247 0.3325138 0.24397031 ENSMUSG00000024247
## ENSMUSG00000028152 0.3784290 0.49260948 ENSMUSG00000028152
## ENSMUSG00000025856 0.6201213 0.49260948 ENSMUSG00000025856
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")
head(nhood_markers)
## GeneID logFC_5 adj.P.Val_5 logFC_6 adj.P.Val_6
## 1 ENSMUSG00000000031 -1.099280300 0.006802748 1.099280300 0.006802748
## 2 ENSMUSG00000000078 -0.201138364 0.363315687 0.201138364 0.363315687
## 3 ENSMUSG00000000088 -0.003216687 0.989599997 0.003216687 0.989599997
## 4 ENSMUSG00000000125 0.047200870 0.823237021 -0.047200870 0.823237021
## 5 ENSMUSG00000000149 0.006896198 0.915659291 -0.006896198 0.915659291
## 6 ENSMUSG00000000184 -0.799865427 0.086422720 0.799865427 0.086422720
or you might use subset.nhoods
to exclude singleton
neighbourhoods, or to subset to the neighbourhoods that show significant
DA.
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)
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...
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.02961019 0.04735079 -2.1933784 0.02846931 0.1423466
## ENSMUSG00000033845 -0.01424245 2.28427474 -0.3275927 0.74327616 1.0000000
## ENSMUSG00000051951 0.00000000 0.00000000 0.0000000 1.00000000 1.0000000
## ENSMUSG00000102343 0.00000000 0.00000000 0.0000000 1.00000000 1.0000000
## ENSMUSG00000025900 0.00000000 0.00000000 0.0000000 1.00000000 1.0000000
## B Nhood.Group
## ENSMUSG00000025902 -8.784328 6
## ENSMUSG00000033845 -11.133300 6
## ENSMUSG00000051951 -11.187000 6
## ENSMUSG00000102343 -11.187000 6
## ENSMUSG00000025900 -11.187000 6
## R version 4.4.2 (2024-10-31)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 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: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] MouseGastrulationData_1.20.0 SpatialExperiment_1.17.0
## [3] MouseThymusAgeing_1.14.0 patchwork_1.3.0
## [5] dplyr_1.1.4 scran_1.35.0
## [7] scater_1.35.0 ggplot2_3.5.1
## [9] scuttle_1.17.0 SingleCellExperiment_1.29.1
## [11] SummarizedExperiment_1.37.0 Biobase_2.67.0
## [13] GenomicRanges_1.59.1 GenomeInfoDb_1.43.2
## [15] IRanges_2.41.1 S4Vectors_0.45.2
## [17] BiocGenerics_0.53.3 generics_0.1.3
## [19] MatrixGenerics_1.19.0 matrixStats_1.4.1
## [21] miloR_2.3.0 edgeR_4.5.0
## [23] limma_3.63.2 BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] RColorBrewer_1.1-3 sys_3.4.3 jsonlite_1.8.9
## [4] magrittr_2.0.3 magick_2.8.5 ggbeeswarm_0.7.2
## [7] farver_2.1.2 rmarkdown_2.29 zlibbioc_1.52.0
## [10] vctrs_0.6.5 memoise_2.0.1 htmltools_0.5.8.1
## [13] S4Arrays_1.7.1 AnnotationHub_3.15.0 curl_6.0.1
## [16] BiocNeighbors_2.1.1 SparseArray_1.7.2 sass_0.4.9
## [19] pracma_2.4.4 bslib_0.8.0 cachem_1.1.0
## [22] buildtools_1.0.0 igraph_2.1.1 mime_0.12
## [25] lifecycle_1.0.4 pkgconfig_2.0.3 rsvd_1.0.5
## [28] Matrix_1.7-1 R6_2.5.1 fastmap_1.2.0
## [31] GenomeInfoDbData_1.2.13 digest_0.6.37 numDeriv_2016.8-1.1
## [34] colorspace_2.1-1 AnnotationDbi_1.69.0 dqrng_0.4.1
## [37] irlba_2.3.5.1 ExperimentHub_2.15.0 RSQLite_2.3.8
## [40] beachmat_2.23.2 labeling_0.4.3 filelock_1.0.3
## [43] fansi_1.0.6 httr_1.4.7 polyclip_1.10-7
## [46] abind_1.4-8 compiler_4.4.2 bit64_4.5.2
## [49] withr_3.0.2 BiocParallel_1.41.0 viridis_0.6.5
## [52] DBI_1.2.3 ggforce_0.4.2 MASS_7.3-61
## [55] rappdirs_0.3.3 DelayedArray_0.33.2 rjson_0.2.23
## [58] bluster_1.17.0 gtools_3.9.5 tools_4.4.2
## [61] vipor_0.4.7 beeswarm_0.4.0 glue_1.8.0
## [64] grid_4.4.2 cluster_2.1.6 gtable_0.3.6
## [67] tidyr_1.3.1 BiocSingular_1.23.0 tidygraph_1.3.1
## [70] ScaledMatrix_1.15.0 metapod_1.15.0 utf8_1.2.4
## [73] XVector_0.47.0 RcppAnnoy_0.0.22 ggrepel_0.9.6
## [76] BiocVersion_3.21.1 pillar_1.9.0 stringr_1.5.1
## [79] BumpyMatrix_1.15.0 splines_4.4.2 tweenr_2.0.3
## [82] BiocFileCache_2.15.0 lattice_0.22-6 FNN_1.1.4.1
## [85] bit_4.5.0 tidyselect_1.2.1 locfit_1.5-9.10
## [88] maketools_1.3.1 Biostrings_2.75.1 knitr_1.49
## [91] gridExtra_2.3 xfun_0.49 graphlayouts_1.2.1
## [94] statmod_1.5.0 stringi_1.8.4 UCSC.utils_1.3.0
## [97] yaml_2.3.10 evaluate_1.0.1 codetools_0.2-20
## [100] ggraph_2.2.1 tibble_3.2.1 BiocManager_1.30.25
## [103] cli_3.6.3 uwot_0.2.2 munsell_0.5.1
## [106] jquerylib_0.1.4 Rcpp_1.0.13-1 dbplyr_2.5.0
## [109] png_0.1-8 parallel_4.4.2 blob_1.2.4
## [112] viridisLite_0.4.2 scales_1.3.0 purrr_1.0.2
## [115] crayon_1.5.3 rlang_1.1.4 cowplot_1.1.3
## [118] KEGGREST_1.47.0