library(MOFA2)
library(tidyverse)
library(pheatmap)

1 Temporal data: Simulate an example data set

To illustrate the MEFISTO method in MOFA2 we simulate a small example data set with 4 different views and one covariates defining a timeline using make_example_data. The simulation is based on 4 factors, two of which vary smoothly along the covariate (with different lengthscales) and two are independent of the covariate.

set.seed(2020)

# set number of samples and time points
N <- 200
time <- seq(0,1,length.out = N)

# generate example data
dd <- make_example_data(sample_cov = time, n_samples = N,
                        n_factors = 4, n_features = 200, n_views = 4,
                        lscales = c(0.5, 0.2, 0, 0))
# input data
data <- dd$data

# covariate matrix with samples in columns
time <- dd$sample_cov
rownames(time) <- "time"

Let’s have a look at the simulated latent temporal processes, which we want to recover:

df <- data.frame(dd$Z, t(time))
df <- gather(df, key = "factor", value = "value", starts_with("simulated_factor"))
ggplot(df, aes(x = time, y = value)) + geom_point() + facet_grid(~factor)

2 MEFISTO framework

Using the MEFISTO framework is very similar to using MOFA2. In addition to the omics data, however, we now additionally specify the time points for each sample. If you are not familiar with the MOFA2 framework, it might be helpful to have a look at MOFA2 tutorials first.

2.1 Create a MOFA object with covariates

To create the MOFA object we need to specify the training data and the covariates for pattern detection and inference of smooth factors. Here, sample_cov is a matrix with samples in columns and one row containing the time points. The sample order must match the order in data columns. Alternatively, a data frame can be provided containing one sample columns with samples names matching the sample names in the data.

First, we start by creating a standard MOFA model.

sm <- create_mofa(data = dd$data)
## Creating MOFA object from a list of matrices (features as rows, sample as columns)...

Now, we can add the additional temporal covariate, that we want to use for training.

sm <- set_covariates(sm, covariates = time)
sm
## Untrained MEFISTO model with the following characteristics: 
##  Number of views: 4 
##  Views names: view_1 view_2 view_3 view_4 
##  Number of features (per view): 200 200 200 200 
##  Number of groups: 1 
##  Groups names: group1 
##  Number of samples (per group): 200 
##  Number of covariates per sample: 1 
## 

We now successfully created a MOFA object that contains 4 views, 1 group and 1 covariate giving the time point for each sample.

2.2 Prepare a MOFA object

Before training, we can specify various options for the model, the training and the data preprocessing. If no options are specified, the model will use the default options. See also get_default_data_options, get_default_model_options and get_default_training_options to have a look at the defaults and change them where required. For illustration, we only use a small number of iterations.

Importantly, to activate the use of the covariate for a functional decomposition (MEFISTO) we now additionally to the standard MOFA options need to specify mefisto_options. For this you can just use the default options (get_default_mefisto_options), unless you want to make use of advanced options such as alignment across groups.

data_opts <- get_default_data_options(sm)

model_opts <- get_default_model_options(sm)
model_opts$num_factors <- 4

train_opts <- get_default_training_options(sm)
train_opts$maxiter <- 100

mefisto_opts <- get_default_mefisto_options(sm)

sm <- prepare_mofa(sm, model_options = model_opts,
                   mefisto_options = mefisto_opts,
                   training_options = train_opts,
                   data_options = data_opts)

2.3 Run MOFA

Now, the MOFA object is ready for training. Using run_mofa we can fit the model, which 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")
sm <- run_mofa(sm, outfile)

2.4 Down-stream analysis

2.4.1 Variance explained per factor

Using plot_variance_explained we can explore which factor is active in which view. plot_factor_cor shows us whether the factors are correlated.

plot_variance_explained(sm)

r <- plot_factor_cor(sm)

2.4.2 Relate factors to the covariate

The MOFA model has learnt scale parameters for each factor, which give us an indication of the smoothness per factor along the covariate (here time) and are between 0 and 1. A scale of 0 means that the factor captures variation independent of time, a value close to 1 tells us that this factor varys very smoothly along time.

get_scales(sm)
##    Factor1    Factor2    Factor3    Factor4 
## 0.00000000 0.01198385 0.99892405 0.99496226

In this example, we find two factors that are non-smooth and two smooth factors. Using plot_factors_vs_cov we can plot the factors along the time line, where we can distinguish smooth and non smooth variation along time.

plot_factors_vs_cov(sm, color_by = "time")

For more customized plots, we can extract the underlying data containing the factor and covariate values for each sample.

df <- plot_factors_vs_cov(sm, color_by = "time",
                    legend = FALSE, return_data = TRUE)
head(df)
##      sample covariate value.covariate  factor value.factor  group   color_by
## 1  sample_1      time      0.00000000 Factor1  -4.98804804 group1 0.00000000
## 2  sample_1      time      0.00000000 Factor2   0.50849870 group1 0.00000000
## 3  sample_1      time      0.00000000 Factor4  -1.62411065 group1 0.00000000
## 4  sample_1      time      0.00000000 Factor3  -2.65808006 group1 0.00000000
## 5 sample_10      time      0.04522613 Factor4  -1.73719393 group1 0.04522613
## 6 sample_10      time      0.04522613 Factor2   0.02898659 group1 0.04522613
##   shape_by
## 1        1
## 2        1
## 3        1
## 4        1
## 5        1
## 6        1

We can compare the above plots to the factors that were simulated above and find that the model recaptured the two smooth as well as two non-smooth patterns in time. Note that factors are invariant to the sign, e.g. Factor 4 is the negative of the simulated factor but we can simply multiply the factors and its weights by -1 to obtain exactly the simulated factor.

2.4.3 Exploration of weights

As with standard MOFA, we can now look deeper into the meaning of these factors by exploring the weights or performing feature set enrichment analysis.

plot_weights(sm, factors = 4, view = 1)

plot_top_weights(sm, factors = 3, view = 2)

In addition, we can take a look at the top feature values per factor along time and see that their patterns are in line with the pattern of the corresponding Factor (here Factor 3).

plot_data_vs_cov(sm, factor=3,
                         features = 2,
                         color_by = "time",
                         dot_size = 1)

2.4.4 Interpolation

Furthermore, we can interpolate or extrapolate a factor to new values. Here, we only show the mean of the prediction, to obtain uncertainties you need to specify the new values before training in get_default_mefisto_options(sm)$new_values.

sm <- interpolate_factors(sm, new_values = seq(0,1.1,0.01))
plot_interpolation_vs_covariate(sm, covariate = "time",
                                factors = "Factor3")

Session Info

sessionInfo()
## R version 4.2.0 RC (2022-04-19 r82224)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.15-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] pheatmap_1.0.12  forcats_0.5.1    stringr_1.4.0    dplyr_1.0.8     
##  [5] purrr_0.3.4      readr_2.1.2      tidyr_1.2.0      tibble_3.1.6    
##  [9] ggplot2_3.3.5    tidyverse_1.3.1  MOFA2_1.6.0      BiocStyle_2.24.0
## 
## loaded via a namespace (and not attached):
##  [1] fs_1.5.2             matrixStats_0.62.0   lubridate_1.8.0     
##  [4] filelock_1.0.2       RColorBrewer_1.1-3   httr_1.4.2          
##  [7] rprojroot_2.0.3      tools_4.2.0          backports_1.4.1     
## [10] bslib_0.3.1          utf8_1.2.2           R6_2.5.1            
## [13] HDF5Array_1.24.0     uwot_0.1.11          DBI_1.1.2           
## [16] BiocGenerics_0.42.0  colorspace_2.0-3     rhdf5filters_1.8.0  
## [19] withr_2.5.0          tidyselect_1.1.2     compiler_4.2.0      
## [22] rvest_1.0.2          cli_3.3.0            basilisk.utils_1.8.0
## [25] xml2_1.3.3           DelayedArray_0.22.0  labeling_0.4.2      
## [28] bookdown_0.26        sass_0.4.1           scales_1.2.0        
## [31] mvtnorm_1.1-3        rappdirs_0.3.3       digest_0.6.29       
## [34] rmarkdown_2.14       basilisk_1.8.0       pkgconfig_2.0.3     
## [37] htmltools_0.5.2      MatrixGenerics_1.8.0 highr_0.9           
## [40] dbplyr_2.1.1         fastmap_1.1.0        rlang_1.0.2         
## [43] readxl_1.4.0         rstudioapi_0.13      farver_2.1.0        
## [46] jquerylib_0.1.4      generics_0.1.2       jsonlite_1.8.0      
## [49] magrittr_2.0.3       Matrix_1.4-1         Rcpp_1.0.8.3        
## [52] munsell_0.5.0        S4Vectors_0.34.0     Rhdf5lib_1.18.0     
## [55] fansi_1.0.3          reticulate_1.24      lifecycle_1.0.1     
## [58] stringi_1.7.6        yaml_2.3.5           rhdf5_2.40.0        
## [61] Rtsne_0.16           plyr_1.8.7           grid_4.2.0          
## [64] parallel_4.2.0       ggrepel_0.9.1        crayon_1.5.1        
## [67] dir.expiry_1.4.0     lattice_0.20-45      haven_2.5.0         
## [70] cowplot_1.1.1        hms_1.1.1            magick_2.7.3        
## [73] knitr_1.38           pillar_1.7.0         reshape2_1.4.4      
## [76] stats4_4.2.0         reprex_2.0.1         glue_1.6.2          
## [79] evaluate_0.15        BiocManager_1.30.17  modelr_0.1.8        
## [82] png_0.1-7            vctrs_0.4.1          tzdb_0.3.0          
## [85] cellranger_1.1.0     gtable_0.3.0         assertthat_0.2.1    
## [88] xfun_0.30            broom_0.8.0          IRanges_2.30.0      
## [91] corrplot_0.92        ellipsis_0.3.2       here_1.0.1