Contents

library(TDbasedUFE)
library(TDbasedUFEadv)
#> 
library(DOSE)
#> DOSE v3.28.0  For help: https://yulab-smu.top/biomedical-knowledge-mining-book/
#> 
#> If you use DOSE in published research, please cite:
#> Guangchuang Yu, Li-Gen Wang, Guang-Rong Yan, Qing-Yu He. DOSE: an R/Bioconductor package for Disease Ontology Semantic and Enrichment analysis. Bioinformatics 2015, 31(4):608-609
library(enrichplot)
library(RTCGA.rnaseq)
#> Loading required package: RTCGA
#> Welcome to the RTCGA (version: 1.32.0). Read more about the project under https://rtcga.github.io/RTCGA/
library(RTCGA.clinical)
library(enrichR)
#> Welcome to enrichR
#> Checking connection ...
#> Enrichr ... Connection is Live!
#> FlyEnrichr ... Connection is Live!
#> WormEnrichr ... Connection is Live!
#> YeastEnrichr ... Connection is Live!
#> FishEnrichr ... Connection is Live!
#> OxEnrichr ... Connection is Live!
library(STRINGdb)

1 Introduction

It might be helpful to demonstrate how to evaluate selected genes by enrichment analysis. Here, we show some of useful tools applied to the output from TDbasedUFEadv In order foe this, we reproduce one example in “How to use TDbasedUFEadv” as follows.

Multi <- list(
  BLCA.rnaseq[seq_len(100), 1 + seq_len(1000)],
  BRCA.rnaseq[seq_len(100), 1 + seq_len(1000)],
  CESC.rnaseq[seq_len(100), 1 + seq_len(1000)],
  COAD.rnaseq[seq_len(100), 1 + seq_len(1000)]
)
Z <- prepareTensorfromList(Multi, 10L)
Z <- aperm(Z, c(2, 1, 3))
Clinical <- list(BLCA.clinical, BRCA.clinical, CESC.clinical, COAD.clinical)
Multi_sample <- list(
  BLCA.rnaseq[seq_len(100), 1, drop = FALSE],
  BRCA.rnaseq[seq_len(100), 1, drop = FALSE],
  CESC.rnaseq[seq_len(100), 1, drop = FALSE],
  COAD.rnaseq[seq_len(100), 1, drop = FALSE]
)
# patient.stage_event.tnm_categories.pathologic_categories.pathologic_m
ID_column_of_Multi_sample <- c(770, 1482, 773, 791)
# patient.bcr_patient_barcode
ID_column_of_Clinical <- c(20, 20, 12, 14)
Z <- PrepareSummarizedExperimentTensor(
  feature = colnames(ACC.rnaseq)[1 + seq_len(1000)],
  sample = array("", 1), value = Z,
  sampleData = prepareCondTCGA(
    Multi_sample, Clinical,
    ID_column_of_Multi_sample, ID_column_of_Clinical
  )
)
HOSVD <- computeHosvd(Z)
#> 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |=======================                                               |  33%
  |                                                                            
  |===============================================                       |  67%
  |                                                                            
  |======================================================================| 100%
cond <- attr(Z, "sampleData")
index <- selectFeatureProj(HOSVD, Multi, cond, de = 1e-3, input_all = 3) # Batch mode

head(tableFeatures(Z, index))
#>       Feature       p value adjusted p value
#> 10    ACTB|60  0.000000e+00     0.000000e+00
#> 11   ACTG1|71  0.000000e+00     0.000000e+00
#> 37  ALDOA|226  0.000000e+00     0.000000e+00
#> 19 ADAM6|8755 5.698305e-299    1.424576e-296
#> 22  AEBP1|165 1.057392e-218    2.114785e-216
#> 9    ACTA2|59 7.862975e-198    1.310496e-195
genes <- unlist(lapply(strsplit(tableFeatures(Z, index)[, 1], "|",
  fixed = TRUE
), "[", 1))
entrez <- unlist(lapply(strsplit(tableFeatures(Z, index)[, 1], "|",
  fixed = TRUE
), "[", 2))

2 Enrichr

Enrichr(Kuleshov et al. 2016) is one of tools that often provides us significant results toward genes selected by TDbasedUFE and TDbasedUFEadv.

setEnrichrSite("Enrichr")
#> Connection changed to https://maayanlab.cloud/Enrichr/
#> Connection is Live!
websiteLive <- TRUE
dbs <- c(
  "GO_Molecular_Function_2015", "GO_Cellular_Component_2015",
  "GO_Biological_Process_2015"
)
enriched <- enrichr(genes, dbs)
#> Uploading data to Enrichr... Done.
#>   Querying GO_Molecular_Function_2015... Done.
#>   Querying GO_Cellular_Component_2015... Done.
#>   Querying GO_Biological_Process_2015... Done.
#> Parsing results... Done.
if (websiteLive) {
  plotEnrich(enriched$GO_Biological_Process_2015,
    showTerms = 20, numChar = 40, y = "Count",
    orderBy = "P.value"
  )
}

Enrichr can provide you huge number of enrichment analyses, many of which have good compatibility with the genes selected by TDbasedUFE as well as TDbasedUFEadv by the experience. Please check Enrichr’s web site to see what kinds of enrichment analyses can be done.

3 STRING

STRING(Szklarczyk et al. 2018) is enrichment analyses based upon protein-protein interaction, which is known to provide often significant results toward genes selected by TDbasedUFE as well as TDbasedUFEadv.

options(timeout = 200)
string_db <- STRINGdb$new(
  version = "11.5",
  species = 9606, score_threshold = 200,
  network_type = "full", input_directory = ""
)
example1_mapped <- string_db$map(data.frame(genes = genes),
  "genes",
  removeUnmappedRows = TRUE
)
#> Warning:  we couldn't map to STRING 1% of your identifiers
hits <- example1_mapped$STRING_id
string_db$plot_network(hits)

4 enrichplot

Although these above can provide us enough number of information to evaluate the genes selected by TDbasedUFE as well as TDbasedUFEadv, one might need all one package for which one does not how to decide which category must be evaluated in enrichment analysis.

In this case, we would recommend Metascape(Zhou et al. 2019) that unfortunately
does not have the ways approached from R. Thus, we recommend RITAN as an alternative. It can list significant ones among multiple categories.

edo <- enrichDGN(entrez)
dotplot(edo, showCategory=30) + ggtitle("dotplot for ORA")

sessionInfo()
#> R version 4.3.1 (2023-06-16)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
#> 
#> 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] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#>  [1] S4Vectors_0.40.0             IRanges_2.36.0              
#>  [3] STRINGdb_2.14.0              enrichR_3.2                 
#>  [5] RTCGA.clinical_20151101.31.0 RTCGA.rnaseq_20151101.31.0  
#>  [7] RTCGA_1.32.0                 enrichplot_1.22.0           
#>  [9] DOSE_3.28.0                  TDbasedUFEadv_1.2.0         
#> [11] TDbasedUFE_1.2.0             BiocStyle_2.30.0            
#> 
#> loaded via a namespace (and not attached):
#>   [1] rTensor_1.4.8                 splines_4.3.1                
#>   [3] later_1.3.1                   ggplotify_0.1.2              
#>   [5] bitops_1.0-7                  filelock_1.0.2               
#>   [7] tibble_3.2.1                  polyclip_1.10-6              
#>   [9] XML_3.99-0.14                 lifecycle_1.0.3              
#>  [11] rstatix_0.7.2                 lattice_0.22-5               
#>  [13] MASS_7.3-60                   backports_1.4.1              
#>  [15] magrittr_2.0.3                sass_0.4.7                   
#>  [17] rmarkdown_2.25                jquerylib_0.1.4              
#>  [19] yaml_2.3.7                    plotrix_3.8-2                
#>  [21] httpuv_1.6.12                 cowplot_1.1.1                
#>  [23] DBI_1.1.3                     RColorBrewer_1.1-3           
#>  [25] abind_1.4-5                   MOFAdata_1.17.1              
#>  [27] zlibbioc_1.48.0               rvest_1.0.3                  
#>  [29] GenomicRanges_1.54.0          purrr_1.0.2                  
#>  [31] ggraph_2.1.0                  BiocGenerics_0.48.0          
#>  [33] RCurl_1.98-1.12               yulab.utils_0.1.0            
#>  [35] hash_2.2.6.3                  WriteXLS_6.4.0               
#>  [37] tweenr_2.0.2                  rappdirs_0.3.3               
#>  [39] GenomeInfoDbData_1.2.11       KMsurv_0.1-5                 
#>  [41] ggrepel_0.9.4                 tidytree_0.4.5               
#>  [43] proto_1.0.0                   codetools_0.2-19             
#>  [45] xml2_1.3.5                    ggforce_0.4.1                
#>  [47] tximportData_1.29.0           tidyselect_1.2.0             
#>  [49] aplot_0.2.2                   farver_2.1.1                 
#>  [51] viridis_0.6.4                 stats4_4.3.1                 
#>  [53] BiocFileCache_2.10.0          jsonlite_1.8.7               
#>  [55] ellipsis_0.3.2                tidygraph_1.2.3              
#>  [57] survival_3.5-7                tools_4.3.1                  
#>  [59] chron_2.3-61                  treeio_1.26.0                
#>  [61] HPO.db_0.99.2                 Rcpp_1.0.11                  
#>  [63] glue_1.6.2                    gridExtra_2.3                
#>  [65] xfun_0.40                     qvalue_2.34.0                
#>  [67] ggthemes_4.2.4                GenomeInfoDb_1.38.0          
#>  [69] dplyr_1.1.3                   withr_2.5.1                  
#>  [71] BiocManager_1.30.22           fastmap_1.1.1                
#>  [73] fansi_1.0.5                   caTools_1.18.2               
#>  [75] digest_0.6.33                 gridGraphics_0.5-1           
#>  [77] R6_2.5.1                      mime_0.12                    
#>  [79] colorspace_2.1-0              GO.db_3.18.0                 
#>  [81] gtools_3.9.4                  RSQLite_2.3.1                
#>  [83] utf8_1.2.4                    tidyr_1.3.0                  
#>  [85] generics_0.1.3                data.table_1.14.8            
#>  [87] graphlayouts_1.0.1            httr_1.4.7                   
#>  [89] scatterpie_0.2.1              sqldf_0.4-11                 
#>  [91] pkgconfig_2.0.3               gtable_0.3.4                 
#>  [93] blob_1.2.4                    XVector_0.42.0               
#>  [95] survMisc_0.5.6                shadowtext_0.1.2             
#>  [97] htmltools_0.5.6.1             carData_3.0-5                
#>  [99] bookdown_0.36                 fgsea_1.28.0                 
#> [101] scales_1.2.1                  Biobase_2.62.0               
#> [103] png_0.1-8                     ggfun_0.1.3                  
#> [105] knitr_1.44                    km.ci_0.5-6                  
#> [107] tzdb_0.4.0                    reshape2_1.4.4               
#> [109] rjson_0.2.21                  nlme_3.1-163                 
#> [111] curl_5.1.0                    cachem_1.0.8                 
#> [113] zoo_1.8-12                    stringr_1.5.0                
#> [115] BiocVersion_3.18.0            KernSmooth_2.23-22           
#> [117] parallel_4.3.1                HDO.db_0.99.1                
#> [119] AnnotationDbi_1.64.0          pillar_1.9.0                 
#> [121] grid_4.3.1                    vctrs_0.6.4                  
#> [123] gplots_3.1.3                  promises_1.2.1               
#> [125] ggpubr_0.6.0                  car_3.1-2                    
#> [127] dbplyr_2.3.4                  xtable_1.8-4                 
#> [129] tximport_1.30.0               evaluate_0.22                
#> [131] readr_2.1.4                   gsubfn_0.7                   
#> [133] cli_3.6.1                     compiler_4.3.1               
#> [135] rlang_1.1.1                   crayon_1.5.2                 
#> [137] ggsignif_0.6.4                labeling_0.4.3               
#> [139] survminer_0.4.9               fs_1.6.3                     
#> [141] plyr_1.8.9                    stringi_1.7.12               
#> [143] viridisLite_0.4.2             BiocParallel_1.36.0          
#> [145] assertthat_0.2.1              MPO.db_0.99.7                
#> [147] munsell_0.5.0                 Biostrings_2.70.0            
#> [149] lazyeval_0.2.2                GOSemSim_2.28.0              
#> [151] Matrix_1.6-1.1                patchwork_1.1.3              
#> [153] hms_1.1.3                     bit64_4.0.5                  
#> [155] ggplot2_3.4.4                 KEGGREST_1.42.0              
#> [157] shiny_1.7.5.1                 interactiveDisplayBase_1.40.0
#> [159] AnnotationHub_3.10.0          igraph_1.5.1                 
#> [161] broom_1.0.5                   memoise_2.0.1                
#> [163] bslib_0.5.1                   ggtree_3.10.0                
#> [165] fastmatch_1.1-4               bit_4.0.5                    
#> [167] ape_5.7-1

Kuleshov, Maxim V., Matthew R. Jones, Andrew D. Rouillard, Nicolas F. Fernandez, Qiaonan Duan, Zichen Wang, Simon Koplev, et al. 2016. “Enrichr: a comprehensive gene set enrichment analysis web server 2016 update.” Nucleic Acids Research 44 (W1): W90–W97. https://doi.org/10.1093/nar/gkw377.

Szklarczyk, Damian, Annika L Gable, David Lyon, Alexander Junge, Stefan Wyder, Jaime Huerta-Cepas, Milan Simonovic, et al. 2018. “STRING v11: protein窶菟rotein association networks with increased coverage, supporting functional discovery in genome-wide experimental datasets.” Nucleic Acids Research 47 (D1): D607–D613. https://doi.org/10.1093/nar/gky1131.

Zhou, Yingyao, Bin Zhou, Lars Pache, Max Chang, Alireza Hadj Khodabakhshi, Olga Tanaseichuk, Christopher Benner, and Sumit K Chanda. 2019. “Metascape provides a biologist-oriented resource for the analysis of systems-level datasets.” Nature Communications 10 (1): 1523. https://doi.org/10.1038/s41467-019-09234-6.