View source: R/da_metagenomeSeq.R
run_metagenomeseq | R Documentation |
Differential expression analysis based on the Zero-inflated Log-Normal mixture model or Zero-inflated Gaussian mixture model using metagenomeSeq.
run_metagenomeseq(
ps,
group,
confounders = character(0),
contrast = NULL,
taxa_rank = "all",
transform = c("identity", "log10", "log10p", "SquareRoot", "CubicRoot", "logit"),
norm = "CSS",
norm_para = list(),
method = c("ZILN", "ZIG"),
p_adjust = c("none", "fdr", "bonferroni", "holm", "hochberg", "hommel", "BH", "BY"),
pvalue_cutoff = 0.05,
...
)
ps |
ps a |
group |
character, the variable to set the group, must be one of the var of the sample metadata. |
confounders |
character vector, the confounding variables to be adjusted.
default |
contrast |
this parameter only used for two groups comparison while there are multiple groups. For more please see the following details. |
taxa_rank |
character to specify taxonomic rank to perform
differential analysis on. Should be one of |
transform |
character, the methods used to transform the microbial
abundance. See
|
norm |
the methods used to normalize the microbial abundance data. See
|
norm_para |
arguments passed to specific normalization methods. |
method |
character, which model used for differential analysis, "ZILN" (Zero-inflated Log-Normal mixture model)" or "ZIG" (Zero-inflated Gaussian mixture model). And the zero-inflated log-normal model is preferred due to the high sensitivity and low FDR. |
p_adjust |
method for multiple test correction, default |
pvalue_cutoff |
numeric, p value cutoff, default 0.05 |
... |
extra arguments passed to the model. more details see
|
metagnomeSeq provides two differential analysis methods, zero-inflated
log-normal mixture model (implemented in
metagenomeSeq::fitFeatureModel()
) and zero-inflated Gaussian mixture
model (implemented in metagenomeSeq::fitZig()
). We recommend
fitFeatureModel over fitZig due to high sensitivity and low FDR. Both
metagenomeSeq::fitFeatureModel()
and metagenomeSeq::fitZig()
require
the abundance profiles before normalization.
For metagenomeSeq::fitZig()
, the output column is the coefficient of
interest, and logFC column in the output of
metagenomeSeq::fitFeatureModel()
is analogous to coefficient. Thus,
logFC is really just the estimate the coefficient of interest in
metagenomeSeq::fitFeatureModel()
. For more details see
these question Difference between fitFeatureModel and fitZIG in metagenomeSeq.
contrast
must be a two length character or NULL
(default). It is only
required to set manually for two groups comparison when there are multiple
groups. The order determines the direction of comparison, the first element
is used to specify the reference group (control). This means that, the first
element is the denominator for the fold change, and the second element is
used as baseline (numerator for fold change). Otherwise, users do required
to concern this paramerter (set as default NULL
), and if there are
two groups, the first level of groups will set as the reference group; if
there are multiple groups, it will perform an ANOVA-like testing to find
markers which difference in any of the groups.
Of note, metagenomeSeq::fitFeatureModel()
is not allows for multiple
groups comparison.
a microbiomeMarker
object.
Yang Cao
Paulson, Joseph N., et al. "Differential abundance analysis for microbial marker-gene surveys." Nature methods 10.12 (2013): 1200-1202.
data(enterotypes_arumugam)
ps <- phyloseq::subset_samples(
enterotypes_arumugam,
Enterotype %in% c("Enterotype 3", "Enterotype 2")
)
run_metagenomeseq(ps, group = "Enterotype")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.