ANCOM Tutorial

knitr::opts_chunk$set(message = FALSE, warning = FALSE, comment = NA,
                      fig.width = 6.25, fig.height = 5)
library(ANCOMBC)
library(tidyverse)

1. Introduction

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.

2. Installation

Download package.

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

Load the package.

library(ANCOMBC)

3. Run ANCOM on a real cross-sectional dataset

3.1 Import example data

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")

# Subset to baseline
pseq = phyloseq::subset_samples(atlas1006, time == 0)

# Re-code the bmi group
meta_data = microbiome::meta(pseq)
meta_data$bmi = recode(meta_data$bmi_group,
                       obese = "obese",
                       severeobese = "obese",
                       morbidobese = "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:
meta_data$bmi = factor(meta_data$bmi, levels = c("obese", "overweight", "lean"))
# You can verify the change by checking:
# levels(meta_data$bmi)

# Create the region variable
meta_data$region = recode(as.character(meta_data$nationality),
                          Scandinavia = "NE", UKIE = "NE", SouthEurope = "SE", 
                          CentralEurope = "CE", EasternEurope = "EE",
                          .missing = "unknown")

phyloseq::sample_data(pseq) = meta_data

# Subset to lean, overweight, and obese subjects
pseq = phyloseq::subset_samples(pseq, bmi %in% c("lean", "overweight", "obese"))
# Discard "EE" as it contains only 1 subject
# Discard subjects with missing values of region
pseq = phyloseq::subset_samples(pseq, ! region %in% c("EE", "unknown"))

print(pseq)
phyloseq-class experiment-level object
otu_table()   OTU Table:         [ 130 taxa and 873 samples ]
sample_data() Sample Data:       [ 873 samples by 12 sample variables ]
tax_table()   Taxonomy Table:    [ 130 taxa by 3 taxonomic ranks ]

3.2 Run ancom function using phyloseq data

set.seed(123)
out = ancom(data = pseq, tax_level = "Family", meta_data = 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, verbose = TRUE)

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 = pseq, tax_level = "Family", meta_data = 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, verbose = TRUE)

3.3 Scatter plot for W statistics

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

3.4 Run ancom function using tse data

tse = mia::makeTreeSummarizedExperimentFromPhyloseq(pseq)

set.seed(123)
out = ancom(data = tse, assay_name = "counts",
            tax_level = "Family", meta_data = 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, verbose = TRUE)

res = out$res

3.5 Run ancom function by directly providing the abundance and metadata

abundance_data = microbiome::abundances(pseq)
aggregate_data = microbiome::abundances(microbiome::aggregate_taxa(pseq, "Family"))
meta_data = microbiome::meta(pseq)

set.seed(123)
out = ancom(data = abundance_data, aggregate_data = aggregate_data, 
            meta_data = meta_data, 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, verbose = TRUE)

res = out$res

4. Run ANCOM on a real longitudinal dataset

4.1 Import example data

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")
print(dietswap)
phyloseq-class experiment-level object
otu_table()   OTU Table:         [ 130 taxa and 222 samples ]
sample_data() Sample Data:       [ 222 samples by 8 sample variables ]
tax_table()   Taxonomy Table:    [ 130 taxa by 3 taxonomic ranks ]

4.2 Run ancom function using phyloseq data

set.seed(123)
out = ancom(data = dietswap, tax_level = "Family", 
            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

4.3 Visualization for W statistics

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

4.4 Run ancom function using tse data

tse = mia::makeTreeSummarizedExperimentFromPhyloseq(dietswap)

set.seed(123)
out = ancom(data = tse, assay_name = "counts", tax_level = "Family",
            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

4.5 Run ancom function by directly providing the abundance and metadata

abundance_data = microbiome::abundances(dietswap)
aggregate_data = microbiome::abundances(microbiome::aggregate_taxa(dietswap, "Family"))
meta_data = microbiome::meta(dietswap)

set.seed(123)
out = ancom(data = abundance_data, aggregate_data = aggregate_data,
            meta_data = meta_data, 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

Session information

sessionInfo()
R version 4.4.2 (2024-10-31)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.1 LTS

Matrix products: default
BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so;  LAPACK version 3.12.0

locale:
 [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
 [3] LC_TIME=en_US.UTF-8        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: Etc/UTC
tzcode source: system (glibc)

attached base packages:
[1] stats     graphics  grDevices utils     datasets  methods   base     

other attached packages:
 [1] lubridate_1.9.3 forcats_1.0.0   stringr_1.5.1   dplyr_1.1.4    
 [5] purrr_1.0.2     readr_2.1.5     tidyr_1.3.1     tibble_3.2.1   
 [9] ggplot2_3.5.1   tidyverse_2.0.0 ANCOMBC_2.9.0   rmarkdown_2.29 

loaded via a namespace (and not attached):
  [1] sys_3.4.3               rstudioapi_0.17.1       jsonlite_1.8.9         
  [4] magrittr_2.0.3          TH.data_1.1-2           farver_2.1.2           
  [7] nloptr_2.1.1            zlibbioc_1.52.0         vctrs_0.6.5            
 [10] multtest_2.63.0         minqa_1.2.8             base64enc_0.1-3        
 [13] htmltools_0.5.8.1       energy_1.7-12           haven_2.5.4            
 [16] Rhdf5lib_1.29.0         cellranger_1.1.0        Formula_1.2-5          
 [19] rhdf5_2.51.0            sass_0.4.9              bslib_0.8.0            
 [22] htmlwidgets_1.6.4       plyr_1.8.9              sandwich_3.1-1         
 [25] rootSolve_1.8.2.4       zoo_1.8-12              cachem_1.1.0           
 [28] buildtools_1.0.0        igraph_2.1.1            lifecycle_1.0.4        
 [31] iterators_1.0.14        pkgconfig_2.0.3         Matrix_1.7-1           
 [34] R6_2.5.1                fastmap_1.2.0           GenomeInfoDbData_1.2.13
 [37] rbibutils_2.3           digest_0.6.37           Exact_3.3              
 [40] numDeriv_2016.8-1.1     colorspace_2.1-1        S4Vectors_0.45.2       
 [43] Hmisc_5.2-0             vegan_2.6-8             labeling_0.4.3         
 [46] fansi_1.0.6             timechange_0.3.0        mgcv_1.9-1             
 [49] httr_1.4.7              compiler_4.4.2          rngtools_1.5.2         
 [52] proxy_0.4-27            bit64_4.5.2             withr_3.0.2            
 [55] doParallel_1.0.17       gsl_2.1-8               htmlTable_2.4.3        
 [58] backports_1.5.0         MASS_7.3-61             biomformat_1.35.0      
 [61] permute_0.9-7           gtools_3.9.5            CVXR_1.0-15            
 [64] gld_2.6.6               tools_4.4.2             foreign_0.8-87         
 [67] ape_5.8                 nnet_7.3-19             glue_1.8.0             
 [70] rhdf5filters_1.19.0     nlme_3.1-166            grid_4.4.2             
 [73] Rtsne_0.17              checkmate_2.3.2         cluster_2.1.6          
 [76] reshape2_1.4.4          ade4_1.7-22             generics_0.1.3         
 [79] microbiome_1.29.0       gtable_0.3.6            tzdb_0.4.0             
 [82] class_7.3-22            data.table_1.16.2       lmom_3.2               
 [85] hms_1.1.3               utf8_1.2.4              XVector_0.47.0         
 [88] BiocGenerics_0.53.3     foreach_1.5.2           pillar_1.9.0           
 [91] splines_4.4.2           lattice_0.22-6          survival_3.7-0         
 [94] gmp_0.7-5               bit_4.5.0               tidyselect_1.2.1       
 [97] maketools_1.3.1         Biostrings_2.75.1       knitr_1.49             
[100] gridExtra_2.3           phyloseq_1.51.0         IRanges_2.41.1         
[103] stats4_4.4.2            xfun_0.49               expm_1.0-0             
[106] Biobase_2.67.0          stringi_1.8.4           UCSC.utils_1.3.0       
[109] yaml_2.3.10             boot_1.3-31             evaluate_1.0.1         
[112] codetools_0.2-20        cli_3.6.3               rpart_4.1.23           
[115] DescTools_0.99.58       Rdpack_2.6.2            munsell_0.5.1          
[118] jquerylib_0.1.4         Rcpp_1.0.13-1           GenomeInfoDb_1.43.1    
[121] readxl_1.4.3            parallel_4.4.2          doRNG_1.8.6            
[124] lme4_1.1-35.5           Rmpfr_1.0-0             mvtnorm_1.3-2          
[127] lmerTest_3.1-3          scales_1.3.0            e1071_1.7-16           
[130] crayon_1.5.3            rlang_1.1.4             multcomp_1.4-26        

References

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, et al. 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.