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'});","}")))
# It appears to be a package compatibility issue between the release version of 
# phyloseq and lme4, a fresh installation of phyloseq might be needed
# See this post: https://github.com/lme4/lme4/issues/743
# remotes::install_github("joey711/phyloseq", force = TRUE)Analysis of Compositions of Microbiomes with Bias Correction 2 (ANCOM-BC2) is a methodology for performing differential abundance (DA) analysis of microbiome count data. This version extends and refines the previously published Analysis of Compositions of Microbiomes with Bias Correction (ANCOM-BC) methodology (Lin and Peddada 2020) in several ways as follows:
Bias correction: ANCOM-BC2 estimates and corrects both the sample-specific (sampling fraction) as well as taxon-specific (sequencing efficiency) biases.
Regularization of variance: Inspired by Significance Analysis of Microarrays (SAM) (Tusher, Tibshirani, and Chu 2001) methodology, a small positive constant is added to the denominator of ANCOM-BC2 test statistic corresponding to each taxon to avoid the significance due to extremely small standard errors, especially for rare taxa. By default, we used the 5-th percentile of the distribution of standard errors for each fixed effect as the regularization factor.
Sensitivity analysis for the pseudo-count addition: Like other differential abundance analysis methods, ANCOM-BC2 applies a log transformation to the observed counts. However, the presence of zero counts poses a challenge, and researchers often consider adding a pseudo-count before the log transformation. However, it has been shown that the choice of pseudo-count can impact the results and lead to an inflated false positive rate (Costea et al. 2014; Paulson, Bravo, and Pop 2014). To address this issue, we conduct a sensitivity analysis to assess the impact of different pseudo-counts on zero counts for each taxon. This involves adding various pseudo-count values (0.1, 0.5, and 1) to the feature table and re-running the ANCOM-BC2 algorithm. The sensitivity score for each taxon is calculated as the proportion of times that the q-value exceeds the specified significance level (alpha). If all q-values consistently show significance or nonsignificance across different pseudo-counts and are consistent with the results obtained without adding pseudo-counts to zero counts (using the complete data), then the taxon is considered not sensitive to the pseudo-count addition. Important: The sensitivity analysis results depend directly on your chosen alpha value (since this threshold determines significance). You should expect sensitivity analysis results to change when you modify the alpha value.
Multi-group comparisons and repeated measurements: The ANCOM-BC2 methodology extends ANCOM-BC for multiple groups and repeated measurements as follows:
Multiple pairwise comparisons: When performning multiple pairwise comparisons, the mixed directional false discover rate (mdFDR) should be taken into account. The mdFDR is the combination of false discovery rate due to multiple testing, multiple pairwise comparisons, and directional tests within each pairwise comparison. For example, suppose we have five taxa and three experimental groups: g1, g2, and g3. Thus, we are performing five tests corresponding to five taxa. For each taxon, we are also conducting three pairwise comparisons (g1 vs. g2, g2 vs. g3, and g1 vs. g3). Within each pairwise comparison, we wish to determine if the abundance has increased or decreased or did not change (direction of the effect size). Errors could occur in each step. The overall false discovery rate is controlled by the mdFDR methodology we adopted from (Guo, Sarkar, and Peddada 2010; Grandhi, Guo, and Peddada 2016).
Multiple pairwise comparisons against a pre-specified group (e.g., Dunnett’s type of test): We use the same set-up as in the multiple pairwise comparisons but use the Dunnett-type modification described in (Grandhi, Guo, and Peddada 2016) to control the mdFDR which is more powerful.
Pattern analysis for ordered groups: In some instances, researchers are interested in discovering abundance patterns of each taxon over the ordered groups, for example, groups based on the health condition of subjects (e.g., lean, overweight, obese). In such cases, in addition to pairwise comparison, one may be interested in identifying taxa whose abundances are increasing or decreasing or have other patterns over the groups. We adopted methodologies from (Jelsema and Peddada 2016) to perform pattern analysis under the ANCOM-BC2 framework.
A clarification regarding Structural zeros: A taxon is considered to have structural zeros in some (>=1) groups if it is completely (or nearly completely) missing in these groups. For instance, suppose there are three groups: g1, g2, and g3. If the counts of taxon A in g1 are 0, but they are nonzero in g2 and g3, then taxon A will be considered to contain structural zeros in g1. In this example, taxon A is declared to be differentially abundant between g1 and g2, g1 and g3, and consequently, it is globally differentially abundant with respect to this group variable. Such taxa are not further analyzed using ANCOM-BC2, but the results are summarized in a separate table.
Download the package.
if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("ANCOMBC")Load the package.
The simulated data contain: one continuous covariate: cont_cov, and one categorical covariate: cat_cov. The categorical covariate has three levels: “1”, “2”, and “3”. Let’s focus on the discussions on the group variable: cat_cov.
The true abundances are generated using the Poisson lognormal (PLN) model based on the mechanism described in the LDM paper (Hu and Satten 2020). The PLN model relates the abundance vector \(Y_i\) with a Gaussian latent vector \(Z_i\) for taxon \(i\) as follows: \[ \text{latent layer: } Z_i \sim N(\mu_i, \Sigma_i) \\ \text{observation layer: } Y_i|Z_i \sim POI(N \exp(Z_i)) \] where \(N\) is the scaling factor. Because of the presence of a latent layer, the PLN model displays a larger variance than the Poisson model (over-dispersion). Also, the covariance (correlation) between abundances has the same sign as the covariance (correlation) between the corresponding latent variables. This property gives enormous flexibility in modeling the variance-covariance structure of microbial abundances since it is easy to specify different variance-covariance matrices in the multivariate Gaussian distribution.
Instead of specifying the variance-covariance matrix, we choose to estimate the variance-covariance matrix from a real dataset: the Quantitative Microbiome Project (QMP) data (Vandeputte et al. 2017). This dataset contains quantitative microbiome count data of 106 samples and 91 OTUs.
data(QMP, package = "ANCOMBC")
set.seed(123)
n = 150
d = ncol(QMP)
diff_prop = 0.1
lfc_cont = 1
lfc_cat2_vs_1 = -2
lfc_cat3_vs_1 = 1
# Generate the true abundances
abn_data = sim_plnm(abn_table = QMP, taxa_are_rows = FALSE, prv_cut = 0.05, 
                    n = n, lib_mean = 1e8, disp = 0.5)
log_abn_data = log(abn_data + 1e-5)
rownames(log_abn_data) = paste0("T", seq_len(d))
colnames(log_abn_data) = paste0("S", seq_len(n))
# Generate the sample and feature meta data
# Sampling fractions are set to differ by batches
smd = data.frame(samp_frac = log(c(runif(n/3, min = 1e-4, max = 1e-3),
                                   runif(n/3, min = 1e-3, max = 1e-2),
                                   runif(n/3, min = 1e-2, max = 1e-1))),
                 cont_cov = rnorm(n),
                 cat_cov = as.factor(rep(seq_len(3), each = n/3)))
rownames(smd) = paste0("S", seq_len(n))
                      
fmd = data.frame(taxon = paste0("T", seq_len(d)),
                 seq_eff = log(runif(d, min = 0.1, max = 1)),
                 lfc_cont = sample(c(0, lfc_cont), 
                                   size = d,
                                   replace = TRUE,
                                   prob = c(1 - diff_prop, diff_prop)),
                 lfc_cat2_vs_1 = sample(c(0, lfc_cat2_vs_1), 
                                        size = d,
                                        replace = TRUE,
                                        prob = c(1 - diff_prop, diff_prop)),
                 lfc_cat3_vs_1 = sample(c(0, lfc_cat3_vs_1), 
                                        size = d,
                                        replace = TRUE,
                                        prob = c(1 - diff_prop, diff_prop))) %>%
    mutate(lfc_cat3_vs_2 = lfc_cat3_vs_1 - lfc_cat2_vs_1)
# Add effect sizes of covariates to the true abundances
smd_dmy = model.matrix(~ 0 + cont_cov + cat_cov, data = smd)
log_abn_data = log_abn_data + outer(fmd$lfc_cont, smd_dmy[, "cont_cov"] )
log_abn_data = log_abn_data + outer(fmd$lfc_cat2_vs_1, smd_dmy[, "cat_cov2"])
log_abn_data = log_abn_data + outer(fmd$lfc_cat3_vs_1, smd_dmy[, "cat_cov3"])
# Add sample- and taxon-specific biases
log_otu_data = t(t(log_abn_data) + smd$samp_frac)
log_otu_data = log_otu_data + fmd$seq_eff
otu_data = round(exp(log_otu_data))set.seed(123)
output = ancombc2(data = otu_data, meta_data = smd, 
                  fix_formula = "cont_cov + cat_cov", rand_formula = NULL,
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "cat_cov", struc_zero = FALSE, neg_lb = FALSE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = FALSE, pairwise = TRUE, 
                  dunnet = FALSE, trend = FALSE,
                  iter_control = list(tol = 1e-5, max_iter = 20, 
                                      verbose = FALSE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = NULL, 
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100), 
                  trend_control = NULL)
res_prim = output$res
res_pair = output$res_pairres_merge1 = res_pair %>%
  dplyr::transmute(taxon, 
                   lfc_est1 = lfc_cat_cov2 * diff_cat_cov2,
                   lfc_est2 = lfc_cat_cov3 * diff_cat_cov3,
                   lfc_est3 = lfc_cat_cov3_cat_cov2 * diff_cat_cov3_cat_cov2) %>%
  dplyr::left_join(fmd %>%
                     dplyr::transmute(taxon, 
                                      lfc_true1 = lfc_cat2_vs_1,
                                      lfc_true2 = lfc_cat3_vs_1,
                                      lfc_true3 = lfc_cat3_vs_2),
                   by = "taxon") %>%
  dplyr::transmute(taxon, 
                   lfc_est1 = case_when(lfc_est1 > 0 ~ 1,
                                        lfc_est1 < 0 ~ -1,
                                        TRUE ~ 0),
                   lfc_est2 = case_when(lfc_est2 > 0 ~ 1,
                                        lfc_est2 < 0 ~ -1,
                                        TRUE ~ 0),
                   lfc_est3 = case_when(lfc_est3 > 0 ~ 1,
                                        lfc_est3 < 0 ~ -1,
                                        TRUE ~ 0),
                   lfc_true1 = case_when(lfc_true1 > 0 ~ 1,
                                         lfc_true1 < 0 ~ -1,
                                         TRUE ~ 0),
                   lfc_true2 = case_when(lfc_true2 > 0 ~ 1,
                                         lfc_true2 < 0 ~ -1,
                                         TRUE ~ 0),
                   lfc_true3 = case_when(lfc_true3 > 0 ~ 1,
                                         lfc_true3 < 0 ~ -1,
                                         TRUE ~ 0))
lfc_est1 = res_merge1$lfc_est1
lfc_true1 = res_merge1$lfc_true1
lfc_est2 = res_merge1$lfc_est2
lfc_true2 = res_merge1$lfc_true2
lfc_est3 = res_merge1$lfc_est3
lfc_true3 = res_merge1$lfc_true3
tp1 = sum(lfc_true1 == 1 & lfc_est1 == 1) +
  sum(lfc_true1 == -1 & lfc_est1 == -1)
fp1 = sum(lfc_true1 == 0 & lfc_est1 != 0) +
  sum(lfc_true1 == 1 & lfc_est1 == -1) +
  sum(lfc_true1 == -1 & lfc_est1 == 1)
fn1 = sum(lfc_true1 != 0 & lfc_est1 == 0)
tp2 = sum(lfc_true2 == 1 & lfc_est2 == 1) +
  sum(lfc_true2 == -1 & lfc_est2 == -1)
fp2 = sum(lfc_true2 == 0 & lfc_est2 != 0) +
  sum(lfc_true2 == 1 & lfc_est2 == -1) +
  sum(lfc_true2 == -1 & lfc_est2 == 1)
fn2 = sum(lfc_true2 != 0 & lfc_est2 == 0)
tp3 = sum(lfc_true3 == 1 & lfc_est3 == 1) +
  sum(lfc_true3 == -1 & lfc_est3 == -1)
fp3 = sum(lfc_true3 == 0 & lfc_est3 != 0) +
  sum(lfc_true3 == 1 & lfc_est3 == -1) +
  sum(lfc_true3 == -1 & lfc_est3 == 1)
fn3 = sum(lfc_true3 != 0 & lfc_est3 == 0)
tp = tp1 + tp2 + tp3
fp = fp1 + fp2 + fp3
fn = fn1 + fn2 + fn3
power1 = tp/(tp + fn)
fdr1 = fp/(tp + fp)
res_merge2 = res_pair %>%
  dplyr::transmute(taxon, 
                   lfc_est1 = lfc_cat_cov2 * diff_robust_cat_cov2,
                   lfc_est2 = lfc_cat_cov3 * diff_robust_cat_cov3,
                   lfc_est3 = lfc_cat_cov3_cat_cov2 * diff_robust_cat_cov3_cat_cov2) %>%
  dplyr::left_join(fmd %>%
                     dplyr::transmute(taxon, 
                                      lfc_true1 = lfc_cat2_vs_1,
                                      lfc_true2 = lfc_cat3_vs_1,
                                      lfc_true3 = lfc_cat3_vs_2),
                   by = "taxon") %>%
  dplyr::transmute(taxon, 
                   lfc_est1 = case_when(lfc_est1 > 0 ~ 1,
                                        lfc_est1 < 0 ~ -1,
                                        TRUE ~ 0),
                   lfc_est2 = case_when(lfc_est2 > 0 ~ 1,
                                        lfc_est2 < 0 ~ -1,
                                        TRUE ~ 0),
                   lfc_est3 = case_when(lfc_est3 > 0 ~ 1,
                                        lfc_est3 < 0 ~ -1,
                                        TRUE ~ 0),
                   lfc_true1 = case_when(lfc_true1 > 0 ~ 1,
                                         lfc_true1 < 0 ~ -1,
                                         TRUE ~ 0),
                   lfc_true2 = case_when(lfc_true2 > 0 ~ 1,
                                         lfc_true2 < 0 ~ -1,
                                         TRUE ~ 0),
                   lfc_true3 = case_when(lfc_true3 > 0 ~ 1,
                                         lfc_true3 < 0 ~ -1,
                                         TRUE ~ 0))
lfc_est1 = res_merge2$lfc_est1
lfc_true1 = res_merge2$lfc_true1
lfc_est2 = res_merge2$lfc_est2
lfc_true2 = res_merge2$lfc_true2
lfc_est3 = res_merge2$lfc_est3
lfc_true3 = res_merge2$lfc_true3
tp1 = sum(lfc_true1 == 1 & lfc_est1 == 1) +
  sum(lfc_true1 == -1 & lfc_est1 == -1)
fp1 = sum(lfc_true1 == 0 & lfc_est1 != 0) +
  sum(lfc_true1 == 1 & lfc_est1 == -1) +
  sum(lfc_true1 == -1 & lfc_est1 == 1)
fn1 = sum(lfc_true1 != 0 & lfc_est1 == 0)
tp2 = sum(lfc_true2 == 1 & lfc_est2 == 1) +
  sum(lfc_true2 == -1 & lfc_est2 == -1)
fp2 = sum(lfc_true2 == 0 & lfc_est2 != 0) +
  sum(lfc_true2 == 1 & lfc_est2 == -1) +
  sum(lfc_true2 == -1 & lfc_est2 == 1)
fn2 = sum(lfc_true2 != 0 & lfc_est2 == 0)
tp3 = sum(lfc_true3 == 1 & lfc_est3 == 1) +
  sum(lfc_true3 == -1 & lfc_est3 == -1)
fp3 = sum(lfc_true3 == 0 & lfc_est3 != 0) +
  sum(lfc_true3 == 1 & lfc_est3 == -1) +
  sum(lfc_true3 == -1 & lfc_est3 == 1)
fn3 = sum(lfc_true3 != 0 & lfc_est3 == 0)
tp = tp1 + tp2 + tp3
fp = fp1 + fp2 + fp3
fn = fn1 + fn2 + fn3
power2 = tp/(tp + fn)
fdr2 = fp/(tp + fp)
tab_summ = data.frame(Comparison = c("Without sensitivity score filter", 
                                     "With sensitivity score filter"),
                      Power = round(c(power1, power2), 2),
                      FDR = round(c(fdr1, fdr2), 2))
tab_summ %>%
    datatable(caption = "Power/FDR Comparison")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 ]bmi labelAfter permutation, the data can be considered null for bmi, meaning that a DA method that effectively controls false positives should produce p-values that are roughly uniformly distributed.
ancombc2 functionThe p-values from ANCOM-BC2, both with and without sensitivity analysis, exhibit a nearly uniform distribution.
Notably, enabling sensitivity analysis significantly reduces the number of significant p-values, thereby decreasing the risk of false positives.
We strongly recommend that users activate the sensitivity analysis and incorporate it into their final assessment of taxon significance.
For researchers aiming to benchmark ANCOM-BC2 against new methods in terms of false positive control, it is essential to account for the sensitivity analysis feature to ensure a fair comparison.
res_perm = output$res %>%
    dplyr::transmute(taxon = taxon,
                     p = p_bmioverweight,
                     ss_pass = passed_ss_bmioverweight) %>%
  rowwise() %>%
  mutate(
      p_update = case_when(
          # For taxa with significant p-values that failed the sensitivity analysis,
          # we assign a random value uniformly drawn from the range [0.05, 1].
          p < 0.05 & ss_pass == FALSE ~ runif(1, min = 0.05, max = 1),
    TRUE ~ p 
  )) %>%
  ungroup()
df_p = res_perm %>%
  dplyr::select(p, p_update) %>%
  pivot_longer(cols = p:p_update, names_to = "type", values_to = "value") 
df_p$type = recode(df_p$type, 
                   `p` = "ANCOM-BC2 (Without Sensitivity Analysis)", 
                   `p_update` = "ANCOM-BC2 (With Sensitivity Analysis)")
num_bins = 30
fig_p = df_p %>%
  ggplot(aes(x = value, fill = type)) +
  geom_histogram(position = "identity", alpha = 0.5, bins = num_bins) +
  geom_hline(yintercept = phyloseq::ntaxa(pseq_perm)/num_bins, linetype = "dashed") +
  theme_minimal() +
  labs(title = "Distribution of P-values",
       x = "P-value",
       y = "Frequency",
       fill = NULL,
       color = NULL) +
    theme_bw() + 
    theme(plot.title = element_text(hjust = 0.5),
          panel.grid.minor.y = element_blank(),
          legend.position = "bottom")
fig_pWe continue using the HITChip Atlas dataset as shown above.
ancombc2 function using the phyloseq objectTo control the FDR arising from multiple testing, we opt for the Holm-Bonferroni method over the Benjamini-Hochberg (BH) procedure, especially when dealing with large sample sizes where statistical power isn’t the primary concern. The Holm-Bonferroni method, accommodating any dependence structure among p-values, is known to be robust against inaccuracies in p-values, an issue often seen in DA analysis. Figures below display only results significant after the Holm-Bonferroni adjustment.
set.seed(123)
# It should be noted that we have set the number of bootstrap samples (B) equal 
# to 10 in the 'trend_control' function for computational expediency. 
# However, it is recommended that users utilize the default value of B, 
# which is 100, or larger values for optimal performance.
output = ancombc2(data = pseq, tax_level = "Family",
                  fix_formula = "age + region + bmi", rand_formula = NULL,
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "bmi", struc_zero = TRUE, neg_lb = TRUE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = TRUE, pairwise = TRUE, dunnet = TRUE, trend = TRUE,
                  iter_control = list(tol = 1e-2, max_iter = 20, 
                                      verbose = TRUE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = lme4::lmerControl(),
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100),
                  trend_control = list(contrast = list(matrix(c(1, 0, -1, 1),
                                                              nrow = 2, 
                                                              byrow = TRUE),
                                                       matrix(c(-1, 0, 1, -1),
                                                              nrow = 2, 
                                                              byrow = TRUE),
                                                       matrix(c(1, 0, 1, -1),
                                                              nrow = 2, 
                                                              byrow = TRUE)),
                                       node = list(2, 2, 1),
                                       solver = "ECOS",
                                       B = 10))Tests for interactions are supported by specifying the interactions in the fix_formula. Please ensure that interactions between variables are included using the * operator, as the : operator is not recognized by the ancombc2 function and will result in an error.
Additionally, when the group variable contains interaction terms, only the main effect will be considered in multi-group comparisons.
set.seed(123)
output_interact = ancombc2(data = pseq, tax_level = "Family",
                          fix_formula = "age * bmi + region", rand_formula = NULL,
                          p_adj_method = "holm", pseudo_sens = TRUE,
                          prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                          group = "bmi", struc_zero = TRUE, neg_lb = TRUE,
                          alpha = 0.05, n_cl = 2, verbose = TRUE)
head(output_interact$res) %>%
    dplyr::select(starts_with("lfc_age")) %>%
    mutate_if(is.numeric, function(x) round(x, 2)) %>%
    datatable(caption = "Example of including interaction terms")The primary output of the ANCOM-BC2 methodology identifies taxa with differential abundance based on the chosen covariate. The results include: 1) log fold changes, 2) standard errors, 3) test statistics, 4) p-values, 5) adjusted p-values, 6) indicators denoting whether the taxon is differentially abundant (TRUE) or not (FALSE), and 7) indicators denoting whether the taxon passed the sensitivity analysis (TRUE) or not (FALSE).
In the subsequent waterfall plot, each bar represents a log fold-change (in natural log) value. Any taxon highlighted in green indicates its successful passage through the sensitivity analysis (diff_robust is TRUE) for pseudo-count addition.
df_age = res_prim %>%
    dplyr::select(taxon, ends_with("age")) 
df_fig_age = df_age %>%
    dplyr::filter(diff_age == 1) %>% 
    dplyr::arrange(desc(lfc_age)) %>%
    dplyr::mutate(direct = ifelse(lfc_age > 0, "Positive LFC", "Negative LFC"),
                  color = ifelse(diff_robust_age, "aquamarine3", "black"))
df_fig_age$taxon = factor(df_fig_age$taxon, levels = df_fig_age$taxon)
df_fig_age$direct = factor(df_fig_age$direct, 
                           levels = c("Positive LFC", "Negative LFC"))
  
fig_age = df_fig_age %>%
    ggplot(aes(x = taxon, y = lfc_age, fill = direct)) + 
    geom_bar(stat = "identity", width = 0.7, color = "black", 
             position = position_dodge(width = 0.4)) +
    geom_errorbar(aes(ymin = lfc_age - se_age, ymax = lfc_age + se_age), 
                  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,
                                     color = df_fig_age$color))
fig_ageIn the subsequent heatmap, each cell represents a log fold-change (in natural log) value. Entries highlighted in green have successfully passed the sensitivity analysis for pseudo-count addition.
df_bmi = res_prim %>%
    dplyr::select(taxon, contains("bmi")) 
df_fig_bmi1 = df_bmi %>%
    dplyr::filter(diff_bmilean == 1 | 
                    diff_bmioverweight == 1) %>%
    dplyr::mutate(lfc1 = ifelse(diff_bmioverweight == 1, 
                                round(lfc_bmioverweight, 2), 0),
                  lfc2 = ifelse(diff_bmilean == 1, 
                                round(lfc_bmilean, 2), 0)) %>%
    tidyr::pivot_longer(cols = lfc1:lfc2, 
                        names_to = "group", values_to = "value") %>%
    dplyr::arrange(taxon)
df_fig_bmi2 = df_bmi %>%
    dplyr::filter(diff_bmilean == 1 | 
                    diff_bmioverweight == 1) %>%
    dplyr::mutate(lfc1 = ifelse(diff_robust_bmioverweight, 
                                "aquamarine3", "black"),
                  lfc2 = ifelse(diff_robust_bmilean, 
                                "aquamarine3", "black")) %>%
    tidyr::pivot_longer(cols = lfc1:lfc2, 
                        names_to = "group", values_to = "color") %>%
    dplyr::arrange(taxon)
df_fig_bmi = df_fig_bmi1 %>%
    dplyr::left_join(df_fig_bmi2, by = c("taxon", "group"))
df_fig_bmi$group = recode(df_fig_bmi$group, 
                          `lfc1` = "Overweight - Obese",
                          `lfc2` = "Lean - Obese")
df_fig_bmi$group = factor(df_fig_bmi$group, 
                          levels = c("Overweight - Obese",
                                     "Lean - Obese"))
  
lo = floor(min(df_fig_bmi$value))
up = ceiling(max(df_fig_bmi$value))
mid = (lo + up)/2
fig_bmi = df_fig_bmi %>%
  ggplot(aes(x = group, 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(group, taxon, label = value, color = color), size = 4) +
  scale_color_identity(guide = "none") +
  labs(x = NULL, y = NULL, title = "Log fold changes as compared to obese subjects") +
  theme_minimal() +
  theme(plot.title = element_text(hjust = 0.5))
fig_bmiThe primary goal of the ANCOM-BC2 global test is to discern taxa that demonstrate differential abundance between a minimum of two groups when analyzing three or more experimental groups.
To illustrate, in our current example, the objective is to pinpoint taxa with differential abundance across the “lean”, “overweight”, and “obese” categories. The results encompass: 1) test statistics, 2) p-values, 3) adjusted p-values, 4) indicators denoting whether the taxon is differentially abundant (TRUE) or not (FALSE), and 5) indicators denoting whether the taxon passed the sensitivity analysis (TRUE) or not (FALSE).
In the subsequent heatmap, each cell represents a log fold-change (in natural log) value. Taxa marked in green have successfully passed the sensitivity analysis for pseudo-count addition.
res_global = output$res_global
df_bmi = res_prim %>%
    dplyr::select(taxon, contains("bmi")) 
df_fig_global = df_bmi %>%
    dplyr::left_join(res_global %>%
                       dplyr::transmute(taxon, 
                                        diff_bmi = diff_abn, 
                                        diff_robust_bmi = diff_robust_abn)) %>%
    dplyr::filter(diff_bmi == 1) %>%
    dplyr::mutate(lfc_overweight = lfc_bmioverweight,
                  lfc_lean = lfc_bmilean,
                  color = ifelse(diff_robust_bmi, "aquamarine3", "black")) %>%
    dplyr::transmute(taxon,
                     `Overweight - Obese` = round(lfc_overweight, 2),
                     `Lean - Obese` = round(lfc_lean, 2), 
                     color = color) %>%
    tidyr::pivot_longer(cols = `Overweight - Obese`:`Lean - Obese`, 
                        names_to = "group", values_to = "value") %>%
    dplyr::arrange(taxon)
df_fig_global$group = factor(df_fig_global$group, 
                             levels = c("Overweight - Obese",
                                        "Lean - Obese"))
  
lo = floor(min(df_fig_global$value))
up = ceiling(max(df_fig_global$value))
mid = (lo + up)/2
fig_global = df_fig_global %>%
  ggplot(aes(x = group, 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(group, 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),
        axis.text.y = element_text(color = df_fig_global %>%
                                       dplyr::distinct(taxon, color) %>%
                                       .$color))
fig_globalThe ANCOM-BC2 methodology for multiple pairwise comparisons is designed to identify taxa that exhibit differential abundance between any two groups within a set of three or more experimental groups, all while maintaining control over the mdFDR.
For instance, in our analysis focusing on the categories “lean”, “overweight”, and “obese”, the output provides: 1) log fold changes, 2) standard errors, 3) test statistics, 4) p-values, 5) adjusted p-values, 6) indicators denoting whether the taxon is differentially abundant (TRUE) or not (FALSE), and 7) indicators denoting whether the taxon passed the sensitivity analysis (TRUE) or not (FALSE).
In the subsequent heatmap, each cell represents a log fold-change (in natural log) value. Entries highlighted in green have successfully passed the sensitivity analysis for pseudo-count addition.
res_pair = output$res_pair
df_fig_pair1 = res_pair %>%
    dplyr::filter(diff_bmioverweight == 1 |
                      diff_bmilean == 1 | 
                      diff_bmilean_bmioverweight == 1) %>%
    dplyr::mutate(lfc1 = ifelse(diff_bmioverweight == 1, 
                                round(lfc_bmioverweight, 2), 0),
                  lfc2 = ifelse(diff_bmilean == 1, 
                                round(lfc_bmilean, 2), 0),
                  lfc3 = ifelse(diff_bmilean_bmioverweight == 1, 
                                round(lfc_bmilean_bmioverweight, 2), 0)) %>%
    tidyr::pivot_longer(cols = lfc1:lfc3, 
                        names_to = "group", values_to = "value") %>%
    dplyr::arrange(taxon)
df_fig_pair2 = res_pair %>%
    dplyr::filter(diff_bmioverweight == 1 |
                      diff_bmilean == 1 | 
                      diff_bmilean_bmioverweight == 1) %>%
    dplyr::mutate(lfc1 = ifelse(diff_robust_bmioverweight, 
                                "aquamarine3", "black"),
                  lfc2 = ifelse(diff_robust_bmilean, 
                                "aquamarine3", "black"),
                  lfc3 = ifelse(diff_robust_bmilean_bmioverweight, 
                                "aquamarine3", "black")) %>%
    tidyr::pivot_longer(cols = lfc1:lfc3, 
                        names_to = "group", values_to = "color") %>%
    dplyr::arrange(taxon)
df_fig_pair = df_fig_pair1 %>%
    dplyr::left_join(df_fig_pair2, by = c("taxon", "group"))
df_fig_pair$group = recode(df_fig_pair$group, 
                          `lfc1` = "Overweight - Obese",
                          `lfc2` = "Lean - Obese",
                          `lfc3` = "Lean - Overweight")
df_fig_pair$group = factor(df_fig_pair$group, 
                          levels = c("Overweight - Obese",
                                     "Lean - Obese", 
                                     "Lean - Overweight"))
lo = floor(min(df_fig_pair$value))
up = ceiling(max(df_fig_pair$value))
mid = (lo + up)/2
fig_pair = df_fig_pair %>%
    ggplot(aes(x = group, 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(group, taxon, label = value, color = color), size = 4) +
    scale_color_identity(guide = FALSE) +
    labs(x = NULL, y = NULL, title = "Log fold changes as compared to obese subjects") +
    theme_minimal() +
    theme(plot.title = element_text(hjust = 0.5))
fig_pairThe Dunnett’s test (Dunnett 1955; Dunnett and Tamhane 1991, 1992) is tailored for comparing multiple experimental groups against a control or reference group. ANCOM-BC2 Dunnett’s type of test applies this framework but also controls the mdFDR. It’s essential to highlight that ANCOM-BC2’s primary results control for multiple testing across taxa but not for multiple comparisons between groups. As such, unlike the ANCOM-BC2 Dunnett’s test, the primary output doesn’t control the mdFDR.
In this illustration, our objective is to pinpoint taxa with differential abundance between “lean”, “overweight”, and the reference group “obese”. The results encompass: 1) log fold changes, 2) standard errors, 3) test statistics, 4) p-values, 5) adjusted p-values, 6) indicators denoting whether the taxon is differentially abundant (TRUE) or not (FALSE), and 7) indicators denoting whether the taxon passed the sensitivity analysis (TRUE) or not (FALSE).
In the subsequent heatmap, each cell represents a log fold-change (in natural log) value. Entries highlighted in green have successfully passed the sensitivity analysis for pseudo-count addition.
res_dunn = output$res_dunn
df_fig_dunn1 = res_dunn %>%
    dplyr::filter(diff_bmioverweight == 1 | 
                      diff_bmilean == 1) %>%
    dplyr::mutate(lfc1 = ifelse(diff_bmioverweight == 1, 
                                round(lfc_bmioverweight, 2), 0),
                  lfc2 = ifelse(diff_bmilean == 1, 
                                round(lfc_bmilean, 2), 0)) %>%
    tidyr::pivot_longer(cols = lfc1:lfc2, 
                        names_to = "group", values_to = "value") %>%
    dplyr::arrange(taxon)
df_fig_dunn2 = res_dunn %>%
    dplyr::filter(diff_bmioverweight == 1 | 
                      diff_bmilean == 1) %>%
    dplyr::mutate(lfc1 = ifelse(diff_robust_bmioverweight, 
                                "aquamarine3", "black"),
                  lfc2 = ifelse(diff_robust_bmilean == 1, 
                                "aquamarine3", "black")) %>%
    tidyr::pivot_longer(cols = lfc1:lfc2, 
                        names_to = "group", values_to = "color") %>%
    dplyr::arrange(taxon)
df_fig_dunn = df_fig_dunn1 %>%
    dplyr::left_join(df_fig_dunn2, by = c("taxon", "group"))
df_fig_dunn$group = recode(df_fig_dunn$group, 
                          `lfc1` = "Overweight - Obese",
                          `lfc2` = "Lean - Obese")
df_fig_dunn$group = factor(df_fig_dunn$group, 
                          levels = c("Overweight - Obese",
                                     "Lean - Obese"))
lo = floor(min(df_fig_dunn$value))
up = ceiling(max(df_fig_dunn$value))
mid = (lo + up)/2
fig_dunn = df_fig_dunn %>%
  ggplot(aes(x = group, 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(group, taxon, label = value, color = color), size = 4) +
    scale_color_identity(guide = FALSE) +
  labs(x = NULL, y = NULL, title = "Log fold changes as compared to obese subjects") +
  theme_minimal() +
  theme(plot.title = element_text(hjust = 0.5))
fig_dunnIn certain experiments, groups inherently follow an order, such as in dose-response studies. In these cases, we’re interested in discerning if microbial abundances align with specific patterns. Such patterns might manifest as monotonically increasing, monotonically decreasing, or an umbrella shape.
ANCOM-BC2 pattern analysis is able to identify potential patterns by testing the contrast: \[Ax \ge 0\] where \(A\) is the contrast matrix and \(x\) is the vector of parameters.
For instance, in this example, we want to identify taxa that are monotonically increasing across “obese”, “overweight”, and “lean”. Note that “obese” is the reference group, and the parameters we can estimate are the differences as compared to the reference group, i.e., \[x = (\text{overweight - obese}, \text{lean - obese})^T\] To test the monotonically increasing trend: \[H_0: \text{obese} = \text{overweight} = \text{lean} \\ H_1: \text{obese} \le \text{overweight} \le \text{lean} \quad \text{with at least one strict inequality}\]
We are essentially testing: \[H_0: 0 = \text{overweight - obese} = \text{lean - obese} \\ H_1: 0 \le \text{overweight - obese} \le \text{lean - obese} \quad \text{with at least one strict inequality}\]
Thus, we shall specify the contrast matrix \(A\) as \[A = \begin{bmatrix} 1 & 0 \\ -1 & 1 \end{bmatrix}\] The first row of \(A\) matrix, \((1, 0)\), indicates that \(\text{"overweight - obese"} \ge 0\), and the second row, \((-1, 1)\) represents \(\text{"lean - obese"} - \text{"overweight - obese"} \ge 0\).
Similarly, to test for the monotonically decreasing trend: \[H_0: \text{obese} = \text{overweight} = \text{lean} \\ H_1: \text{obese} \ge \text{overweight} \ge \text{lean} \quad \text{with at least one strict inequality}\] We shall specify the contrast matrix \(A\) as \[A = \begin{bmatrix} -1 & 0 \\ 1 & -1 \end{bmatrix}\] Lastly, to test for an umbrella trend: \[H_0: \text{obese} = \text{overweight} = \text{lean} \\ H_1: \text{obese} \le \text{overweight} \ge \text{lean} \quad \text{with at least one strict inequality}\] We shall specify the contrast matrix \(A\) as \[A = \begin{bmatrix} 1 & 0 \\ 1 & -1 \end{bmatrix}\]
For testing monotonic trend (increasing or decreasing), one should specify the node parameter in trend_control as the last position of x. In this example, the vector of parameters \(x\) is of length 2, thus, the last position is 2. For testing umbrella shape, such as in the above umbrella shape example, one should set node as the position of the turning point of x. In this example, the turning position is overweight, thus, node = 1.
We will test the monotonically increasing and decreasing trends, as well as the umbrella trend in this example. The result contains: 1) log fold changes, 2) standard errors, 3) test statistics, 4) p-values, 5) adjusted p-values, 6) indicators denoting whether the taxon is differentially abundant (TRUE) or not (FALSE), and 7) indicators denoting whether the taxon passed the sensitivity analysis (TRUE) or not (FALSE).
In the subsequent heatmap, each cell represents a log fold-change (in natural log) value. Taxa marked in green have successfully passed the sensitivity analysis for pseudo-count addition.
res_trend = output$res_trend
df_fig_trend = res_trend %>%
    dplyr::filter(diff_abn == 1) %>%
    dplyr::mutate(lfc1 = round(lfc_bmioverweight, 2),
                  lfc2 = round(lfc_bmilean, 2),
                  color = ifelse(diff_robust_abn, "aquamarine3", "black")) %>%
    tidyr::pivot_longer(cols = lfc1:lfc2, 
                        names_to = "group", values_to = "value") %>%
    dplyr::arrange(taxon)
df_fig_trend$group = recode(df_fig_trend$group, 
                          `lfc1` = "Overweight - Obese",
                          `lfc2` = "Lean - Obese")
df_fig_trend$group = factor(df_fig_trend$group, 
                          levels = c("Overweight - Obese", 
                                     "Lean - Obese"))
lo = floor(min(df_fig_trend$value))
up = ceiling(max(df_fig_trend$value))
mid = (lo + up)/2
fig_trend = df_fig_trend %>%
    ggplot(aes(x = group, 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(group, taxon, 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),
          axis.text.y = element_text(color = df_fig_trend %>%
                                         dplyr::distinct(taxon, color) %>%
                                         .$color))
fig_trendancombc2 function using the tse objectOne can also run the ancombc2 function using the TreeSummarizedExperiment object.
tse = mia::makeTreeSummarizedExperimentFromPhyloseq(atlas1006)
tse = tse[, tse$time == 0]
tse$bmi = recode(tse$bmi_group,
                 obese = "obese",
                 severeobese = "obese",
                 morbidobese = "obese")
tse = tse[, tse$bmi %in% c("lean", "overweight", "obese")]
tse$bmi = factor(tse$bmi, levels = c("obese", "overweight", "lean"))
tse$region = recode(as.character(tse$nationality),
                    Scandinavia = "NE", UKIE = "NE", SouthEurope = "SE", 
                    CentralEurope = "CE", EasternEurope = "EE",
                    .missing = "unknown")
tse = tse[, ! tse$region %in% c("EE", "unknown")]
set.seed(123)
output = ancombc2(data = tse, assay_name = "counts", tax_level = "Family",
                  fix_formula = "age + region + bmi", rand_formula = NULL,
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "bmi", struc_zero = TRUE, neg_lb = TRUE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = TRUE, pairwise = TRUE, dunnet = TRUE, trend = TRUE,
                  iter_control = list(tol = 1e-2, max_iter = 20, 
                                      verbose = TRUE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = lme4::lmerControl(),
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100),
                  trend_control = list(contrast = list(matrix(c(1, 0, -1, 1),
                                                              nrow = 2, 
                                                              byrow = TRUE),
                                                       matrix(c(-1, 0, 1, -1),
                                                              nrow = 2, 
                                                              byrow = TRUE),
                                                       matrix(c(1, 0, 1, -1),
                                                              nrow = 2, 
                                                              byrow = TRUE)),
                                       node = list(2, 2, 1),
                                       solver = "ECOS",
                                       B = 10))ancombc2 function by directly providing the abundance and metadataPlease ensure that you provide the following: 1) The abundance data at its lowest possible taxonomic level. 2) The aggregated data at the desired taxonomic level; if no aggregation is performed, it can be the same as the original abundance data. 3) The sample metadata.
abundance_data = microbiome::abundances(pseq)
aggregate_data = microbiome::abundances(microbiome::aggregate_taxa(pseq, "Family"))
meta_data = microbiome::meta(pseq)
set.seed(123)
output = ancombc2(data = abundance_data, 
                  aggregate_data = aggregate_data,
                  meta_data = meta_data,
                  fix_formula = "age + region + bmi", rand_formula = NULL,
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "bmi", struc_zero = TRUE, neg_lb = TRUE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = TRUE, pairwise = TRUE, dunnet = TRUE, trend = TRUE,
                  iter_control = list(tol = 1e-2, max_iter = 20, 
                                      verbose = TRUE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = lme4::lmerControl(),
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100),
                  trend_control = list(contrast = list(matrix(c(1, 0, -1, 1),
                                                              nrow = 2, 
                                                              byrow = TRUE),
                                                       matrix(c(-1, 0, 1, -1),
                                                              nrow = 2, 
                                                              byrow = TRUE),
                                                       matrix(c(1, 0, 1, -1),
                                                              nrow = 2, 
                                                              byrow = TRUE)),
                                       node = list(2, 2, 1),
                                       solver = "ECOS",
                                       B = 10))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.
ancombc2 function using the phyloseq objectIn 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”
Procedures of ANCOM-BC2 global test, pairwise directional test, Dunnett’s type of test, and trend test are the same as those for the cross-sectional data shown above.
set.seed(123)
# It should be noted that we have set the number of bootstrap samples (B) equal 
# to 10 in the 'trend_control' function for computational expediency. 
# However, it is recommended that users utilize the default value of B, 
# which is 100, or larger values for optimal performance.
output = ancombc2(data = dietswap, tax_level = "Family",
                  fix_formula = "nationality + timepoint + group",
                  rand_formula = "(timepoint | subject)",
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "group", struc_zero = TRUE, neg_lb = TRUE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = TRUE, pairwise = TRUE, dunnet = TRUE, trend = TRUE,
                  iter_control = list(tol = 1e-2, max_iter = 20, 
                                      verbose = TRUE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = lme4::lmerControl(),
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100),
                  trend_control = list(contrast = list(matrix(c(1, 0, -1, 1),
                                                              nrow = 2, 
                                                              byrow = TRUE)),
                                       node = list(2),
                                       solver = "ECOS",
                                       B = 10))Again, tests for interactions are supported by specifying the interactions in the fix_formula. Please ensure that interactions between variables are included using the * operator, as the : operator is not recognized by the ancombc2 function and will result in an error.
Additionally, when the group variable contains interaction terms, only the main effect will be considered in multi-group comparisons.
set.seed(123)
output = ancombc2(data = dietswap, tax_level = "Family",
                  fix_formula = "nationality + timepoint * group",
                  rand_formula = "(timepoint | subject)",
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "group", struc_zero = TRUE, neg_lb = TRUE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = TRUE, pairwise = TRUE, dunnet = TRUE, trend = TRUE,
                  iter_control = list(tol = 1e-2, max_iter = 20, 
                                      verbose = TRUE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = lme4::lmerControl(),
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100),
                  trend_control = list(contrast = list(matrix(c(1, 0, -1, 1),
                                                              nrow = 2, 
                                                              byrow = TRUE)),
                                       node = list(2),
                                       solver = "ECOS",
                                       B = 10))ancombc2 function using the tse objecttse = mia::makeTreeSummarizedExperimentFromPhyloseq(dietswap)
set.seed(123)
output = ancombc2(data = tse, assay_name = "counts", tax_level = "Family",
                  fix_formula = "nationality + timepoint + group",
                  rand_formula = "(timepoint | subject)",
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "group", struc_zero = TRUE, neg_lb = TRUE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = TRUE, pairwise = TRUE, dunnet = TRUE, trend = TRUE,
                  iter_control = list(tol = 1e-2, max_iter = 20, 
                                      verbose = TRUE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = lme4::lmerControl(),
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100),
                  trend_control = list(contrast = list(matrix(c(1, 0, -1, 1),
                                                              nrow = 2, 
                                                              byrow = TRUE)),
                                       node = list(2),
                                       solver = "ECOS",
                                       B = 10))ancombc2 function by directly providing the abundance and metadataabundance_data = microbiome::abundances(dietswap)
aggregate_data = microbiome::abundances(microbiome::aggregate_taxa(dietswap, "Family"))
meta_data = microbiome::meta(dietswap)
set.seed(123)
output = ancombc2(data = abundance_data, 
                  aggregate_data = aggregate_data,
                  meta_data = meta_data,
                  fix_formula = "nationality + timepoint + group",
                  rand_formula = "(timepoint | subject)",
                  p_adj_method = "holm", pseudo_sens = TRUE,
                  prv_cut = 0.10, lib_cut = 1000, s0_perc = 0.05,
                  group = "group", struc_zero = TRUE, neg_lb = TRUE,
                  alpha = 0.05, n_cl = 2, verbose = TRUE,
                  global = TRUE, pairwise = TRUE, dunnet = TRUE, trend = TRUE,
                  iter_control = list(tol = 1e-2, max_iter = 20, 
                                      verbose = TRUE),
                  em_control = list(tol = 1e-5, max_iter = 100),
                  lme_control = lme4::lmerControl(),
                  mdfdr_control = list(fwer_ctrl_method = "holm", B = 100),
                  trend_control = list(contrast = list(matrix(c(1, 0, -1, 1),
                                                              nrow = 2, 
                                                              byrow = TRUE)),
                                       node = list(2),
                                       solver = "ECOS",
                                       B = 10))It is important to acknowledge that the estimation of sampling fractions in ANCOM-BC2 is limited to an additive constant. This means that only the difference between bias-corrected log abundances is meaningful, rather than the absolute values themselves.
Moreover, within each taxon, the bias-corrected log abundances are centered across samples. ANCOM-BC2 operates on the assumption that while these taxon-specific biases differ between taxa, they stay consistent within a given taxon across various samples. This assumption enables ANCOM-BC2 to accommodate intra-taxon variations through the centering of the bias-corrected log abundances.
bias_correct_log_table = output$bias_correct_log_table
# By default, ANCOM-BC2 does not add pseudo-counts to zero counts, which can 
# result in NAs in the bias-corrected log abundances. Users have the option to 
# either leave the NAs as they are or replace them with zeros. 
# This replacement is equivalent to adding pseudo-counts of ones to the zero counts. 
bias_correct_log_table[is.na(bias_correct_log_table)] = 0
# Show the first 6 samples
round(bias_correct_log_table[, 1:6], 2) %>% 
  datatable(caption = "Bias-corrected log abundances")R version 4.5.0 RC (2025-04-04 r88126)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.2 LTS
Matrix products: default
BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.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.2   rngtools_1.5.2  foreach_1.5.2   DT_0.33        
 [5] lubridate_1.9.4 forcats_1.0.0   stringr_1.5.1   dplyr_1.1.4    
 [9] purrr_1.0.4     readr_2.1.5     tidyr_1.3.1     tibble_3.2.1   
[13] ggplot2_3.5.2   tidyverse_2.0.0 ANCOMBC_2.10.1 
loaded via a namespace (and not attached):
  [1] RColorBrewer_1.1-3      rstudioapi_0.17.1       jsonlite_2.0.0         
  [4] magrittr_2.0.3          TH.data_1.1-3           farver_2.1.2           
  [7] nloptr_2.2.1            rmarkdown_2.29          fs_1.6.6               
 [10] vctrs_0.6.5             multtest_2.64.0         minqa_1.2.8            
 [13] base64enc_0.1-3         htmltools_0.5.8.1       energy_1.7-12          
 [16] haven_2.5.4             Rhdf5lib_1.30.0         cellranger_1.1.0       
 [19] Formula_1.2-5           rhdf5_2.52.0            sass_0.4.10            
 [22] bslib_0.9.0             htmlwidgets_1.6.4       plyr_1.8.9             
 [25] sandwich_3.1-1          rootSolve_1.8.2.4       zoo_1.8-14             
 [28] cachem_1.1.0            igraph_2.1.4            lifecycle_1.0.4        
 [31] iterators_1.0.14        pkgconfig_2.0.3         Matrix_1.7-3           
 [34] R6_2.6.1                fastmap_1.2.0           GenomeInfoDbData_1.2.14
 [37] rbibutils_2.3           digest_0.6.37           Exact_3.3              
 [40] numDeriv_2016.8-1.1     colorspace_2.1-1        S4Vectors_0.46.0       
 [43] crosstalk_1.2.1         Hmisc_5.2-3             vegan_2.6-10           
 [46] labeling_0.4.3          timechange_0.3.0        mgcv_1.9-3             
 [49] httr_1.4.7              compiler_4.5.0          proxy_0.4-27           
 [52] bit64_4.6.0-1           withr_3.0.2             doParallel_1.0.17      
 [55] gsl_2.1-8               htmlTable_2.4.3         backports_1.5.0        
 [58] MASS_7.3-65             biomformat_1.36.0       permute_0.9-7          
 [61] gtools_3.9.5            CVXR_1.0-15             gld_2.6.7              
 [64] tools_4.5.0             foreign_0.8-90          ape_5.8-1              
 [67] nnet_7.3-20             glue_1.8.0              rhdf5filters_1.20.0    
 [70] nlme_3.1-168            grid_4.5.0              Rtsne_0.17             
 [73] checkmate_2.3.2         cluster_2.1.8.1         reshape2_1.4.4         
 [76] ade4_1.7-23             generics_0.1.4          microbiome_1.30.0      
 [79] gtable_0.3.6            tzdb_0.5.0              class_7.3-23           
 [82] data.table_1.17.2       lmom_3.2                hms_1.1.3              
 [85] XVector_0.48.0          BiocGenerics_0.54.0     pillar_1.10.2          
 [88] splines_4.5.0           lattice_0.22-7          survival_3.8-3         
 [91] gmp_0.7-5               bit_4.6.0               tidyselect_1.2.1       
 [94] Biostrings_2.76.0       knitr_1.50              reformulas_0.4.1       
 [97] gridExtra_2.3           phyloseq_1.52.0         IRanges_2.42.0         
[100] stats4_4.5.0            xfun_0.52               expm_1.0-0             
[103] Biobase_2.68.0          stringi_1.8.7           UCSC.utils_1.4.0       
[106] yaml_2.3.10             boot_1.3-31             evaluate_1.0.3         
[109] codetools_0.2-20        cli_3.6.5               rpart_4.1.24           
[112] DescTools_0.99.60       Rdpack_2.6.4            jquerylib_0.1.4        
[115] dichromat_2.0-0.1       Rcpp_1.0.14             GenomeInfoDb_1.44.0    
[118] readxl_1.4.5            parallel_4.5.0          lme4_1.1-37            
[121] Rmpfr_1.1-0             mvtnorm_1.3-3           lmerTest_3.1-3         
[124] scales_1.4.0            e1071_1.7-16            crayon_1.5.3           
[127] rlang_1.1.6             multcomp_1.4-28        Costea, Paul I, Georg Zeller, Shinichi Sunagawa, and Peer Bork. 2014. “A Fair Comparison.” Nature Methods 11 (4): 359–59.
Dunnett, Charles W. 1955. “A Multiple Comparison Procedure for Comparing Several Treatments with a Control.” Journal of the American Statistical Association 50 (272): 1096–1121.
Dunnett, Charles W, and Ajit C Tamhane. 1991. “Step-down Multiple Tests for Comparing Treatments with a Control in Unbalanced One-Way Layouts.” Statistics in Medicine 10 (6): 939–47.
———. 1992. “A Step-up Multiple Test Procedure.” Journal of the American Statistical Association 87 (417): 162–70.
Grandhi, Anjana, Wenge Guo, and Shyamal D Peddada. 2016. “A Multiple Testing Procedure for Multi-Dimensional Pairwise Comparisons with Application to Gene Expression Studies.” BMC Bioinformatics 17 (1): 1–12.
Guo, Wenge, Sanat K Sarkar, and Shyamal D Peddada. 2010. “Controlling False Discoveries in Multidimensional Directional Decisions, with Applications to Gene Expression Data on Ordered Categories.” Biometrics 66 (2): 485–92.
Hu, Yi-Juan, and Glen A Satten. 2020. “Testing Hypotheses About the Microbiome Using the Linear Decomposition Model (Ldm).” Bioinformatics 36 (14): 4106–15.
Jelsema, Casey M, and Shyamal D Peddada. 2016. “CLME: An R Package for Linear Mixed Effects Models Under Inequality Constraints.” Journal of Statistical Software 75.
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.
Paulson, Joseph N, Héctor Corrada Bravo, and Mihai Pop. 2014. “Reply to:" A Fair Comparison".” Nature Methods 11 (4): 359–60.
Tusher, Virginia Goss, Robert Tibshirani, and Gilbert Chu. 2001. “Significance Analysis of Microarrays Applied to the Ionizing Radiation Response.” Proceedings of the National Academy of Sciences 98 (9): 5116–21.
Vandeputte, Doris, Gunter Kathagen, Kevin D’hoe, Sara Vieira-Silva, Mireia Valles-Colomer, João Sabino, Jun Wang, et al. 2017. “Quantitative Microbiome Profiling Links Gut Community Variation to Microbial Load.” Nature 551 (7681): 507–11.