Contents

1 Overview

augere.gsea implements augere pipeline functions for gene set enrichment analyses (GSEA). Each analysis pipeline creates a self-contained Rmarkdown report that documents each step in the DE analysis, after which it compiles the report and returns the analysis results to the user in the current R session. We can then inspect the report to determine exactly how the analysis was performed. We can also modify the report to customize the analysis beyond the options provided by the pipeline function.

2 Quick start

To install this package, follow the usual instructions for Bioconductor:

install.packages("BiocManager") # if not already available
BiocManager::install("augere.gsea")

Let’s demonstrate on a dataset derived from the airway package. First, we perform a differential expression analysis with edgeR:

library(augere.de)
se <- loadExampleDataset()

# Testing for DE between 'trt' and 'untrt' in the 'dex' grouping factor.
output.de <- tempfile()
res.de <- runVoom(se, group="dex", comparisons=c("trt", "untrt"), output.dir=output.de)
res.de$results[[1]]
## DataFrame with 63677 rows and 6 columns
##                      LogFC   AveExpr           t    PValue       FDR         B
##                  <numeric> <numeric>   <numeric> <numeric> <numeric> <numeric>
## ENSG00000000003 -0.3816595   5.02415 -2.85280352 0.0187438 0.0933503  -3.98117
## ENSG00000000005         NA        NA          NA        NA        NA        NA
## ENSG00000000419  0.1895912   4.60047  2.31869869 0.0452045 0.1637096  -4.81906
## ENSG00000000457  0.0242607   3.47124  0.24152618 0.8144866 0.9082962  -6.90369
## ENSG00000000460 -0.0015551   1.40713 -0.00511864 0.9960261 0.9977997  -6.43023
## ...                    ...       ...         ...       ...       ...       ...
## ENSG00000273489         NA        NA          NA        NA        NA        NA
## ENSG00000273490         NA        NA          NA        NA        NA        NA
## ENSG00000273491         NA        NA          NA        NA        NA        NA
## ENSG00000273492         NA        NA          NA        NA        NA        NA
## ENSG00000273493         NA        NA          NA        NA        NA        NA

We then choose the collection of gene sets to test. Any gene sets can be used as long as they share the same gene identifiers as those in the row names of the DE table. Here, we’ll use the first few GO terms from org.Hs.eg.db:

library(org.Hs.eg.db)
sets <- select(org.Hs.eg.db, keys=keys(org.Hs.eg.db, "GO"), keytype="GO", columns="ENSEMBL")
sets <- split(sets$ENSEMBL, sets$GO)
sets <- head(sets, 100) # to cut down on vignette runtime.
length(sets)
## [1] 100

We use the pre-computed DE results to check each gene set for significant enrichment of significant genes. This is done by calling the runPrecomputed() function, which performs a variety of different tests to test for enrichment:

library(augere.gsea)
output.gsea <- tempfile()

# Set the seed just in case any of the methods require randomization.
# runPrecomputed() also hard-codes some set.seed() calls in the generated
# Rmarkdown report to ensure that the report's results are reproducible.
set.seed(100000)

res.gsea <- runPrecomputed(
    res.de$results[[1]],
    sets,
    signif.field="FDR",
    signif.threshold=0.05,
    rank.field="t",
    output.dir=output.gsea
) 
res.gsea
## $hypergeometric
## DataFrame with 100 rows and 4 columns
##             NumGenes    NumSig    PValue       FDR
##            <integer> <integer> <numeric> <numeric>
## GO:0000009         2         0  1.000000         1
## GO:0000012        13         1  0.869174         1
## GO:0000014         8         2  0.326251         1
## GO:0000015         4         0  1.000000         1
## GO:0000016         0         0        NA        NA
## ...              ...       ...       ...       ...
## GO:0000244        37         1  0.996952         1
## GO:0000245        34         2  0.966982         1
## GO:0000246         1         0  1.000000         1
## GO:0000247         1         0  1.000000         1
## GO:0000248         2         0  1.000000         1
## 
## $goseq
## DataFrame with 100 rows and 4 columns
##             NumGenes    NumSig    PValue       FDR
##            <integer> <integer> <numeric> <numeric>
## GO:0000009         2         0  1.000000         1
## GO:0000012        13         1  0.858395         1
## GO:0000014         8         2  0.377812         1
## GO:0000015         4         0  1.000000         1
## GO:0000016         0         0        NA        NA
## ...              ...       ...       ...       ...
## GO:0000244        37         1  0.999196         1
## GO:0000245        34         2  0.988993         1
## GO:0000246         1         0  1.000000         1
## GO:0000247         1         0  1.000000         1
## GO:0000248         2         0  1.000000         1
## 
## $fgsea
## DataFrame with 100 rows and 5 columns
##             NumGenes        ES       NES    PValue       FDR
##            <integer> <numeric> <numeric> <numeric> <numeric>
## GO:0000009         2  0.408141  0.760599  0.696304  0.998002
## GO:0000012        13  0.477022  1.068831  0.414585  0.998002
## GO:0000014         8  0.562902  1.180735  0.281718  0.998002
## GO:0000015         4  0.452687  0.875743  0.643357  0.998002
## GO:0000016         0        NA        NA        NA        NA
## ...              ...       ...       ...       ...       ...
## GO:0000244        37 0.2017021  0.495753  0.998002  0.998002
## GO:0000245        34 0.2638262  0.644800  0.979021  0.998002
## GO:0000246         1 0.6996445  1.373488  0.333666  0.998002
## GO:0000247         1 0.5862449  1.150871  0.425574  0.998002
## GO:0000248         2 0.0854533  0.159248  0.974026  0.998002
## 
## $cameraPR
## DataFrame with 100 rows and 3 columns
##             NumGenes    PValue       FDR
##            <integer> <numeric> <numeric>
## GO:0000009         2 0.4856659  0.913107
## GO:0000012        13 0.0576934  0.520878
## GO:0000014         8 0.7417102  0.966919
## GO:0000015         4 0.5560069  0.913107
## GO:0000016        NA        NA        NA
## ...              ...       ...       ...
## GO:0000244        37  0.957846  0.985176
## GO:0000245        34  0.921241  0.985176
## GO:0000246         1  0.244644  0.901829
## GO:0000247         1  0.382609  0.905196
## GO:0000248         2  0.975020  0.985176

Like its counterparts from augere.de, runPrecomputed() also generates an Rmarkdown report containing (almost) all of the steps required to reproduce the analysis. Let’s have a look at it:

fname <- file.path(output.gsea, "report.Rmd")
cat(head(readLines(fname), 50), sep="\n")
## ---
## title: GSEA with precomputed statistics
## author:
##   - biocbuild
## date: "`r Sys.Date()`"
## output: BiocStyle::html_document
## ---
## 
## ```{r, echo=FALSE}
## knitr::opts_chunk$set(error=FALSE, warning=FALSE, message=FALSE)
## ```
## 
## # Loading statistics
## 
## We construct the data frame of precomputed statistics for each gene.
## This is typically from a differential expression analysis.
## 
## ```{r}
## tab <- local({ # augere.core input (1)
##     stop("insert commands to generate 'tab' here")
## })
## tab
## ```
## 
## We remove rows with `NA` values, e.g., because they correspond to low-abundance genes.
## 
## ```{r}
## discard <- lapply(tab[,c("FDR", "t"),drop=FALSE], function(x) is.na(x))
## discard <- Reduce(`|`, discard)
## tab <- tab[!discard,,drop=FALSE]
## nrow(tab)
## ```
## 
## We also extract some statistics for use in all methods:
## 
## ```{r}
## is.sig <- tab[,"FDR"] <= 0.05
## summary(is.sig)
## rank.stat <- tab[,"t"]
## summary(rank.stat)
## ```
## 
## # Loading gene sets
## 
## We set up the collection of gene sets. 
## 
## ```{r}
## sets <- local({ # augere.core input (2)
##     stop("insert commands to generate 'sets' here")
## })

The function will also save the results to file inside a results subdirectory, which can be loaded back into the session using the readResult() function:

reloaded <- readResult(file.path(output.gsea, "results", "hypergeometric"))

# The result itself:
reloaded$x
## DataFrame with 100 rows and 4 columns
##             NumGenes    NumSig    PValue       FDR
##            <integer> <integer> <numeric> <numeric>
## GO:0000009         2         0  1.000000         1
## GO:0000012        13         1  0.869174         1
## GO:0000014         8         2  0.326251         1
## GO:0000015         4         0  1.000000         1
## GO:0000016         0         0        NA        NA
## ...              ...       ...       ...       ...
## GO:0000244        37         1  0.996952         1
## GO:0000245        34         2  0.966982         1
## GO:0000246         1         0  1.000000         1
## GO:0000247         1         0  1.000000         1
## GO:0000248         2         0  1.000000         1
# Along with some metadata:
str(reloaded$metadata)
## List of 4
##  $ title                          : chr "Hypergeometric test results"
##  $ authors                        :List of 1
##   ..$ : chr "biocbuild"
##  $ type                           : chr "precomputed_gene_set_enrichment"
##  $ precomputed_gene_set_enrichment:List of 2
##   ..$ method     : chr "hypergeometric"
##   ..$ alternative: chr "mixed"

3 Considering directionality

By default, the alternative hypothesis in runPrecomputed() is non-directional, i.e., we test each set for over-representation of any significant genes. However, we can also consider the sign of the differential expression of the genes. For example, if we are only interested in significant genes that are upregulated in trt versus untrt, we could do:

output.up <- tempfile()
set.seed(100001)
res.up <- runPrecomputed(
    res.de$results[[1]],
    sets,
    alternative="up",
    signif.field="FDR",
    signif.threshold=0.05,
    rank.field="t",
    sign.field="LogFC", # supplying the sign for each gene. 
    output.dir=output.up,
    save.results=FALSE # skipping the save, to reduce run-time.
)
res.up$hypergeometric
## DataFrame with 100 rows and 4 columns
##             NumGenes    NumSig    PValue       FDR
##            <integer> <integer> <numeric> <numeric>
## GO:0000009         2         0  1.000000         1
## GO:0000012        13         0  1.000000         1
## GO:0000014         8         1  0.500057         1
## GO:0000015         4         0  1.000000         1
## GO:0000016         0         0        NA        NA
## ...              ...       ...       ...       ...
## GO:0000244        37         1  0.959616         1
## GO:0000245        34         2  0.785997         1
## GO:0000246         1         0  1.000000         1
## GO:0000247         1         0  1.000000         1
## GO:0000248         2         0  1.000000         1

Alternatively, we can test for sets that are overrepresented in genes changing in either direction. This adds an extra Direction field to the result data frames, specifying the net direction in which the genes are changing within the set. The “either” approach is subtly different from the default of “mixed”; the latter doesn’t care about direction at all, whereas the former tests both directions and reports the more significant one.

output.either <- tempfile()
set.seed(100002)
res.either <- runPrecomputed(
    res.de$results[[1]],
    sets,
    alternative="either",
    signif.field="FDR",
    signif.threshold=0.05,
    rank.field="t",
    sign.field="LogFC",
    output.dir=output.either,
    save.results=FALSE
)
res.either$hypergeometric
## DataFrame with 100 rows and 5 columns
##             NumGenes    NumSig    PValue       FDR   Direction
##            <integer> <integer> <numeric> <numeric> <character>
## GO:0000009         2         0  2.000000         1        down
## GO:0000012        13         1  1.127366         1        down
## GO:0000014         8         2  0.799382         1        down
## GO:0000015         4         0  2.000000         1        down
## GO:0000016         0         0        NA        NA          NA
## ...              ...       ...       ...       ...         ...
## GO:0000244        37         1   1.91923         1          up
## GO:0000245        34         2   1.57199         1          up
## GO:0000246         1         0   2.00000         1        down
## GO:0000247         1         0   2.00000         1        down
## GO:0000248         2         0   2.00000         1        down

In some situations, the test statistic from the DE analysis is unsigned because it is a squared representation of a signed test statistic. For example, the quasi-likelihood framework of edgeR reports F-statistics that are always non-negative. These are converted back to signed statistics by taking the square root and multiplying by the sign of the log-fold change, effectively mimicking a t-statistic.

output.de2 <- tempfile()
res.de2 <- runEdgeR(se, group="dex", comparisons=c("trt", "untrt"), output.dir=output.de2)
res.de2$results[[1]]
## DataFrame with 63677 rows and 5 columns
##                      LogFC   AveExpr         F    PValue       FDR
##                  <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000000003 -0.3862784   5.05487 6.7849013 0.0287436  0.130627
## ENSG00000000005         NA        NA        NA        NA        NA
## ENSG00000000419  0.1967653   4.60947 6.0508824 0.0364214  0.150396
## ENSG00000000457  0.0254907   3.48214 0.0637723 0.8063598  0.911506
## ENSG00000000460 -0.1225528   1.48030 0.2075723 0.6595808  0.827011
## ...                    ...       ...       ...       ...       ...
## ENSG00000273489         NA        NA        NA        NA        NA
## ENSG00000273490         NA        NA        NA        NA        NA
## ENSG00000273491         NA        NA        NA        NA        NA
## ENSG00000273492         NA        NA        NA        NA        NA
## ENSG00000273493         NA        NA        NA        NA        NA
output.up2 <- tempfile()
set.seed(100003)
res.up2 <- runPrecomputed(
    res.de2$results[[1]],
    sets,
    signif.field="FDR",
    signif.threshold=0.05,
    rank.field="F",
    rank.sqrt=TRUE, # square-rooting the F statistic to mimic a t-test.
    sign.field="LogFC", # using this to restore the sign to the test statistic.
    alternative="up",
    output.dir=output.up2,
    save.results=FALSE,
)
res.up2$hypergeometric
## DataFrame with 100 rows and 4 columns
##             NumGenes    NumSig    PValue       FDR
##            <integer> <integer> <numeric> <numeric>
## GO:0000009         2         0  1.000000         1
## GO:0000012        13         0  1.000000         1
## GO:0000014         8         1  0.436034         1
## GO:0000015         4         0  1.000000         1
## GO:0000016         0         0        NA        NA
## ...              ...       ...       ...       ...
## GO:0000244        37         0         1         1
## GO:0000245        34         0         1         1
## GO:0000246         1         0         1         1
## GO:0000247         1         0         1         1
## GO:0000248         2         0         1         1

Of course, these considerations are unnecessary for the default alternative="mixed", where the sign is ignored anyway.

4 Differential gene set tests

runContrast() performs differential gene set tests that compare the expression values across all genes in the set. This means that we need to supply the counts themselves along with the design and contrast:

output.dgst <- tempfile()
set.seed(200000)
res.dgst <- runContrast(
    se,
    sets,
    group="dex",
    comparison=c("trt", "untrt"),
    save.results=FALSE,
    output.dir=output.dgst
) 
res.dgst
## $fry
## DataFrame with 100 rows and 6 columns
##             NumGenes   Direction      PValue       FDR PValue.Mixed   FDR.Mixed
##            <integer> <character>   <numeric> <numeric>    <numeric>   <numeric>
## GO:0000009         2          Up 0.408345874 0.5801702  7.33654e-01 8.33697e-01
## GO:0000012        13        Down 0.000532368 0.0111803  4.52691e-07 1.96822e-06
## GO:0000014         8          Up 0.412698365 0.5801702  6.43704e-04 1.56197e-03
## GO:0000015         4        Down 0.961220866 0.9777157  5.17536e-01 6.38933e-01
## GO:0000016         0          Up         NaN       NaN  0.00000e+00 0.00000e+00
## ...              ...         ...         ...       ...          ...         ...
## GO:0000244        37          Up   0.6045212 0.6980781   0.00191857 0.004170815
## GO:0000245        34          Up   0.7146045 0.7788387   0.00011890 0.000312894
## GO:0000246         1          Up   0.1101072 0.2427364   0.11010723 0.183512058
## GO:0000247         1        Down   0.0075379 0.0415309   0.00753790 0.015383468
## GO:0000248         2        Down   0.9848608 0.9848608   1.00000000 1.000000000
## 
## $camera
## DataFrame with 100 rows and 4 columns
##             NumGenes   Direction    PValue       FDR
##            <integer> <character> <numeric> <numeric>
## GO:0000009         2          Up 0.8273552  0.932769
## GO:0000012        13        Down 0.0161556  0.261181
## GO:0000014         8          Up 0.6146707  0.895889
## GO:0000015         4        Down 0.9328063  0.972927
## GO:0000016         0          Up       NaN       NaN
## ...              ...         ...       ...       ...
## GO:0000244        37          Up  0.579835  0.895889
## GO:0000245        34          Up  0.875072  0.932769
## GO:0000246         1          Up  0.278953  0.733098
## GO:0000247         1        Down  0.392637  0.863806
## GO:0000248         2          Up  0.991322  0.996229

This returns the same type of output as runPrecomputed() but using different tests, i.e., ROAST, fry, CAMERA, ROMER. These tests typically provide more accurate type I error control than their pre-computed counterparts, as the former can use the expression values to account for the correlations between genes. This avoids the assumption of independent genes that is commonly made by the pre-computed tests.

Some of the runContrast() tests (namely, ROAST and fry()) are “self-contained” tests, where the null hypothesis is that none of the genes in the set are DE. This differs from the other “competitive” tests (CAMERA, ROMER and those in runPrecomputed()), where the null hypothesis is that there are no more/stronger DE genes in the set compared to those outside of the set. These two categories of tests answer different questions and must be chosen accordingly. To illustrate, consider a simple experiment where we treat a cell line with a drug, causing many pathways to be activated.

runContrast() effectively does a voom() analysis before applying the relevant tests. This means that many of the options in augere.de’s runVoom() can also be used here. For example, if we want to supply a custom design matrix:

output.dgst2 <- tempfile()
set.seed(200001)
res.dgst2 <- runContrast(
    se,
    sets,
    design=~dex,
    contrast="dexuntrt",
    output.dir=output.dgst2
) 

5 More output options

As shown above, we can set save.results=FALSE to instruct the various pipelines to not write results to disk. This is more efficient if we only want the results in the current R session.

Setting dry.run=TRUE will instruct each pipeline function to only generate the Rmarkdown report but not evaluate it. This can be helpful if further modifications to the report are intended.

output.dry <- tempfile()
set.seed(300000)
runPrecomputed(
    res.de$results[[1]],
    sets,
    signif.field="FDR",
    signif.threshold=0.05,
    rank.field="t",
    output.dir=output.dry,
    dry.run=TRUE
)
## NULL
list.files(output.dry)
## [1] "report.Rmd"

The generated report is not quite self-contained as runPrecomputed() doesn’t know how the result table was constructed. In the absence of this information, the report contains a placeholder stop() to remind the user to insert the relevant commands:

## We construct the data frame of precomputed statistics for each gene.
## This is typically from a differential expression analysis.
## 
## ```{r}
## tab <- local({ # augere.core input (1)
##     stop("insert commands to generate 'tab' here")
## })
## tab
## ```
## 
## We remove rows with `NA` values, e.g., because they correspond to low-abundance genes.

If we know how the object was generated, we can pass this information to runPrecomputed():

set.seed(300001)
runPrecomputed(
    wrapInput(res.de$results[[1]], commands=sprintf(
        # Pulling the result out of the directory saved by runVoom().
        "augere.core::readResult(%s)$x", deparse(file.path(output.de, "results", "de-1"))
    )),
    wrapInput(sets, commands="
library(org.Hs.eg.db)
sets <- select(org.Hs.eg.db, keys=keys(org.Hs.eg.db, 'GO'), keytype='GO', columns='ENSEMBL')
sets <- split(sets$ENSEMBL, sets$GO)
head(sets, 100)"), 
    signif.field="FDR",
    signif.threshold=0.05,
    rank.field="t",
    output.dir=output.dry,
    dry.run=TRUE
)
## NULL

This allows runPrecomputed() to insert the relevant commands into the report for full reproducibility:

## We construct the data frame of precomputed statistics for each gene.
## This is typically from a differential expression analysis.
## 
## ```{r}
## tab <- local({ # augere.core input (1)
##     augere.core::readResult("/tmp/RtmpRJN1zU/file1f86fa23a006af/results/de-1")$x
## })
## tab
## ```
## 
## We remove rows with `NA` values, e.g., because they correspond to low-abundance genes.

We can also attach gene set-specific annotation to the result tables, if the sets are represented as S4Vectors::List objects with annotations in their mcols().

library(GO.db)
annotated.sets <- List(sets)
mcols(annotated.sets)$name <- mapIds(
    GO.db,
    keys=names(annotated.sets),
    keytype="GOID",
    column="TERM"
)

output.anno.dir <- tempfile()
set.seed(300002)
res.anno <- runPrecomputed(
    res.de$results[[1]],
    annotated.sets,
    signif.field="FDR",
    signif.threshold=0.05,
    rank.field="t",
    save.results=FALSE,
    annotation="name",
    output.dir=output.anno.dir
)
# The desired columns of the 'mcols' are now included: 
res.anno$hypergeometric
## DataFrame with 100 rows and 5 columns
##                              name  NumGenes    NumSig    PValue       FDR
##                       <character> <integer> <integer> <numeric> <numeric>
## GO:0000009 alpha-1,6-mannosyltr..         2         0  1.000000         1
## GO:0000012 single strand break ..        13         1  0.869174         1
## GO:0000014 single-stranded DNA ..         8         2  0.326251         1
## GO:0000015 phosphopyruvate hydr..         4         0  1.000000         1
## GO:0000016       lactase activity         0         0        NA        NA
## ...                           ...       ...       ...       ...       ...
## GO:0000244 spliceosomal tri-snR..        37         1  0.996952         1
## GO:0000245 spliceosomal complex..        34         2  0.966982         1
## GO:0000246 Delta24(24-1) sterol..         1         0  1.000000         1
## GO:0000247 C-8 sterol isomerase..         1         0  1.000000         1
## GO:0000248 C-5 sterol desaturas..         2         0  1.000000         1

Session information

sessionInfo()
## R version 4.6.0 RC (2026-04-17 r89917)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.24-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.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] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] GO.db_3.23.1         augere.gsea_0.99.1   org.Hs.eg.db_3.23.1 
##  [4] AnnotationDbi_1.75.0 IRanges_2.47.0       S4Vectors_0.51.1    
##  [7] Biobase_2.73.1       BiocGenerics_0.59.0  generics_0.1.4      
## [10] limma_3.69.0         augere.de_0.99.2     BiocStyle_2.41.0    
## 
## loaded via a namespace (and not attached):
##   [1] DBI_1.3.0                   bitops_1.0-9               
##   [3] httr2_1.2.2                 BiasedUrn_2.0.12           
##   [5] biomaRt_2.69.0              rlang_1.2.0                
##   [7] magrittr_2.0.5              otel_0.2.0                 
##   [9] matrixStats_1.5.0           compiler_4.6.0             
##  [11] RSQLite_2.4.6               mgcv_1.9-4                 
##  [13] GenomicFeatures_1.65.0      png_0.1-9                  
##  [15] vctrs_0.7.3                 txdbmaker_1.9.0            
##  [17] stringr_1.6.0               pkgconfig_2.0.3            
##  [19] crayon_1.5.3                fastmap_1.2.0              
##  [21] dbplyr_2.5.2                magick_2.9.1               
##  [23] XVector_0.53.0              Rsamtools_2.29.0           
##  [25] rmarkdown_2.31              UCSC.utils_1.9.0           
##  [27] tinytex_0.59                bit_4.6.0                  
##  [29] xfun_0.57                   cachem_1.1.0               
##  [31] cigarillo_1.3.0             GenomeInfoDb_1.49.0        
##  [33] jsonlite_2.0.0              progress_1.2.3             
##  [35] blob_1.3.0                  rhdf5filters_1.25.0        
##  [37] DelayedArray_0.39.1         Rhdf5lib_2.1.0             
##  [39] BiocParallel_1.47.0         prettyunits_1.2.0          
##  [41] parallel_4.6.0              R6_2.6.1                   
##  [43] RColorBrewer_1.1-3          stringi_1.8.7              
##  [45] bslib_0.10.0                rtracklayer_1.73.0         
##  [47] GenomicRanges_1.65.0        jquerylib_0.1.4            
##  [49] Rcpp_1.1.1-1.1              Seqinfo_1.3.0              
##  [51] bookdown_0.46               SummarizedExperiment_1.43.0
##  [53] knitr_1.51                  BiocBaseUtils_1.15.0       
##  [55] splines_4.6.0               Matrix_1.7-5               
##  [57] tidyselect_1.2.1            dichromat_2.0-0.1          
##  [59] abind_1.4-8                 yaml_2.3.12                
##  [61] codetools_0.2-20            curl_7.1.0                 
##  [63] lattice_0.22-9              tibble_3.3.1               
##  [65] S7_0.2.2                    KEGGREST_1.53.0            
##  [67] evaluate_1.0.5              goseq_1.65.0               
##  [69] BiocFileCache_3.3.0         alabaster.schemas_1.13.0   
##  [71] Biostrings_2.81.1           pillar_1.11.1              
##  [73] BiocManager_1.30.27         filelock_1.0.3             
##  [75] MatrixGenerics_1.25.0       RCurl_1.98-1.18            
##  [77] ggplot2_4.0.3               hms_1.1.4                  
##  [79] scales_1.4.0                alabaster.base_1.13.0      
##  [81] glue_1.8.1                  alabaster.ranges_1.13.0    
##  [83] augere.core_0.99.3          alabaster.matrix_1.13.0    
##  [85] tools_4.6.0                 BiocIO_1.23.3              
##  [87] data.table_1.18.4           fgsea_1.39.0               
##  [89] locfit_1.5-9.12             GenomicAlignments_1.49.0   
##  [91] geneLenDataBase_1.49.0      XML_3.99-0.23              
##  [93] fastmatch_1.1-8             cowplot_1.2.0              
##  [95] rhdf5_2.57.0                grid_4.6.0                 
##  [97] edgeR_4.11.0                nlme_3.1-169               
##  [99] HDF5Array_1.41.0            restfulr_0.0.16            
## [101] cli_3.6.6                   rappdirs_0.3.4             
## [103] S4Arrays_1.13.0             dplyr_1.2.1                
## [105] gtable_0.3.6                alabaster.se_1.13.0        
## [107] sass_0.4.10                 digest_0.6.39              
## [109] SparseArray_1.13.2          farver_2.1.2               
## [111] rjson_0.2.23                memoise_2.0.1              
## [113] htmltools_0.5.9             lifecycle_1.0.5            
## [115] h5mread_1.5.0               httr_1.4.8                 
## [117] statmod_1.5.1               bit64_4.8.0