/
evaluation.Rmd
146 lines (105 loc) · 6.17 KB
/
evaluation.Rmd
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
---
title: "Evaluating integration results by CIDER"
output:
rmarkdown::html_vignette:
toc: TRUE
number_sections: true
vignette: >
%\VignetteIndexEntry{Evaluating integration results by CIDER}
%\VignetteEngine{knitr::rmarkdown}
%\VignetteEncoding{UTF-8}
---
```{r, include = FALSE}
knitr::opts_chunk$set(
collapse = TRUE,
comment = "#>", warning = FALSE, message = FALSE
)
```
# Introduction
Integration and batch correction methods have become a **popular** component in the bioinformatic workflows for scRNA-Seq data analysis, whilst the integration results (mostly corrected PCs or less commonly corrected read counts) are **rarely** validated or evaluated with an objective metric.
To assess the correctness of integration (i.e. whether cells belonging to the same population are gathered and ones belonging to different populations stay separate after integration), the existing evaluation metrics require the existence of ground truth for cell population annotations.
CIDER provides a **ground-truth-free** approach to evaluate the integration results. This vignette focuses how showing the process using the example data of dendritic cells.
# Set up
Apart from **CIDER**, the following packages also need to be loaded:
```{r setup}
library(CIDER)
library(Seurat)
library(cowplot)
library(ggplot2)
```
# Load dendritic data
The example data can be downloaded from https://figshare.com/s/d5474749ca8c711cc205. This dataset contains 26593 genes and 564 cells from two batches.
```{r}
load("../data/dendritic.rda")
dim(dendritic)
```
```{r}
table(dendritic$Batch)
```
# Perform integration
First an integration method$^1$ is applied on the dendritic data. You can apply other integration methods to the your data, as long as the correct PCs are stored in your Seurat object, i.e. `Reductions(seu.integrated, "pca")` or `seu.integrated@reductions$pca`.
```{r integration}
seu.list <- SplitObject(dendritic, split.by = "Batch")
for (i in 1:length(seu.list)) {
seu.list[[i]] <- NormalizeData(seu.list[[i]], verbose = FALSE)
seu.list[[i]] <- FindVariableFeatures(seu.list[[i]],
selection.method = "vst",
nfeatures = 1000, verbose = FALSE)
}
seu.anchors <- FindIntegrationAnchors(object.list = seu.list,
dims = 1:15, verbose = FALSE)
seu.integrated <- IntegrateData(anchorset = seu.anchors,
dims = 1:15, verbose = FALSE)
DefaultAssay(seu.integrated) <- "integrated"
seu.integrated <- ScaleData(seu.integrated, verbose = FALSE)
seu.integrated <- RunPCA(seu.integrated, verbose = FALSE)
seu.integrated <- RunTSNE(seu.integrated, reduction = "pca", dims = 1:5)
seu.integrated@reductions$pca
```
Clear the intermediate outcome.
```{r}
rm(seu.list, seu.anchors)
gc()
```
# Evaluate by CIDER
## Calculate similarity and p values (essential)
CIDER evaluates integration results in three steps:
1. Clustering based on the corrected PCs (`hdbscan.seurat`). This step uses HDBSCAN, which is a density-based clustering algorithm$^2$. The clustering results are stored in `seu.integrated$dbscan_cluster`. Clusters are further divided into batch-specific clusters by concatenating dbscan_cluster and batch, stored in `seu.integrated$initial_cluster`.
1. Compute IDER-based similarity matrix (`getIDEr`) among the batch-specific initial clusters. If multiple CPUs are availble, you can set `use.parallel = TRUE` and `n.cores` to the number of available cores to speed it up.
1. Assign the similarity and estimate empirical p values (`estimateProb`) for the correctness of integration. High similarity values and low p values indicate that the cell are similar to the surrounding cells and likely integrated correctly.
```{r}
seu.integrated <- hdbscan.seurat(seu.integrated)
ider <- getIDEr(seu.integrated, use.parallel = FALSE, verbose = FALSE)
seu.integrated <- estimateProb(seu.integrated, ider)
```
## Evaluation scores
The evaluation scores can be viewed by the `scatterPlot` as below. As shown cells with dbscan_cluster of 2 and 3 have low regional similarity and high empirical p values, suggesting that they can be incorrectly integrated.
```{r, fig.height=1.8, fig.width=7}
p1 <- scatterPlot(seu.integrated, "tsne", "dbscan_cluster")
p2 <- scatterPlot(seu.integrated, "tsne", colour.by = "similarity") + labs(fill = "Similarity")
p3 <- scatterPlot(seu.integrated, "tsne", colour.by = "pvalue") + labs(fill = "Prob of \nrejection")
plot_grid(p1,p2,p3, ncol = 3)
```
## The IDER-based similarity matrix
To have more insight, we can view the IDER-based similarity matrix by functions `plotNetwork` or `plotHeatmap`. Both of them require the input of a Seurat object and the output of `getIDEr`. In this example, 1_Batch1 and 1_Batch2 as well as 4_Batch1 and 4_Batch2 have high similarity.
`plotNetwork` generates a graph where vertexes are initial clusters and edge widths are similarity values. The parameter `weight.factor` controls the scale of edge widths; larger `weight.factor` will give bolder edges proportionally.
```{r, fig.height=5, fig.width=5}
plotNetwork(seu.integrated, ider, weight.factor = 3)
```
`plotHeatmap` generates a heatmap where each cell is coloured and labeled by the similarity values.
```{r, fig.height=5, fig.width=5}
plotHeatmap(seu.integrated, ider)
```
# Refering to ground-truth annotation
So far the evaluation have completed and CIDER has not used the ground truth at all!
Let's peep at the ground truth before the closure of this vignette. As shown in the figure below, the clusters having low IDER-based similarity and high p values actually have at least two popuplations (CD1C and CD141), verifying that CIDER spots the wrongly integrated cells.
```{r, fig.height=3, fig.width=5}
scatterPlot(seu.integrated, "tsne", colour.by = "Group") + labs(fill = "Group\n (ground truth)")
```
# Technical
```{r sessionInfo}
sessionInfo()
```
# References
1. Stuart and Butler et al. Comprehensive Integration of Single-Cell Data. Cell (2019).
2. Campello, Ricardo JGB, Davoud Moulavi, and Jörg Sander. “Density-based clustering based on hierarchical density estimates.” Pacific-Asia conference on knowledge discovery and data mining. Springer, Berlin, Heidelberg, 2013.