miaSim 1.10.0
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.
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)
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))
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.
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)
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)
Stochastic logistic model is used to determine dead and alive counts in community.
tse_logistic <- simulateStochasticLogistic(n_species = 5)
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.
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)
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)
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.
sessionInfo()
## R version 4.4.0 beta (2024-04-15 r86425)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.19-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_US.UTF-8 LC_COLLATE=en_US.UTF-8
## [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.12.0 ggraph_2.2.1
## [3] ggplot2_3.5.1 mia_1.12.0
## [5] MultiAssayExperiment_1.30.0 miaSim_1.10.0
## [7] TreeSummarizedExperiment_2.12.0 Biostrings_2.72.0
## [9] XVector_0.44.0 SingleCellExperiment_1.26.0
## [11] SummarizedExperiment_1.34.0 Biobase_2.64.0
## [13] GenomicRanges_1.56.0 GenomeInfoDb_1.40.0
## [15] IRanges_2.38.0 S4Vectors_0.42.0
## [17] BiocGenerics_0.50.0 MatrixGenerics_1.16.0
## [19] matrixStats_1.3.0 BiocStyle_2.32.0
##
## loaded via a namespace (and not attached):
## [1] jsonlite_1.8.8 magrittr_2.0.3
## [3] magick_2.8.3 ggbeeswarm_0.7.2
## [5] farver_2.1.1 rmarkdown_2.26
## [7] fs_1.6.4 zlibbioc_1.50.0
## [9] vctrs_0.6.5 memoise_2.0.1
## [11] DelayedMatrixStats_1.26.0 ggtree_3.12.0
## [13] tinytex_0.50 htmltools_0.5.8.1
## [15] S4Arrays_1.4.0 BiocNeighbors_1.22.0
## [17] deSolve_1.40 gridGraphics_0.5-1
## [19] SparseArray_1.4.0 sass_0.4.9
## [21] pracma_2.4.4 bslib_0.7.0
## [23] plyr_1.8.9 DECIPHER_3.0.0
## [25] cachem_1.0.8 igraph_2.0.3
## [27] lifecycle_1.0.4 pkgconfig_2.0.3
## [29] rsvd_1.0.5 Matrix_1.7-0
## [31] R6_2.5.1 fastmap_1.1.1
## [33] GenomeInfoDbData_1.2.12 digest_0.6.35
## [35] aplot_0.2.2 colorspace_2.1-0
## [37] ggnewscale_0.4.10 patchwork_1.2.0
## [39] scater_1.32.0 irlba_2.3.5.1
## [41] vegan_2.6-4 beachmat_2.20.0
## [43] labeling_0.4.3 fansi_1.0.6
## [45] httr_1.4.7 polyclip_1.10-6
## [47] abind_1.4-5 mgcv_1.9-1
## [49] compiler_4.4.0 withr_3.0.0
## [51] BiocParallel_1.38.0 viridis_0.6.5
## [53] DBI_1.2.2 highr_0.10
## [55] ggforce_0.4.2 MASS_7.3-60.2
## [57] poweRlaw_0.80.0 DelayedArray_0.30.0
## [59] bluster_1.14.0 permute_0.9-7
## [61] tools_4.4.0 vipor_0.4.7
## [63] beeswarm_0.4.0 ape_5.8
## [65] glue_1.7.0 nlme_3.1-164
## [67] grid_4.4.0 cluster_2.1.6
## [69] reshape2_1.4.4 generics_0.1.3
## [71] gtable_0.3.5 tidyr_1.3.1
## [73] BiocSingular_1.20.0 tidygraph_1.3.1
## [75] ScaledMatrix_1.12.0 utf8_1.2.4
## [77] ggrepel_0.9.5 pillar_1.9.0
## [79] stringr_1.5.1 yulab.utils_0.1.4
## [81] splines_4.4.0 dplyr_1.1.4
## [83] tweenr_2.0.3 treeio_1.28.0
## [85] lattice_0.22-6 tidyselect_1.2.1
## [87] DirichletMultinomial_1.46.0 scuttle_1.14.0
## [89] knitr_1.46 gridExtra_2.3
## [91] bookdown_0.39 xfun_0.43
## [93] graphlayouts_1.1.1 stringi_1.8.3
## [95] UCSC.utils_1.0.0 lazyeval_0.2.2
## [97] ggfun_0.1.4 yaml_2.3.8
## [99] evaluate_0.23 codetools_0.2-20
## [101] tibble_3.2.1 BiocManager_1.30.22
## [103] ggplotify_0.1.2 cli_3.6.2
## [105] munsell_0.5.1 jquerylib_0.1.4
## [107] Rcpp_1.0.12 parallel_4.4.0
## [109] sparseMatrixStats_1.16.0 decontam_1.24.0
## [111] viridisLite_0.4.2 tidytree_0.4.6
## [113] scales_1.3.0 purrr_1.0.2
## [115] crayon_1.5.2 rlang_1.1.3