This vignette contains a detailed tutorial on how to train a MOFA model using R. A concise template script can be found here. Many more examples on application of MOFA to various multi-omics data sets can be found here.
library(data.table)
library(MOFA2)
MOFA (and factor analysis models in general) are useful to uncover variation in complex data sets that contain multiple sources of heterogeneity. This requires a relatively large sample size (at least ~15 samples). In addition, MOFA needs the multi-modal measurements to be derived from the same samples. It is fine if you have samples that are missing some data modality, but there has to be a significant degree of matched measurements.
Proper normalisation of the data is critical. The model can handle three types of data: continuous (modelled with a gaussian likelihood), small counts (modelled with a Poisson likelihood) and binary measurements (modelled with a bernoulli likelihood). Non-gaussian likelihoods give non-optimal results, we recommend the user to apply data transformations to obtain continuous measurements. For example, for count-based data such as RNA-seq or ATAC-seq we recommend size factor normalisation + variance stabilisation (i.e. a log transformation).
It is strongly recommended that you select highly variable features (HVGs) per assay before fitting the model. This ensures a faster training and a more robust inference procedure. Also, for data modalities that have very different dimensionalities we suggest a stronger feature selection fort he bigger views, with the aim of reducing the feature imbalance between data modalities.
To create a MOFA object you need to specify three dimensions: samples, features and view(s). Optionally, a group can also be specified for each sample (no group structure by default). MOFA objects can be created from a wide range of input formats, including:
A list of matrices, where each entry corresponds to one view. Samples are stored in columns and features in rows.
Let’s simulate some data to start with
data <- make_example_data(
n_views = 2,
n_samples = 200,
n_features = 1000,
n_factors = 10
)[[1]]
lapply(data,dim)
## $view_1
## [1] 1000 200
##
## $view_2
## [1] 1000 200
Create the MOFA object:
MOFAobject <- create_mofa(data)
Plot the data overview
plot_data_overview(MOFAobject)
In case you are using the multi-group functionality, the groups can be specified
using the groups
argument as a vector with the group ID for each sample.
Keep in mind that the multi-group functionality is a rather advanced option that we
discourage for beginners. For more details on how the multi-group inference works,
read the FAQ section and
check this vignette.
N = ncol(data[[1]])
groups = c(rep("A",N/2), rep("B",N/2))
MOFAobject <- create_mofa(data, groups=groups)
Plot the data overview
plot_data_overview(MOFAobject)
A long data.frame with columns sample
, feature
, view
, group
(optional), value
might be the best format for complex data sets with multiple omics and potentially multiple groups of data. Also, there is no need to add rows that correspond to missing data:
filepath <- system.file("extdata", "test_data.RData", package = "MOFA2")
load(filepath)
head(dt)
## sample feature view value
## 1: sample_0_group_1 feature_0_view_0 view_0 2.08
## 2: sample_1_group_1 feature_0_view_0 view_0 0.01
## 3: sample_2_group_1 feature_0_view_0 view_0 -0.11
## 4: sample_3_group_1 feature_0_view_0 view_0 -0.82
## 5: sample_4_group_1 feature_0_view_0 view_0 -1.13
## 6: sample_5_group_1 feature_0_view_0 view_0 -0.25
Create the MOFA object
MOFAobject <- create_mofa(dt)
## Creating MOFA object from a data.frame...
print(MOFAobject)
## Untrained MOFA model with the following characteristics:
## Number of views: 2
## Views names: view_0 view_1
## Number of features (per view): 1000 1000
## Number of groups: 1
## Groups names: single_group
## Number of samples (per group): 100
##
Plot data overview
plot_data_overview(MOFAobject)
FALSE
FALSE
data_opts <- get_default_data_options(MOFAobject)
head(data_opts)
## $scale_views
## [1] FALSE
##
## $scale_groups
## [1] FALSE
##
## $center_groups
## [1] TRUE
##
## $use_float32
## [1] TRUE
##
## $views
## [1] "view_0" "view_1"
##
## $groups
## [1] "single_group"
FALSE
.TRUE
.TRUE
if using multiple groups.TRUE
if using multiple views.Only change the default model options if you are familiar with the underlying mathematical model.
model_opts <- get_default_model_options(MOFAobject)
model_opts$num_factors <- 10
head(model_opts)
## $likelihoods
## view_0 view_1
## "gaussian" "gaussian"
##
## $num_factors
## [1] 10
##
## $spikeslab_factors
## [1] FALSE
##
## $spikeslab_weights
## [1] FALSE
##
## $ard_factors
## [1] FALSE
##
## $ard_weights
## [1] TRUE
train_opts <- get_default_training_options(MOFAobject)
head(train_opts)
## $maxiter
## [1] 1000
##
## $convergence_mode
## [1] "fast"
##
## $drop_factor_threshold
## [1] -1
##
## $verbose
## [1] FALSE
##
## $startELBO
## [1] 1
##
## $freqELBO
## [1] 5
Prepare the MOFA object
MOFAobject <- prepare_mofa(
object = MOFAobject,
data_options = data_opts,
model_options = model_opts,
training_options = train_opts
)
Train the MOFA model. Remember that in this step the MOFA2
R package connets with the
mofapy2
Python package using reticulate
. This is the source of most problems when
running MOFA. See our FAQ section
if you have issues. The output is saved in the file specified as outfile
. If none is specified, the output is saved in a temporary location.
outfile = file.path(tempdir(),"model.hdf5")
MOFAobject.trained <- run_mofa(MOFAobject, outfile)
## Warning: Output file /tmp/RtmpRtXFUO/model.hdf5 already exists, it will be replaced
## Connecting to the mofapy2 python package using reticulate (use_basilisk = FALSE)...
## Please make sure to manually specify the right python binary when loading R with reticulate::use_python(..., force=TRUE) or the right conda environment with reticulate::use_condaenv(..., force=TRUE)
## If you prefer to let us automatically install a conda environment with 'mofapy2' installed using the 'basilisk' package, please use the argument 'use_basilisk = TRUE'
If everything is successful, you should observe an output analogous to the following:
######################################
## Training the model with seed 1 ##
######################################
Iteration 1: time=0.03, ELBO=-52650.68, deltaELBO=837116.802 (94.082647669%), Factors=10
(...)
Iteration 9: time=0.04, ELBO=-50114.43, deltaELBO=23.907 (0.002686924%), Factors=10
#######################
## Training finished ##
#######################
Saving model in `/var/folders/.../model.hdf5.../tmp/RtmpRtXFUO/model.hdf5.
This finishes the tutorial on how to train a MOFA object from R. To continue with the downstream analysis, follow this tutorial
Session Info
sessionInfo()
## R version 4.2.1 (2022-06-23)
## 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
## [3] LC_TIME=en_GB 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
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] data.table_1.14.4 pheatmap_1.0.12 forcats_0.5.2 stringr_1.4.1
## [5] dplyr_1.0.10 purrr_0.3.5 readr_2.1.3 tidyr_1.2.1
## [9] tibble_3.1.8 ggplot2_3.3.6 tidyverse_1.3.2 MOFA2_1.8.0
## [13] BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] googledrive_2.0.0 Rtsne_0.16 colorspace_2.0-3
## [4] ggsignif_0.6.4 ellipsis_0.3.2 rprojroot_2.0.3
## [7] fs_1.5.2 ggpubr_0.4.0 farver_2.1.1
## [10] ggrepel_0.9.1 fansi_1.0.3 mvtnorm_1.1-3
## [13] lubridate_1.8.0 xml2_1.3.3 codetools_0.2-18
## [16] splines_4.2.1 cachem_1.0.6 knitr_1.40
## [19] jsonlite_1.8.3 broom_1.0.1 dbplyr_2.2.1
## [22] png_0.1-7 uwot_0.1.14 HDF5Array_1.26.0
## [25] BiocManager_1.30.19 compiler_4.2.1 httr_1.4.4
## [28] basilisk_1.10.0 backports_1.4.1 assertthat_0.2.1
## [31] Matrix_1.5-1 fastmap_1.1.0 gargle_1.2.1
## [34] cli_3.4.1 htmltools_0.5.3 tools_4.2.1
## [37] gtable_0.3.1 glue_1.6.2 reshape2_1.4.4
## [40] Rcpp_1.0.9 carData_3.0-5 cellranger_1.1.0
## [43] jquerylib_0.1.4 vctrs_0.5.0 rhdf5filters_1.10.0
## [46] nlme_3.1-160 xfun_0.34 rvest_1.0.3
## [49] irlba_2.3.5.1 lifecycle_1.0.3 rstatix_0.7.0
## [52] googlesheets4_1.0.1 scales_1.2.1 basilisk.utils_1.10.0
## [55] hms_1.1.2 MatrixGenerics_1.10.0 parallel_4.2.1
## [58] rhdf5_2.42.0 RColorBrewer_1.1-3 yaml_2.3.6
## [61] reticulate_1.26 sass_0.4.2 reshape_0.8.9
## [64] stringi_1.7.8 highr_0.9 S4Vectors_0.36.0
## [67] corrplot_0.92 BiocGenerics_0.44.0 filelock_1.0.2
## [70] rlang_1.0.6 pkgconfig_2.0.3 matrixStats_0.62.0
## [73] evaluate_0.17 lattice_0.20-45 Rhdf5lib_1.20.0
## [76] labeling_0.4.2 cowplot_1.1.1 tidyselect_1.2.0
## [79] RcppAnnoy_0.0.20 here_1.0.1 GGally_2.1.2
## [82] plyr_1.8.7 magrittr_2.0.3 bookdown_0.29
## [85] R6_2.5.1 IRanges_2.32.0 magick_2.7.3
## [88] generics_0.1.3 DelayedArray_0.24.0 DBI_1.1.3
## [91] mgcv_1.8-41 pillar_1.8.1 haven_2.5.1
## [94] withr_2.5.0 abind_1.4-5 dir.expiry_1.6.0
## [97] modelr_0.1.9 crayon_1.5.2 car_3.1-1
## [100] utf8_1.2.2 tzdb_0.3.0 rmarkdown_2.17
## [103] grid_4.2.1 readxl_1.4.1 reprex_2.0.2
## [106] digest_0.6.30 stats4_4.2.1 munsell_0.5.0
## [109] bslib_0.4.0