1 vissE
This package implements the vissE algorithm to summarise results of gene-set analyses. Usually, the results of a gene-set enrichment analysis (e.g using limma::fry, singscore or GSEA) consist of a long list of gene-sets. Biologists then have to search through these lists to determines emerging themes to explain the altered biological processes. This task can be labour intensive therefore we need solutions to summarise large sets of results from such analyses.
This package provides an approach to provide summaries of results from gene-set enrichment analyses. It exploits the relatedness between gene-sets and the inherent hierarchical structure that may exist in pathway databases and gene ontologies to cluster results. For each cluster of gene-sets vissE identifies, it performs text-mining to automate characterisation of biological functions and processes represented by the cluster.
An additional power of vissE is to perform a novel type of gene-set enrichment analysis based on the network of similarity between gene-sets. Given a list of genes (e.g. from a DE analysis), vissE can characterise said list by first identifying all other gene-sets that are similar to it, following up with clustering the resulting gene-sets and finally performing text-mining to reveal emerging themes.
In addition to these analyses, it provides visualisations to assist the users in understanding the results of their experiment. This document will demonstrate these functions across the two use-cases. The vissE package can be downloaded as follows:
2 Summarising the results of a gene-set enrichment analysis
Often, the results of a gene-set enrichment analysis (be it an over representation analysis of a functional class scoring approach) is a list of gene-sets that are accompanied by their statistics and p-values or false discovery rates (FDR). These results are mostly scanned through by biologists who then extract relevant themes pertaining to the experiment of interest. The approach here, vissE, will allow automated extraction of themes.
The example below can be used with the results of any enrichment analysis. The data below is simulated to demonstrate the workflow.
library(msigdb)
library(GSEABase)
#load the MSigDB from the msigdb package
msigdb_hs = getMsigdb()
#append KEGG gene-sets
msigdb_hs = appendKEGG(msigdb_hs)
#select h, c2, and c5 collections (recommended)
msigdb_hs = subsetCollection(msigdb_hs, c('h', 'c2', 'c5'))
#randomly sample gene-sets to simulate the results of an enrichment analysis
set.seed(360)
geneset_res = sample(sapply(msigdb_hs, setName), 2500)
#create a GeneSetCollection using the gene-set analysis results
geneset_gsc = msigdb_hs[geneset_res]
geneset_gsc
#> GeneSetCollection
#>   names: GOMF_GLUCOSE_6-PHOSPHATE:INORGANIC_PHOSPHATE_ANTIPORTER_ACTIVITY, GOCC_MULTIVESICULAR_BODY,_INTERNAL_VESICLE, ..., REACTOME_LYSINE_CATABOLISM (2500 total)
#>   unique identifiers: SLC37A4, SLC37A1, ..., IGIP (20039 total)
#>   types in collection:
#>     geneIdType: SymbolIdentifier (1 total)
#>     collectionType: BroadCollection (1 total)A vissE analysis involves 3 steps:
- Compute gene-set overlaps and the gene-set overlap network
- Identify clusters of gene-sets based on their overlap
- Characterise clusters using text mining
- (Optional) Visualise gene-level statistics
2.1 Compute gene-set overlap
The default approach to computing overlaps is using the Jaccard index. Overlap is computed based on the gene overlap between gene-sets. Alternatively, the overlap coefficient can be used. The latter can be used to highlight hierarchical overlaps (such as those present in the gene ontology).
library(vissE)
#compute gene-set overlap
gs_ovlap = computeMsigOverlap(geneset_gsc, thresh = 0.25)
#create an overlap network
gs_ovnet = computeMsigNetwork(gs_ovlap, msigdb_hs)
#plot the network
set.seed(36) #set seed for reproducible layout
plotMsigNetwork(gs_ovnet)The overlap network plot above is annotated using the MSigDB category. If gene-set statistics are available, they can be projected onto the network too. Gene-set statistics can be passed onto the plotting function as a named vector.
#simulate gene-set statistics
geneset_stats = rnorm(2500)
names(geneset_stats) = geneset_res
head(geneset_stats)
#> GOMF_GLUCOSE_6-PHOSPHATE:INORGANIC_PHOSPHATE_ANTIPORTER_ACTIVITY 
#>                                                      -0.02780652 
#>                       GOCC_MULTIVESICULAR_BODY,_INTERNAL_VESICLE 
#>                                                      -0.54903451 
#>          HP_APLASIA_HYPOPLASIA_OF_THE_ABDOMINAL_WALL_MUSCULATURE 
#>                                                      -0.05165514 
#>              GOBP_REGULATION_OF_ANTIBACTERIAL_PEPTIDE_PRODUCTION 
#>                                                      -0.12148471 
#>       GOBP_POSITIVE_REGULATION_OF_NITRIC-OXIDE_SYNTHASE_ACTIVITY 
#>                                                      -1.32827497 
#>                                            GOCC_PCNA-P21_COMPLEX 
#>                                                       0.74960550
#plot the network and overlay gene-set statistics
set.seed(36) #set seed for reproducible layout
plotMsigNetwork(gs_ovnet, genesetStat = geneset_stats)2.2 Identify clusters of gene-sets
Related gene-sets likely represent related processes. The next step is to identify clusters of gene-sets so that they can be assessed for biological themes. The specific clustering approach can be selected by the user though we recommend graph clustering approaches to use the information provided in the overlap graph. We recommend using the igraph::cluster_walktrap() algorithm as it works well with dense graphs. Many other algorithms are implemented in the igraph package and these can be used instead of the walktrap algorithm.
library(igraph)
#identify clusters
grps = cluster_walktrap(gs_ovnet)
#extract clustering results
grps = groups(grps)
#sort by cluster size
grps = grps[order(sapply(grps, length), decreasing = TRUE)]
#plot the top 12 clusters
set.seed(36) #set seed for reproducible layout
plotMsigNetwork(gs_ovnet, markGroups = grps[1:6], genesetStat = geneset_stats)2.3 Characterise gene-set clusters
Gene-set clusters identified can be assessed for their biological similarities using text-mining approaches. Here, we perform a frequency analysis (adjusted for using the inverse document frequency) on the gene-set names or their short descriptions to assess recurring biological themes in clusters. These results are then presented as word clouds.
#compute and plot the results of text-mining
#using gene-set Names
plotMsigWordcloud(msigdb_hs, grps[1:6], type = 'Name')2.4 Visualise gene-level statistics for gene-set clusters
Gene-level statistics for each gene-set cluster can be visualised to better understand the genes contributing to significance of gene-sets. Gene-level statistics can be passed onto the plotting function as a named vector. A jitter is applied on the x-axis (due to its discrete nature).
library(ggplot2)
#simulate gene-set statistics
set.seed(36)
genes = unique(unlist(geneIds(geneset_gsc)))
gene_stats = rnorm(length(genes))
names(gene_stats) = genes
head(gene_stats)
#>    SLC37A4    SLC37A1    SLC37A2       CD63       EGFR    LAPTM4B 
#>  0.3117314  0.8498291  0.7055331  1.6999284 -1.3455710 -0.5698134
#plot the gene-level statistics
plotGeneStats(gene_stats, msigdb_hs, grps[1:6]) +
  geom_hline(yintercept = 0, colour = 2, lty = 2)2.5 Visualise protein-protein interactions (PPI) in each cluster
An alternative line of evidence for a common functional role of genes are the protein-protein interactions between them. Genes involved in a biological process are likely to interact with each other to achieve the desired function. We can therefore investigate protein-protein interactions within each cluster and thus assess evidence of a common process. In vissE, this can be done by inducing the protein-protein interaction of all genes in a gene-set cluster. Furthermore, the individual nodes in the network can be mapped onto properties such as the gene-level statistic. Networks can then be filtered based on the gene-level statistic, the confidence value of each interaction and the frequency of each gene in the cluster (i.e., how many gene-sets it belongs to).
We will retrieve the PPI from the msigdb R/Bioconductor package. Setting inferred to TRUE will allow PPIs inferred from across organisms to be used in the analysis.
#load PPI from the msigdb package
ppi = getIMEX('hs', inferred = TRUE)
#create the PPI plot
set.seed(36)
plotMsigPPI(
  ppi,
  msigdb_hs,
  grps[1:6],
  geneStat = gene_stats,
  threshStatistic = 1.5,
  threshConfidence = 0.5
)2.6 Combine results to interpret results
Results of a vissE analysis are best presented and interpreted as paneled plots that combine all of the above plots. This allows for collective interpretation of the gene-set clusters.
library(patchwork)
#create independent plots
set.seed(36) #set seed for reproducible layout
p1 = plotMsigWordcloud(msigdb_hs, grps[1:6], type = 'Name')
p2 = plotMsigNetwork(gs_ovnet, markGroups = grps[1:6], genesetStat = geneset_stats)
p3 = plotGeneStats(gene_stats, msigdb_hs, grps[1:6]) +
  geom_hline(yintercept = 0, colour = 2, lty = 2)
p4 = plotMsigPPI(
  ppi,
  msigdb_hs,
  grps[1:6],
  geneStat = gene_stats,
  threshStatistic = 1.5,
  threshConfidence = 0.5
)
#combine using functions from ggpubr
p1 + p2 + p3 + p4 + plot_layout(2, 2)3 Session information
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] patchwork_1.1.1      ggplot2_3.3.5        igraph_1.2.7        
#>  [4] vissE_1.2.2          GSEABase_1.56.0      graph_1.72.0        
#>  [7] annotate_1.72.0      XML_3.99-0.8         AnnotationDbi_1.56.0
#> [10] IRanges_2.28.0       S4Vectors_0.32.0     Biobase_2.54.0      
#> [13] BiocGenerics_0.40.0  msigdb_1.1.5        
#> 
#> loaded via a namespace (and not attached):
#>   [1] colorspace_2.0-2              ellipsis_0.3.2               
#>   [3] qdapRegex_0.7.2               XVector_0.34.0               
#>   [5] farver_2.1.0                  graphlayouts_0.7.1           
#>   [7] scico_1.2.0                   ggwordcloud_0.5.0            
#>   [9] ggrepel_0.9.1                 bit64_4.0.5                  
#>  [11] interactiveDisplayBase_1.32.0 fansi_0.5.0                  
#>  [13] textclean_0.9.3               xml2_1.3.2                   
#>  [15] cachem_1.0.6                  knitr_1.36                   
#>  [17] polyclip_1.10-0               jsonlite_1.7.2               
#>  [19] dbplyr_2.1.1                  png_0.1-7                    
#>  [21] ggforce_0.3.3                 shiny_1.7.1                  
#>  [23] BiocManager_1.30.16           compiler_4.1.1               
#>  [25] httr_1.4.2                    Matrix_1.3-4                 
#>  [27] assertthat_0.2.1              fastmap_1.1.0                
#>  [29] later_1.3.0                   tweenr_1.0.2                 
#>  [31] htmltools_0.5.2               tools_4.1.1                  
#>  [33] NLP_0.2-1                     gtable_0.3.0                 
#>  [35] glue_1.4.2                    GenomeInfoDbData_1.2.7       
#>  [37] reshape2_1.4.4                dplyr_1.0.7                  
#>  [39] rappdirs_0.3.3                Rcpp_1.0.7                   
#>  [41] slam_0.1-48                   lexicon_1.2.1                
#>  [43] jquerylib_0.1.4               vctrs_0.3.8                  
#>  [45] Biostrings_2.62.0             ExperimentHub_2.2.0          
#>  [47] ggraph_2.0.5                  xfun_0.27                    
#>  [49] stringr_1.4.0                 textstem_0.1.4               
#>  [51] mime_0.12                     lifecycle_1.0.1              
#>  [53] AnnotationHub_3.2.0           org.Hs.eg.db_3.14.0          
#>  [55] zlibbioc_1.40.0               MASS_7.3-54                  
#>  [57] scales_1.1.1                  tidygraph_1.2.0              
#>  [59] BiocStyle_2.22.0              promises_1.2.0.1             
#>  [61] parallel_4.1.1                RColorBrewer_1.1-2           
#>  [63] yaml_2.2.1                    curl_4.3.2                   
#>  [65] memoise_2.0.0                 gridExtra_2.3                
#>  [67] sass_0.4.0                    stringi_1.7.5                
#>  [69] RSQLite_2.2.8                 BiocVersion_3.14.0           
#>  [71] highr_0.9                     filelock_1.0.2               
#>  [73] sylly.en_0.1-3                GenomeInfoDb_1.30.0          
#>  [75] rlang_0.4.12                  pkgconfig_2.0.3              
#>  [77] bitops_1.0-7                  lattice_0.20-45              
#>  [79] evaluate_0.14                 koRpus.lang.en_0.1-4         
#>  [81] purrr_0.3.4                   labeling_0.4.2               
#>  [83] bit_4.0.4                     tidyselect_1.1.1             
#>  [85] plyr_1.8.6                    magrittr_2.0.1               
#>  [87] R6_2.5.1                      generics_0.1.1               
#>  [89] DBI_1.1.1                     pillar_1.6.4                 
#>  [91] withr_2.4.2                   prettydoc_0.4.1              
#>  [93] KEGGREST_1.34.0               RCurl_1.98-1.5               
#>  [95] textshape_1.7.3               tibble_3.1.5                 
#>  [97] crayon_1.4.1                  utf8_1.2.2                   
#>  [99] sylly_0.1-6                   BiocFileCache_2.2.0          
#> [101] rmarkdown_2.11                koRpus_0.13-8                
#> [103] viridis_0.6.2                 syuzhet_1.0.6                
#> [105] grid_4.1.1                    data.table_1.14.2            
#> [107] blob_1.2.2                    digest_0.6.28                
#> [109] xtable_1.8-4                  tm_0.7-8                     
#> [111] tidyr_1.1.4                   httpuv_1.6.3                 
#> [113] munsell_0.5.0                 viridisLite_0.4.0            
#> [115] bslib_0.3.1