dreamlet() evaluates precision-weighted linear (mixed)
models on each gene that passes standard filters. The linear mixed model
used by dream() can be a little fragile for small sample
sizes and correlated covariates. dreamlet() runs
variancePartition::dream() in the backend for each cell
cluster. dream() reports model failures for each cell
cluster and dreamlet() reports these failures to the user.
dreamlet() returns all successful model
fits to be used for downstream analysis.
Due to a recent bug in the dependency
Matrix package, all random effects models may fail for
technical reasons.
This case can be detected and resolve as follows:
library(lme4)
# Fit simple mixed model
lmer(Reaction ~ (1 | Subject), sleepstudy)
# Error in initializePtr() :
# function 'chm_factor_ldetL2' not provided by package 'Matrix'
This error indicates incompatible installs of Matrix and
lme4.
This can be solved with
install.packages("lme4", type = "source")
followed by restarting R.
The most common issue is when dreamlet() analysis
succeeds for most genes, but a handful of genes fail in each cell
cluster. These genes can fail if the iterative process of fitting the
linear mixed model does not converge, or if the estimated covariance
matrix that is supposed be positive definite has an eigen-value that is
negative or too close to zero due to rounding errors in floating point
arithmetic.
In these cases, dreamlet() stores a summary of these
failures for all cell clusters that is accessible with
details(). Specific failure messages for each cell cluster
and gene can be extracted using seeErrors()
Here we demonstrate how dreamlet() handles model
failures:
library(dreamlet)
library(muscat)
library(SingleCellExperiment)
data(example_sce)
# create pseudobulk for each sample and cell cluster
pb <- aggregateToPseudoBulk(example_sce,
assay = "counts",
cluster_id = "cluster_id",
sample_id = "sample_id",
verbose = FALSE
)
# voom-style normalization for each cell cluster
res.proc <- processAssays(
pb[1:300, ],
~group_id
)
# Redundant formula
# This example is an extreme example of redundancy
# but more subtle cases often show up in real data
form <- ~ group_id + (1 | group_id)
# fit dreamlet model
res.dl <- dreamlet(res.proc, form)
## B cells...7.9 secs
## CD14+ Monocytes...10 secs
## CD4 T cells...9 secs
## CD8 T cells...4.4 secs
## FCGR3A+ Monocytes...11 secs
##
## Of 1,062 models fit across all assays, 96.2% failed
# summary of models
res.dl
## class: dreamletResult
## assays(5): B cells CD14+ Monocytes CD4 T cells CD8 T cells FCGR3A+ Monocytes
## Genes:
## min: 3
## max: 11
## details(7): assay n_retain ... n_errors error_initial
## coefNames(2): (Intercept) group_idstim
##
## Of 1,062 models fit across all assays, 96.2% failed
# summary of models for each cell cluster
details(res.dl)
## assay n_retain formula formDropsTerms n_genes n_errors error_initial
## 1 B cells 4 ~group_id + (1 | group_id) FALSE 201 190 FALSE
## 2 CD14+ Monocytes 4 ~group_id + (1 | group_id) FALSE 269 263 FALSE
## 3 CD4 T cells 4 ~group_id + (1 | group_id) FALSE 216 207 FALSE
## 4 CD8 T cells 4 ~group_id + (1 | group_id) FALSE 118 115 FALSE
## 5 FCGR3A+ Monocytes 4 ~group_id + (1 | group_id) FALSE 258 247 FALSE
# Extract errors as a tibble
seeErrors(res.dl)
## # A tibble: 1,022 × 3
## assay feature errorText
## <chr> <chr> <chr>
## 1 B cells ISG15 "Error in lmerTest:::as_lmerModLT(model, devfun, tol = tol)…
## 2 B cells AURKAIP1 "Error in lmerTest:::as_lmerModLT(model, devfun, tol = tol)…
# Extract gene-level errors for B cells
clustID <- "B cells"
attr(res.dl[[clustID]])[1:2]
## ISG15
## "Error in lmerTest:::as_lmerModLT(model, devfun, tol = tol): (converted from warning)
## Model may not have converged with 1 eigenvalue close to zero: 3.2e-11\n"
##
## AURKAIP1
## "Error in lmerTest:::as_lmerModLT(model, devfun, tol = tol): (converted from warning)
## Model may not have converged with 1 eigenvalue close to zero: -2.9e-11\n"
Some model failures affect every gene in an assay.
dreamlet() tests an initial model fit before analysis of
each gene. If the initial fit fails for a cell cluster, then
details(res.dl)$error_initial will be
TRUE.
## R version 4.3.0 (2023-04-21)
## Platform: x86_64-apple-darwin22.4.0 (64-bit)
## Running under: macOS 14.2.1
##
## Matrix products: default
## BLAS: /Users/gabrielhoffman/prog/R-4.3.0/lib/libRblas.dylib
## LAPACK: /usr/local/Cellar/r/4.3.0_1/lib/R/lib/libRlapack.dylib; LAPACK version 3.11.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 datasets utils methods base
##
## loaded via a namespace (and not attached):
## [1] vctrs_0.6.3 cli_3.6.1 knitr_1.43 rlang_1.1.1
## [5] xfun_0.39 processx_3.8.1 purrr_1.0.2 styler_1.10.1
## [9] jsonlite_1.8.5 glue_1.6.2 clipr_0.8.0 htmltools_0.5.5
## [13] ps_1.7.5 sass_0.4.6 fansi_1.0.4 rmarkdown_2.22
## [17] R.cache_0.16.0 jquerylib_0.1.4 evaluate_0.21 tibble_3.2.1
## [21] fastmap_1.1.1 yaml_2.3.7 lifecycle_1.0.3 compiler_4.3.0
## [25] fs_1.6.2 pkgconfig_2.0.3 rstudioapi_0.15.0 R.oo_1.25.0
## [29] R.utils_2.12.2 digest_0.6.33 R6_2.5.1 reprex_2.0.2
## [33] utf8_1.2.3 pillar_1.9.0 callr_3.7.3 magrittr_2.0.3
## [37] bslib_0.4.2 R.methodsS3_1.8.2 tools_4.3.0 withr_2.5.0
## [41] cachem_1.0.8