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.

Analysis of Composition of Microbiomes (ANCOM) (Mandal et al. 2015) is also a 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. Larger the value of W, the more likely the taxon is differentially abundant. For more details, please refer to the ANCOM paper.

Sparse Estimation of Correlations among Microbiomes (SECOM) is a methodology which aims to detect both linear and nonlinear relationships between a pair of taxa within an ecosystem (e.g. gut) or across ecosystems (e.g. gut and tongue). SECOM corrects both sample-specific and taxon-specific biases, obtains a consistent estimator for the correlation matrix of microbial absolute abundances, while maintaining the underlying true sparsity.

2. Installation

Download package.

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

Load the package.

library(ANCOMBC)

3. Example Data

Cross-sectional data

The HITChip Atlas data set (Lahti et al. 2014) is available via the microbiome R package (Lahti et al. 2017) in phyloseq (McMurdie and Holmes 2013) format.

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

# Re-code the bmi group
sample_data(pseq)$bmi_group = recode(sample_data(pseq)$bmi_group,
                                     lean = "lean",
                                     overweight = "overweight",
                                     obese = "obese",
                                     severeobese = "obese",
                                     morbidobese = "obese")
# Subset to lean and obese subjects
pseq = subset_samples(pseq, bmi_group %in% c("lean", "obese"))

# Create the region variable
sample_data(pseq)$region = recode(sample_data(pseq)$nationality,
                                  Scandinavia = "NE",
                                  UKIE = "NE",
                                  SouthEurope = "SE",
                                  CentralEurope = "CE",
                                  EasternEurope = "EE")

# Discard "EE" as it contains only 1 subject
pseq = subset_samples(pseq, region != "EE")
# Genus level data
genus_data = aggregate_taxa(pseq, "Genus")
# Family level data
family_data = aggregate_taxa(pseq, "Family")
# Phylum level data
phylum_data = aggregate_taxa(pseq, "Phylum")

print(genus_data)
phyloseq-class experiment-level object
otu_table()   OTU Table:         [ 130 taxa and 722 samples ]
sample_data() Sample Data:       [ 722 samples by 11 sample variables ]
tax_table()   Taxonomy Table:    [ 130 taxa by 3 taxonomic ranks ]
print(family_data)
phyloseq-class experiment-level object
otu_table()   OTU Table:         [ 22 taxa and 722 samples ]
sample_data() Sample Data:       [ 722 samples by 11 sample variables ]
tax_table()   Taxonomy Table:    [ 22 taxa by 3 taxonomic ranks ]
print(phylum_data)
phyloseq-class experiment-level object
otu_table()   OTU Table:         [ 8 taxa and 722 samples ]
sample_data() Sample Data:       [ 722 samples by 11 sample variables ]
tax_table()   Taxonomy Table:    [ 8 taxa by 2 taxonomic ranks ]

Longitudinal data

A two-week diet swap study between western (USA) and traditional (rural Africa) diets (Lahti et al. 2014). The data set is available via the microbiome R package (Lahti et al. 2017) in phyloseq (McMurdie and Holmes 2013) format.

data(dietswap)

# Aggregate to family level
family_data2 = aggregate_taxa(dietswap, "Family")

print(family_data2)
phyloseq-class experiment-level object
otu_table()   OTU Table:         [ 22 taxa and 222 samples ]
sample_data() Sample Data:       [ 222 samples by 8 sample variables ]
tax_table()   Taxonomy Table:    [ 22 taxa by 3 taxonomic ranks ]

4. ANCOMBC Implementation

4.1 Run ancombc function

out = ancombc(phyloseq = family_data, formula = "age + region + bmi_group", 
              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)

res = out$res
res_global = out$res_global

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("Age", "NE - CE", "SE - CE", "US - CE", "Obese - Lean")
colnames(tab_lfc) = col_name
tab_lfc %>% 
  datatable(caption = "Log Fold Changes from the Primary Result") %>%
  formatRound(col_name, digits = 2)

SE

tab_se = res$se
colnames(tab_se) = col_name
tab_se %>% 
  datatable(caption = "SEs from the Primary Result") %>%
  formatRound(col_name, 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, 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, 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, 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")

Visualization for age

df_lfc = data.frame(res$lfc * res$diff_abn, check.names = FALSE) %>% 
  rownames_to_column("taxon_id")
df_se = data.frame(res$se * res$diff_abn, check.names = FALSE) %>% 
  rownames_to_column("taxon_id")
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 = "Waterfall Plot 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 %>% 
  dplyr::left_join(df_se, by = "taxon_id") %>%
  dplyr::transmute(taxon_id, bmi = bmi_groupobese, bmiSE = bmi_groupobeseSE) %>%
  dplyr::filter(bmi != 0) %>% 
  dplyr::arrange(desc(bmi)) %>%
  dplyr::mutate(direct = ifelse(bmi > 0, "Positive LFC", "Negative LFC"))
df_fig_bmi$taxon_id = factor(df_fig_bmi$taxon_id, levels = df_fig_bmi$taxon_id)
df_fig_bmi$direct = factor(df_fig_bmi$direct, 
                        levels = c("Positive LFC", "Negative LFC"))
  
p_bmi = ggplot(data = df_fig_bmi, 
           aes(x = taxon_id, y = bmi, fill = direct, color = direct)) + 
  geom_bar(stat = "identity", width = 0.7, 
           position = position_dodge(width = 0.4)) +
  geom_errorbar(aes(ymin = bmi - bmiSE, ymax = bmi + bmiSE), width = 0.2,
                position = position_dodge(0.05), color = "black") + 
  labs(x = NULL, y = "Log fold change", 
       title = "Waterfall Plot of BMI (Obese - Lean)") + 
  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_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[, "W", drop = FALSE]
colnames(tab_w) = "Region"
tab_w %>% datatable(caption = "Test Statistics 
                    from the Global Test Result") %>%
      formatRound(c("Region"), digits = 2)

P-values

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

Adjusted p-values

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

Differentially abundant taxa

tab_diff = res_global[, "diff_abn", drop = FALSE]
colnames(tab_diff) = "Region"
tab_diff %>% datatable(caption = "Differentially Abundant Taxa 
                       from the Global Test Result")

Bias-corrected abundances

Step 1: estimate sample-specific sampling fractions (in log scale). Note that for each sample, if it contains missing values for any variable specified in the formula, the corresponding sampling fraction estimate for this sample will be NA since the sampling fraction is not estimable with the presence of missing values.

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

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(abundances(family_data) + 1) 
# Adjust the log observed abundances
log_obs_abn_adj = t(t(log_obs_abn) - samp_frac)
# Show the first 6 samples
round(log_obs_abn_adj[, 1:6], 2) %>% 
  datatable(caption = "Bias-corrected log observed abundances")

Visualization

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

df_sig = as.data.frame(t(log_obs_abn_adj[sig_taxa, ])) %>%
  tibble::rownames_to_column("sample") %>%
  dplyr::left_join(meta(family_data) %>%
                     select(sample, region),
                   by = "sample") %>%
  dplyr::filter(!is.na(region)) %>%
  tidyr::pivot_longer(cols = -one_of("sample", "region"), 
                      names_to = "taxon", values_to = "value")

df_heat = df_sig %>%
  dplyr::group_by(region, taxon) %>%
  dplyr::summarise_if(is.numeric, mean, na.rm = TRUE) %>%
  dplyr::mutate(value = round(value, 2)) %>%
  dplyr::arrange(region)
df_heat$taxon = factor(df_heat$taxon, levels = 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 = "white") +
  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 = "Heat Map of Region") +
  theme_minimal() +
  theme(plot.title = element_text(hjust = 0.5))
p_heat

5. ANCOM Implementation

5.1 Cross-sectional data

Run ancom function

set.seed(123)
out1 = ancom(phyloseq = family_data,  p_adj_method = "holm",
             prv_cut = 0.10, lib_cut = 1000, main_var = "bmi_group", 
             adj_formula = "age + region", rand_formula = NULL,
             lme_control = NULL, struc_zero = TRUE, neg_lb = TRUE,
             alpha = 0.05, n_cl = 2)

res1 = out1$res

# Similarly, if the main variable of interest is continuous, such as age, the
# ancom model can be specified as
# out1 = ancom(phyloseq = family_data,  p_adj_method = "holm",
#              prv_cut = 0.10, lib_cut = 0, main_var = "age",
#              adj_formula = "region + bmi_group", rand_formula = NULL,
#              lme_control = NULL, struc_zero = FALSE, neg_lb = FALSE,
#              alpha = 0.05, n_cl = 2)

Visualization for W statistics

q_val = out1$q_data
beta_val = out1$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(out1$zero_ind), 
                nrow(feature_table), 
                sum(apply(out1$zero_ind, 1, sum) == 0))
# Cutoff values for declaring differentially abundant taxa
cut_off = 0.7 * (n_taxa - 1)

df_fig_w = res1 %>%
  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_id = factor(df_fig_w$taxon_id, levels = df_fig_w$taxon_id)
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_id, 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

Compare ANCOM-BC and ANCOM results

df_compare = res$diff_abn %>%
  tibble::rownames_to_column("taxon_id") %>%
  dplyr::left_join(res1, by = "taxon_id") %>%
  dplyr::transmute(ancombc = bmi_groupobese,
                   ancom = detected_0.7)
vdc = vennCounts(df_compare)
class(vdc) = "matrix"
df_vdc = as.data.frame(vdc)[-1, ] %>%
  dplyr::mutate(x = c(1, -1, 0),
                y = c(0, 0, 0)) %>%
  dplyr::transmute(x, y, label = Counts)

df_fig_venn = data.frame(x = c(-0.866, 0.866),
                     y = c(0, 0),
                     labels = c("ANCOM-BC", "ANCOM"))
df_fig_venn$labels = factor(df_fig_venn$labels, 
                            levels = c("ANCOM-BC", "ANCOM"))
p_venn = df_fig_venn %>%
  ggplot(aes(x0 = x, y0 = y, r = 1.5, fill = labels)) +
  geom_circle(alpha = .3, size = 1, color = "grey") +
  geom_text(data = df_vdc, aes(x, y, label = label),
            inherit.aes = FALSE, size = 8)+
  labs(fill = NULL) +
  coord_fixed() +
  theme_void()
p_venn

5.2 Longitudinal data

Run ancom function

set.seed(123)
out2 = ancom(phyloseq = family_data2,  p_adj_method = "holm",
             prv_cut = 0.10, lib_cut = 1000, main_var = "group", 
             adj_formula = "nationality + timepoint + bmi_group", 
             rand_formula = "~ timepoint | subject",
             lme_control = list(maxIter = 100, msMaxIter = 100, opt = "optim"), 
             struc_zero = TRUE, neg_lb = TRUE, alpha = 0.05, n_cl = 2)

res2 = out2$res

Visualization for W statistics

q_val = out2$q_data
beta_val = out2$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(out2$zero_ind), 
                nrow(feature_table), 
                sum(apply(out2$zero_ind, 1, sum) == 0))
# Cutoff values for declaring differentially abundant taxa
cut_off = 0.7 * (n_taxa - 1)

df_fig_w = res2 %>%
  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_id = factor(df_fig_w$taxon_id, levels = df_fig_w$taxon_id)
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_id, 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

6. SECOM Implementation

6.1 Run secom functions

set.seed(123)
# Linear relationships
res_linear = secom_linear(pseqs = list(c(genus_data, phylum_data)), pseudo = 0, 
                          prv_cut = 0.5, lib_cut = 1000, corr_cut = 0.5, 
                          wins_quant = c(0.05, 0.95), method = "pearson", 
                          soft = FALSE, thresh_len = 20, n_cv = 10, 
                          thresh_hard = 0.3, max_p = 0.005, n_cl = 2)

# Nonlinear relationships
res_dist = secom_dist(pseqs = list(c(genus_data, phylum_data)), pseudo = 0, 
                      prv_cut = 0.5, lib_cut = 1000, corr_cut = 0.5, 
                      wins_quant = c(0.05, 0.95), R = 1000, 
                      thresh_hard = 0.3, max_p = 0.005, n_cl = 2)

6.2 Visualizations

Pearson correlation with thresholding

corrplot(res_linear$corr_th, method = "color", addgrid.col = "grey",
         order = "alphabet", type = "lower", col.lim = c(-1, 1), 
         diag = FALSE, addCoef.col = "black", number.cex = 0.8,
         col = COL2("BrBG"), title = "Pearson (Thresholding)", 
         mar = c(1, 0, 1, 0))

Pearson correlation with p-value filtering

corrplot(res_linear$corr_fl, method = "color", addgrid.col = "grey",
         order = "alphabet", type = "lower", col.lim = c(-1, 1), 
         diag = FALSE, addCoef.col = "black", number.cex = 0.8,
         col = COL2("BrBG"), title = "Pearson (Thresholding)", 
         mar = c(1, 0, 1, 0))

Distance correlation with p-value filtering

corrplot(res_dist$dcorr_fl, method = "color", addgrid.col = "grey",
         order = "alphabet", type = "lower", col.lim = c(-1, 1), 
         diag = FALSE, addCoef.col = "black", number.cex = 0.8,
         col = COL2("BrBG"), title = "Pearson (Thresholding)", 
         mar = c(1, 0, 1, 0))

Session information

sessionInfo()
R version 4.2.0 (2022-04-22)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 20.04.4 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so

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       

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

other attached packages:
 [1] DT_0.23           ggforce_0.3.3     limma_3.52.1      corrplot_0.92    
 [5] qwraps2_0.5.2     magrittr_2.0.3    microbiome_1.18.0 phyloseq_1.40.0  
 [9] forcats_0.5.1     stringr_1.4.0     dplyr_1.0.9       purrr_0.3.4      
[13] readr_2.1.2       tidyr_1.2.0       tibble_3.1.7      ggplot2_3.3.6    
[17] tidyverse_1.3.1   ANCOMBC_1.6.1    

loaded via a namespace (and not attached):
  [1] readxl_1.4.0           backports_1.4.1        Hmisc_4.7-0           
  [4] plyr_1.8.7             igraph_1.3.1           splines_4.2.0         
  [7] crosstalk_1.2.0        GenomeInfoDb_1.32.2    digest_0.6.29         
 [10] foreach_1.5.2          htmltools_0.5.2        fansi_1.0.3           
 [13] checkmate_2.1.0        cluster_2.1.3          doParallel_1.0.17     
 [16] tzdb_0.3.0             Biostrings_2.64.0      modelr_0.1.8          
 [19] jpeg_0.1-9             colorspace_2.0-3       rvest_1.0.2           
 [22] haven_2.5.0            rbibutils_2.2.8        xfun_0.31             
 [25] crayon_1.5.1           RCurl_1.98-1.7         jsonlite_1.8.0        
 [28] Exact_3.1              survival_3.3-1         iterators_1.0.14      
 [31] ape_5.6-2              glue_1.6.2             polyclip_1.10-0       
 [34] gtable_0.3.0           zlibbioc_1.42.0        XVector_0.36.0        
 [37] Rhdf5lib_1.18.2        BiocGenerics_0.42.0    scales_1.2.0          
 [40] mvtnorm_1.1-3          DBI_1.1.2              rngtools_1.5.2        
 [43] Rcpp_1.0.8.3           htmlTable_2.4.0        foreign_0.8-82        
 [46] proxy_0.4-27           Formula_1.2-4          stats4_4.2.0          
 [49] htmlwidgets_1.5.4      httr_1.4.3             RColorBrewer_1.1-3    
 [52] ellipsis_0.3.2         farver_2.1.0           pkgconfig_2.0.3       
 [55] nnet_7.3-17            sass_0.4.1             dbplyr_2.2.0          
 [58] utf8_1.2.2             labeling_0.4.2         tidyselect_1.1.2      
 [61] rlang_1.0.2            reshape2_1.4.4         munsell_0.5.0         
 [64] cellranger_1.1.0       tools_4.2.0            cli_3.3.0             
 [67] generics_0.1.2         ade4_1.7-19            broom_0.8.0           
 [70] evaluate_0.15          biomformat_1.24.0      fastmap_1.1.0         
 [73] yaml_2.3.5             knitr_1.39             fs_1.5.2              
 [76] rootSolve_1.8.2.3      nlme_3.1-157           doRNG_1.8.2           
 [79] xml2_1.3.3             compiler_4.2.0         rstudioapi_0.13       
 [82] png_0.1-7              e1071_1.7-11           reprex_2.0.1          
 [85] tweenr_1.0.2           bslib_0.3.1            DescTools_0.99.45     
 [88] stringi_1.7.6          gsl_2.1-7.1            highr_0.9             
 [91] lattice_0.20-45        Matrix_1.4-1           nloptr_2.0.3          
 [94] vegan_2.6-2            permute_0.9-7          multtest_2.52.0       
 [97] vctrs_0.4.1            pillar_1.7.0           lifecycle_1.0.1       
[100] rhdf5filters_1.8.0     Rdpack_2.3.1           jquerylib_0.1.4       
[103] data.table_1.14.2      bitops_1.0-7           lmom_2.9              
[106] R6_2.5.1               latticeExtra_0.6-29    gridExtra_2.3         
[109] IRanges_2.30.0         gld_2.6.4              codetools_0.2-18      
[112] boot_1.3-28            energy_1.7-10          MASS_7.3-57           
[115] assertthat_0.2.1       rhdf5_2.40.0           withr_2.5.0           
[118] S4Vectors_0.34.0       GenomeInfoDbData_1.2.8 mgcv_1.8-40           
[121] expm_0.999-6           parallel_4.2.0         hms_1.1.1             
[124] grid_4.2.0             rpart_4.1.16           class_7.3-20          
[127] rmarkdown_2.14         Rtsne_0.16             Biobase_2.56.0        
[130] lubridate_1.8.0        base64enc_0.1-3       

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.

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.