Evaluating single-cell clustering with ClustAssess

In this vignette we will illustrate several of the tools available in ClustAssess on a small single-cell RNA-seq dataset.

library(Seurat)
#> Attaching SeuratObject
library(ClustAssess)

# we will use the pbmc_small single-cell RNA-seq dataset available via Seurat
DimPlot(pbmc_small, group.by='letter.idents')

Proportion of ambiguously clustered pairs

The proportion of ambiguously clustered pairs (PAC) uses consensus clustering to infer the optimal number of clusters for the data, by observing how variably pairs of elements cluster together. The lower the PAC, the more stable the clustering. PAC was originally presented in https://doi.org/10.1038/srep06207.

# retrieve scaled data for PAC calculation
pbmc.data = GetAssayData(pbmc_small, assay='RNA', slot='scale.data')

# perform consensus clustering
cc.res = consensus_cluster(t(pbmc.data), 
                           k_max=30, 
                           n_reps=100,
                           p_sample=0.8,
                           p_feature=0.8,
                           verbose=TRUE)
#> Calculating consensus clustering

# assess the PAC convergence for a few values of k - each curve should 
# have converged to some value
k.plot = c(4,6,8,10)
pac_convergence(cc.res, k.plot)


# visualize the final PAC across k - there seems to be a local maximum at k=7, 
# indicating that 7 clusters leads to a less stable clustering of the data than 
# nearby values of k
pac_landscape(cc.res)

Element-centric clustering comparison

We compare the similarity of clustering results between methods (in this case, between Louvain community detection and k-means) using element-centric similarity, which quantifies the clustering similarity per cell. Higher ECS indicates that an observation is clustered similarly across methods. ECS was introduced in https://doi.org/10.1038/s41598-019-44892-y.

# first, cluster with Louvain algorithm
pbmc_small = FindClusters(pbmc_small, resolution=0.8, verbose=FALSE)
DimPlot(pbmc_small, group.by='seurat_clusters')


# also cluster with PCA+k-means
pbmc_pca = Embeddings(pbmc_small, 'pca')
pbmc_small@meta.data$kmeans_clusters = kmeans(pbmc_pca, 
                                              centers=3, 
                                              nstart=10, 
                                              iter.max=1e3)$cluster
DimPlot(pbmc_small, group.by='kmeans_clusters')


# where are the clustering results more similar?
pbmc_small@meta.data$ecs = element_sim_elscore(pbmc_small@meta.data$seurat_clusters, 
                                               pbmc_small@meta.data$kmeans_clusters)
FeaturePlot(pbmc_small, 'ecs')

mean(pbmc_small@meta.data$ecs)
#> [1] 0.6632927

Jaccard similarity of cluster markers

As a common step in computational single-cell RNA-seq analyses, discriminative marker genes are identified for each cluster. These markers are then used to infer the cell type of the respective cluster. Here, we compare the markers obtained for each clustering method to ask: how similarly would each cell be interpreted across clustering methods? We compare the markers per cell using the Jaccard similarity (defined as the size of the intersect divided by the size of the overlap) of the marker gene lists. The higher the JSI, the more similar are the marker genes for that cell.

# first, we calculate the markers on the Louvain clustering
Idents(pbmc_small) = pbmc_small@meta.data$seurat_clusters
louvain.markers = FindAllMarkers(pbmc_small, 
                                 logfc.threshold=1, 
                                 min.pct=0.0, 
                                 test.use='roc', 
                                 verbose=FALSE)

# then we get the markers on the k-means clustering
Idents(pbmc_small) = pbmc_small@meta.data$kmeans_clusters
kmeans.markers = FindAllMarkers(pbmc_small, 
                                logfc.threshold=1, 
                                min.pct=0.0, 
                                test.use='roc', 
                                verbose=FALSE)

# next, compare the top 10 markers per cell
pbmc_small@meta.data$marker.gene.jsi = marker_overlap(louvain.markers, 
                                          kmeans.markers, 
                                          pbmc_small@meta.data$seurat_clusters, 
                                          pbmc_small@meta.data$kmeans_clusters, 
                                          n=10, 
                                          rank_by='power')
#> Warning: Unknown or uninitialised column: `avg_logFC`.

#> Warning: Unknown or uninitialised column: `avg_logFC`.

# which cells have the same markers, regardless of clustering?
FeaturePlot(pbmc_small, 'marker.gene.jsi')

mean(pbmc_small@meta.data$marker.gene.jsi)
#> [1] 0.6352273

Element-wise consistency

How consistently are cells clustered by k-means? We will rerun k-means clustering 20 times to investigate.

clustering.list = list()
n.reps = 20
for (i in 1:n.reps){
  # we set nstart=1 and a fairly high iter.max - this should mean that
  # the algorithm converges, and that the variability in final clusterings
  # depends mainly on the random initial cluster assignments
  km.res = kmeans(pbmc_pca, centers=3, nstart=1, iter.max=1e3)$cluster
  clustering.list[[i]] = km.res
}

# now, we calculate the element-wise consistency (aka frustration), which
# performs pair-wise comparisons between all 20 clusterings; the 
# consistency is the average per-cell ECS across all comparisons. The higher
# the consistency, the more consistently is that cell clustered across
# random seeds.
pbmc_small@meta.data$consistency = element_consistency(clustering.list)

# which cells are clustered more consistently?
FeaturePlot(pbmc_small, 'consistency')

mean(pbmc_small@meta.data$consistency)
#> [1] 0.6712922

Session info

sessionInfo()
#> R version 4.0.3 (2020-10-10)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Debian GNU/Linux 10 (buster)
#> 
#> Matrix products: default
#> BLAS:   /usr/lib/x86_64-linux-gnu/blas/libblas.so.3.8.0
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.8.0
#> 
#> locale:
#>  [1] LC_CTYPE=en_GB.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_GB.UTF-8        LC_COLLATE=C              
#>  [5] LC_MONETARY=en_GB.UTF-8    LC_MESSAGES=en_GB.UTF-8   
#>  [7] LC_PAPER=en_GB.UTF-8       LC_NAME=C                 
#>  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
#> [11] LC_MEASUREMENT=en_GB.UTF-8 LC_IDENTIFICATION=C       
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] ClustAssess_0.3.0  SeuratObject_4.0.4 Seurat_4.1.0      
#> 
#> loaded via a namespace (and not attached):
#>   [1] Rtsne_0.15            colorspace_2.0-2      deldir_1.0-6         
#>   [4] ellipsis_0.3.2        ggridges_0.5.3        spatstat.data_2.1-2  
#>   [7] farver_2.1.0          leiden_0.3.9          listenv_0.8.0        
#>  [10] ggrepel_0.9.1         fansi_1.0.0           codetools_0.2-18     
#>  [13] splines_4.0.3         knitr_1.37            polyclip_1.10-0      
#>  [16] jsonlite_1.7.2        ica_1.0-2             cluster_2.1.2        
#>  [19] png_0.1-7             uwot_0.1.11           shiny_1.7.1          
#>  [22] sctransform_0.3.3     spatstat.sparse_2.1-0 compiler_4.0.3       
#>  [25] httr_1.4.2            assertthat_0.2.1      Matrix_1.4-0         
#>  [28] fastmap_1.1.0         lazyeval_0.2.2        later_1.3.0          
#>  [31] prettyunits_1.1.1     htmltools_0.5.2       tools_4.0.3          
#>  [34] igraph_1.2.11         gtable_0.3.0          glue_1.6.0           
#>  [37] RANN_2.6.1            reshape2_1.4.4        dplyr_1.0.7          
#>  [40] Rcpp_1.0.8            scattermore_0.7       jquerylib_0.1.4      
#>  [43] vctrs_0.3.8           nlme_3.1-155          iterators_1.0.13     
#>  [46] lmtest_0.9-39         fastcluster_1.2.3     xfun_0.29            
#>  [49] stringr_1.4.0         globals_0.14.0        mime_0.12            
#>  [52] miniUI_0.1.1.1        lifecycle_1.0.1       irlba_2.3.5          
#>  [55] goftest_1.2-3         future_1.23.0         MASS_7.3-55          
#>  [58] zoo_1.8-9             scales_1.1.1          spatstat.core_2.3-2  
#>  [61] hms_1.1.1             promises_1.2.0.1      spatstat.utils_2.3-0 
#>  [64] parallel_4.0.3        RColorBrewer_1.1-2    yaml_2.2.1           
#>  [67] reticulate_1.22       pbapply_1.5-0         gridExtra_2.3        
#>  [70] ggplot2_3.3.5         sass_0.4.0            rpart_4.1-15         
#>  [73] stringi_1.7.6         highr_0.9             foreach_1.5.1        
#>  [76] rlang_0.4.12          pkgconfig_2.0.3       matrixStats_0.61.0   
#>  [79] evaluate_0.14         lattice_0.20-45       ROCR_1.0-11          
#>  [82] purrr_0.3.4           tensor_1.5            labeling_0.4.2       
#>  [85] patchwork_1.1.1       htmlwidgets_1.5.4     cowplot_1.1.1        
#>  [88] tidyselect_1.1.1      parallelly_1.30.0     RcppAnnoy_0.0.19     
#>  [91] plyr_1.8.6            magrittr_2.0.1        R6_2.5.1             
#>  [94] generics_0.1.1        DBI_1.1.2             pillar_1.6.4         
#>  [97] mgcv_1.8-38           fitdistrplus_1.1-6    survival_3.2-13      
#> [100] abind_1.4-5           tibble_3.1.6          future.apply_1.8.1   
#> [103] crayon_1.4.2          KernSmooth_2.23-20    utf8_1.2.2           
#> [106] spatstat.geom_2.3-1   plotly_4.10.0         rmarkdown_2.11       
#> [109] progress_1.2.2        grid_4.0.3            data.table_1.14.2    
#> [112] digest_0.6.29         xtable_1.8-4          tidyr_1.1.4          
#> [115] httpuv_1.6.5          munsell_0.5.0         viridisLite_0.4.0    
#> [118] bslib_0.3.1