DA_metagenomeSeq {benchdamic} | R Documentation |
Fast run for the metagenomeSeq's differential abundance detection method.
DA_metagenomeSeq( object, pseudo_count = FALSE, design = NULL, coef = 2, norm = c("TMM", "TMMwsp", "RLE", "upperquartile", "posupperquartile", "none", "ratio", "poscounts", "iterate", "TSS", "CSSmedian", "CSSdefault"), verbose = TRUE )
object |
phyloseq object. |
pseudo_count |
add 1 to all counts if TRUE (default
|
design |
The model for the count distribution. Can be the variable name, or a character similar to "~ 1 + group", or a formula, or a 'model.matrix' object. |
coef |
integer or character index vector indicating which coefficients of the linear model are to be tested equal to zero. |
norm |
name of the normalization method used to compute the
scaling factors to use in the differential abundance analysis. If |
verbose |
an optional logical value. If |
A list object containing the matrix of p-values 'pValMat', the matrix of summary statistics for each tag 'statInfo', and a suggested 'name' of the final object considering the parameters passed to the function.
fitZig
for the Zero-Inflated Gaussian
regression model estimation and MRfulltable
for results extraction.
set.seed(1) # Create a very simple phyloseq object counts <- matrix(rnbinom(n = 60, size = 3, prob = 0.5), nrow = 10, ncol = 6) metadata <- data.frame("Sample" = c("S1", "S2", "S3", "S4", "S5", "S6"), "group" = as.factor(c("A", "A", "A", "B", "B", "B"))) ps <- phyloseq::phyloseq(phyloseq::otu_table(counts, taxa_are_rows = TRUE), phyloseq::sample_data(metadata)) # Calculate the CSSdefault scaling factors ps_NF <- norm_CSS(object = ps, method = "default") # The phyloseq object now contains the scaling factors: scaleFacts <- phyloseq::sample_data(ps_NF)[, "NF.CSSdefault"] head(scaleFacts) # Differential abundance DA_metagenomeSeq(object = ps_NF, pseudo_count = FALSE, design = ~ group, coef = 2, norm = "CSSdefault")