edgerWrp | R Documentation |
edgerWrp
is a wrapper using functions from the edgeR
package (Robinson et al. 2010, Bioinformatics; McCarthy et al. 2012,
Nucleic Acids Research) to fit models and perform a moderated test for
each entity.
edgerWrp(
count,
lib_size = NULL,
option = c("glm", "glmQL"),
design,
contrast = NULL,
normalize = TRUE,
normalize_method = "TMM",
...
)
count |
A matrix with features (e.g., genes or microbes) in rows and samples in columns. |
lib_size |
A numeric vector with library sizes for each sample. If
|
option |
Either |
design |
A numeric design matrix, e.g. created by
|
contrast |
A numeric vector specifying one contrast of
the linear model coefficients to be tested. Its length
must equal the number of columns of |
normalize |
A logical scalar, specifying whether normalization factors
should be calculated (using |
normalize_method |
Normalization method to be used. Please refer to
|
... |
More arguments to pass to |
The function performs the following steps:
Create a DGEList
object. If lib_size
is
given, set the library sizes to these values, otherwise use the column sums
of the count matrix.
If normalize
is TRUE
, estimate normalization factors
using calcNormFactors
.
Estimate dispersions with estimateDisp
.
Depending on the value of option
, apply either the LRT or
QLF edgeR workflows (i.e., either glmFit
+
glmLRT
or glmQLFit
+
glmQLFTest
), testing for the specified contrast.
The output of glmQLFTest
or
glmLRT
depending on the specified option
.
Ruizhu Huang
suppressPackageStartupMessages({
library(TreeSummarizedExperiment)
})
## Read example data
x <- readRDS(system.file("extdata/da_sim_100_30_18de.rds",
package = "treeclimbR"))
## Run differential abundance analysis
out <- edgerWrp(count = assay(x), option = "glm",
design = model.matrix(~ group, data = colData(x)),
contrast = c(0, 1))
## The output is an edgeR DGELRT object
class(out)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.