mfa
is an R package for fitting a Bayesian mixture of
factor analysers to infer developmental trajectories with bifurcations
from single-cell gene expression data. It is able to jointly infer
pseudotimes, branching, and genes differentially regulated across
branches using a generative, Bayesian hierarchical model. Inference is
performed using fast Gibbs sampling.
mfa
can be installed in one of two ways:
We first create some synthetic data for 100 cells and 40 genes
calling the mfa
function create_synthetic
.
This returns a list with gene expression, pseudotime, branch allocation,
and various parameter estimates:
## List of 7
## $ X : num [1:100, 1:40] 13.35 5.37 10.26 15.11 14.43 ...
## ..- attr(*, "dimnames")=List of 2
## .. ..$ : chr [1:100] "cell1" "cell2" "cell3" "cell4" ...
## .. ..$ : chr [1:40] "feature1" "feature2" "feature3" "feature4" ...
## $ branch : int [1:100] 0 1 1 0 0 0 0 0 1 1 ...
## $ pst : num [1:100] 0.864 0.118 0.207 0.705 0.456 ...
## $ k : num [1:40, 1:2] 8.91 6.86 7.62 8.71 5.43 ...
## $ phi : num [1:40, 1:2] 7.39 9.45 8.6 9.5 9.8 ...
## $ delta : num [1:40, 1:2] 0.118 0.392 0.395 0.448 0.277 ...
## $ p_transient: num 0
## NULL
We can then PCA and put into a tidy format:
df_synth <- as_data_frame(prcomp(synth$X)$x[,1:2]) %>%
mutate(pseudotime = synth$pst,
branch = factor(synth$branch))
and have a look at a PCA representation, coloured by both pseudotime and branch allocation:
mfa
The input to mfa
is either an ExpressionSet
(e.g. from using the package Scater)
or a cell-by-gene expression matrix. If an ExpressionSet
is
provided then the values in the exprs
slot are used for
gene expression.
We invoke mfa
with a call to the mfa(...)
function. Depending on the size of the dataset and number of MCMC
iterations used, this may take some time:
## MFA fit with
## 100 cells and 40 genes
## ( 2000 iterations )
Particular care must be paid to the initialisation of the
pseudotimes: by default they are initialised to the first principal
component, though if the researcher suspects (based on plotting marker
genes) that the trajectory corresponds to a different PC, this can be
set using the pc_initialise
argument.
As in any MCMC analysis, basic care is needed to make sure the samples have converged to something resembling the stationary distribution (see e.g. Cowles and Carlin (1996) for a full discussion).
For a quick summary of these, mfa
provides two
functions: plot_mfa_trace
and
plot_mfa_autocorr
for quick plotting of the trace and
autocorrelation of the posterior log-likelihood:
We can extract posterior mean estimates along with credible intervals
using the summary
function:
## # A tibble: 6 × 5
## pseudotime branch branch_certainty pseudotime_lower pseudotime_upper
## <dbl> <fct> <dbl> <dbl> <dbl>
## 1 -1.09 2 1 -1.24 -0.950
## 2 1.19 1 1 0.953 1.58
## 3 0.752 1 1 0.427 1.03
## 4 -0.641 2 1 -0.735 -0.497
## 5 0.138 2 1 0.00421 0.228
## 6 -0.702 2 1 -0.811 -0.572
This has six entries:
pseudotime
The MAP pseudotime estimatebranch
The MAP branch estimatebranch_certainty
The proportion of MCMC traces (after
burn-in) for which the cell was assigned to the MAP branchpseudotime_lower
and pseudotime_upper
: the
lower and upper 95% highest-probability-density posterior credible
intervalsWe can compare the inferred pseudotimes to the true values:
qplot(synth$pst, ms$pseudotime, color = factor(synth$branch)) +
xlab('True pseudotime') + ylab('Inferred pseudotime') +
scale_color_discrete(name = 'True\nbranch')
And we can equivalently plot the PCA representation coloured by MAP branch:
A unique part of this model is that through an ARD-like prior
structure on the loading matrices we can automatically infer which genes
are involved in the bifurcation process. For a quick-and-dirty look we
can use the plot_chi
function, where larger values of
inverse-chi imply the gene is associated with the bifurcation:
To calculate the MAP values for chi we can call the
calculate_chi
function, which returns a
data_frame
with the feature names and values:
## # A tibble: 6 × 2
## feature chi_map
## <chr> <dbl>
## 1 feature1 0.433
## 2 feature2 0.850
## 3 feature3 1.63
## 4 feature4 1.16
## 5 feature5 0.387
## 6 feature6 0.577
mfa
classA call to mfa(...)
returns an mfa
object
that contains all the information about the dataset and the MCMC
inference performed. Note that it does not contain a copy of
the original data. We can see the structure by calling str
on an mfa
object:
## List of 10
## $ traces :List of 10
## $ iter : num 2000
## $ thin : num 1
## $ burn : num 1000
## $ b : num 2
## $ collapse : logi FALSE
## $ N : int 100
## $ G : int 40
## $ feature_names: chr [1:40] "feature1" "feature2" "feature3" "feature4" ...
## $ cell_names : chr [1:100] "cell1" "cell2" "cell3" "cell4" ...
## - attr(*, "class")= chr "mfa"
This contains the following slots:
traces
- the raw MCMC traces (discussed in following
section)iter
- the number of MCMC iterationsthin
- the thinning of the MCMC chainburn
- the number of MCMC iterations thrown away as
burn-inb
- the number of branches modelledcollapse
- whether collapsed Gibbs sampling was
implementedN
- the number of cellsG
- the number of features (e.g. genes)feature_names
- the names of the features
(e.g. genes)cell_names
- the names of the cellsMCMC traces can be accessed through the traces
slot of
an mfa
object. This gives a list with an element for each
variable, along with the log-likelihood:
## [1] "tau_trace" "gamma_trace" "pst_trace"
## [4] "theta_trace" "lambda_theta_trace" "chi_trace"
## [7] "eta_trace" "k_trace" "c_trace"
## [10] "lp_trace"
For non-branch-specific variables this is simply a matrix. For example, for the variable τ is just an interation-by-gene matrix:
## num [1:1000, 1:40] 4.94 5.48 7.59 6.05 6.9 ...
## - attr(*, "dimnames")=List of 2
## ..$ : NULL
## ..$ : chr [1:40] "tau[1]" "tau[2]" "tau[3]" "tau[4]" ...
We can easily get the posterior mean by calling
colMeans
. More fancy posterior density estimation can be
perfomed using the MCMCglmm
package, such as
posterior.mode(...)
for MAP estimation (though in practice
this is often similar to posterior mean). We can estimate posterior
intervals using the HPDInterval(...)
function from the
coda
package (note that traces must be converted to
coda
objects before calling either of these).
Some variables are branch dependent, meaning the traces returned are
arrays (or tensors in fashionable speak) that have dimension
iteration x gene x branch
. An example is the k variable:
## num [1:1000, 1:40, 1:2] -1.101 -1.033 -1.135 -1.184 -0.984 ...
To get posterior means (or modes, or intervals) we then need to use
the apply
function to iterate over the branches. To find
the posterior means of k
, we then call
## num [1:40, 1:2] -1.11 -1.12 -1.11 -1.06 -1.17 ...
This returns a gene-by-branch matrix of posterior estimates.
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 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: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] dplyr_1.1.4 ggplot2_3.5.1 mfa_1.29.0 BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] gtable_0.3.6 tensorA_0.36.2.1 xfun_0.48
## [4] bslib_0.8.0 Biobase_2.67.0 GGally_2.2.1
## [7] lattice_0.22-6 vctrs_0.6.5 tools_4.4.1
## [10] generics_0.1.3 parallel_4.4.1 tibble_3.2.1
## [13] fansi_1.0.6 highr_0.11 pkgconfig_2.0.3
## [16] Matrix_1.7-1 RColorBrewer_1.1-3 ggmcmc_1.5.1.1
## [19] lifecycle_1.0.4 cubature_2.1.1 compiler_4.4.1
## [22] farver_2.1.2 MatrixModels_0.5-3 mcmc_0.9-8
## [25] munsell_0.5.1 codetools_0.2-20 SparseM_1.84-2
## [28] quantreg_5.99 htmltools_0.5.8.1 sys_3.4.3
## [31] buildtools_1.0.0 sass_0.4.9 yaml_2.3.10
## [34] tidyr_1.3.1 pillar_1.9.0 jquerylib_0.1.4
## [37] MASS_7.3-61 cachem_1.1.0 nlme_3.1-166
## [40] ggstats_0.7.0 tidyselect_1.2.1 digest_0.6.37
## [43] purrr_1.0.2 maketools_1.3.1 labeling_0.4.3
## [46] splines_4.4.1 fastmap_1.2.0 grid_4.4.1
## [49] colorspace_2.1-1 cli_3.6.3 magrittr_2.0.3
## [52] survival_3.7-0 utf8_1.2.4 ape_5.8
## [55] corpcor_1.6.10 withr_3.0.2 scales_1.3.0
## [58] MCMCglmm_2.36 rmarkdown_2.28 coda_0.19-4.1
## [61] evaluate_1.0.1 knitr_1.48 rlang_1.1.4
## [64] MCMCpack_1.7-1 Rcpp_1.0.13 glue_1.8.0
## [67] BiocManager_1.30.25 BiocGenerics_0.53.0 jsonlite_1.8.9
## [70] R6_2.5.1 plyr_1.8.9