Contents

1 Introduction

miaSim implements tools for microbiome data simulation based on varying ecological modeling assumptions. These can be used to simulate species abundance matrices, including time series. Detailed function documentation is available at the function reference

The miaSim package supports the R/Bioconductor framework based on the TreeSummarizedExperiment data container. For more information on operating with this data format in microbial ecology, see the online tutorial.

1.1 Installation

The stable Bioconductor release version can be installed as follows.

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
if (!requireNamespace("miaSim", quietly = TRUE))    
    BiocManager::install("miaSim")

The experimental Bioconductor devel version can be installed as follows.

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
# The following initializes usage of Bioc devel
BiocManager::install(version='devel')
BiocManager::install("miaSim")

Load the library

library(miaSim)

1.2 Examples

1.2.1 Generating species interaction matrices

Some of the models rely on interaction matrices that represents interaction heterogeneity between species. The interaction matrix can be generated with different distributional assumptions.

Generate interactions from normal distribution:

A_normal <- powerlawA(n_species = 4, alpha = 3)

Generate interactions from uniform distribution:

A_uniform <- randomA(n_species = 10,
                 diagonal = -0.4,
                     connectance = 0.5,
             interactions = runif(n = 10^2, min = -0.8, max = 0.8))

1.2.2 Ricker model

Ricker model is a discrete version of the gLV:

tse_ricker <- simulateRicker(n_species=4, A = A_normal, t_end=100, norm = FALSE)

The number of species specified in the interaction matrix must be the same as the species used in the models.

1.2.3 Hubbell model

Hubbell Neutral simulation model characterizes diversity and relative abundance of species in ecological communities assuming migration, births and deaths but no interactions. Losses become replaced by migration or birth.

tse_hubbell <- simulateHubbell(n_species = 8,
                               M = 10,
                   carrying_capacity = 1000,
                               k_events = 50,
                   migration_p = 0.02,
                   t_end = 100)

One can also simulate parameters for the Hubbell model.

params_hubbell <- simulateHubbellRates(x0 = c(0,5,10),
    migration_p = 0.1, metacommunity_probability = NULL, k_events = 1, 
    growth_rates = NULL, norm = FALSE, t_end=1000)

1.2.4 Self-Organised Instability (SOI)

The Self-Organised Instability (SOI) model generates time series for communities and accelerates stochastic simulation.

tse_soi <- simulateSOI(n_species = 4, carrying_capacity = 1000,
                       A = A_normal, k_events=5,
               x0 = NULL,t_end = 150, norm = TRUE)

1.2.5 Stochastic logistic model

Stochastic logistic model is used to determine dead and alive counts in community.

tse_logistic <- simulateStochasticLogistic(n_species = 5)

1.2.6 Consumer-resource model

The consumer resource model requires the randomE function. This returns a matrix containing the production rates and consumption rates of each species. The resulting matrix is used as a determination of resource consumption efficiency.

# Consumer-resource model as a TreeSE object
tse_crm <- simulateConsumerResource(n_species = 2,
                                    n_resources = 4,
                                    E = randomE(n_species = 2, n_resources = 4))

You could visualize the simulated dynamics using tools from the miaTime package.

1.2.7 Generalized Lotka-Volterra (gLV)

The generalized Lotka-Volterra simulation model generates time-series assuming microbial population dynamics and interaction.

tse_glv <- simulateGLV(n_species = 4,
                       A = A_normal,
               t_start = 0, 
                       t_store = 1000,
               stochastic = FALSE,
               norm = FALSE)

1.3 Data containers

The simulated data sets are returned as TreeSummarizedExperiment objects. This provides access to a broad range of tools for microbiome analysis that support this format (see microbiome.github.io). More examples on the object manipulation and analysis can be found at OMA Online Manual.

For instance, to plot population density we can use the miaViz package:

library(miaViz)
p1 <- plotAbundanceDensity(tse_hubbell, assay.type = "counts")
p2 <- plotSeries(tse_hubbell, x = "time")
print(p1+p2)

1.4 Case studies

Source code for replicating the published case studies using the miaSim package (Gao et al. 2023) is available in the phyloseq folder (based on the phyloseq data container). Some of the original case studies have now been replicated also with TreeSummarized, see the TreeSummarizedExperiment folder.

2 Session info

sessionInfo()
## R version 4.3.1 (2023-06-16)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
## 
## 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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] miaViz_1.10.0                   ggraph_2.1.0                   
##  [3] ggplot2_3.4.4                   mia_1.10.0                     
##  [5] MultiAssayExperiment_1.28.0     miaSim_1.8.0                   
##  [7] TreeSummarizedExperiment_2.10.0 Biostrings_2.70.0              
##  [9] XVector_0.42.0                  SingleCellExperiment_1.24.0    
## [11] SummarizedExperiment_1.32.0     Biobase_2.62.0                 
## [13] GenomicRanges_1.54.0            GenomeInfoDb_1.38.0            
## [15] IRanges_2.36.0                  S4Vectors_0.40.0               
## [17] BiocGenerics_0.48.0             MatrixGenerics_1.14.0          
## [19] matrixStats_1.0.0               BiocStyle_2.30.0               
## 
## loaded via a namespace (and not attached):
##   [1] jsonlite_1.8.7              magrittr_2.0.3             
##   [3] magick_2.8.1                ggbeeswarm_0.7.2           
##   [5] farver_2.1.1                rmarkdown_2.25             
##   [7] fs_1.6.3                    zlibbioc_1.48.0            
##   [9] vctrs_0.6.4                 memoise_2.0.1              
##  [11] DelayedMatrixStats_1.24.0   RCurl_1.98-1.12            
##  [13] ggtree_3.10.0               htmltools_0.5.6.1          
##  [15] S4Arrays_1.2.0              BiocNeighbors_1.20.0       
##  [17] deSolve_1.38                gridGraphics_0.5-1         
##  [19] SparseArray_1.2.0           sass_0.4.7                 
##  [21] pracma_2.4.2                bslib_0.5.1                
##  [23] plyr_1.8.9                  DECIPHER_2.30.0            
##  [25] cachem_1.0.8                igraph_1.5.1               
##  [27] lifecycle_1.0.3             pkgconfig_2.0.3            
##  [29] rsvd_1.0.5                  Matrix_1.6-1.1             
##  [31] R6_2.5.1                    fastmap_1.1.1              
##  [33] GenomeInfoDbData_1.2.11     aplot_0.2.2                
##  [35] digest_0.6.33               colorspace_2.1-0           
##  [37] ggnewscale_0.4.9            patchwork_1.1.3            
##  [39] scater_1.30.0               irlba_2.3.5.1              
##  [41] RSQLite_2.3.1               vegan_2.6-4                
##  [43] beachmat_2.18.0             labeling_0.4.3             
##  [45] fansi_1.0.5                 polyclip_1.10-6            
##  [47] abind_1.4-5                 mgcv_1.9-0                 
##  [49] compiler_4.3.1              bit64_4.0.5                
##  [51] withr_2.5.1                 BiocParallel_1.36.0        
##  [53] viridis_0.6.4               DBI_1.1.3                  
##  [55] ggforce_0.4.1               MASS_7.3-60                
##  [57] poweRlaw_0.70.6             DelayedArray_0.28.0        
##  [59] bluster_1.12.0              permute_0.9-7              
##  [61] tools_4.3.1                 vipor_0.4.5                
##  [63] beeswarm_0.4.0              ape_5.7-1                  
##  [65] glue_1.6.2                  nlme_3.1-163               
##  [67] grid_4.3.1                  cluster_2.1.4              
##  [69] reshape2_1.4.4              generics_0.1.3             
##  [71] gtable_0.3.4                tidyr_1.3.0                
##  [73] BiocSingular_1.18.0         tidygraph_1.2.3            
##  [75] ScaledMatrix_1.10.0         utf8_1.2.4                 
##  [77] ggrepel_0.9.4               pillar_1.9.0               
##  [79] stringr_1.5.0               yulab.utils_0.1.0          
##  [81] splines_4.3.1               dplyr_1.1.3                
##  [83] tweenr_2.0.2                treeio_1.26.0              
##  [85] lattice_0.22-5              bit_4.0.5                  
##  [87] tidyselect_1.2.0            DirichletMultinomial_1.44.0
##  [89] scuttle_1.12.0              knitr_1.44                 
##  [91] gridExtra_2.3               bookdown_0.36              
##  [93] xfun_0.40                   graphlayouts_1.0.1         
##  [95] stringi_1.7.12              ggfun_0.1.3                
##  [97] lazyeval_0.2.2              yaml_2.3.7                 
##  [99] evaluate_0.22               codetools_0.2-19           
## [101] tibble_3.2.1                BiocManager_1.30.22        
## [103] ggplotify_0.1.2             cli_3.6.1                  
## [105] munsell_0.5.0               jquerylib_0.1.4            
## [107] Rcpp_1.0.11                 parallel_4.3.1             
## [109] blob_1.2.4                  sparseMatrixStats_1.14.0   
## [111] bitops_1.0-7                decontam_1.22.0            
## [113] viridisLite_0.4.2           tidytree_0.4.5             
## [115] scales_1.2.1                purrr_1.0.2                
## [117] crayon_1.5.2                rlang_1.1.1