ANCOM-BC Tutorial

Huang Lin\(^1\)

\(^1\)NIEHS, Research Triangle Park, NC 27709, USA

October 24, 2023

knitr::opts_chunk$set(message = FALSE, warning = FALSE, comment = NA, 
                      fig.width = 6.25, fig.height = 5)
library(ANCOMBC)
library(tidyverse)
library(DT)
options(DT.options = list(
  initComplete = JS("function(settings, json) {",
  "$(this.api().table().header()).css({'background-color': 
  '#000', 'color': '#fff'});","}")))

1. Introduction

Analysis of Compositions of Microbiomes with Bias Correction (ANCOM-BC) (Lin and Peddada 2020) is a methodology of differential abundance (DA) analysis for microbial absolute abundances. ANCOM-BC estimates the unknown sampling fractions, corrects the bias induced by their differences through a log linear regression model including the estimated sampling fraction as an offset terms, and identifies taxa that are differentially abundant according to the variable of interest. For more details, please refer to the ANCOM-BC paper.

2. Installation

Download package.

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

Load the package.

library(ANCOMBC)

3. 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:

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

4 ANCOM-BC Implementation

4.1 Run ancombc function

out = ancombc(data = tse, assay_name = "counts", 
              tax_level = "Family", phyloseq = NULL, 
              formula = "age + region + bmi", 
              p_adj_method = "holm", prv_cut = 0.10, lib_cut = 1000, 
              group = "bmi", struc_zero = TRUE, neg_lb = TRUE, tol = 1e-5, 
              max_iter = 100, conserve = TRUE, alpha = 0.05, global = TRUE,
              n_cl = 1, verbose = TRUE)

res = out$res
res_global = out$res_global

# ancombc also supports importing data in phyloseq format
# tse_alt = agglomerateByRank(tse, "Family")
# pseq = makePhyloseqFromTreeSummarizedExperiment(tse_alt)
# out = ancombc(data = NULL, assay_name = NULL,
#               tax_level = "Family", phyloseq = pseq,
#               formula = "age + region + bmi",
#               p_adj_method = "holm", prv_cut = 0.10, lib_cut = 1000,
#               group = "region", struc_zero = TRUE, neg_lb = TRUE, tol = 1e-5,
#               max_iter = 100, conserve = TRUE, alpha = 0.05, global = TRUE,
#               n_cl = 1, verbose = TRUE)

4.2 ANCOMBC primary result

Result from the ANCOM-BC log-linear model to determine taxa that are differentially abundant according to the covariate of interest. It contains: 1) log fold changes; 2) standard errors; 3) test statistics; 4) p-values; 5) adjusted p-values; 6) indicators whether the taxon is differentially abundant (TRUE) or not (FALSE).

LFC

tab_lfc = res$lfc
col_name = c("Taxon", "Intercept", "Age", "NE - CE", "SE - CE", 
             "US - CE", "Overweight - Obese", "Lean - Obese")
colnames(tab_lfc) = col_name
tab_lfc %>% 
  datatable(caption = "Log Fold Changes from the Primary Result") %>%
  formatRound(col_name[-1], digits = 2)

SE

tab_se = res$se
colnames(tab_se) = col_name
tab_se %>% 
  datatable(caption = "SEs from the Primary Result") %>%
  formatRound(col_name[-1], digits = 2)

Test statistic

tab_w = res$W
colnames(tab_w) = col_name
tab_w %>% 
  datatable(caption = "Test Statistics from the Primary Result") %>%
  formatRound(col_name[-1], digits = 2)

P-values

tab_p = res$p_val
colnames(tab_p) = col_name
tab_p %>% 
  datatable(caption = "P-values from the Primary Result") %>%
  formatRound(col_name[-1], digits = 2)

Adjusted p-values

tab_q = res$q
colnames(tab_q) = col_name
tab_q %>% 
  datatable(caption = "Adjusted p-values from the Primary Result") %>%
  formatRound(col_name[-1], digits = 2)

Differentially abundant taxa

tab_diff = res$diff_abn
colnames(tab_diff) = col_name
tab_diff %>% 
  datatable(caption = "Differentially Abundant Taxa from the Primary Result")

Bias-corrected abundances

To obtain bias-corrected abundances, the following steps can be taken:

Step 1: Calculate the estimated sample-specific sampling fractions, in log scale.

Step 2: Correct the log observed abundances by subtracting the estimated sampling fraction from the log observed abundances of each sample.

It is important to note that we can only estimate sampling fractions up to an additive constant, meaning that only the difference between bias-corrected abundances is meaningful. Additionally, taxon-specific biases are not taken into account in the calculation of bias-corrected abundances, as it is assumed that these biases vary across taxa but remain constant across samples within a taxon.

samp_frac = out$samp_frac
# Replace NA with 0
samp_frac[is.na(samp_frac)] = 0 
# Add pesudo-count (1) to avoid taking the log of 0
log_obs_abn = log(out$feature_table + 1)
# Adjust the log observed abundances
log_corr_abn = t(t(log_obs_abn) - samp_frac)
# Show the first 6 samples
round(log_corr_abn[, 1:6], 2) %>% 
  datatable(caption = "Bias-corrected log observed abundances")

Visualization for age

df_lfc = data.frame(res$lfc[, -1] * res$diff_abn[, -1], check.names = FALSE) %>%
    mutate(taxon_id = res$diff_abn$taxon) %>%
    dplyr::select(taxon_id, everything())
df_se = data.frame(res$se[, -1] * res$diff_abn[, -1], check.names = FALSE) %>% 
  mutate(taxon_id = res$diff_abn$taxon) %>%
    dplyr::select(taxon_id, everything())
colnames(df_se)[-1] = paste0(colnames(df_se)[-1], "SE")

df_fig_age = df_lfc %>% 
  dplyr::left_join(df_se, by = "taxon_id") %>%
  dplyr::transmute(taxon_id, age, ageSE) %>%
  dplyr::filter(age != 0) %>% 
  dplyr::arrange(desc(age)) %>%
  dplyr::mutate(direct = ifelse(age > 0, "Positive LFC", "Negative LFC"))
df_fig_age$taxon_id = factor(df_fig_age$taxon_id, levels = df_fig_age$taxon_id)
df_fig_age$direct = factor(df_fig_age$direct, 
                        levels = c("Positive LFC", "Negative LFC"))
  
p_age = ggplot(data = df_fig_age, 
           aes(x = taxon_id, y = age, fill = direct, color = direct)) + 
  geom_bar(stat = "identity", width = 0.7, 
           position = position_dodge(width = 0.4)) +
  geom_errorbar(aes(ymin = age - ageSE, ymax = age + ageSE), width = 0.2,
                position = position_dodge(0.05), color = "black") + 
  labs(x = NULL, y = "Log fold change", 
       title = "Log fold changes as one unit increase of age") + 
  scale_fill_discrete(name = NULL) +
  scale_color_discrete(name = NULL) +
  theme_bw() + 
  theme(plot.title = element_text(hjust = 0.5),
        panel.grid.minor.y = element_blank(),
        axis.text.x = element_text(angle = 60, hjust = 1))
p_age

Visualization for BMI

df_fig_bmi = df_lfc %>% 
  filter(bmioverweight != 0 | bmilean != 0) %>%
  transmute(taxon_id, 
            `Overweight vs. Obese` = round(bmioverweight, 2),
            `Lean vs. Obese` = round(bmilean, 2)) %>%
  pivot_longer(cols = `Overweight vs. Obese`:`Lean vs. Obese`, 
               names_to = "group", values_to = "value") %>%
  arrange(taxon_id)
lo = floor(min(df_fig_bmi$value))
up = ceiling(max(df_fig_bmi$value))
mid = (lo + up)/2
p_bmi = df_fig_bmi %>%
  ggplot(aes(x = group, y = taxon_id, fill = value)) + 
  geom_tile(color = "black") +
  scale_fill_gradient2(low = "blue", high = "red", mid = "white", 
                       na.value = "white", midpoint = mid, limit = c(lo, up),
                       name = NULL) +
  geom_text(aes(group, taxon_id, label = value), color = "black", size = 4) +
  labs(x = NULL, y = NULL, title = "Log fold changes as compared to obese subjects") +
  theme_minimal() +
  theme(plot.title = element_text(hjust = 0.5))
p_bmi

4.3 ANCOMBC global test result

Result from the ANCOM-BC global test to determine taxa that are differentially abundant between at least two groups across three or more different groups. In this example, we want to identify taxa that are differentially abundant between at least two regions across CE, NE, SE, and US. The result contains: 1) test statistics; 2) p-values; 3) adjusted p-values; 4) indicators whether the taxon is differentially abundant (TRUE) or not (FALSE).

Test statistics

tab_w = res_global[, c("taxon", "W")]
tab_w %>% datatable(caption = "Test Statistics 
                    from the Global Test Result") %>%
      formatRound(c("W"), digits = 2)

P-values

tab_p = res_global[, c("taxon", "p_val")]
tab_p %>% datatable(caption = "P-values 
                    from the Global Test Result") %>%
      formatRound(c("p_val"), digits = 2)

Adjusted p-values

tab_q = res_global[, c("taxon", "q_val")]
tab_q %>% datatable(caption = "Adjusted p-values 
                    from the Global Test Result") %>%
      formatRound(c("q_val"), digits = 2)

Differentially abundant taxa

tab_diff = res_global[, c("taxon", "diff_abn")]
tab_diff %>% datatable(caption = "Differentially Abundant Taxa 
                       from the Global Test Result")

Visualization

sig_taxa = res_global %>%
  dplyr::filter(diff_abn == TRUE) %>%
  .$taxon

df_bmi = tab_lfc %>%
    dplyr::select(Taxon, `Overweight - Obese`, `Lean - Obese`) %>%
    filter(Taxon %in% sig_taxa)

df_heat = df_bmi %>%
    pivot_longer(cols = -one_of("Taxon"),
                 names_to = "region", values_to = "value") %>%
    mutate(value = round(value, 2))
df_heat$Taxon = factor(df_heat$Taxon, levels = sort(sig_taxa))

lo = floor(min(df_heat$value))
up = ceiling(max(df_heat$value))
mid = (lo + up)/2
p_heat = df_heat %>%
  ggplot(aes(x = region, y = Taxon, fill = value)) + 
  geom_tile(color = "black") +
  scale_fill_gradient2(low = "blue", high = "red", mid = "white", 
                       na.value = "white", midpoint = mid, limit = c(lo, up),
                       name = NULL) +
  geom_text(aes(region, Taxon, label = value), color = "black", size = 4) +
  labs(x = NULL, y = NULL, 
       title = "Log fold changes for globally significant taxa") +
  theme_minimal() +
  theme(plot.title = element_text(hjust = 0.5))
p_heat

Session information

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] stats     graphics  grDevices utils     datasets  methods   base     

other attached packages:
 [1] doRNG_1.8.6                     rngtools_1.5.2                 
 [3] foreach_1.5.2                   DT_0.30                        
 [5] TreeSummarizedExperiment_2.10.0 GenomicRanges_1.54.0           
 [7] SummarizedExperiment_1.32.0     SingleCellExperiment_1.24.0    
 [9] IRanges_2.36.0                  S4Vectors_0.40.0               
[11] phyloseq_1.46.0                 lubridate_1.9.3                
[13] forcats_1.0.0                   stringr_1.5.0                  
[15] dplyr_1.1.3                     purrr_1.0.2                    
[17] readr_2.1.4                     tidyr_1.3.0                    
[19] tibble_3.2.1                    ggplot2_3.4.4                  
[21] tidyverse_2.0.0                 ANCOMBC_2.4.0                  

loaded via a namespace (and not attached):
  [1] splines_4.3.1               bitops_1.0-7               
  [3] cellranger_1.1.0            rpart_4.1.21               
  [5] DirichletMultinomial_1.44.0 lifecycle_1.0.3            
  [7] Rdpack_2.5                  doParallel_1.0.17          
  [9] lattice_0.22-5              MASS_7.3-60                
 [11] crosstalk_1.2.0             MultiAssayExperiment_1.28.0
 [13] backports_1.4.1             magrittr_2.0.3             
 [15] Hmisc_5.1-1                 sass_0.4.7                 
 [17] rmarkdown_2.25              jquerylib_0.1.4            
 [19] yaml_2.3.7                  gld_2.6.6                  
 [21] DBI_1.1.3                   minqa_1.2.6                
 [23] ade4_1.7-22                 multcomp_1.4-25            
 [25] abind_1.4-5                 zlibbioc_1.48.0            
 [27] expm_0.999-7                BiocGenerics_0.48.0        
 [29] RCurl_1.98-1.12             TH.data_1.1-2              
 [31] yulab.utils_0.1.0           nnet_7.3-19                
 [33] sandwich_3.0-2              GenomeInfoDbData_1.2.11    
 [35] ggrepel_0.9.4               irlba_2.3.5.1              
 [37] tidytree_0.4.5              vegan_2.6-4                
 [39] permute_0.9-7               DelayedMatrixStats_1.24.0  
 [41] codetools_0.2-19            DelayedArray_0.28.0        
 [43] scuttle_1.12.0              energy_1.7-11              
 [45] tidyselect_1.2.0            farver_2.1.1               
 [47] lme4_1.1-34                 gmp_0.7-2                  
 [49] ScaledMatrix_1.10.0         viridis_0.6.4              
 [51] matrixStats_1.0.0           stats4_4.3.1               
 [53] base64enc_0.1-3             jsonlite_1.8.7             
 [55] multtest_2.58.0             BiocNeighbors_1.20.0       
 [57] e1071_1.7-13                ellipsis_0.3.2             
 [59] decontam_1.22.0             mia_1.10.0                 
 [61] Formula_1.2-5               survival_3.5-7             
 [63] scater_1.30.0               iterators_1.0.14           
 [65] tools_4.3.1                 treeio_1.26.0              
 [67] DescTools_0.99.50           Rcpp_1.0.11                
 [69] glue_1.6.2                  gridExtra_2.3              
 [71] SparseArray_1.2.0           xfun_0.40                  
 [73] mgcv_1.9-0                  MatrixGenerics_1.14.0      
 [75] GenomeInfoDb_1.38.0         withr_2.5.1                
 [77] numDeriv_2016.8-1.1         fastmap_1.1.1              
 [79] rhdf5filters_1.14.0         boot_1.3-28.1              
 [81] bluster_1.12.0              fansi_1.0.5                
 [83] digest_0.6.33               rsvd_1.0.5                 
 [85] timechange_0.2.0            R6_2.5.1                   
 [87] colorspace_2.1-0            gtools_3.9.4               
 [89] RSQLite_2.3.1               utf8_1.2.4                 
 [91] generics_0.1.3              data.table_1.14.8          
 [93] DECIPHER_2.30.0             class_7.3-22               
 [95] CVXR_1.0-11                 httr_1.4.7                 
 [97] htmlwidgets_1.6.2           S4Arrays_1.2.0             
 [99] pkgconfig_2.0.3             gtable_0.3.4               
[101] Exact_3.2                   Rmpfr_0.9-3                
[103] blob_1.2.4                  XVector_0.42.0             
[105] htmltools_0.5.6.1           biomformat_1.30.0          
[107] scales_1.2.1                Biobase_2.62.0             
[109] lmom_3.0                    knitr_1.44                 
[111] rstudioapi_0.15.0           tzdb_0.4.0                 
[113] reshape2_1.4.4              checkmate_2.2.0            
[115] nlme_3.1-163                nloptr_2.0.3               
[117] rhdf5_2.46.0                proxy_0.4-27               
[119] cachem_1.0.8                zoo_1.8-12                 
[121] rootSolve_1.8.2.4           parallel_4.3.1             
[123] vipor_0.4.5                 foreign_0.8-85             
[125] pillar_1.9.0                grid_4.3.1                 
[127] vctrs_0.6.4                 BiocSingular_1.18.0        
[129] beachmat_2.18.0             cluster_2.1.4              
[131] beeswarm_0.4.0              htmlTable_2.4.1            
[133] evaluate_0.22               mvtnorm_1.2-3              
[135] cli_3.6.1                   compiler_4.3.1             
[137] rlang_1.1.1                 crayon_1.5.2               
[139] labeling_0.4.3              plyr_1.8.9                 
[141] fs_1.6.3                    ggbeeswarm_0.7.2           
[143] stringi_1.7.12              viridisLite_0.4.2          
[145] BiocParallel_1.36.0         lmerTest_3.1-3             
[147] munsell_0.5.0               Biostrings_2.70.0          
[149] gsl_2.1-8                   lazyeval_0.2.2             
[151] Matrix_1.6-1.1              hms_1.1.3                  
[153] sparseMatrixStats_1.14.0    bit64_4.0.5                
[155] Rhdf5lib_1.24.0             rbibutils_2.2.15           
[157] igraph_1.5.1                memoise_2.0.1              
[159] bslib_0.5.1                 bit_4.0.5                  
[161] readxl_1.4.3                ape_5.7-1                  

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, and others. 2017. “Tools for Microbiome Analysis in R.” Version 1: 10013.

Lin, Huang, and Shyamal Das Peddada. 2020. “Analysis of Compositions of Microbiomes with Bias Correction.” Nature Communications 11 (1): 1–11.

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.