1 Introduction

The rapid advancement of new technologies for omics data acquisition has resulted in the production of massive amounts of data that combine different omics types, such as gene expression, methylation and copy number variation data. These technologies enable the generation of large-scale, high-resolution and multidimensional data that can reveal the molecular mechanisms and interactions underlying biological phenomena. For instance, The Cancer Genome Atlas (TCGA), a landmark cancer genomics program, molecularly characterized over 20,000 primary cancers and matched normal samples spanning 33 cancer types. These data offer a rich source of information, but they are still underutilized because the new techniques and statistical models for data analysis have not kept up with the pace of data acquisition technologies. Indeed, most of the conventional approaches for omics data modeling depend on the comparison of only one data type across different groups of samples. Moreover, the analysis and integration of omics data pose significant challenges in terms of data interpretation. Therefore, new technologies for omics data acquisition require the development of novel methods and models for data processing and mining that can cope with the volume, variety and complexity of omics data.

In this work, we present gINTomics, our new Multi Omics integration R package, approaching Multi-Omics integration from a new perspective, with the goal to assess the impact of different omics on the final outcome of regulatory networks, that is gene expression. Therefore, for each gene/miRNA, we try to determine the association between its expression and that of its regulators, taking into account also genomic modifications such as copy number variations (CNV) and methylation.

2 Installation

To install this package:

if (!require("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("gINTomics")

#devtools::install_github("angelovelle96/gINTomics")

3 How to use gINTomics

3.1 Structure of the package

gINTomics is designed to perform both two-Omics and Multi Omics integrations. There are four different functions for the two-omics integration. The run_cnv_integration function can be used to integrate Gene or miRNA Expression data with Copy Number Variation data. The run_met_integration is designed for the integration between Gene Expression data with methylation data. When both gene CNV and methylation data are available, the run_genomic_integration function can be used to integrate them together with gene expression. Finally, the run_tf_integration function can be used for the integration between the Expression of a target gene or miRNA and every kind of regulator, such as Transcription Factors or miRNAs expression. The package also gives the possibility to perform Multi Omics integration. The run_multiomics function takes as input a MultiAssayExperiment, that can be generated with the help of the create_multiassay function, and will perform all the possible integration available for the omics provided by the user. All the integration functions exploit different statistical models depending on the input data type, the available models are negative binomial and linear regression. When the number of covariates is too high, a random forest model will select only the most importants which will be included in the integration models. In order to make the results more interpretable, gINTomics provides a comprehensive and interactive shiny app for the graphical representation of the results, the shiny app can be called with the run_shiny function, which takes as input the output of the run_multiomics function. Moreover additional functions are available to build single plots to visualize the results of the integration outside of the shiny app. In the following sections, we use an example MultiAssayExperiment of ovarian cancer to show how to use gINTomics with the different available integrations. The object contains all the available input data types: Gene expression data, miRNA expression data, gene methylation data, gene Copy Number Variations and miRNA Copy Number Variations.

3.2 Generate a MultiAssayExperiment

The package contains a pre built MultiAssayExperiment, anyway in this section we will divide it in signle omics matrices and see how to build a proper MultiAssayExperiment from zero.

# loading packages
library(gINTomics)
library(MultiAssayExperiment)
library(shiny)
data("mmultiassay_ov")
mmultiassay_ov
## A MultiAssayExperiment object of 5 listed
##  experiments with user-defined names and respective classes.
##  Containing an ExperimentList class object of length 5:
##  [1] gene_exp: RangedSummarizedExperiment with 1000 rows and 20 columns
##  [2] methylation: SummarizedExperiment with 1000 rows and 20 columns
##  [3] cnv_data: SummarizedExperiment with 1000 rows and 20 columns
##  [4] miRNA_exp: SummarizedExperiment with 1000 rows and 20 columns
##  [5] miRNA_cnv_data: SummarizedExperiment with 1000 rows and 20 columns
## Functionality:
##  experiments() - obtain the ExperimentList instance
##  colData() - the primary/phenotype DataFrame
##  sampleMap() - the sample coordination DataFrame
##  `$`, `[`, `[[` - extract colData columns, subset, or experiment
##  *Format() - convert into a long or wide DataFrame
##  assays() - convert ExperimentList to a SimpleList of matrices
##  exportClass() - save data to flat files

Here we will extract the single omics from the MultiAssayExperiment

## Here we just select part of the data o speed up the process
tmp <- lapply(experiments(mmultiassay_ov), function(x) x[1:400,])
mmultiassay_ov <- MultiAssayExperiment(experiments = tmp)

gene_exp_matrix <- as.matrix(assay(mmultiassay_ov[["gene_exp"]]))
miRNA_exp_matrix <- as.matrix(assay(mmultiassay_ov[["miRNA_exp"]]))
meth_matrix <- as.matrix(assay(mmultiassay_ov[["methylation"]]))
gene_cnv_matrix <- as.matrix(assay(mmultiassay_ov[["cnv_data"]]))
miRNA_cnv_matrix <- as.matrix(assay(mmultiassay_ov[["miRNA_cnv_data"]]))

Now let’s build a proper MultiAssayExperiment starting from single matrices

new_multiassay <- create_multiassay(methylation = meth_matrix, 
                                    gene_exp = gene_exp_matrix,
                                    cnv_data = gene_cnv_matrix,
                                    miRNA_exp = miRNA_exp_matrix,
                                    miRNA_cnv_data = miRNA_cnv_matrix)

new_multiassay
## A MultiAssayExperiment object of 5 listed
##  experiments with user-defined names and respective classes.
##  Containing an ExperimentList class object of length 5:
##  [1] gene_exp: matrix with 400 rows and 20 columns
##  [2] methylation: matrix with 400 rows and 20 columns
##  [3] cnv_data: matrix with 400 rows and 20 columns
##  [4] miRNA_exp: matrix with 400 rows and 20 columns
##  [5] miRNA_cnv_data: matrix with 400 rows and 20 columns
## Functionality:
##  experiments() - obtain the ExperimentList instance
##  colData() - the primary/phenotype DataFrame
##  sampleMap() - the sample coordination DataFrame
##  `$`, `[`, `[[` - extract colData columns, subset, or experiment
##  *Format() - convert into a long or wide DataFrame
##  assays() - convert ExperimentList to a SimpleList of matrices
##  exportClass() - save data to flat files

3.3 Run Genomic integration

In this section we will see how to perform a gene_expression~CNV+met integration. The input data should be provided as matrices or data.frames. Rows represent samples, while each column represents the different response variables/covariates of the models. By default all the integration functions will assume that expression data are obtained through sequencing, you should set the sequencing_data argument to FALSE if they are not. Expression data could be both normalized or not normalized, the function is able to normalize data by setting the normalize argument to TRUE (default). If you want you can specify custom interaction through the interactions argument otherwise the function will automatically look for the genes of expression in cnv_data. NOTE: if you customize genomics interactions, linking more covariates to a single response, the output may be not more compatible with the shiny visualization.

gene_genomic_integration <- run_genomic_integration(expression = t(gene_exp_matrix),
                                                    cnv_data = t(gene_cnv_matrix),
                                                    methylation = t(meth_matrix))
summary(gene_genomic_integration)
##           Length Class      Mode
## coef_data  3     data.frame list
## pval_data  3     data.frame list
## fdr_data   3     data.frame list
## residuals 20     data.frame list
## data       3     -none-     list

3.4 Run CNV integration

In this section we will see how to perform an expression~CNV integration. The input data (for both expression and CNV) should be provided as matrices or data.frames. Rows represent samples, while each column represents the different response variables/covariates of the models. By default all the integration functions will assume that expression data are obtained through sequencing, you should set the sequencing_data argument to FALSE if they are not. Expression data could be both normalized or not normalized, the function is able to normalize data by setting the normalize argument to TRUE (default). If you want you can specify custom interaction through the interactions argument otherwise the function will automatically look for the genes of expression in cnv_data.

gene_cnv_integration <- run_cnv_integration(expression = t(gene_exp_matrix),
                                            cnv_data = t(gene_cnv_matrix))
summary(gene_cnv_integration)
##           Length Class      Mode
## coef_data  2     data.frame list
## pval_data  2     data.frame list
## fdr_data   2     data.frame list
## residuals 20     data.frame list
## data       3     -none-     list

3.5 Run methylation integration

In this section we will see how to perform an expression~methylation integration. The input data should be the same of run_cnv_integration, but the covariates matrix will contain methylation data instead of Copy Number Variations data. Expression data could be both normalized or not normalized, the function is able to normalize data by setting the normalize argument to TRUE (default). If you want you can specify custom interaction through the interactions argument otherwise the function will automatically look for the genes of expression in methylation.

gene_met_integration <- run_met_integration(expression = t(gene_exp_matrix),
                                            methylation = t(meth_matrix))
summary(gene_met_integration)
##           Length Class      Mode
## coef_data  2     data.frame list
## pval_data  2     data.frame list
## fdr_data   2     data.frame list
## residuals 20     data.frame list
## data       3     -none-     list

3.6 Run TF-target integration

In this section we will see how to perform an expression~TF integration. In this case you can use as input data a single gene expression matrix if both your TF and targets are genes and are contained in the gene expression matrix. If you want the package to automatically download the interactions between TFs and targets you have to set the argument type to “tf”. Otherwise you can also specify your custom interactions providing them with the interactions argument. You can handle data normalization as in the previous functions through the normalize argument. This function is designed to integrate TF expression data but it can handle every type of numerical data representing a gene expression regulator. So you can pass the regulators matrix to the tf_expression argument and specify your custom interactions. If expression data are not obtained through sequencing, remember to set sequencing data to FALSE.

tf_target_integration <- run_tf_integration(expression = t(gene_exp_matrix),
                                            type = "tf")
## Warning: One or more parsing issues, call `problems()` on your data frame for details,
## e.g.:
##   dat <- vroom(...)
##   problems(dat)
summary(tf_target_integration)
##           Length Class      Mode
## coef_data  2     data.frame list
## pval_data  2     data.frame list
## fdr_data   2     data.frame list
## residuals 20     data.frame list
## data       3     -none-     list

3.7 Run miRNA-target integration

In this section we will see how to perform an expression~miRNA integration. Gene expression data will be provided through the expression argument while miRNA expression data through the tf_expression argument. Input matrices should follow the same rules of the previous functions. If you want the package to automatically download the interactions between miRNA and targets you have to set the argument type to “miRNA_target”. Otherwise you can also specify your custom interactions providing them with the interactions argument. You can handle data normalization as in the previous functions through the normalize argument.

miRNA_target_integration <- run_tf_integration(expression = t(gene_exp_matrix),
                                               tf_expression = t(miRNA_exp_matrix),
                                               type = "miRNA_target")
summary(miRNA_target_integration)
##           Length Class      Mode
## coef_data  4     data.frame list
## pval_data  4     data.frame list
## fdr_data   4     data.frame list
## residuals 20     data.frame list
## data       3     -none-     list

3.8 Run TF-miRNA integration

In this section we will see how to perform an miRNA~TF integration. miRNA expression data will be provided through the expression argument while gene expression data through the tf_expression argument. Input matrices should follow the same rules of the previous functions. If you want the package to automatically download the interactions between TF and miRNA you have to set the argument type to “tf_miRNA”. Otherwise you can also specify your custom interactions providing them with the interactions argument. You can handle data normalization as in the previous functions through the normalize argument.

tf_miRNA_integration <- run_tf_integration(expression = t(miRNA_exp_matrix),
                                               tf_expression = t(gene_exp_matrix),
                                               type = "tf_miRNA")
summary(tf_miRNA_integration)
##           Length Class      Mode
## coef_data  7     data.frame list
## pval_data  7     data.frame list
## fdr_data   7     data.frame list
## residuals 20     data.frame list
## data       3     -none-     list

3.9 Run complete Multi-Omics integration

Finally we will see how to perform a complete integration using all the available omics. In order to run this function you need to generate a MultiAssayExperiment as showed at the beginning of this vignette. The function will automatically use all the omics contained in the MultiAssayExperiment to perform all the possible integrations showed before. Moreover, if genomic data are available (CNV and/or Methylation), the first step will be the genomic integration and all the following integrations that contain gene expression data as response variable will use instead the residuals of the genomic model in order to filter out the effect of CNV and/or methylation. This framework is used also for miRNA, but in this case only CNV data are supported.

## Here we run the model
multiomics_integration <- run_multiomics(data = new_multiassay)
summary(multiomics_integration)
##                  Length Class  Mode
## gene_genomic_res 5      -none- list
## mirna_cnv_res    5      -none- list
## tf_res           5      -none- list
## tf_mirna_res     5      -none- list
## mirna_target_res 5      -none- list

4 Visualization

Now let’s see how you can visualize the results of your integration models.

4.1 Shiny app

gINTomics provides an interactive environment, based on Shiny, that allows the user to easily visualize output results from integration models and to save them for downstream tasks and reports creation. Once multiomics integration has been performed users can provide the results to the run_shiny function. NOTE: Only the output of the run_multiomics function is compatible with run_shiny

run_shiny(multiomics_integration)

4.2 plots

4.2.1 network plot

Network plot shows the significant interactions between transcriptional regulators (TFs and miRNAs, if present) and their targets genes/miRNA. Nodes and edges are selected ordering them by the most high coefficient values (absolute value) and by default the top 200 interactions are showed. You can change the number of interactions with the num_interactions argument. Here you can see the code to generate a network from a multiomics integration

data_table <- extract_model_res(multiomics_integration)
##   2135 genes were dropped because they have exons located on both strands
##   of the same reference sequence or on more than one reference sequence,
##   so cannot be represented by a single genomic range.
##   Use 'single.strand.genes.only=FALSE' to get all the genes in a
##   GRangesList object, or use suppressMessages() to suppress this message.
## 'select()' returned 1:many mapping between keys and columns
data_table <- data_table[data_table$cov!="(Intercept)",]
# plot_network(data_table, num_interactions = 200)

4.2.2 Venn Diagram

The Venn Diagram is designed for the genomic integration. It can help to identify genes which are significantly regulated by both CNV and methylation. Here you can see the code to generate a Venn Diagram from a multiomics integration

# plot_venn(data_table)

4.2.3 Volcano plot

The Volcano Plot shows the distribution of integration coefficients for every integration type associated with a genomic class (cnv, met, cnv-mirna). For each integration coefficient, on the y axis you have the -log10 of Pvalue/FDR and on the x axis the value of the coefficient. Here you can see the code to generate a Volcano plot for CNV and one for methylation from a multiomics integration

# plot_volcano(data_table, omics = "gene_genomic_res", cnv_met = "cnv")
# plot_volcano(data_table, omics = "gene_genomic_res", cnv_met = "met")

4.2.4 Ridgeline plot

The ridgeline plot is designed to compare different distributions, it has been integrated in the package with the aim to compare the distribution of significant and non significant coefficients returned by our integration models. For each distribution, on the y axis you have the frequencies and on the x axis the values of the coefficients. Here you can see the code to generate a Ridgeline plot for CNV and one for methylation from a multiomics integration

# plot_ridge(data_table, omics = "gene_genomic_res", cnv_met = "cnv")
# plot_ridge(data_table, omics = "gene_genomic_res", cnv_met = "met")

4.2.5 Chromosome distribution plot

This barplot highlights the distribution of significant and non significant covariates among chromosomes. This kind of visualization will identify chromosomes in which the type of regulation under analysis is particularly active. Here you can see the code to generate a chromosome distribution plot specif for the genomic integration, starting from the results of the multiomics integration.

# plot_chr_distribution(data_table = data_table,
#                     omics = "gene_genomic_res")

4.2.6 TF distribution plot

This barplot highlights the number of significant tragets for each TFs. Here you can see the code to generate a TF distribution plot starting from the results of the multiomics integration.

# plot_tf_distribution(data_table = data_table)

4.2.7 Enrichment plot

The enrichment plot shows the enrichment results obtained with enrichGO and enrichKEGG (clusterProfiler). The genomic enrichment is performed providing the list of genes significantly regulated by methylation or CNV, while the transcriptional one with the list of genes significantly regulated by each transcription factor (we run an enrichment for each TF that significantly regulates at least 12 targets) Here you can see the code to run a genomic enrichment starting from the results of the multiomics integration and to visualize the top results with a DotPlot

#gen_enr <- run_genomic_enrich(multiomics_integration, qvalueCutoff = 1, pvalueCutoff = 0.05, pAdjustMethod = "none")
#dot_plotly(gen_enr$cnv[[1]]$go)

5 Session info

Here is the output of sessionInfo() on the system on which this document was compiled.

sessionInfo()
## R version 4.4.0 beta (2024-04-15 r86425)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 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
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [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    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] shiny_1.8.1.1               MultiAssayExperiment_1.30.0
##  [3] SummarizedExperiment_1.34.0 Biobase_2.64.0             
##  [5] GenomicRanges_1.56.0        GenomeInfoDb_1.40.0        
##  [7] IRanges_2.38.0              S4Vectors_0.42.0           
##  [9] BiocGenerics_0.50.0         MatrixGenerics_1.16.0      
## [11] matrixStats_1.3.0           gINTomics_1.0.0            
## [13] BiocStyle_2.32.0           
## 
## loaded via a namespace (and not attached):
##   [1] fs_1.6.4                                 
##   [2] bitops_1.0-7                             
##   [3] enrichplot_1.24.0                        
##   [4] fontawesome_0.5.2                        
##   [5] lubridate_1.9.3                          
##   [6] HDO.db_0.99.1                            
##   [7] httr_1.4.7                               
##   [8] RColorBrewer_1.1-3                       
##   [9] doParallel_1.0.17                        
##  [10] tools_4.4.0                              
##  [11] backports_1.4.1                          
##  [12] utf8_1.2.4                               
##  [13] R6_2.5.1                                 
##  [14] DT_0.33                                  
##  [15] lazyeval_0.2.2                           
##  [16] GetoptLong_1.0.5                         
##  [17] withr_3.0.0                              
##  [18] graphite_1.50.0                          
##  [19] prettyunits_1.2.0                        
##  [20] gridExtra_2.3                            
##  [21] cli_3.6.2                                
##  [22] scatterpie_0.2.2                         
##  [23] sass_0.4.9                               
##  [24] randomForest_4.7-1.1                     
##  [25] readr_2.1.5                              
##  [26] ggridges_0.5.6                           
##  [27] Rsamtools_2.20.0                         
##  [28] systemfonts_1.0.6                        
##  [29] yulab.utils_0.1.4                        
##  [30] gson_0.1.0                               
##  [31] DOSE_3.30.0                              
##  [32] svglite_2.1.3                            
##  [33] InteractiveComplexHeatmap_1.12.0         
##  [34] limma_3.60.0                             
##  [35] readxl_1.4.3                             
##  [36] rstudioapi_0.16.0                        
##  [37] RSQLite_2.3.6                            
##  [38] visNetwork_2.1.2                         
##  [39] generics_0.1.3                           
##  [40] gridGraphics_0.5-1                       
##  [41] shape_1.4.6.1                            
##  [42] BiocIO_1.14.0                            
##  [43] vroom_1.6.5                              
##  [44] gtools_3.9.5                             
##  [45] dplyr_1.1.4                              
##  [46] GO.db_3.19.1                             
##  [47] Matrix_1.7-0                             
##  [48] fansi_1.0.6                              
##  [49] logger_0.3.0                             
##  [50] abind_1.4-5                              
##  [51] lifecycle_1.0.4                          
##  [52] yaml_2.3.8                               
##  [53] edgeR_4.2.0                              
##  [54] qvalue_2.36.0                            
##  [55] SparseArray_1.4.0                        
##  [56] BiocFileCache_2.12.0                     
##  [57] grid_4.4.0                               
##  [58] blob_1.2.4                               
##  [59] promises_1.3.0                           
##  [60] shinydashboard_0.7.2                     
##  [61] crayon_1.5.2                             
##  [62] lattice_0.22-6                           
##  [63] cowplot_1.1.3                            
##  [64] GenomicFeatures_1.56.0                   
##  [65] chromote_0.2.0                           
##  [66] KEGGREST_1.44.0                          
##  [67] pillar_1.9.0                             
##  [68] knitr_1.46                               
##  [69] ComplexHeatmap_2.20.0                    
##  [70] fgsea_1.30.0                             
##  [71] rjson_0.2.21                             
##  [72] clisymbols_1.2.0                         
##  [73] codetools_0.2-20                         
##  [74] fastmatch_1.1-4                          
##  [75] glue_1.7.0                               
##  [76] ggvenn_0.1.10                            
##  [77] ggfun_0.1.4                              
##  [78] data.table_1.15.4                        
##  [79] vctrs_0.6.5                              
##  [80] png_0.1-8                                
##  [81] treeio_1.28.0                            
##  [82] org.Mm.eg.db_3.19.1                      
##  [83] cellranger_1.1.0                         
##  [84] gtable_0.3.5                             
##  [85] cachem_1.0.8                             
##  [86] OmnipathR_3.12.0                         
##  [87] xfun_0.43                                
##  [88] S4Arrays_1.4.0                           
##  [89] mime_0.12                                
##  [90] tidygraph_1.3.1                          
##  [91] iterators_1.0.14                         
##  [92] statmod_1.5.0                            
##  [93] nlme_3.1-164                             
##  [94] ggtree_3.12.0                            
##  [95] bit64_4.0.5                              
##  [96] progress_1.2.3                           
##  [97] filelock_1.0.3                           
##  [98] bslib_0.7.0                              
##  [99] colorspace_2.1-0                         
## [100] DBI_1.2.2                                
## [101] tidyselect_1.2.1                         
## [102] processx_3.8.4                           
## [103] bit_4.0.5                                
## [104] compiler_4.4.0                           
## [105] curl_5.2.1                               
## [106] shiny.gosling_1.0.0                      
## [107] rvest_1.0.4                              
## [108] httr2_1.0.1                              
## [109] graph_1.82.0                             
## [110] xml2_1.3.6                               
## [111] plotly_4.10.4                            
## [112] TxDb.Mmusculus.UCSC.mm10.knownGene_3.10.0
## [113] DelayedArray_0.30.0                      
## [114] bookdown_0.39                            
## [115] shadowtext_0.1.3                         
## [116] rtracklayer_1.64.0                       
## [117] checkmate_2.3.1                          
## [118] scales_1.3.0                             
## [119] callr_3.7.6                              
## [120] rappdirs_0.3.3                           
## [121] stringr_1.5.1                            
## [122] digest_0.6.35                            
## [123] rmarkdown_2.26                           
## [124] XVector_0.44.0                           
## [125] htmltools_0.5.8.1                        
## [126] pkgconfig_2.0.3                          
## [127] highr_0.10                               
## [128] dbplyr_2.5.0                             
## [129] fastmap_1.1.1                            
## [130] rlang_1.1.3                              
## [131] GlobalOptions_0.1.2                      
## [132] htmlwidgets_1.6.4                        
## [133] UCSC.utils_1.0.0                         
## [134] farver_2.1.1                             
## [135] jquerylib_0.1.4                          
## [136] jsonlite_1.8.8                           
## [137] BiocParallel_1.38.0                      
## [138] GOSemSim_2.30.0                          
## [139] RCurl_1.98-1.14                          
## [140] magrittr_2.0.3                           
## [141] kableExtra_1.4.0                         
## [142] GenomeInfoDbData_1.2.12                  
## [143] ggplotify_0.1.2                          
## [144] patchwork_1.2.0                          
## [145] munsell_0.5.1                            
## [146] Rcpp_1.0.12                              
## [147] ape_5.8                                  
## [148] viridis_0.6.5                            
## [149] stringi_1.8.3                            
## [150] ggraph_2.2.1                             
## [151] zlibbioc_1.50.0                          
## [152] MASS_7.3-60.2                            
## [153] org.Hs.eg.db_3.19.1                      
## [154] plyr_1.8.9                               
## [155] parallel_4.4.0                           
## [156] ggrepel_0.9.5                            
## [157] Biostrings_2.72.0                        
## [158] graphlayouts_1.1.1                       
## [159] splines_4.4.0                            
## [160] hms_1.1.3                                
## [161] circlize_0.4.16                          
## [162] locfit_1.5-9.9                           
## [163] ps_1.7.6                                 
## [164] igraph_2.0.3                             
## [165] reshape2_1.4.4                           
## [166] biomaRt_2.60.0                           
## [167] XML_3.99-0.16.1                          
## [168] evaluate_0.23                            
## [169] BiocManager_1.30.22                      
## [170] selectr_0.4-2                            
## [171] tzdb_0.4.0                               
## [172] foreach_1.5.2                            
## [173] tweenr_2.0.3                             
## [174] httpuv_1.6.15                            
## [175] tidyr_1.3.1                              
## [176] purrr_1.0.2                              
## [177] polyclip_1.10-6                          
## [178] clue_0.3-65                              
## [179] ggplot2_3.5.1                            
## [180] BiocBaseUtils_1.6.0                      
## [181] ReactomePA_1.48.0                        
## [182] ggforce_0.4.2                            
## [183] xtable_1.8-4                             
## [184] restfulr_0.0.15                          
## [185] reactome.db_1.88.0                       
## [186] tidytree_0.4.6                           
## [187] later_1.3.2                              
## [188] viridisLite_0.4.2                        
## [189] TxDb.Hsapiens.UCSC.hg38.knownGene_3.18.0 
## [190] tibble_3.2.1                             
## [191] clusterProfiler_4.12.0                   
## [192] websocket_1.4.1                          
## [193] aplot_0.2.2                              
## [194] memoise_2.0.1                            
## [195] AnnotationDbi_1.66.0                     
## [196] GenomicAlignments_1.40.0                 
## [197] cluster_2.1.6                            
## [198] timechange_0.3.0