dream analysis
Differential expression testing with linear mixed models for repeated measures
Developed by Gabriel Hoffman
Run on 2024-06-12 13:01:04
Source:vignettes/dream.Rmd
dream.Rmd
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:
- flexible modeling of repeated measures gene expression data
- precision weights to model measurement error in RNA-seq
counts
- ability to model multiple random effects
- random effects estimated separately for each gene
- hypothesis testing for fixed effects in linear mixed models
- small sample size hypothesis test
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.
-
voomWithDreamWeights()
replacesvoom()
to estimate precision weights -
dream()
replaceslmFit()
to estimate regression coefficients.
-
variancePartition::eBayes()
replaceslimma::eBayes()
to apply empircal Bayes shrinkage on linear mixed models. -
variancePartition::topTable()
replaceslimma::topTable()
to give seamless access to results fromdream()
.
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
# For the hypothesis testing, by default,
# dream() uses the KR method for <= 20 samples,
# otherwise it uses the Satterthwaite approximation
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 32.66731 1.934507e-21 1.934507e-18 38.78853
## ENST00000257181.9 gene=PRPF38A 1.380549 4.398270 24.57867 1.493752e-18 7.468760e-16 32.48845
## ENST00000525790.1 gene=TDRKH 1.508341 3.184931 19.32367 3.712012e-16 1.025502e-13 27.06758
## z.std
## ENST00000283033.5 gene=TXNDC11 9.508490
## ENST00000257181.9 gene=PRPF38A 8.790139
## ENST00000525790.1 gene=TDRKH 8.147605
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.350356 0.0000224728 0.0224728 0.8858852
## ENST00000593466.1 gene=DDA1 -1.7265710 3.901579 -3.664795 0.0013569332 0.6784666 -1.5402278
## ENST00000200676.3 gene=CETP 1.4777422 3.723438 3.873141 0.0040308941 0.9994857 -1.8241804
## z.std
## ENST00000355624.3 gene=RAB11FIP2 -4.238790
## ENST00000593466.1 gene=DDA1 -3.203659
## ENST00000200676.3 gene=CETP 2.875734
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
## ENST00000418210.2 gene=TMEM64 -1.0343236 4.715367 -7.000080 3.300023e-07 3.054320e-06 8.484489
## ENST00000456159.1 gene=MET -0.9830788 2.458926 -6.062227 6.108639e-07 3.054320e-06 5.836340
## ENST00000555834.1 gene=RPS6KL1 -0.8954795 5.272063 -5.450341 3.957534e-06 1.319178e-05 3.997352
## z.std
## ENST00000418210.2 gene=TMEM64 -5.105453
## ENST00000456159.1 gene=MET -4.987750
## ENST00000555834.1 gene=RPS6KL1 -4.613600
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 823.4480 3.966577e-23
## ENST00000555834.1 gene=RPS6KL1 5.662699 5.719196 5.272063 791.4321 6.389688e-23
## ENST00000589123.1 gene=NFIC 6.545195 6.181023 5.855335 651.9138 6.549307e-22
## adj.P.Val F.std
## ENST00000418210.2 gene=TMEM64 3.194844e-22 51.58155
## ENST00000555834.1 gene=RPS6KL1 3.194844e-22 51.10477
## ENST00000589123.1 gene=NFIC 2.183102e-21 48.77751
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.
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.
- Specify parameters with the BPPARAM argument.
# 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.3.0 (2023-04-21)
## Platform: x86_64-apple-darwin22.4.0 (64-bit)
## Running under: macOS 14.5
##
## Matrix products: default
## BLAS: /Users/gabrielhoffman/prog/R-4.3.0/lib/libRblas.dylib
## LAPACK: /usr/local/Cellar/r/4.4.0_1/lib/R/lib/libRlapack.dylib; LAPACK version 3.12.0
##
## locale:
## [1] en_US.UTF-8/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
##
## time zone: America/New_York
## tzcode source: internal
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] knitr_1.43 edgeR_3.42.4 variancePartition_1.33.15
## [4] BiocParallel_1.34.2 limma_3.56.2 ggplot2_3.4.4
## [7] pander_0.6.5
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.0 farver_2.1.1 dplyr_1.1.2 bitops_1.0-7
## [5] fastmap_1.1.1 digest_0.6.33 lifecycle_1.0.3 statmod_1.5.0
## [9] magrittr_2.0.3 compiler_4.3.0 rlang_1.1.1 sass_0.4.9
## [13] tools_4.3.0 utf8_1.2.3 yaml_2.3.7 labeling_0.4.2
## [17] htmlwidgets_1.6.2 plyr_1.8.8 KernSmooth_2.23-21 withr_2.5.0
## [21] purrr_1.0.2 numDeriv_2016.8-1.1 BiocGenerics_0.46.0 desc_1.4.2
## [25] grid_4.3.0 aod_1.3.2 fansi_1.0.4 caTools_1.18.2
## [29] colorspace_2.1-0 scales_1.2.1 gtools_3.9.4 iterators_1.0.14
## [33] MASS_7.3-60 cli_3.6.1 mvtnorm_1.2-2 rmarkdown_2.22
## [37] ragg_1.2.5 generics_0.1.3 reshape2_1.4.4 minqa_1.2.5
## [41] cachem_1.0.8 stringr_1.5.0 splines_4.3.0 parallel_4.3.0
## [45] matrixStats_1.0.0 vctrs_0.6.3 boot_1.3-28.1 Matrix_1.6-5
## [49] jsonlite_1.8.5 pbkrtest_0.5.2 systemfonts_1.0.4 locfit_1.5-9.7
## [53] jquerylib_0.1.4 tidyr_1.3.1 snow_0.4-4 glue_1.6.2
## [57] pkgdown_2.0.9 nloptr_2.0.3 codetools_0.2-19 stringi_1.7.12
## [61] gtable_0.3.3 EnvStats_2.7.0 lme4_1.1-35.1 lmerTest_3.1-3
## [65] munsell_0.5.0 tibble_3.2.1 remaCor_0.0.18 pillar_1.9.0
## [69] htmltools_0.5.8.1 gplots_3.1.3 R6_2.5.1 textshaping_0.3.6
## [73] Rdpack_2.4 rprojroot_2.0.3 evaluate_0.21 lattice_0.21-8
## [77] Biobase_2.60.0 highr_0.10 rbibutils_2.2.13 backports_1.4.1
## [81] RhpcBLASctl_0.23-42 memoise_2.0.1 broom_1.0.5 fANCOVA_0.6-1
## [85] corpcor_1.6.10 bslib_0.7.0 Rcpp_1.0.11 nlme_3.1-162
## [89] xfun_0.39 fs_1.6.2 pkgconfig_2.0.3
<>