Results from the univariate regressions performed using can be combined in a post-processing step to perform multivariate hypothesis testing. In this example, we fit on transcript-level counts and then perform multivariate hypothesis testing by combining transcripts at the gene-level. This is done with the function.

Import transcript-level counts

Read in transcript counts from the package.

library(readr)
library(tximport)
library(tximportData)

# specify directory
path <- system.file("extdata", package = "tximportData")

# read sample meta-data
samples <- read.table(file.path(path, "samples.txt"), header = TRUE)
samples.ext <- read.table(file.path(path, "samples_extended.txt"), header = TRUE, sep = "\t")

# read assignment of transcripts to genes
# remove genes on the PAR, since these are present twice
tx2gene <- read_csv(file.path(path, "tx2gene.gencode.v27.csv"))
tx2gene <- tx2gene[grep("PAR_Y", tx2gene$GENEID, invert = TRUE), ]

# read transcript-level quatifictions
files <- file.path(path, "salmon", samples$run, "quant.sf.gz")
txi <- tximport(files, type = "salmon", txOut = TRUE)

# Create metadata simulating two conditions
sampleTable <- data.frame(condition = factor(rep(c("A", "B"), each = 3)))
rownames(sampleTable) <- paste0("Sample", 1:6)

Standard dream analysis

Perform standard analysis at the transcript-level

library(variancePartition)
library(edgeR)

# Prepare transcript-level reads
dge <- DGEList(txi$counts)
design <- model.matrix(~condition, data = sampleTable)
isexpr <- filterByExpr(dge, design)
dge <- dge[isexpr, ]
dge <- calcNormFactors(dge)

# Estimate precision weights
vobj <- voomWithDreamWeights(dge, ~condition, sampleTable)

# Fit regression model one transcript at a time
fit <- dream(vobj, ~condition, sampleTable)
fit <- eBayes(fit)

Multivariate analysis

Combine the transcript-level results at the gene-level. The mapping between transcript and gene is stored in as a list.

# Prepare transcript to gene mapping
# keep only transcripts present in vobj
# then convert to list with key GENEID and values TXNAMEs
keep <- tx2gene$TXNAME %in% rownames(vobj)
tx2gene.lst <- unstack(tx2gene[keep, ])

# Run multivariate test on entries in each feature set
# Default method is "FE.empirical", but use "FE" here to reduce runtime
res <- mvTest(fit, vobj, tx2gene.lst, coef = "conditionB", method = "FE")

# truncate gene names since they have version numbers
# ENST00000498289.5 -> ENST00000498289
res$ID.short <- gsub("\\..+", "", res$ID)

Gene set analysis

Perform gene set analysis using on the gene-level test statistics.

# must have zenith > v1.0.2
library(zenith)
library(GSEABase)

gs <- get_MSigDB("C1", to = "ENSEMBL")

df_gsa <- zenithPR_gsa(res$stat, res$ID.short, gs, inter.gene.cor = .05)

head(df_gsa)
##                 NGenes Correlation     delta       se      p.less    p.greater       PValue
## M14982_chr7p13      25        0.05  7.906128 2.082903 0.999926118 7.388155e-05 0.0001477631
## M5824_chr11p13      30        0.05 -6.028272 2.007216 0.001337443 9.986626e-01 0.0026748867
## M7314_chr4p14       25        0.05 -5.077180 2.084132 0.007428521 9.925715e-01 0.0148570424
## M3783_chr2q37       71        0.05  3.749809 1.768579 0.982999226 1.700077e-02 0.0340015489
## M6742_chr2q36       20        0.05 -4.344151 2.194043 0.023861861 9.761381e-01 0.0477237227
## M14543_chr18q22     17        0.05  4.325680 2.286441 0.970737601 2.926240e-02 0.0585247982
##                 Direction        FDR         Geneset     coef
## M14982_chr7p13         Up 0.03664525  M14982_chr7p13 zenithPR
## M5824_chr11p13       Down 0.33168596  M5824_chr11p13 zenithPR
## M7314_chr4p14        Down 0.99873228   M7314_chr4p14 zenithPR
## M3783_chr2q37          Up 0.99873228   M3783_chr2q37 zenithPR
## M6742_chr2q36        Down 0.99873228   M6742_chr2q36 zenithPR
## M14543_chr18q22        Up 0.99873228 M14543_chr18q22 zenithPR

Session info

## R Under development (unstable) (2024-10-21 r87258)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 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
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C               LC_TIME=en_GB             
##  [4] LC_COLLATE=C               LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                  LC_ADDRESS=C              
## [10] LC_TELEPHONE=C             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   base     
## 
## other attached packages:
##  [1] org.Hs.eg.db_3.20.0      msigdbr_7.5.1            GSEABase_1.69.0         
##  [4] graph_1.85.0             annotate_1.85.0          XML_3.99-0.17           
##  [7] AnnotationDbi_1.69.0     IRanges_2.41.0           S4Vectors_0.45.1        
## [10] Biobase_2.67.0           BiocGenerics_0.53.2      generics_0.1.3          
## [13] zenith_1.9.0             tximportData_1.35.0      tximport_1.35.0         
## [16] readr_2.1.5              edgeR_4.5.0              pander_0.6.5            
## [19] variancePartition_1.37.1 BiocParallel_1.41.0      limma_3.63.2            
## [22] ggplot2_3.5.1            knitr_1.49              
## 
## loaded via a namespace (and not attached):
##   [1] jsonlite_1.8.9              magrittr_2.0.3              farver_2.1.2               
##   [4] nloptr_2.1.1                rmarkdown_2.29              zlibbioc_1.53.0            
##   [7] vctrs_0.6.5                 memoise_2.0.1               minqa_1.2.8                
##  [10] RCurl_1.98-1.16             progress_1.2.3              S4Arrays_1.7.1             
##  [13] htmltools_0.5.8.1           curl_6.0.0                  broom_1.0.7                
##  [16] SparseArray_1.7.1           sass_0.4.9                  KernSmooth_2.23-24         
##  [19] bslib_0.8.0                 pbkrtest_0.5.3              plyr_1.8.9                 
##  [22] cachem_1.1.0                lifecycle_1.0.4             iterators_1.0.14           
##  [25] pkgconfig_2.0.3             Matrix_1.7-1                R6_2.5.1                   
##  [28] fastmap_1.2.0               GenomeInfoDbData_1.2.13     rbibutils_2.3              
##  [31] MatrixGenerics_1.19.0       digest_0.6.37               numDeriv_2016.8-1.1        
##  [34] colorspace_2.1-1            GenomicRanges_1.59.0        RSQLite_2.3.7              
##  [37] filelock_1.0.3              labeling_0.4.3              RcppZiggurat_0.1.6         
##  [40] fansi_1.0.6                 abind_1.4-8                 httr_1.4.7                 
##  [43] compiler_4.5.0              bit64_4.5.2                 aod_1.3.3                  
##  [46] withr_3.0.2                 backports_1.5.0             DBI_1.2.3                  
##  [49] gplots_3.2.0                MASS_7.3-61                 DelayedArray_0.33.1        
##  [52] corpcor_1.6.10              gtools_3.9.5                caTools_1.18.3             
##  [55] tools_4.5.0                 remaCor_0.0.18              glue_1.8.0                 
##  [58] nlme_3.1-166                grid_4.5.0                  reshape2_1.4.4             
##  [61] snow_0.4-4                  gtable_0.3.6                tzdb_0.4.0                 
##  [64] tidyr_1.3.1                 hms_1.1.3                   utf8_1.2.4                 
##  [67] XVector_0.47.0              pillar_1.9.0                stringr_1.5.1              
##  [70] babelgene_22.9              vroom_1.6.5                 splines_4.5.0              
##  [73] dplyr_1.1.4                 BiocFileCache_2.15.0        lattice_0.22-6             
##  [76] bit_4.5.0                   tidyselect_1.2.1            locfit_1.5-9.10            
##  [79] Biostrings_2.75.1           SummarizedExperiment_1.37.0 RhpcBLASctl_0.23-42        
##  [82] xfun_0.49                   statmod_1.5.0               matrixStats_1.4.1          
##  [85] KEGGgraph_1.67.0            stringi_1.8.4               UCSC.utils_1.3.0           
##  [88] yaml_2.3.10                 boot_1.3-31                 evaluate_1.0.1             
##  [91] codetools_0.2-20            archive_1.1.10              tibble_3.2.1               
##  [94] Rgraphviz_2.51.0            cli_3.6.3                   RcppParallel_5.1.9         
##  [97] xtable_1.8-4                Rdpack_2.6.1                munsell_0.5.1              
## [100] jquerylib_0.1.4             Rcpp_1.0.13-1               GenomeInfoDb_1.43.0        
## [103] EnvStats_3.0.0              dbplyr_2.5.0                png_0.1-8                  
## [106] Rfast_2.1.0                 parallel_4.5.0              blob_1.2.4                 
## [109] prettyunits_1.2.0           bitops_1.0-9                lme4_1.1-35.5              
## [112] mvtnorm_1.3-2               lmerTest_3.1-3              scales_1.3.0               
## [115] purrr_1.0.2                 crayon_1.5.3                fANCOVA_0.6-1              
## [118] rlang_1.1.4                 EnrichmentBrowser_2.37.0    KEGGREST_1.47.0

<>

References