Simplify Functional Enrichment Results

Zuguang Gu (z.gu@dkfz.de)

2023-11-06

The simplifyEnrichment package clusters functional terms into groups by clustering the similarity matrix of the terms with a new proposed method “binary cut” which recursively applies partition around medoids (PAM) with two groups on the similarity matrix and in each iteration step, a score is assigned to decide whether the group of gene sets that corresponds to the current sub-matrix should be split or not. For more details of the method, please refer to the simplifyEnrichment paper.

Simplify GO enrichment results

The major use case for simplifyEnrichment is for simplying the GO enrichment results by clustering the corresponding semantic similarity matrix of the significant GO terms. To demonstrate the usage, we first generate a list of random GO IDs from the Biological Process (BP) ontology category:

library(simplifyEnrichment)
set.seed(888)
go_id = random_GO(500)

simplifyEnrichment starts with the GO similarity matrix. Users can use their own similarity matrices or use the GO_similarity() function to calculate the semantic similarity matrix. The GO_similarity() function is simply a wrapper on GOSemSim::termSim(). The function accepts a vector of GO IDs. Note the GO terms should only belong to one same ontology (i.e., BP, CC or MF).

mat = GO_similarity(go_id)

By default, GO_similarity() uses Rel method in GOSemSim::termSim(). Other methods to calculate GO similarities can be set by measure argument, e.g.:

GO_similarity(go_id, measure = "Wang")

With the similarity matrix mat, users can directly apply simplifyGO() function to perform the clustering as well as visualizing the results.

df = simplifyGO(mat)

On the right side of the heatmap there are the word cloud annotations which summarize the functions with keywords in every GO cluster. Additionally, enrichment is done on keywords compared to GO background vocabulary and the significance corresponds to the font size of the keywords.

Note there is no word cloud for the cluster that is merged from small clusters (size < 5).

The returned variable df is a data frame with GO IDs and the cluster labels:

head(df)
##           id cluster
## 1 GO:0086066       1
## 2 GO:0090461       2
## 3 GO:0032912       3
## 4 GO:0090220       4
## 5 GO:0032495       5
## 6 GO:0070585       4

The size of GO clusters can be retrieved by:

sort(table(df$cluster))
## 
##  15  16  17  18  19   8   9  14  10  12  13   6  11   2   7   4   5   1   3 
##   1   1   1   1   1   2   2   3   5   5   5   8   8  13  36  69  84 120 135

Or split the data frame by the cluster labels:

split(df, df$cluster)

plot argument can be set to FALSE in simplifyGO(), so that no plot is generated and only the data frame is returned.

If the aim is only to cluster GO terms, binary_cut() or cluster_terms() functions can be directly applied:

binary_cut(mat)
##   [1]  1  2  3  4  5  4  4  4  1  3  3  5  1  5  1  4  1  4  5  4  5  1  4  3  6
##  [26]  3  7  4  3  1  1  3  8  7  3  3  5  1  4  4  5  2  4  9  5  1  1  7  3  7
##  [51]  5  3 10  7  1  7  3 10  3  3  7  1  5  5  3  3  1  3  1  3  3  4 11  3  1
##  [76]  3  4  7  3  3  1  5  4  6  3  3  1  1  4  5  7  7  3  4  5  3  1  1  4  6
## [101]  4  3  4  4  1  4  3  5  3  7  3  1  3  5  3  1  1  1  9  4  4 12  5  1  1
## [126]  1  1  5  4  1  4  1  3  5  5  5  1  5  3  5  5  5  1  3  5  5 13  1  2  4
## [151] 12  3  7  1  3  5  8  1  4  5  1  5  1  1  5  3  3  4  1  1  3  3  3  1  4
## [176]  4  1  3  4  3  3  5  5  4  4  3  1  3  5  3  5  3  3  3  1  1  1  5  4  3
## [201]  2  5  4  4  2  3  1  1  3  3  3  2  3  3  3  4  5  3  4  3  4  6  6  4  7
## [226]  4  5  5  1  1  3  1 14 12  7  7  5  5  3  3  7  3  1  4  1  1  5 14  3  5
## [251]  5  1 15  1  3 11  4  1  5 13  3  1  7  3  1  5  1  3  5  5  3  6  1  5  3
## [276] 11  1  2  7  4  5 13  1  1 10  3  5  3  3  3  7  4  5  3  3  4  1  1  4  7
## [301]  3  2  4  5  3  3  5  3  1  3  2  5  2  3  3  1  1  5  5  4 11  1  1  5  3
## [326]  3  1  4  5  1  1  1  1  7  3  2  3  2  5  3  5  1 13  3  3  4  7  1  4  6
## [351]  4  7  1  3  1  4  7  3  5  3  1  2 11  4  3  3  1  1  1  5  1  1  1  1  3
## [376]  5  4  7  3 10  1  3  1  7  3  1  3 11  3 16  4  6  3  1  7  7 14  1  5  4
## [401]  5  3  1  3  3  5  3  3  7  4  5  1  7  5  7  3  1 12  3  7 10  1  1  4  3
## [426]  1  3 13  1  3  3  1  1  5  3  3  1  1  1  1 17  5  1  3  7  3  1  5  1  5
## [451]  7  3  4  5  4  4  1  1  3  5  7 11  4  4  5  3  4  4  1  1 18  5 12  5  3
## [476]  3  3  3  5  5  3  3  3  4  1 11  5  1  3  3  1  5  3  1  3 19  7  5  3  1

or

cluster_terms(mat, method = "binary_cut")

binary_cut() and cluster_terms() basically generate the same clusterings, but the labels of clusters might differ.

Simplify general functional enrichment results

Semantic measures can be used for the similarity of GO terms. However, there are still a lot of ontologies (e.g. MsigDB gene sets) that are only represented as a list of genes where the similarity between gene sets are mainly measured by gene overlap. simplifyEnrichment provides the term_similarity() and other related functions (term_similarity_from_enrichResult(), term_similarity_from_KEGG(), term_similarity_from_Reactome(), term_similarity_from_MSigDB() and term_similarity_from_gmt()) which calculate the similarity of terms by the gene overlapping, with methods of Jaccard coefficient, Dice coefficient, overlap coefficient and kappa coefficient.

The similarity can be calculated by providing:

  1. A list of gene sets where each gene set contains a vector of genes.
  2. A enrichResult object which is normally from the ‘clusterProfiler’, ‘DOSE’, ‘meshes’ or ‘ReactomePA’ package.
  3. A list of KEGG/Reactome/MsigDB IDs. The gene set names can also be provided for MsigDB ontologies.
  4. A gmt file and the corresponding gene set IDs.

Once you have the similarity matrix, you can send it to simplifyEnrichment() function. But note, as we benchmarked in the manuscript, the clustering on the gene overlap similarity performs much worse than on the semantic similarity.

Comparing clustering methods

In the simplifyEnrichment package, there are also functions that compare clustering results from different methods. Here we still use previously generated variable mat which is the similarity matrix from the 500 random GO terms. Simply running compare_clustering_methods() function performs all supported methods (in all_clustering_methods()) excluding mclust, because mclust usually takes very long time to run. The function generates a figure with three panels:

  1. A heatmap of the similarity matrix with different clusterings as row annotations.
  2. A heatmap of the pair-wise concordance of the clustering from every two methods.
  3. Barplots of the difference scores for each method, the number of clusters (total clusters and the clusters with size >= 5) and the mean similarity of the terms that are in the same clusters (block mean).

In the barplots, the three metrics are defined as follows:

  1. Different score: This is the difference between the similarity values for the terms that belong to the same clusters and different clusters. For a similarity matrix \(M\), for term \(i\) and term \(j\) where \(i \ne j\), the similarity value \(x_{i,j}\) is saved to the vector \(\mathbf{x_1}\) only when term \(i\) and \(j\) are in a same cluster. \(x_{i,j}\) is saved to the vector \(\mathbf{x_2}\) when term \(i\) and \(j\) are not in the same cluster. The difference score measures the distribution difference between \(\mathbf{x_1}\) and \(\mathbf{x_2}\), calculated as the Kolmogorov-Smirnov statistic between the two distributions.
  2. Number of clusters: For each clustering, there are two numbers: the number of total clusters and the number of clusters with size >= 5 (only the big clusters).
  3. Block mean: Mean similarity values of the diagonal blocks in the similarity heatmap. Using the same convention as for the difference score, the block mean is the mean value of \(\mathbf{x_1}\).
compare_clustering_methods(mat)

If plot_type argument is set to heatmap. There are heatmaps for the similarity matrix under different clusterings methods. The last panel is a table with the number of clusters.

compare_clustering_methods(mat, plot_type = "heatmap")

Please note, the clustering methods might have randomness, which means, different runs of compare_clustering_methods() may generate different clusterings (slightly different). Thus, if users want to compare the plots between compare_clustering_methods(mat) and compare_clustering_methods(mat, plot_type = "heatmap"), they should set the same random seed before executing the function.

set.seed(123)
compare_clustering_methods(mat)
set.seed(123)
compare_clustering_methods(mat, plot_type = "heatmap")

compare_clustering_methods() is simply a wrapper on cmp_make_clusters() and cmp_make_plot() functions where the former function performs clustering with different methods and the latter visualizes the results. To compare different plots, users can also use the following code without specifying the random seed.

clt = cmp_make_clusters(mat) # just a list of cluster labels
cmp_make_plot(mat, clt)
cmp_make_plot(mat, clt, plot_type = "heatmap")

Register new clustering methods

New clustering methods can be added by register_clustering_methods(), removed by remove_clustering_methods() and reset to the default methods by reset_clustering_methods(). All the supported methods can be retrieved by all_clustering_methods(). compare_clustering_methods() runs all the clustering methods in all_clustering_methods().

The new clustering methods should be as user-defined functions and sent to register_clustering_methods() as named arguments, e.g.:

register_clustering_methods(
    method1 = function(mat, ...) ...,
    method2 = function(mat, ...) ...,
    ...
)

The functions should accept at least one argument which is the input matrix (mat in above example). The second optional argument should always be ... so that parameters for the clustering function can be passed by control argument from cluster_terms() or simplifyGO(). If users forget to add ..., it is added internally.

Please note, the user-defined function should automatically identify the optimized number of clusters. The function should return a vector of cluster labels. Internally it is converted to numeric labels.

Examples

There are following examples which we did for the benchmarking in the manuscript:

Apply to multiple lists of GO IDs

It is always very common that users have multiple lists of GO enrichment results (e.g. from multiple groups of genes) and they want to compare the significant terms between different lists, e.g. to see which biological functions are more specific in a certain list. There is a function simplifyGOFromMultipleLists() in the package which helps this type of analysis.

The input data for simplifyGOFromMultipleLists() (with the argument lt) can have three types of formats:

If the GO enrichment results is directly from upstream analysis, e.g. the package clusterProfiler or other similar packages, the results are most probably represented as a list of data frames, thus, we first demonstrate the usage on a list of data frames.

The function functional_enrichment() in cola package applies functional enrichment on different groups of signature genes from consensus clustering. The function internally uses clusterProfiler and returns a list of data frames:

# perform functional enrichment on the signatures genes from cola anlaysis 
library(cola)
data(golub_cola) 
res = golub_cola["ATC:skmeans"]

library(hu6800.db)
x = hu6800ENTREZID
mapped_probes = mappedkeys(x)
id_mapping = unlist(as.list(x[mapped_probes]))

lt = functional_enrichment(res, k = 3, id_mapping = id_mapping)
## - 2058/4116 significant genes are taken from 3-group comparisons
## - on k-means group 1/4, 531 genes
##   - 478/531 (90%) genes left after id mapping
##   - gene set enrichment, GO:BP
## - on k-means group 2/4, 811 genes
##   - 640/811 (78.9%) genes left after id mapping
##   - gene set enrichment, GO:BP
## - on k-means group 3/4, 315 genes
##   - 276/315 (87.6%) genes left after id mapping
##   - gene set enrichment, GO:BP
## - on k-means group 4/4, 401 genes
##   - 374/401 (93.3%) genes left after id mapping
##   - gene set enrichment, GO:BP
names(lt)
## [1] "BP_km1" "BP_km2" "BP_km3" "BP_km4"
head(lt[[1]][, 1:7])
##                    ID                Description GeneRatio   BgRatio
## GO:0006974 GO:0006974        DNA damage response    79/471 884/18870
## GO:0006281 GO:0006281                 DNA repair    63/471 601/18870
## GO:1903047 GO:1903047 mitotic cell cycle process    69/471 774/18870
## GO:0000278 GO:0000278         mitotic cell cycle    74/471 928/18870
## GO:0006260 GO:0006260            DNA replication    38/471 278/18870
## GO:0051276 GO:0051276    chromosome organization    57/471 626/18870
##                  pvalue     p.adjust       qvalue
## GO:0006974 1.875493e-23 8.807316e-20 6.680704e-20
## GO:0006281 2.435061e-22 5.717522e-19 4.336971e-19
## GO:1903047 1.934970e-20 3.028872e-17 2.297522e-17
## GO:0000278 4.277686e-19 5.022004e-16 3.809392e-16
## GO:0006260 1.196027e-17 1.123308e-14 8.520746e-15
## GO:0051276 2.126220e-17 1.664121e-14 1.262303e-14

By default, simplifyGOFromMultipleLists() automatically identifies the columns that contain GO IDs and adjusted p-values, so here we directly send lt to simplifyGOFromMultipleLists(). We additionally set padj_cutoff to 0.001 because under the default cutoff 0.01, there are too many GO IDs and to save the running time, we set a more strict cutoff.

simplifyGOFromMultipleLists(lt, padj_cutoff = 0.001)

Next we demonstrate two other data types for simplifyGOFromMultipleLists(). Both usages are straightforward. The first is a list of numeric vectors:

lt2 = lapply(lt, function(x) structure(x$p.adjust, names = x$ID))
simplifyGOFromMultipleLists(lt2, padj_cutoff = 0.001)

And the second is a list of character vectors of GO IDs:

lt3 = lapply(lt, function(x) x$ID[x$p.adjust < 0.001])
simplifyGOFromMultipleLists(lt3)

The process of this analysis is as follows. Let’s assume there are \(n\) GO lists, we first construct a global matrix where columns correspond to the \(n\) GO lists and rows correspond to the “union” of all GO IDs in the \(n\) lists. The value for the ith GO ID and in the jth list are taken from the corresponding numeric vector in lt. If the jth vector in lt does not contain the ith GO ID, the value defined by default argument is taken there (e.g. in most cases the numeric values are adjusted p-values, thus default is set to 1). Let’s call this matrix as \(M_0\).

Next step is to filter \(M_0\) so that we only take a subset of GO IDs of interest. We define a proper function via argument filter to remove GO IDs that are not important for the analysis. Function for filter is applied to every row in \(M_0\) and filter function needs to return a logical value to decide whether to keep or remove the current GO ID. For example, if the values in lt are adjusted p-values, the filter function can be set as function(x) any(x < padj_cutoff) so that the GO ID is kept as long as it is signfiicant in at least one list. After the filtering, let’s call the filtered matrix \(M_1\).

GO IDs in \(M_1\) (row names of \(M_1\)) are used for clustering. A heatmap of \(M_1\) is attached to the left of the GO similarity heatmap so that the group-specific (or list-specific) patterns can be easily observed and to corresponded to GO functions.

Argument heatmap_param controls several parameters for heatmap \(M_1\):

Session Info

sessionInfo()
## R Under development (unstable) (2023-10-22 r85388)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
## 
## Random number generation:
##  RNG:     L'Ecuyer-CMRG 
##  Normal:  Inversion 
##  Sample:  Rejection 
##  
## 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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    grid      stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] hu6800.db_3.13.0          org.Hs.eg.db_3.18.0      
##  [3] AnnotationDbi_1.65.2      IRanges_2.37.0           
##  [5] S4Vectors_0.41.1          Biobase_2.63.0           
##  [7] cola_2.9.0                simplifyEnrichment_1.13.1
##  [9] BiocGenerics_0.49.1       knitr_1.45               
## 
## loaded via a namespace (and not attached):
##   [1] splines_4.4.0                 later_1.3.1                  
##   [3] ggplotify_0.1.2               bitops_1.0-7                 
##   [5] filelock_1.0.2                tibble_3.2.1                 
##   [7] polyclip_1.10-6               XML_3.99-0.15                
##   [9] lifecycle_1.0.3               doParallel_1.0.17            
##  [11] NLP_0.2-1                     lattice_0.22-5               
##  [13] prabclus_2.3-3                MASS_7.3-60.1                
##  [15] magrittr_2.0.3                sass_0.4.7                   
##  [17] rmarkdown_2.25                jquerylib_0.1.4              
##  [19] yaml_2.3.7                    httpuv_1.6.12                
##  [21] doRNG_1.8.6                   flexmix_2.3-19               
##  [23] cowplot_1.1.1                 DBI_1.1.3                    
##  [25] RColorBrewer_1.1-3            eulerr_7.0.0                 
##  [27] zlibbioc_1.49.0               expm_0.999-7                 
##  [29] purrr_1.0.2                   ggraph_2.1.0                 
##  [31] RCurl_1.98-1.13               yulab.utils_0.1.0            
##  [33] nnet_7.3-19                   tweenr_2.0.2                 
##  [35] rappdirs_0.3.3                circlize_0.4.15              
##  [37] GenomeInfoDbData_1.2.11       enrichplot_1.23.0            
##  [39] ggrepel_0.9.4                 tm_0.7-11                    
##  [41] irlba_2.3.5.1                 tidytree_0.4.5               
##  [43] genefilter_1.85.0             annotate_1.81.0              
##  [45] brew_1.0-8                    commonmark_1.9.0             
##  [47] codetools_0.2-19              ggforce_0.4.1                
##  [49] DOSE_3.29.0                   xml2_1.3.5                   
##  [51] tidyselect_1.2.0              shape_1.4.6                  
##  [53] aplot_0.2.2                   farver_2.1.1                 
##  [55] viridis_0.6.4                 matrixStats_1.0.0            
##  [57] BiocFileCache_2.11.1          dynamicTreeCut_1.63-1        
##  [59] jsonlite_1.8.7                GetoptLong_1.0.5             
##  [61] tidygraph_1.2.3               ellipsis_0.3.2               
##  [63] survival_3.5-7                iterators_1.0.14             
##  [65] foreach_1.5.2                 dbscan_1.1-11                
##  [67] tools_4.4.0                   treeio_1.27.0                
##  [69] HPO.db_0.99.2                 Rcpp_1.0.11                  
##  [71] glue_1.6.2                    gridExtra_2.3                
##  [73] xfun_0.41                     qvalue_2.35.0                
##  [75] MatrixGenerics_1.15.0         GenomeInfoDb_1.39.0          
##  [77] dplyr_1.1.3                   withr_2.5.2                  
##  [79] BiocManager_1.30.22           fastmap_1.1.1                
##  [81] fansi_1.0.5                   digest_0.6.33                
##  [83] gridGraphics_0.5-1            R6_2.5.1                     
##  [85] mime_0.12                     microbenchmark_1.4.10        
##  [87] colorspace_2.1-0              GO.db_3.18.0                 
##  [89] Cairo_1.6-1                   markdown_1.11                
##  [91] RSQLite_2.3.3                 diptest_0.76-0               
##  [93] tidyr_1.3.0                   utf8_1.2.4                   
##  [95] generics_0.1.3                data.table_1.14.8            
##  [97] robustbase_0.99-0             class_7.3-22                 
##  [99] graphlayouts_1.0.2            httr_1.4.7                   
## [101] scatterpie_0.2.1              MCL_1.0                      
## [103] pkgconfig_2.0.3               gtable_0.3.4                 
## [105] modeltools_0.2-23             blob_1.2.4                   
## [107] ComplexHeatmap_2.19.0         impute_1.77.0                
## [109] XVector_0.43.0                shadowtext_0.1.2             
## [111] clusterProfiler_4.11.0        htmltools_0.5.7              
## [113] fgsea_1.29.0                  clue_0.3-65                  
## [115] scales_1.2.1                  png_0.1-8                    
## [117] ggfun_0.1.3                   reshape2_1.4.4               
## [119] rjson_0.2.21                  nlme_3.1-163                 
## [121] curl_5.1.0                    cachem_1.0.8                 
## [123] GlobalOptions_0.1.2           stringr_1.5.0                
## [125] BiocVersion_3.19.1            parallel_4.4.0               
## [127] HDO.db_0.99.1                 pillar_1.9.0                 
## [129] proxyC_0.3.4                  vctrs_0.6.4                  
## [131] slam_0.1-50                   promises_1.2.1               
## [133] dbplyr_2.4.0                  xtable_1.8-4                 
## [135] cluster_2.1.4                 evaluate_0.23                
## [137] magick_2.8.1                  cli_3.6.1                    
## [139] compiler_4.4.0                rlang_1.1.2                  
## [141] crayon_1.5.2                  rngtools_1.5.2               
## [143] labeling_0.4.3                mclust_6.0.0                 
## [145] fs_1.6.3                      skmeans_0.2-16               
## [147] plyr_1.8.9                    stringi_1.7.12               
## [149] viridisLite_0.4.2             BiocParallel_1.37.0          
## [151] MPO.db_0.99.7                 munsell_0.5.0                
## [153] Biostrings_2.71.1             lazyeval_0.2.2               
## [155] GOSemSim_2.29.0               Matrix_1.6-1.1               
## [157] patchwork_1.1.3               bit64_4.0.5                  
## [159] ggplot2_3.4.4                 KEGGREST_1.43.0              
## [161] fpc_2.2-10                    shiny_1.7.5.1                
## [163] interactiveDisplayBase_1.41.0 highr_0.10                   
## [165] apcluster_1.4.11              AnnotationHub_3.11.0         
## [167] kernlab_0.9-32                gridtext_0.1.5               
## [169] igraph_1.5.1                  memoise_2.0.1                
## [171] RcppParallel_5.1.7            bslib_0.5.1                  
## [173] ggtree_3.11.0                 fastmatch_1.1-4              
## [175] DEoptimR_1.1-3                bit_4.0.5                    
## [177] gson_0.1.0                    ape_5.7-1