zenith performs gene set analysis on the result of differential expression using linear (mixed) modeling with dream by considering the correlation between gene expression traits. This package implements the camera method from the limma package proposed by Wu and Smyth (2012). zenith() is a simple extension of camera() to be compatible with linear (mixed) models implemented in dream().

Standard workflow

# Load packages
library(zenith)
library(edgeR)
library(variancePartition)
library(tweeDEseqCountData)
library(kableExtra)

# Load RNA-seq data from LCL's
data(pickrell)
geneCounts = exprs(pickrell.eset)
df_metadata = pData(pickrell.eset)

# Filter genes
# Note this is low coverage data, so just use as code example
dsgn = model.matrix(~ gender, df_metadata)
keep = filterByExpr(geneCounts, dsgn, min.count=5)

# Compute library size normalization
dge = DGEList(counts = geneCounts[keep,])
dge = calcNormFactors(dge)

# Estimate precision weights using voom
vobj = voomWithDreamWeights(dge, ~ gender, df_metadata)

# Apply dream analysis
fit = dream(vobj, ~ gender, df_metadata)
fit = eBayes(fit)

# Load get_MSigDB database, Hallmark genes
# use gene 'SYMBOL', or 'ENSEMBL' id
# use get_GeneOntology() to load Gene Ontology
msdb.gs = get_MSigDB("H", to="ENSEMBL")
   
# Run zenith analysis, and specific which coefficient to evaluate
res.gsa = zenith_gsa(fit, msdb.gs, 'gendermale', progressbar=FALSE )

# Show top gene sets: head(res.gsa)
kable_styling(kable(head(res.gsa), row.names=FALSE))
coef Geneset NGenes Correlation delta se p.less p.greater PValue Direction FDR
gendermale HALLMARK_TNFA_SIGNALING_VIA_NFKB 118 0.01 -0.9999389 0.1621467 0.0000000 1.0000000 0.0000000 Down 0.0000021
gendermale HALLMARK_CHOLESTEROL_HOMEOSTASIS 37 0.01 -1.0915363 0.2296170 0.0000055 0.9999945 0.0000110 Down 0.0002582
gendermale HALLMARK_INFLAMMATORY_RESPONSE 93 0.01 -0.7819360 0.1722474 0.0000120 0.9999880 0.0000241 Down 0.0003769
gendermale HALLMARK_IL2_STAT5_SIGNALING 105 0.01 -0.6493291 0.1672385 0.0001194 0.9998806 0.0002389 Down 0.0028066
gendermale HALLMARK_ESTROGEN_RESPONSE_LATE 81 0.01 -0.5424138 0.1789531 0.0017312 0.9982688 0.0034624 Down 0.0279921
gendermale HALLMARK_EPITHELIAL_MESENCHYMAL_TRANSITION 53 0.01 -0.6136470 0.2031832 0.0017867 0.9982133 0.0035735 Down 0.0279921
# for each cell type select 3 genesets with largest t-statistic
# and 1 geneset with the lowest
# Grey boxes indicate the gene set could not be evaluted because
#    to few genes were represented
plotZenithResults(res.gsa)

Session Info

## R Under development (unstable) (2025-03-13 r87965)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.21-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] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
##  [1] kableExtra_1.4.0          tweeDEseqCountData_1.45.0
##  [3] Biobase_2.67.0            BiocGenerics_0.53.6      
##  [5] generics_0.1.3            variancePartition_1.37.2 
##  [7] BiocParallel_1.41.5       ggplot2_3.5.1            
##  [9] edgeR_4.5.10              zenith_1.9.3             
## [11] limma_3.63.11             knitr_1.50               
## 
## loaded via a namespace (and not attached):
##   [1] rstudioapi_0.17.1           jsonlite_2.0.0             
##   [3] magrittr_2.0.3              farver_2.1.2               
##   [5] nloptr_2.2.1                rmarkdown_2.29             
##   [7] vctrs_0.6.5                 memoise_2.0.1              
##   [9] minqa_1.2.8                 RCurl_1.98-1.17            
##  [11] htmltools_0.5.8.1           S4Arrays_1.7.3             
##  [13] progress_1.2.3              broom_1.0.8                
##  [15] SparseArray_1.7.7           sass_0.4.9                 
##  [17] KernSmooth_2.23-26          bslib_0.9.0                
##  [19] pbkrtest_0.5.3              plyr_1.8.9                 
##  [21] cachem_1.1.0                lifecycle_1.0.4            
##  [23] iterators_1.0.14            pkgconfig_2.0.3            
##  [25] Matrix_1.7-3                R6_2.6.1                   
##  [27] fastmap_1.2.0               GenomeInfoDbData_1.2.14    
##  [29] rbibutils_2.3               MatrixGenerics_1.19.1      
##  [31] digest_0.6.37               numDeriv_2016.8-1.1        
##  [33] colorspace_2.1-1            AnnotationDbi_1.69.0       
##  [35] S4Vectors_0.45.4            GenomicRanges_1.59.1       
##  [37] RSQLite_2.3.9               labeling_0.4.3             
##  [39] httr_1.4.7                  abind_1.4-8                
##  [41] compiler_4.6.0              withr_3.0.2                
##  [43] bit64_4.6.0-1               aod_1.3.3                  
##  [45] backports_1.5.0             DBI_1.2.3                  
##  [47] gplots_3.2.0                MASS_7.3-65                
##  [49] DelayedArray_0.33.6         corpcor_1.6.10             
##  [51] gtools_3.9.5                caTools_1.18.3             
##  [53] tools_4.6.0                 msigdbr_10.0.1             
##  [55] remaCor_0.0.18              glue_1.8.0                 
##  [57] nlme_3.1-168                grid_4.6.0                 
##  [59] reshape2_1.4.4              gtable_0.3.6               
##  [61] tidyr_1.3.1                 hms_1.1.3                  
##  [63] xml2_1.3.8                  XVector_0.47.2             
##  [65] pillar_1.10.1               stringr_1.5.1              
##  [67] babelgene_22.9              splines_4.6.0              
##  [69] dplyr_1.1.4                 lattice_0.22-7             
##  [71] bit_4.6.0                   annotate_1.85.0            
##  [73] tidyselect_1.2.1            locfit_1.5-9.12            
##  [75] Biostrings_2.75.4           reformulas_0.4.0           
##  [77] IRanges_2.41.3              SummarizedExperiment_1.37.0
##  [79] svglite_2.1.3               RhpcBLASctl_0.23-42        
##  [81] stats4_4.6.0                xfun_0.52                  
##  [83] statmod_1.5.0               matrixStats_1.5.0          
##  [85] KEGGgraph_1.67.0            stringi_1.8.7              
##  [87] UCSC.utils_1.3.1            yaml_2.3.10                
##  [89] boot_1.3-31                 evaluate_1.0.3             
##  [91] codetools_0.2-20            tibble_3.2.1               
##  [93] Rgraphviz_2.51.9            graph_1.85.3               
##  [95] cli_3.6.4                   RcppParallel_5.1.10        
##  [97] systemfonts_1.2.1           xtable_1.8-4               
##  [99] Rdpack_2.6.3                munsell_0.5.1              
## [101] jquerylib_0.1.4             Rcpp_1.0.14                
## [103] GenomeInfoDb_1.43.4         zigg_0.0.2                 
## [105] EnvStats_3.0.0              png_0.1-8                  
## [107] XML_3.99-0.18               parallel_4.6.0             
## [109] Rfast_2.1.5.1               assertthat_0.2.1           
## [111] blob_1.2.4                  prettyunits_1.2.0          
## [113] bitops_1.0-9                lme4_1.1-37                
## [115] viridisLite_0.4.2           mvtnorm_1.3-3              
## [117] GSEABase_1.69.1             lmerTest_3.1-3             
## [119] scales_1.3.0                purrr_1.0.4                
## [121] crayon_1.5.3                fANCOVA_0.6-1              
## [123] rlang_1.1.5                 EnrichmentBrowser_2.37.0   
## [125] KEGGREST_1.47.0