Differential expression for repeated measures (dream) uses a linear model model to increase power and decrease false positives for RNA-seq datasets with multiple measurements per individual. The analysis fits seamlessly into the widely used workflow of limma/voom (Law et al. 2014). Dream uses a linear model model to increase power and decrease false positives for RNA-seq datasets with repeated measurements. Dream achieves this by combining multiple statistical concepts into a single statistical model. The model includes:

Dream also includes multi-threaded analysis across thousands of genes on a multi-core machine.

Standard RNA-seq processing

This tutorial assumes that the reader is familiar with the limma/voom workflow for RNA-seq. Process raw count data using limma/voom.

library('variancePartition')
library('edgeR')
library('BiocParallel')
data(varPartDEdata)

# filter genes by number of counts
isexpr = rowSums(cpm(countMatrix)>0.1) >= 5

# Standard usage of limma/voom
dge = DGEList( countMatrix[isexpr,] )
dge = calcNormFactors( dge )

# make this vignette faster by analyzing a subset of genes
dge = dge[1:1000,]

Limma Analysis

Limma has a built-in approach for analyzing repeated measures data using duplicateCorrelation(). The model can handle a single random effect, and forces the magnitude of the random effect to be the same across all genes.

# apply duplicateCorrelation is two rounds
design = model.matrix( ~ Disease, metadata)
vobj_tmp = voom( dge, design, plot=FALSE)
dupcor <- duplicateCorrelation(vobj_tmp,design,block=metadata$Individual)

# run voom considering the duplicateCorrelation results
# in order to compute more accurate precision weights
# Otherwise, use the results from the first voom run
vobj = voom( dge, design, plot=FALSE, block=metadata$Individual, correlation=dupcor$consensus)

# Estimate linear mixed model with a single variance component
# Fit the model for each gene, 
dupcor <- duplicateCorrelation(vobj, design, block=metadata$Individual)

# But this step uses only the genome-wide average for the random effect
fitDupCor <- lmFit(vobj, design, block=metadata$Individual, correlation=dupcor$consensus)

# Fit Empirical Bayes for moderated t-statistics
fitDupCor <- eBayes( fitDupCor )

Dream Analysis

The dream method replaces 4 core functions of limma with a linear mixed model.

  1. voomWithDreamWeights() replaces voom() to estimate precision weights
  2. dream() replaces lmFit() to estimate regression coefficients.
  3. variancePartition::eBayes() replaces limma::eBayes() to apply empircal Bayes shrinkage on linear mixed models.
  4. variancePartition::topTable() replaces limma::topTable() to give seamless access to results from dream().

For models with only fixed effects, variancePartition::eBayes(), and variancePartition::topTable() work seamlessly and give results equivalent to the limma functions with the same name. From the user perspective, the dream() workflow is the same as limma since the statistical differences are handled behind the scenes.

# Specify parallel processing parameters
# this is used implicitly by dream() to run in parallel
param = SnowParam(4, "SOCK", progressbar=TRUE)

# The variable to be tested must be a fixed effect
form <- ~ Disease + (1|Individual) 

# estimate weights using linear mixed model of dream
vobjDream = voomWithDreamWeights( dge, form, metadata, BPPARAM=param )

# Fit the dream model on each gene
# By default, uses the Satterthwaite approximation for the hypothesis test
fitmm = dream( vobjDream, form, metadata )
fitmm = eBayes(fitmm)
# Examine design matrix
head(fitmm$design, 3)
##           (Intercept) Disease1
## sample_01           1        0
## sample_02           1        0
## sample_03           1        0
# Get results of hypothesis test on coefficients of interest
topTable( fitmm, coef='Disease1', number=3 )
##                                   logFC  AveExpr        t      P.Value    adj.P.Val        B
## ENST00000283033.5 gene=TXNDC11 1.556233 3.567624 31.02068 3.099042e-21 3.099042e-18 38.15636
## ENST00000257181.9 gene=PRPF38A 1.380549 4.398270 21.59222 1.683871e-17 8.419354e-15 29.93271
## ENST00000525790.1 gene=TDRKH   1.508341 3.184931 19.51845 1.767605e-16 5.892017e-14 27.78202
##                                   z.std
## ENST00000283033.5 gene=TXNDC11 9.459336
## ENST00000257181.9 gene=PRPF38A 8.513757
## ENST00000525790.1 gene=TDRKH   8.236879

Since dream uses an estimated degrees of freedom value for each hypothsis test, the degrees of freedom is different for each gene here. Therefore, the t-statistics are not directly comparable since they have different degrees of freedom. In order to be able to compare test statistics, we report z.std which is the p-value transformed into a signed z-score. This can be used for downstream analysis.

Note that if a random effect is not specified, dream() automatically uses lmFit(), but the user must run eBayes() afterward.

Advanced hypothesis testing

Using contrasts to compare coefficients

You can also perform a hypothesis test of the difference between two or more coefficients by using a contrast matrix. The contrasts are evaluated at the time of the model fit and the results can be extracted with topTable(). This behaves like makeContrasts() and contrasts.fit() in limma.

Multiple contrasts can be evaluated at the same time, in order to save computation time. Make sure to inspect your contrast matrix to confirm it is testing what you intend.

form <- ~ 0 + DiseaseSubtype + Sex + (1|Individual) 

L = makeContrastsDream( form, metadata, 
  contrasts = c(compare2_1 = "DiseaseSubtype2 - DiseaseSubtype1", 
                compare1_0 = "DiseaseSubtype1 - DiseaseSubtype0"))

# Visualize contrast matrix
plotContrasts(L) 

# fit dream model with contrasts
fit = dream( vobjDream, form, metadata, L)
fit = eBayes(fit)

# get names of available coefficients and contrasts for testing
colnames(fit)
## [1] "compare2_1"      "compare1_0"      "DiseaseSubtype0" "DiseaseSubtype1" "DiseaseSubtype2"
## [6] "SexM"
# extract results from first contrast
topTable( fit, coef="compare2_1", number=3 )
##                                       logFC  AveExpr         t      P.Value  adj.P.Val          B
## ENST00000355624.3 gene=RAB11FIP2 -0.9493146 5.260280 -5.114589 3.698015e-05 0.03698015  0.1470281
## ENST00000593466.1 gene=DDA1      -1.7265709 3.901579 -3.700902 1.205688e-03 0.60284415 -1.2751588
## ENST00000200676.3 gene=CETP       1.4777422 3.723438  3.910146 3.397163e-03 0.99951347 -1.4383267
##                                      z.std
## ENST00000355624.3 gene=RAB11FIP2 -4.125580
## ENST00000593466.1 gene=DDA1      -3.237531
## ENST00000200676.3 gene=CETP       2.929309

Comparing multiple coefficients

So far contrasts have only involved the difference between two coefficients. But contrasts can also compare any linear combination of coefficients. Here, consider comparing DiseaseSubtype0 to the mean of DiseaseSubtype1 and DiseaseSubtype2. Note you can also customize the name of the contrast.

L2 = makeContrastsDream( form, metadata, contrasts = 
  c(Test1 = "DiseaseSubtype0 - (DiseaseSubtype1 + DiseaseSubtype2)/2"))

plotContrasts(L2)

# fit dream model to evaluate contrasts
fit = dream( vobjDream[1:10,], form, metadata, L=L2)
fit = eBayes(fit)

topTable(fit, coef="Test1", number=3)
##                                    logFC  AveExpr         t      P.Value    adj.P.Val        B
## ENST00000456159.1 gene=MET    -0.9830788 2.458926 -6.975125 3.529900e-08 3.529900e-07 8.662257
## ENST00000570099.1 gene=YPEL3  -0.9512918 2.063331 -6.648837 6.060468e-07 2.172997e-06 7.230068
## ENST00000418210.2 gene=TMEM64 -1.0343236 4.715367 -6.655755 6.518991e-07 2.172997e-06 7.694537
##                                   z.std
## ENST00000456159.1 gene=MET    -5.512889
## ENST00000570099.1 gene=YPEL3  -4.989280
## ENST00000418210.2 gene=TMEM64 -4.975171

Joint hypothesis test of multiple coefficients

Joint hypothesis testing of multiple coefficients at the same time can be performed by using an F-test. Just like in limma, the results can be extracted using topTable()

# extract results from first contrast
topTable( fit, coef=c("DiseaseSubtype2", "DiseaseSubtype1"), number=3 )
##                                DiseaseSubtype2 DiseaseSubtype1  AveExpr        F      P.Value
## ENST00000418210.2 gene=TMEM64         5.301001        5.211674 4.715367 744.4317 5.216147e-23
## ENST00000555834.1 gene=RPS6KL1        5.662699        5.719196 5.272063 657.9520 2.377852e-22
## ENST00000589123.1 gene=NFIC           6.545195        6.181023 5.855335 519.2464 4.314164e-21
##                                   adj.P.Val    F.std
## ENST00000418210.2 gene=TMEM64  5.216147e-22 51.30770
## ENST00000555834.1 gene=RPS6KL1 1.188926e-21 49.79067
## ENST00000589123.1 gene=NFIC    1.438055e-20 46.89238

Since dream uses an estimated degrees of freedom value for each hypothsis test, the degrees of freedom is different for each gene here. Therefore, the F-statistics are not directly comparable since they have different degrees of freedom. In order to be able to compare test statistics, we report F.std which is the p-value transformed into an F-statistic with \(df_1\) the number of coefficients tested and \(df_2=\infty\). This can be used for downstream analysis.

Small-sample method

For small datasets, the Kenward-Roger method can be more powerful. But it is substantially more computationally intensive.

fitmmKR = dream( vobjDream, form, metadata, ddf="Kenward-Roger")
fitmmKR = eBayes(fitmmKR)

variancePartition plot

Dream and variancePartition share the same underlying linear mixed model framework. A variancePartition analysis can indicate important variables that should be included as fixed or random effects in the dream analysis.

# Note: this could be run with either vobj from voom()
# or vobjDream from voomWithDreamWeights()
# The resuylts are similar
form = ~ (1|Individual) + (1|Disease)
vp = fitExtractVarPartModel( vobj, form, metadata)

plotVarPart( sortCols(vp))

Comparing p-values

Here we compare p-values from from dream() and duplicateCorrelation. In order to understand the empircal difference between them, we can plot the \(-\log_{10}\) p-values from each method.

# Compare p-values and make plot
p1 = topTable(fitDupCor, coef="Disease1", number=Inf, sort.by="none")$P.Value
p2 = topTable(fitmm, number=Inf, sort.by="none")$P.Value

plotCompareP( p1, p2, vp$Individual, dupcor$consensus)

The duplicateCorrelation method estimates a single variance term genome-wide even though the donor contribution of a particular gene can vary substantially from the genome-wide trend. Using a single value genome-wide for the within-donor variance can reduce power and increase the false positive rate in a particular, reproducible way. Let \(\tau^2_g\) be the value of the donor component for gene \(g\) and \(\bar{\tau}^2\) be the genome-wide mean. For genes where \(\tau^2_g>\bar{\tau}^2\), using \(\bar{\tau}^2\) under-corrects for the donor component so that it increases the false positive rate compared to using \(\tau^2_g\). Conversely, for genes where \(\tau^2_g<\bar{\tau}^2\), using \(\bar{\tau}^2\) over-corrects for the donor component so that it decreases power. Increasing sample size does not overcome this issue. The dream method overcomes this issue by using a \(\tau^2_g\).

Here, the \(-\log_{10}\) p-values from both methods are plotted and colored by the donor contribution estiamted by variancePartition. The green value indicates \(\bar{\tau}^2\), while red and blue indicate higher and lower values, respectively. When only one variance component is used and the contrast matrix is simple, the effect of using dream versus duplicateCorrelation is determined by the comparison of \(\tau^2_g\) to \(\bar{\tau}^2\):

dream can increase the \(-\log_{10}\) p-value for genes with a lower donor component (i.e. \(\tau^2_g<\bar{\tau}^2\)) and decrease \(-\log_{10}\) p-value for genes with a higher donor component (i.e. \(\tau^2_g>\bar{\tau}^2\))

Note that using more variance components or a more complicated contrast matrix can make the relationship more complicated.

Parallel processing

variancePartition functions including dream(), fitExtractVarPartModel() and fitVarPartModel() can take advange of multicore machines to speed up analysis. It uses the BiocParallel package to manage the parallelization.

# Request 4 cores, and enable the progress bar
# This is the ideal for Linux, OS X and Windows
param = SnowParam(4, "SOCK", progressbar=TRUE)
fitmm = dream( vobjDream, form, metadata, BPPARAM=param)

By default BPPARAM = SerialParam().

Session info

## R version 4.2.2 (2022-10-31)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
## 
## 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       
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] edgeR_3.40.2             pander_0.6.5             variancePartition_1.28.9
## [4] BiocParallel_1.32.6      limma_3.54.2             ggplot2_3.4.1           
## [7] knitr_1.42              
## 
## loaded via a namespace (and not attached):
##  [1] Biobase_2.58.0          sass_0.4.5              tidyr_1.3.0             jsonlite_1.8.4         
##  [5] splines_4.2.2           foreach_1.5.2           RhpcBLASctl_0.23-42     gtools_3.9.4           
##  [9] bslib_0.4.2             Rdpack_2.4              statmod_1.5.0           highr_0.10             
## [13] yaml_2.3.7              progress_1.2.2          numDeriv_2016.8-1.1     pillar_1.8.1           
## [17] backports_1.4.1         lattice_0.20-45         glue_1.6.2              RUnit_0.4.32           
## [21] digest_0.6.31           rbibutils_2.2.13        minqa_1.2.5             colorspace_2.1-0       
## [25] htmltools_0.5.4         Matrix_1.5-3            plyr_1.8.8              pkgconfig_2.0.3        
## [29] broom_1.0.4             purrr_1.0.1             mvtnorm_1.1-3           snow_0.4-4             
## [33] scales_1.2.1            lme4_1.1-32             tibble_3.2.0            farver_2.1.1           
## [37] generics_0.1.3          ellipsis_0.3.2          cachem_1.0.7            withr_2.5.0            
## [41] BiocGenerics_0.44.0     remaCor_0.0.11          pbkrtest_0.5.2          cli_3.6.0              
## [45] magrittr_2.0.3          crayon_1.5.2            evaluate_0.20           fansi_1.0.4            
## [49] doParallel_1.0.17       nlme_3.1-162            MASS_7.3-58.3           gplots_3.1.3           
## [53] tools_4.2.2             prettyunits_1.1.1       hms_1.1.2               lifecycle_1.0.3        
## [57] stringr_1.5.0           locfit_1.5-9.7          munsell_0.5.0           aod_1.3.2              
## [61] compiler_4.2.2          jquerylib_0.1.4         clusterGeneration_1.3.7 caTools_1.18.2         
## [65] rlang_1.1.0             grid_4.2.2              nloptr_2.0.3            iterators_1.0.14       
## [69] labeling_0.4.2          bitops_1.0-7            rmarkdown_2.20          boot_1.3-28.1          
## [73] lmerTest_3.1-3          gtable_0.3.2            codetools_0.2-19        reshape2_1.4.4         
## [77] R6_2.5.1                dplyr_1.1.0             fastmap_1.1.1           utf8_1.2.3             
## [81] KernSmooth_2.23-20      stringi_1.7.12          parallel_4.2.2          Rcpp_1.0.10            
## [85] vctrs_0.6.0             tidyselect_1.2.0        xfun_0.37

<>

References

Law, C. W., Y. Chen, W. Shi, and G. K. Smyth. 2014. “Voom: precision weights unlock linear model analysis tools for RNA-seq read counts.” Genome Biology 15 (2): R29. https://doi.org/10.1186/gb-2014-15-2-r29.