View source: R/differentialTest.R
differentialTest | R Documentation |
Identify differentially-abundant and differentially-variable taxa
differentialTest(
formula,
phi.formula,
formula_null,
phi.formula_null,
data,
link = "logit",
phi.link = "logit",
test,
boot = FALSE,
B = 1000,
sample_data = NULL,
taxa_are_rows = TRUE,
filter_discriminant = TRUE,
fdr_cutoff = 0.05,
fdr = "fdr",
full_output = FALSE,
inits = NULL,
inits_null = NULL,
try_only = NULL,
verbose = FALSE,
robust = FALSE,
...
)
formula |
an object of class |
phi.formula |
an object of class |
formula_null |
Formula for mean under null, without response |
phi.formula_null |
Formula for overdispersion under null, without response |
data |
a data frame containing the OTU table, or |
link |
link function for abundance covariates, defaults to |
phi.link |
link function for dispersion covariates, defaults to |
test |
Character. Hypothesis testing procedure to use. One of |
boot |
Boolean. Defaults to |
B |
Optional integer. Number of bootstrap iterations. Ignored if |
sample_data |
Data frame or matrix. Defaults to |
taxa_are_rows |
Boolean. Optional. If |
filter_discriminant |
Boolean. Defaults to |
fdr_cutoff |
Integer. Defaults to |
fdr |
Character. Defaults to |
full_output |
Boolean. Opetional. Defaults to |
inits |
Optional initializations for model fit using |
inits_null |
Optional initializations for model fit using |
try_only |
Optional numeric. Will try only the |
verbose |
Boolean. Defaults to |
robust |
Should robust standard errors be used? If not, model-based standard errors are used. Logical, defaults to |
... |
Optional additional arguments for |
See package vignette for details and example usage. Make sure the number of columns in all of the initializations are correct! inits
probably shouldn't match inits_null
. To use a contrast matrix, see contrastsTest
.
An object of class differentialTest
. List with elements p
containing the p-values, p_fdr
containing the p-values after false discovery rate control, significant_taxa
containing the taxa names of the statistically significant taxa, significant_models
containing a list of the model fits for the significant taxa, all_models
containing a list of the model fits for all taxa, restrictions_DA
containing a list of covariates that were tested for differential abundance, restrictions_DV
containing a list of covariates that were tested for differential variability, discriminant_taxa_DA
containing the taxa for which at least one covariate associated with the abundance was perfectly discriminant, discriminant_taxa_DV
containing the taxa for which at least one covariate associated with the dispersion was perfectly discriminant, data
containing the data used to fit the models. If full_output = TRUE
, it will also include full_output
, a list of all model output from bbdml
.
# data frame example
data(soil_phylum_small_sample)
data(soil_phylum_small_otu)
da_analysis <- differentialTest(formula = ~ DayAmdmt,
phi.formula = ~ DayAmdmt,
formula_null = ~ 1,
phi.formula_null = ~ DayAmdmt,
test = "Wald", boot = FALSE,
data = soil_phylum_small_otu,
sample_data = soil_phylum_small_sample,
fdr_cutoff = 0.05,
try_only = 1:5)
# phyloseq example (only run if you have phyloseq installed)
## Not run:
data_phylo <- phyloseq::phyloseq(phyloseq::sample_data(soil_phylum_small_sample),
phyloseq::otu_table(soil_phylum_small_otu, taxa_are_rows = TRUE))
da_analysis <- differentialTest(formula = ~ DayAmdmt,
phi.formula = ~ DayAmdmt,
formula_null = ~ 1,
phi.formula_null = ~ DayAmdmt,
test = "Wald", boot = FALSE,
data = data_phylo,
fdr_cutoff = 0.05,
try_only = 1:5)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.