This package implements the CBEA approach for performing set-based enrichment analysis for microbiome relative abundance data. A preprint of the package can be found on bioXriv. In summary, CBEA (Competitive Balances for taxonomic Enrichment Analysis) provides an estimate of the activity of a set by transforming an input taxa-by-sample data matrix into a corresponding set-by-sample data matrix. The resulting output can be used for additional downstream analyses such as differential abundance, classification, clustering, etc. using set-based features instead of the original units.
The transformation that CBEA applies is based on the isometric log
ratio transformation:
$$
CBEA_{i,\mathbb{S}} =
\sqrt{\frac{|\mathbb{S}||\mathbb{S_c}|}{|\mathbb{S}| + |\mathbb{S_c}|}}
\ln \frac{g(X_{i,j | j\in \mathbb{S}})}{g(X_{i,j | j \notin
\mathbb{S}})}
$$ Where π is the set of
interest, πC is
itβs complement, g() is the
geometric mean operation, and X is the original data matrix where
i is the index representing
samples and j is the index
representing variables (or taxa).
The inference procedure is performed through estimating the null distribution of the test statistic. This can be done either via permutations or a parametric fit of a distributional form on the permuted scores. Users can also adjust for variance inflation due to inter-taxa correlation. Please refer to the main manuscript for any additional details.
CBEA
CBEA is an
R
package available via the Bioconductor repository for packages.
It requires installing the R
open source statistical
programming language, which can be accessed on any operating system from
CRAN. After which you can
install CBEA by
using the following commands in your R
session:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("CBEA")
## Check that you have a valid Bioconductor installation
BiocManager::valid()
If there are any issues with the installation procedure or package
features, the best place would be to file an issue at the GitHub repository. For any
additional support you can use the Bioconductor support site
and use the CBEA
tag and check the older posts.
Please note that if you want to receive help you should adhere to the posting
guidelines. It is particularly critical that you provide a small
reproducible example and your session information so package developers
can track down the source of the error.
First, we can load some pre-packaged data sets contained in CBEA. Here
weβre loading the data from the Human Microbiome Project (HMP) in
TreeSummarizedExperiment
data container from the TreeSummarizedExperiment.
This package does not support phyloseq
from the phyloseq
package but users can leverage the mia package
to convert between the data types.
In addition, users can also input raw matrices or data frames,
however those require additional arguments. The
taxa_are_rows
argument requires users specify whether the
data frame/matrix has taxa abundances as rows (or as columns). The
id_col
argument requires users to specify (for data frames
only) a vector of names of row metadata that will be excluded from the
analysis.
data("hmp_gingival")
abun <- hmp_gingival$data
metab_sets <- hmp_gingival$set
abun # this is a TreeSummarizedExperiment object
#> class: TreeSummarizedExperiment
#> dim: 5378 369
#> metadata(4): experimentData phylogeneticTree experimentData
#> phylogeneticTree
#> assays(1): 16SrRNA
#> rownames(5378): OTU_97.10005 OTU_97.10006 ... OTU_97.9991 OTU_97.9995
#> rowData names(7): CONSENSUS_LINEAGE SUPERKINGDOM ... FAMILY GENUS
#> colnames(369): 700014427 700014521 ... 700111587 700111758
#> colData names(7): RSID VISITNO ... HMP_BODY_SUBSITE SRS_SAMPLE_ID
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> rowLinks: NULL
#> rowTree: NULL
#> colLinks: NULL
#> colTree: NULL
CBEA accepts any type of sets, as long as it is in the
BiocSet
format where the elements in the sets can be
matched to taxa names in the data set. The main function will check if
these names match.
metab_sets
#> class: BiocSet
#>
#> es_element():
#> # A tibble: 4,828 Γ 1
#> element
#> <chr>
#> 1 OTU_97.10005
#> 2 OTU_97.10053
#> 3 OTU_97.10090
#> # βΉ 4,825 more rows
#>
#> es_set():
#> # A tibble: 3 Γ 1
#> set
#> <chr>
#> 1 F Anaerobic
#> 2 Aerobic
#> 3 Anaerobic
#>
#> es_elementset() <active>:
#> # A tibble: 4,828 Γ 2
#> element set
#> <chr> <chr>
#> 1 OTU_97.10005 F Anaerobic
#> 2 OTU_97.10053 F Anaerobic
#> 3 OTU_97.10090 F Anaerobic
#> # βΉ 4,825 more rows
For more information on BiocSet
, please refer to the
documentation from BiocSet.
However, simply speaking, BiocSet
acts similar to a list of
three data frames and can be used in conjunction with dplyr/tidyr.
After specifying the inputs, cbea
is the main function
to apply the method. If there are zeros in the abundance data, the
cbea
will add a pseudocount to avoid issues with the
log-ratio transformation (but will throw a warning). If a different
zero-handling approach is desired, users should pre-process the
abundance data with the appropriate method. For parametric fits,
cbea
relies on the fitdistrplus
and mixtools
packages to estimate the parameters of the null. Specific arguments to
control this fitting procedure can be provided as a named list in the
control
argument.
Applying cbea
is one command:
results <- cbea(abun, set = metab_sets, abund_values = "16SrRNA",
output = "cdf", distr = "mnorm", adj = TRUE, thresh = 0.05, n_perm = 10)
#> Warning in .local(obj, set, output, distr, adj, n_perm, parametric, thresh, : Taxonomic count table contains zeros,
#> which would invalidate the log-ratio transform.
#> Adding a pseudocount of 1e-5...
#> number of iterations= 644
#> number of iterations= 1253
#> number of iterations= 468
#> number of iterations= 847
#> number of iterations= 612
#> number of iterations= 105
results
#> CBEA output of class 'CBEAout' with 369 samples and 3 sets
#> Fit type: Parametric with 2-component Gaussian Mixture Distribution
#> Number of permutations: 10
#> Output type: CDF values (cdf)
Some important arguments to control the behaviour of CBEA.
output
: This controls what type of output is being
returned. CBEA usually estimates a parametric null and users can specify
what they want in return. If users want to perform downstream analysis
with set-level features, they can return CDF values or z-scores of each
raw score computed against that distribution (options cdf
or zscore
). Alternatively, users can just return the raw
scores themselves (no distribution fitting will be performed) using
raw
as the option. Users can also use this distribution to
estimate unadjusted p-values (option pval
) to see whether a
set is enriched at each sample. These unadjusted p-values can be
converted based on a threshold (based on thresh
which is
default to be set at 0.05) into a dummy variable indicating enrichment
(option sig
). Note: CDF values and
Z-scores are not available for non-parametric null estimations.parametric
: This is a logical argument to specify
whether a the null distribution will be specified via parametric fit or
via non-parametric permutation testing. If parametric
is
TRUE
, users need to specify distr
and
adj
. If parametric
is FALSE
,
users need to increase n_perm
.distr
: The form of the distribution if parametric fit
is desired. As of now only supports norm
,
mnorm
.adj
: Whether the distribution should be adjusted for
variance inflation. This procedure is done by combining the mean
estimate from scores computed from permuted data set and the variance
estimate from raw scores (computed on the unpermuted data set).The output object is of class CBEAout
, which is an S3
object. The underlying data structure is a list of lists, where the
outer lists represent different aspects of the output. For example
R
represent the final scores while diagnostic
represent certain goodness-of-fit statistics.
Within each aspect, there is a list of size equivalent to the total
number of sets evaluated. For example, the results
object
is of size 3 representing the evaluated sets.
str(results$R)
#> List of 3
#> $ Aerobic : num [1:369] 0.9694 0.0286 0.9953 0.8842 0.5617 ...
#> $ Anaerobic : num [1:369] 0.000461 0.58112 0.203075 0.034167 0.088761 ...
#> $ F Anaerobic: num [1:369] 0.997 0.89 0.303 0.96 0.863 ...
Users can use tidy
and glance
following the
broom to
process CBEAout
into nice objects. The tidy
function returns a tibble
of scores (samples by set). The
glance
function returns some diagnostics. There are two
options for the glance
function:
fit_comparison
allows users to compare the l-moments of the
data, the permuted data, and the final fitted distribution;
fit_diagnostic
shows goodness-of-fit statistics of the
distribution fitting procedure itself, with log-likelihoods and
Anderson-Darling (column βadβ) statistics.
tidy(results)
#> # A tibble: 369 Γ 4
#> sample_ids Aerobic Anaerobic F_Anaerobic
#> <chr> <dbl> <dbl> <dbl>
#> 1 700014427 0.969 0.000461 0.997
#> 2 700014521 0.0286 0.581 0.890
#> 3 700014603 0.995 0.203 0.303
#> 4 700014749 0.884 0.0342 0.960
#> 5 700014791 0.562 0.0888 0.863
#> 6 700014917 0.720 0.134 0.898
#> 7 700014989 0.0249 0.992 0.0462
#> 8 700015076 0.685 0.0698 0.911
#> 9 700015149 0.402 0.000344 1.00
#> 10 700015215 0.0976 0.729 0.605
#> # βΉ 359 more rows
glance(results, "fit_comparison")
#> # A tibble: 9 Γ 7
#> set_ids final_param distr l_location l_scale l_skewness l_kurtosis
#> <chr> <named list> <chr> <dbl> <dbl> <dbl> <dbl>
#> 1 Aerobic <dbl [6]> data 6.98 15.2 0.0300 0.114
#> 2 Aerobic <dbl [6]> perm -0.0557 2.58 -0.0000123 0.131
#> 3 Aerobic <dbl [6]> fit 0.0768 13.0 -0.00401 0.127
#> 4 Anaerobic <dbl [6]> data -6.42 20.0 0.0360 0.0917
#> 5 Anaerobic <dbl [6]> perm 0.758 2.54 0.0146 0.122
#> 6 Anaerobic <dbl [6]> fit 0.986 15.0 0.0104 0.132
#> 7 F_Anaerobic <dbl [6]> data 5.88 16.5 0.0481 0.0952
#> 8 F_Anaerobic <dbl [6]> perm -0.465 2.71 0.000482 0.125
#> 9 F_Anaerobic <dbl [6]> fit -0.591 15.1 0.000309 0.123
glance(results, "fit_diagnostic")
#> # A tibble: 6 Γ 5
#> set_ids final_param loglik ad type
#> <chr> <named list> <dbl> <dbl> <chr>
#> 1 Aerobic <dbl [6]> -10853. 0.229 permuted_distr
#> 2 Aerobic <dbl [6]> -1735. 0.141 unperm_distr
#> 3 Anaerobic <dbl [6]> -10791. 0.327 permuted_distr
#> 4 Anaerobic <dbl [6]> -1832. 0.185 unperm_distr
#> 5 F_Anaerobic <dbl [6]> -11036. 0.321 permuted_distr
#> 6 F_Anaerobic <dbl [6]> -1762. 0.277 unperm_distr
CBEA
has in-built capacity to perform calculations
paralelled across the total number of sets. The engine for
parallelization is BiocParallel.
If NULL
, SerialParam
backend will be used.
BiocParallel::registered()
#> $MulticoreParam
#> class: MulticoreParam
#> bpisup: FALSE; bpnworkers: 2; bptasks: 0; bpjobname: BPJOB
#> bplog: FALSE; bpthreshold: INFO; bpstopOnError: TRUE
#> bpRNGseed: ; bptimeout: NA; bpprogressbar: FALSE
#> bpexportglobals: TRUE; bpexportvariables: FALSE; bpforceGC: FALSE
#> bpfallback: TRUE
#> bplogdir: NA
#> bpresultdir: NA
#> cluster type: FORK
#>
#> $SnowParam
#> class: SnowParam
#> bpisup: FALSE; bpnworkers: 2; bptasks: 0; bpjobname: BPJOB
#> bplog: FALSE; bpthreshold: INFO; bpstopOnError: TRUE
#> bpRNGseed: ; bptimeout: NA; bpprogressbar: FALSE
#> bpexportglobals: TRUE; bpexportvariables: TRUE; bpforceGC: FALSE
#> bpfallback: TRUE
#> bplogdir: NA
#> bpresultdir: NA
#> cluster type: SOCK
#>
#> $SerialParam
#> class: SerialParam
#> bpisup: FALSE; bpnworkers: 1; bptasks: 0; bpjobname: BPJOB
#> bplog: FALSE; bpthreshold: INFO; bpstopOnError: TRUE
#> bpRNGseed: ; bptimeout: NA; bpprogressbar: FALSE
#> bpexportglobals: FALSE; bpexportvariables: FALSE; bpforceGC: FALSE
#> bpfallback: FALSE
#> bplogdir: NA
#> bpresultdir: NA
CBEA
We hope that CBEA will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
## Citation info
citation("CBEA")
#> To cite package 'CBEA' in publications use:
#>
#> Nguyen Q (2024). _CBEA: R package for performing CBEA approach_.
#> doi:10.18129/B9.bioc.CBEA <https://doi.org/10.18129/B9.bioc.CBEA>,
#> https://github.com/qpmnguyen/CBEA - R package version 1.7.0,
#> <http://www.bioconductor.org/packages/CBEA>.
#>
#> Nguyen Q (2024). "CBEA: Competitive balances for taxonomic enrichment
#> analysis." _bioRxiv_. doi:10.1101/TODO
#> <https://doi.org/10.1101/TODO>,
#> <https://www.biorxiv.org/content/10.1101/TODO>.
#>
#> To see these entries in BibTeX format, use 'print(<citation>,
#> bibtex=TRUE)', 'toBibtex(.)', or set
#> 'options(citation.bibtex.max=999)'.
The CBEA package (Nguyen, 2024) was made possible thanks to:
This package was developed using biocthis.
Code for creating the vignette
## Create the vignette
library("rmarkdown")
system.time(render("basic_usage.Rmd", "BiocStyle::html_document"))
## Extract the R code
library("knitr")
knit("basic_usage.Rmd", tangle = TRUE)
Date the vignette was generated.
#> [1] "2024-12-19 03:25:42 UTC"
Wallclock time spent generating the vignette.
#> Time difference of 14.891 secs
R
session information.
#> β Session info βββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββ
#> setting value
#> version R version 4.4.2 (2024-10-31)
#> os Ubuntu 24.04.1 LTS
#> system x86_64, linux-gnu
#> ui X11
#> language (EN)
#> collate C
#> ctype en_US.UTF-8
#> tz Etc/UTC
#> date 2024-12-19
#> pandoc 3.2.1 @ /usr/local/bin/ (via rmarkdown)
#>
#> β Packages βββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββ
#> package * version date (UTC) lib source
#> abind 1.4-8 2024-09-12 [2] RSPM (R 4.4.0)
#> AnnotationDbi 1.69.0 2024-11-29 [2] https://bioc.r-universe.dev (R 4.4.2)
#> ape 5.8-1 2024-12-16 [2] RSPM (R 4.4.0)
#> backports 1.5.0 2024-05-23 [2] RSPM (R 4.4.0)
#> bibtex 0.5.1 2023-01-26 [2] RSPM (R 4.4.0)
#> Biobase 2.67.0 2024-11-30 [2] https://bioc.r-universe.dev (R 4.4.2)
#> BiocGenerics 0.53.3 2024-12-15 [2] https://bioc.r-universe.dev (R 4.4.2)
#> BiocIO 1.17.1 2024-11-22 [2] https://bioc.r-universe.dev (R 4.4.2)
#> BiocManager 1.30.25 2024-08-28 [2] RSPM (R 4.4.0)
#> BiocParallel 1.41.0 2024-11-29 [2] https://bioc.r-universe.dev (R 4.4.2)
#> BiocSet * 1.21.0 2024-11-19 [2] https://bioc.r-universe.dev (R 4.4.2)
#> BiocStyle * 2.35.0 2024-12-19 [2] https://bioc.r-universe.dev (R 4.4.2)
#> Biostrings 2.75.3 2024-12-16 [2] https://bioc.r-universe.dev (R 4.4.2)
#> bit 4.5.0.1 2024-12-03 [2] RSPM (R 4.4.0)
#> bit64 4.5.2 2024-09-22 [2] RSPM (R 4.4.0)
#> blob 1.2.4 2023-03-17 [2] RSPM (R 4.4.0)
#> bslib 0.8.0 2024-07-29 [2] RSPM (R 4.4.0)
#> buildtools 1.0.0 2024-12-16 [3] local (/pkg)
#> cachem 1.1.0 2024-05-16 [2] RSPM (R 4.4.0)
#> CBEA * 1.7.0 2024-12-19 [1] https://bioc.r-universe.dev (R 4.4.2)
#> cli 3.6.3 2024-06-21 [2] RSPM (R 4.4.0)
#> codetools 0.2-20 2024-03-31 [2] RSPM (R 4.4.0)
#> colorspace 2.1-1 2024-07-26 [2] RSPM (R 4.4.0)
#> crayon 1.5.3 2024-06-20 [2] RSPM (R 4.4.0)
#> data.table 1.16.4 2024-12-06 [2] RSPM (R 4.4.0)
#> DBI 1.2.3 2024-06-02 [2] RSPM (R 4.4.0)
#> DelayedArray 0.33.3 2024-12-03 [2] https://bioc.r-universe.dev (R 4.4.2)
#> digest 0.6.37 2024-08-19 [2] RSPM (R 4.4.0)
#> dplyr * 1.1.4 2023-11-17 [2] RSPM (R 4.4.0)
#> evaluate 1.0.1 2024-10-10 [2] RSPM (R 4.4.0)
#> fastmap 1.2.0 2024-05-15 [2] RSPM (R 4.4.0)
#> fitdistrplus 1.2-1 2024-07-12 [2] RSPM (R 4.4.0)
#> forcats * 1.0.0 2023-01-29 [2] RSPM (R 4.4.0)
#> fs 1.6.5 2024-10-30 [2] RSPM (R 4.4.0)
#> generics 0.1.3 2022-07-05 [2] RSPM (R 4.4.0)
#> GenomeInfoDb 1.43.2 2024-11-28 [2] https://bioc.r-universe.dev (R 4.4.2)
#> GenomeInfoDbData 1.2.13 2024-12-19 [2] Bioconductor
#> GenomicRanges 1.59.1 2024-12-19 [2] https://bioc.r-universe.dev (R 4.4.2)
#> ggplot2 * 3.5.1 2024-04-23 [2] RSPM (R 4.4.0)
#> glue 1.8.0 2024-09-30 [2] RSPM (R 4.4.0)
#> goftest 1.2-3 2021-10-07 [2] RSPM (R 4.4.0)
#> gtable 0.3.6 2024-10-25 [2] RSPM (R 4.4.0)
#> hms 1.1.3 2023-03-21 [2] RSPM (R 4.4.0)
#> htmltools 0.5.8.1 2024-04-04 [2] RSPM (R 4.4.0)
#> htmlwidgets 1.6.4 2023-12-06 [2] RSPM (R 4.4.0)
#> httr 1.4.7 2023-08-15 [2] RSPM (R 4.4.0)
#> IRanges 2.41.2 2024-12-03 [2] https://bioc.r-universe.dev (R 4.4.2)
#> jquerylib 0.1.4 2021-04-26 [2] RSPM (R 4.4.0)
#> jsonlite 1.8.9 2024-09-20 [2] RSPM (R 4.4.0)
#> KEGGREST 1.47.0 2024-11-29 [2] https://bioc.r-universe.dev (R 4.4.2)
#> kernlab 0.9-33 2024-08-13 [2] RSPM (R 4.4.0)
#> knitr 1.49 2024-11-08 [2] RSPM (R 4.4.0)
#> lattice 0.22-6 2024-03-20 [2] RSPM (R 4.4.0)
#> lazyeval 0.2.2 2019-03-15 [2] RSPM (R 4.4.0)
#> lifecycle 1.0.4 2023-11-07 [2] RSPM (R 4.4.0)
#> lmom 3.2 2024-09-30 [2] RSPM (R 4.4.0)
#> lubridate * 1.9.4 2024-12-08 [2] RSPM (R 4.4.0)
#> magrittr 2.0.3 2022-03-30 [2] RSPM (R 4.4.0)
#> maketools 1.3.1 2024-10-04 [3] RSPM (R 4.4.0)
#> MASS 7.3-61 2024-06-13 [2] RSPM (R 4.4.0)
#> Matrix 1.7-1 2024-10-18 [2] RSPM (R 4.4.0)
#> MatrixGenerics 1.19.0 2024-12-06 [2] https://bioc.r-universe.dev (R 4.4.2)
#> matrixStats 1.4.1 2024-09-08 [2] RSPM (R 4.4.0)
#> memoise 2.0.1 2021-11-26 [2] RSPM (R 4.4.0)
#> mixtools 2.0.0 2022-12-05 [2] RSPM (R 4.4.0)
#> munsell 0.5.1 2024-04-01 [2] RSPM (R 4.4.0)
#> nlme 3.1-166 2024-08-14 [2] RSPM (R 4.4.0)
#> ontologyIndex 2.12 2024-02-27 [2] RSPM (R 4.4.0)
#> pillar 1.10.0 2024-12-17 [2] RSPM (R 4.4.0)
#> pkgconfig 2.0.3 2019-09-22 [2] RSPM (R 4.4.0)
#> plotly 4.10.4 2024-01-13 [2] RSPM (R 4.4.0)
#> plyr 1.8.9 2023-10-02 [2] RSPM (R 4.4.0)
#> png 0.1-8 2022-11-29 [2] RSPM (R 4.4.0)
#> purrr * 1.0.2 2023-08-10 [2] RSPM (R 4.4.0)
#> R6 2.5.1 2021-08-19 [2] RSPM (R 4.4.0)
#> Rcpp 1.0.13-1 2024-11-02 [2] RSPM (R 4.4.0)
#> readr * 2.1.5 2024-01-10 [2] RSPM (R 4.4.0)
#> RefManageR * 1.4.0 2022-09-30 [2] RSPM (R 4.4.0)
#> rlang 1.1.4 2024-06-04 [2] RSPM (R 4.4.0)
#> rmarkdown 2.29 2024-11-04 [2] RSPM (R 4.4.0)
#> RSQLite 2.3.9 2024-12-03 [2] RSPM (R 4.4.0)
#> S4Arrays 1.7.1 2024-12-18 [2] https://bioc.r-universe.dev (R 4.4.2)
#> S4Vectors 0.45.2 2024-12-16 [2] https://bioc.r-universe.dev (R 4.4.2)
#> sass 0.4.9 2024-03-15 [2] RSPM (R 4.4.0)
#> scales 1.3.0 2023-11-28 [2] RSPM (R 4.4.0)
#> segmented 2.1-3 2024-10-25 [2] RSPM (R 4.4.0)
#> sessioninfo * 1.2.2 2021-12-06 [2] RSPM (R 4.4.0)
#> SingleCellExperiment 1.29.1 2024-12-09 [2] https://bioc.r-universe.dev (R 4.4.2)
#> SparseArray 1.7.2 2024-12-15 [2] https://bioc.r-universe.dev (R 4.4.2)
#> stringi 1.8.4 2024-05-06 [2] RSPM (R 4.4.0)
#> stringr * 1.5.1 2023-11-14 [2] RSPM (R 4.4.0)
#> SummarizedExperiment 1.37.0 2024-11-21 [2] https://bioc.r-universe.dev (R 4.4.2)
#> survival 3.8-3 2024-12-17 [2] RSPM (R 4.4.0)
#> sys 3.4.3 2024-10-04 [2] RSPM (R 4.4.0)
#> tibble * 3.2.1 2023-03-20 [2] RSPM (R 4.4.0)
#> tidyr * 1.3.1 2024-01-24 [2] RSPM (R 4.4.0)
#> tidyselect 1.2.1 2024-03-11 [2] RSPM (R 4.4.0)
#> tidytree 0.4.6 2023-12-12 [2] RSPM (R 4.4.0)
#> tidyverse * 2.0.0 2023-02-22 [2] RSPM (R 4.4.0)
#> timechange 0.3.0 2024-01-18 [2] RSPM (R 4.4.0)
#> treeio 1.31.0 2024-11-19 [2] https://bioc.r-universe.dev (R 4.4.2)
#> TreeSummarizedExperiment 2.15.0 2024-11-26 [2] https://bioc.r-universe.dev (R 4.4.2)
#> tzdb 0.4.0 2023-05-12 [2] RSPM (R 4.4.0)
#> UCSC.utils 1.3.0 2024-11-30 [2] https://bioc.r-universe.dev (R 4.4.2)
#> utf8 1.2.4 2023-10-22 [2] RSPM (R 4.4.0)
#> vctrs 0.6.5 2023-12-01 [2] RSPM (R 4.4.0)
#> viridisLite 0.4.2 2023-05-02 [2] RSPM (R 4.4.0)
#> withr 3.0.2 2024-10-28 [2] RSPM (R 4.4.0)
#> xfun 0.49 2024-10-31 [2] RSPM (R 4.4.0)
#> xml2 1.3.6 2023-12-04 [2] RSPM (R 4.4.0)
#> XVector 0.47.0 2024-11-21 [2] https://bioc.r-universe.dev (R 4.4.2)
#> yaml 2.3.10 2024-07-26 [2] RSPM (R 4.4.0)
#> yulab.utils 0.1.8 2024-11-07 [2] RSPM (R 4.4.0)
#> zlibbioc 1.52.0 2024-10-29 [2] Bioconductor 3.20 (R 4.4.2)
#>
#> [1] /tmp/RtmpyceNDl/Rinst2b3079392d32
#> [2] /github/workspace/pkglib
#> [3] /usr/local/lib/R/site-library
#> [4] /usr/lib/R/site-library
#> [5] /usr/lib/R/library
#>
#> ββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββ
This vignette was generated using BiocStyle (OleΕ, 2024) with knitr (Xie, 2024) and rmarkdown (Allaire, Xie, Dervieux et al., 2024) running behind the scenes.
Citations made with RefManageR (McLean, 2017).
[1] J. Allaire, Y. Xie, C. Dervieux, et al. rmarkdown: Dynamic Documents for R. R package version 2.29. 2024. URL: https://github.com/rstudio/rmarkdown.
[2] T. Benaglia, D. Chauveau, D. R. Hunter, et al. βmixtools: An R Package for Analyzing Finite Mixture Modelsβ. In: Journal of Statistical Software 32.6 (2009), pp.Β 1β29. URL: https://www.jstatsoft.org/v32/i06/.
[3] M. L. Delignette-Muller and C. Dutang. βfitdistrplus: An R Package for Fitting Distributionsβ. In: Journal of Statistical Software 64.4 (2015), pp.Β 1β34. DOI: 10.18637/jss.v064.i04.
[4] M. W. McLean. βRefManageR: Import and Manage BibTeX and BibLaTeX References in Rβ. In: The Journal of Open Source Software (2017). DOI: 10.21105/joss.00338.
[5] P. J. McMurdie and S. Holmes. βphyloseq: An R package for reproducible interactive analysis and graphics of microbiome census dataβ. In: PLoS ONE 8.4 (2013), p.Β e61217. URL: http://dx.plos.org/10.1371/journal.pone.0061217.
[6] M. Morgan, J. Wang, V. Obenchain, et al. BiocParallel: Bioconductor facilities for parallel evaluation. R package version 1.41.0. 2024. URL: https://github.com/Bioconductor/BiocParallel.
[7] Q. Nguyen. CBEA: R package for performing CBEA approach. https://github.com/qpmnguyen/CBEA - R package version 1.7.0. 2024. DOI: 10.18129/B9.bioc.CBEA. URL: http://www.bioconductor.org/packages/CBEA.
[8] A. OleΕ. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.35.0. 2024. URL: https://github.com/Bioconductor/BiocStyle.
[9] R Core Team. R: A Language and Environment for Statistical Computing. R Foundation for Statistical Computing. Vienna, Austria, 2024. URL: https://www.R-project.org/.
[10] D. Robinson, A. Hayes, and S. Couch. broom: Convert Statistical Objects into Tidy Tibbles. R package version 1.0.7, https://github.com/tidymodels/broom. 2024. URL: https://broom.tidymodels.org/.
[11] H. Wickham. βtestthat: Get Started with Testingβ. In: The R Journal 3 (2011), pp.Β 5β10. URL: https://journal.r-project.org/archive/2011-1/RJournal_2011-1_Wickham.pdf.
[12] H. Wickham, M. Averick, J. Bryan, et al. βWelcome to the tidyverseβ. In: Journal of Open Source Software 4.43 (2019), p.Β 1686. DOI: 10.21105/joss.01686.
[13] H. Wickham, W. Chang, R. Flight, et al. sessioninfo: R Session Information. R package version 1.2.2, https://r-lib.github.io/sessioninfo/. 2021. URL: https://github.com/r-lib/sessioninfo#readme.
[14] Y. Xie. knitr: A General-Purpose Package for Dynamic Report Generation in R. R package version 1.49. 2024. URL: https://yihui.org/knitr/.