Skip to contents

This function calculates beta diversity measures between time points, performs linear modeling, and generates results.

Usage

generate_beta_change_test_pair(
  data.obj,
  dist.obj = NULL,
  time.var = NULL,
  subject.var,
  group.var,
  adj.vars,
  change.base = NULL,
  dist.name = c("BC", "Jaccard", "UniFrac", "GUniFrac", "WUniFrac", "JS")
)

Arguments

data.obj

A list object in a format specific to MicrobiomeStat, which can include components such as feature.tab (matrix), feature.ann (matrix), meta.dat (data.frame), tree, and feature.agg.list (list). The data.obj can be converted from other formats using several functions from the MicrobiomeStat package, including: 'mStat_convert_DGEList_to_data_obj', 'mStat_convert_DESeqDataSet_to_data_obj', 'mStat_convert_phyloseq_to_data_obj', 'mStat_convert_SummarizedExperiment_to_data_obj', 'mStat_import_qiime2_as_data_obj', 'mStat_import_mothur_as_data_obj', 'mStat_import_dada2_as_data_obj', and 'mStat_import_biom_as_data_obj'. Alternatively, users can construct their own data.obj. Note that not all components of data.obj may be required for all functions in the MicrobiomeStat package.

dist.obj

Distance matrix between samples, usually calculated using mStat_calculate_beta_diversity function. If NULL, beta diversity will be automatically computed from data.obj using mStat_calculate_beta_diversity.

time.var

The name of the column in metadata containing the time variable. This should be a column with time points for each sample. Required to identify pairs of samples for the same subject across time.

subject.var

The name of the column in metadata containing the subject IDs. This should uniquely identify each subject in the study. Required to identify samples that belong to the same subject.

group.var

The name of the column in metadata containing the grouping variable to use in linear modeling. This grouping variable will be used as a predictor in the linear models for beta diversity change. Optional.

adj.vars

Vector of names of additional variables in metadata to be included as covariates in the linear models. Can be empty.

change.base

The baseline time point value in the time variable to be used as the reference for calculating beta diversity change. Required if time.var contains multiple time points. If NULL, the first time point will be used as change.base automatically.

dist.name

A character vector specifying which beta diversity indices to calculate. Supported indices are "BC" (Bray-Curtis), "Jaccard", "UniFrac" (unweighted UniFrac), "GUniFrac" (generalized UniFrac), "WUniFrac" (weighted UniFrac), and "JS" (Jensen-Shannon divergence). If a name is provided but the corresponding object does not exist within dist.obj, it will be computed internally. If the specific index is not supported, an error message will be returned.

Value

A named list containing linear modeling results for each beta diversity metric. Each list element corresponds to one of the distance metrics specified in dist.name. It contains a coefficient table from fitting a linear model with the beta diversity change as response and the group_var and adj_vars as predictors. If group_var has multiple levels, ANOVA results are also included after the coefficients. Column names include: Term, Estimate, Std.Error, Statistic, P.Value

Details

For each subject, it calculates the change in beta diversity between baseline time point and later time point. This within-subject change is used as the outcome in linear models.

Adjustment for covariates is supported. Both ANOVA and coefficient tables are provided if a grouping variable is specified.

Examples

if (FALSE) { # \dontrun{

# Load packages
library(vegan)

# Load data
data(peerj32.obj)
generate_beta_change_test_pair(
  data.obj = peerj32.obj,
  dist.obj = NULL,
  time.var = "time",
  subject.var = "subject",
  group.var = "group",
  adj.vars = NULL,
  change.base = "1",
  dist.name = c('BC', 'Jaccard')
)
generate_beta_change_test_pair(
  data.obj = peerj32.obj,
  dist.obj = NULL,
  time.var = "time",
  subject.var = "subject",
  group.var = "group",
  adj.vars = "sex",
  change.base = "1",
  dist.name = c('BC', 'Jaccard')
)

data("subset_pairs.obj")
generate_beta_change_test_pair(
  data.obj = subset_pairs.obj,
  dist.obj = NULL,
  time.var = "Antibiotic",
  subject.var = "MouseID",
  group.var = "Sex",
  adj.vars = NULL,
  change.base = "Baseline",
  dist.name = c('BC', 'Jaccard')
)
} # }