knitr::opts_chunk$set(message = FALSE, warning = FALSE, comment = NA,
fig.width = 6.25, fig.height = 5)
library(ANCOMBC)
library(tidyverse)
Analysis of Composition of Microbiomes (ANCOM) (Mandal et al. 2015) is a differential abundance (DA) analysis for microbial absolute abundances. It accounts for the compositionality of microbiome data by performing the additive log ratio (ALR) transformation. ANCOM employs a heuristic strategy to declare taxa that are significantly differentially abundant. For a given taxon, the output W statistic represents the number ALR transformed models where the taxon is differentially abundant with regard to the variable of interest. The larger the value of W, the more likely the taxon is differentially abundant. For more details, please refer to the ANCOM paper.
Download package.
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("ANCOMBC")
Load the package.
The HITChip Atlas dataset contains genus-level microbiota profiling with HITChip for 1006 western adults with no reported health complications, reported in (Lahti et al. 2014). The dataset is available via the microbiome R package (Lahti et al. 2017) in phyloseq (McMurdie and Holmes 2013) format. In this tutorial, we consider the following covariates:
Continuous covariates: “age”
Categorical covariates: “region”, “bmi”
The group variable of interest: “bmi”
Three groups: “lean”, “overweight”, “obese”
The reference group: “obese”
data(atlas1006, package = "microbiome")
tse = mia::makeTreeSummarizedExperimentFromPhyloseq(atlas1006)
# Subset to baseline
tse = tse[, tse$time == 0]
# Re-code the bmi group
tse$bmi = recode(tse$bmi_group,
obese = "obese",
severeobese = "obese",
morbidobese = "obese")
# Subset to lean, overweight, and obese subjects
tse = tse[, tse$bmi %in% c("lean", "overweight", "obese")]
# Note that by default, levels of a categorical variable in R are sorted
# alphabetically. In this case, the reference level for `bmi` will be
# `lean`. To manually change the reference level, for instance, setting `obese`
# as the reference level, use:
tse$bmi = factor(tse$bmi, levels = c("obese", "overweight", "lean"))
# You can verify the change by checking:
# levels(sample_data(tse)$bmi)
# Create the region variable
tse$region = recode(as.character(tse$nationality),
Scandinavia = "NE", UKIE = "NE", SouthEurope = "SE",
CentralEurope = "CE", EasternEurope = "EE",
.missing = "unknown")
# Discard "EE" as it contains only 1 subject
# Discard subjects with missing values of region
tse = tse[, ! tse$region %in% c("EE", "unknown")]
print(tse)
class: TreeSummarizedExperiment
dim: 130 873
metadata(0):
assays(1): counts
rownames(130): Actinomycetaceae Aerococcus ... Xanthomonadaceae
Yersinia et rel.
rowData names(3): Phylum Family Genus
colnames(873): Sample-1 Sample-2 ... Sample-1005 Sample-1006
colData names(12): age sex ... bmi region
reducedDimNames(0):
mainExpName: NULL
altExpNames(0):
rowLinks: NULL
rowTree: NULL
colLinks: NULL
colTree: NULL
set.seed(123)
out = ancom(data = tse, assay_name = "counts",
tax_level = "Family", phyloseq = NULL,
p_adj_method = "holm", prv_cut = 0.10, lib_cut = 1000,
main_var = "bmi", adj_formula = "age + region",
rand_formula = NULL, lme_control = NULL, struc_zero = TRUE,
neg_lb = TRUE, alpha = 0.05, n_cl = 2)
res = out$res
# Similarly, if the main variable of interest is continuous, such as age, the
# ancom model can be specified as
# out = ancom(data = tse, assay_name = "counts",
# tax_level = "Family", phyloseq = NULL,
# p_adj_method = "holm", prv_cut = 0.10, lib_cut = 1000,
# main_var = "age", adj_formula = "bmi + region",
# rand_formula = NULL, lme_control = NULL, struc_zero = FALSE,
# neg_lb = FALSE, alpha = 0.05, n_cl = 2)
# ancom also supports importing data in phyloseq format
# tse_alt = agglomerateByRank(tse, "Family")
# pseq = makePhyloseqFromTreeSummarizedExperiment(tse_alt)
# out = ancom(data = NULL, assay_name = NULL,
# tax_level = "Family", phyloseq = pseq,
# p_adj_method = "holm", prv_cut = 0.10, lib_cut = 1000,
# main_var = "bmi", adj_formula = "age + region",
# rand_formula = NULL, lme_control = NULL, struc_zero = TRUE,
# neg_lb = TRUE, alpha = 0.05, n_cl = 2)
q_val = out$q_data
beta_val = out$beta_data
# Only consider the effect sizes with the corresponding q-value less than alpha
beta_val = beta_val * (q_val < 0.05)
# Choose the maximum of beta's as the effect size
beta_pos = apply(abs(beta_val), 2, which.max)
beta_max = vapply(seq_along(beta_pos), function(i)
beta_val[beta_pos[i], i], FUN.VALUE = double(1))
# Number of taxa except structural zeros
n_taxa = ifelse(is.null(out$zero_ind),
nrow(tse),
sum(apply(out$zero_ind[, -1], 1, sum) == 0))
# Cutoff values for declaring differentially abundant taxa
cut_off = 0.7 * (n_taxa - 1)
df_fig_w = res %>%
dplyr::mutate(beta = beta_max,
direct = case_when(
detected_0.7 == TRUE & beta > 0 ~ "Positive",
detected_0.7 == TRUE & beta <= 0 ~ "Negative",
TRUE ~ "Not Significant"
)) %>%
dplyr::arrange(W)
df_fig_w$taxon = factor(df_fig_w$taxon, levels = df_fig_w$taxon)
df_fig_w$W = replace(df_fig_w$W, is.infinite(df_fig_w$W), n_taxa - 1)
df_fig_w$direct = factor(df_fig_w$direct,
levels = c("Negative", "Positive", "Not Significant"))
p_w = df_fig_w %>%
ggplot(aes(x = taxon, y = W, color = direct)) +
geom_point(size = 2, alpha = 0.6) +
labs(x = "Taxon", y = "W") +
scale_color_discrete(name = NULL) +
geom_hline(yintercept = cut_off, linetype = "dotted",
color = "blue", size = 1.5) +
geom_text(aes(x = 2, y = cut_off + 0.5, label = "W[0.7]"),
size = 5, vjust = -0.5, hjust = 0, color = "orange", parse = TRUE) +
theme_bw() +
theme(axis.text.x = element_blank(),
axis.ticks.x = element_blank(),
panel.grid.major = element_blank())
p_w
A two-week diet swap study between western (USA) and traditional (rural Africa) diets (Lahti et al. 2014). The dataset is available via the microbiome R package (Lahti et al. 2017) in phyloseq (McMurdie and Holmes 2013) format. In this tutorial, we consider the following fixed effects:
Continuous covariates: “timepoint”
Categorical covariates: “nationality”
The group variable of interest: “group”
Three groups: “DI”, “ED”, “HE”
The reference group: “DI”
and the following random effects:
A random intercept
A random slope: “timepoint”
data(dietswap, package = "microbiome")
tse = mia::makeTreeSummarizedExperimentFromPhyloseq(dietswap)
print(tse)
class: TreeSummarizedExperiment
dim: 130 222
metadata(0):
assays(1): counts
rownames(130): Actinomycetaceae Aerococcus ... Xanthomonadaceae
Yersinia et rel.
rowData names(3): Phylum Family Genus
colnames(222): Sample-1 Sample-2 ... Sample-221 Sample-222
colData names(8): subject sex ... timepoint.within.group bmi_group
reducedDimNames(0):
mainExpName: NULL
altExpNames(0):
rowLinks: NULL
rowTree: NULL
colLinks: NULL
colTree: NULL
set.seed(123)
out = ancom(data = tse, assay_name = "counts",
tax_level = "Family", phyloseq = NULL,
p_adj_method = "holm", prv_cut = 0.10, lib_cut = 1000,
main_var = "group",
adj_formula = "nationality + timepoint",
rand_formula = "(timepoint | subject)",
lme_control = lme4::lmerControl(),
struc_zero = TRUE, neg_lb = TRUE, alpha = 0.05, n_cl = 2)
res = out$res
q_val = out$q_data
beta_val = out$beta_data
# Only consider the effect sizes with the corresponding q-value less than alpha
beta_val = beta_val * (q_val < 0.05)
# Choose the maximum of beta's as the effect size
beta_pos = apply(abs(beta_val), 2, which.max)
beta_max = vapply(seq_along(beta_pos), function(i) beta_val[beta_pos[i], i],
FUN.VALUE = double(1))
# Number of taxa except structural zeros
n_taxa = ifelse(is.null(out$zero_ind),
nrow(tse),
sum(apply(out$zero_ind[, -1], 1, sum) == 0))
# Cutoff values for declaring differentially abundant taxa
cut_off = 0.7 * (n_taxa - 1)
df_fig_w = res %>%
dplyr::mutate(beta = beta_max,
direct = case_when(
detected_0.7 == TRUE & beta > 0 ~ "Positive",
detected_0.7 == TRUE & beta <= 0 ~ "Negative",
TRUE ~ "Not Significant"
)) %>%
dplyr::arrange(W)
df_fig_w$taxon = factor(df_fig_w$taxon, levels = df_fig_w$taxon)
df_fig_w$W = replace(df_fig_w$W, is.infinite(df_fig_w$W), n_taxa - 1)
df_fig_w$direct = factor(df_fig_w$direct,
levels = c("Negative", "Positive", "Not Significant"))
p_w = df_fig_w %>%
ggplot(aes(x = taxon, y = W, color = direct)) +
geom_point(size = 2, alpha = 0.6) +
labs(x = "Taxon", y = "W") +
scale_color_discrete(name = NULL) +
geom_hline(yintercept = cut_off, linetype = "dotted",
color = "blue", size = 1.5) +
geom_text(aes(x = 2, y = cut_off + 0.5, label = "W[0.7]"),
size = 5, vjust = -0.5, hjust = 0, color = "orange", parse = TRUE) +
theme_bw() +
theme(axis.text.x = element_blank(),
axis.ticks.x = element_blank(),
panel.grid.major = element_blank())
p_w
R version 4.4.0 RC (2024-04-16 r86468)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 22.04.4 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.20-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] stats graphics grDevices utils datasets methods base
other attached packages:
[1] phyloseq_1.49.0 lubridate_1.9.3 forcats_1.0.0 stringr_1.5.1
[5] dplyr_1.1.4 purrr_1.0.2 readr_2.1.5 tidyr_1.3.1
[9] tibble_3.2.1 ggplot2_3.5.1 tidyverse_2.0.0 ANCOMBC_2.7.0
loaded via a namespace (and not attached):
[1] splines_4.4.0 cellranger_1.1.0
[3] rpart_4.1.23 DirichletMultinomial_1.47.0
[5] lifecycle_1.0.4 Rdpack_2.6
[7] doParallel_1.0.17 lattice_0.22-6
[9] MASS_7.3-60.2 MultiAssayExperiment_1.31.0
[11] backports_1.4.1 magrittr_2.0.3
[13] Hmisc_5.1-2 sass_0.4.9
[15] rmarkdown_2.26 jquerylib_0.1.4
[17] yaml_2.3.8 doRNG_1.8.6
[19] gld_2.6.6 DBI_1.2.2
[21] minqa_1.2.6 ade4_1.7-22
[23] multcomp_1.4-25 abind_1.4-5
[25] zlibbioc_1.51.0 expm_0.999-9
[27] GenomicRanges_1.57.0 BiocGenerics_0.51.0
[29] yulab.utils_0.1.4 nnet_7.3-19
[31] TH.data_1.1-2 sandwich_3.1-0
[33] GenomeInfoDbData_1.2.12 IRanges_2.39.0
[35] S4Vectors_0.43.0 ggrepel_0.9.5
[37] irlba_2.3.5.1 tidytree_0.4.6
[39] vegan_2.6-4 permute_0.9-7
[41] DelayedMatrixStats_1.27.0 codetools_0.2-20
[43] DelayedArray_0.31.0 scuttle_1.15.0
[45] energy_1.7-11 tidyselect_1.2.1
[47] farver_2.1.1 UCSC.utils_1.1.0
[49] lme4_1.1-35.3 gmp_0.7-4
[51] ScaledMatrix_1.13.0 viridis_0.6.5
[53] matrixStats_1.3.0 stats4_4.4.0
[55] base64enc_0.1-3 jsonlite_1.8.8
[57] multtest_2.61.0 BiocNeighbors_1.23.0
[59] e1071_1.7-14 decontam_1.25.0
[61] mia_1.13.0 Formula_1.2-5
[63] survival_3.6-4 scater_1.33.0
[65] iterators_1.0.14 foreach_1.5.2
[67] tools_4.4.0 treeio_1.29.0
[69] DescTools_0.99.54 Rcpp_1.0.12
[71] glue_1.7.0 gridExtra_2.3
[73] SparseArray_1.5.0 xfun_0.43
[75] mgcv_1.9-1 MatrixGenerics_1.17.0
[77] GenomeInfoDb_1.41.0 TreeSummarizedExperiment_2.13.0
[79] withr_3.0.0 numDeriv_2016.8-1.1
[81] fastmap_1.1.1 rhdf5filters_1.17.0
[83] boot_1.3-30 bluster_1.15.0
[85] fansi_1.0.6 digest_0.6.35
[87] rsvd_1.0.5 timechange_0.3.0
[89] R6_2.5.1 colorspace_2.1-0
[91] gtools_3.9.5 utf8_1.2.4
[93] generics_0.1.3 data.table_1.15.4
[95] DECIPHER_3.1.0 class_7.3-22
[97] CVXR_1.0-12 httr_1.4.7
[99] htmlwidgets_1.6.4 S4Arrays_1.5.0
[101] pkgconfig_2.0.3 gtable_0.3.5
[103] Exact_3.2 Rmpfr_0.9-5
[105] SingleCellExperiment_1.27.0 XVector_0.45.0
[107] htmltools_0.5.8.1 biomformat_1.33.0
[109] scales_1.3.0 Biobase_2.65.0
[111] lmom_3.0 knitr_1.46
[113] rstudioapi_0.16.0 tzdb_0.4.0
[115] reshape2_1.4.4 checkmate_2.3.1
[117] nlme_3.1-164 nloptr_2.0.3
[119] rhdf5_2.49.0 proxy_0.4-27
[121] cachem_1.0.8 zoo_1.8-12
[123] rootSolve_1.8.2.4 parallel_4.4.0
[125] vipor_0.4.7 foreign_0.8-86
[127] pillar_1.9.0 grid_4.4.0
[129] vctrs_0.6.5 BiocSingular_1.21.0
[131] beachmat_2.21.0 cluster_2.1.6
[133] beeswarm_0.4.0 htmlTable_2.4.2
[135] evaluate_0.23 mvtnorm_1.2-4
[137] cli_3.6.2 compiler_4.4.0
[139] rlang_1.1.3 crayon_1.5.2
[141] rngtools_1.5.2 labeling_0.4.3
[143] plyr_1.8.9 fs_1.6.4
[145] ggbeeswarm_0.7.2 stringi_1.8.3
[147] viridisLite_0.4.2 BiocParallel_1.39.0
[149] lmerTest_3.1-3 munsell_0.5.1
[151] Biostrings_2.73.0 gsl_2.1-8
[153] lazyeval_0.2.2 Matrix_1.7-0
[155] hms_1.1.3 sparseMatrixStats_1.17.0
[157] bit64_4.0.5 Rhdf5lib_1.27.0
[159] highr_0.10 SummarizedExperiment_1.35.0
[161] rbibutils_2.2.16 igraph_2.0.3
[163] memoise_2.0.1 bslib_0.7.0
[165] bit_4.0.5 readxl_1.4.3
[167] ape_5.8
Lahti, Leo, Jarkko Salojärvi, Anne Salonen, Marten Scheffer, and Willem M De Vos. 2014. “Tipping Elements in the Human Intestinal Ecosystem.” Nature Communications 5 (1): 1–10.
Lahti, Leo, Sudarshan Shetty, T Blake, J Salojarvi, and others. 2017. “Tools for Microbiome Analysis in R.” Version 1: 10013.
Mandal, Siddhartha, Will Van Treuren, Richard A White, Merete Eggesbø, Rob Knight, and Shyamal D Peddada. 2015. “Analysis of Composition of Microbiomes: A Novel Method for Studying Microbial Composition.” Microbial Ecology in Health and Disease 26 (1): 27663.
McMurdie, Paul J, and Susan Holmes. 2013. “Phyloseq: An R Package for Reproducible Interactive Analysis and Graphics of Microbiome Census Data.” PloS One 8 (4): e61217.