performMirnaDE()
and performGeneDE()
are two functions provided by MIRit
to conduct miRNA and gene differential expression analysis, respectively.
In particular, these functions allow the user to compute differential
expression through different methods, namely edgeR
(Quasi-Likelihood
framework), DESeq2
, limma-voom
and limma
. Data deriving from NGS
experiments and microarray technology are all suitable for these functions.
For precise indications about how to use these functions, please refer to
the details section.
Usage
performMirnaDE(
mirnaObj,
group,
contrast,
design,
method = "edgeR",
logFC = 0,
pCutoff = 0.05,
pAdjustment = "fdr",
filterByExpr.args = list(),
calcNormFactors.args = list(),
estimateDisp.args = list(robust = TRUE),
glmQLFit.args = list(),
glmQLFTest.args = list(),
DESeq.args = list(),
useVoomWithQualityWeights = TRUE,
voom.args = list(),
lmFit.args = list(),
eBayes.args = list(),
useArrayWeights = TRUE,
useWsva = FALSE,
wsva.args = list(),
arrayWeights.args = list(),
useDuplicateCorrelation = FALSE,
correlationBlockVariable = NULL,
duplicateCorrelation.args = list()
)
performGeneDE(
mirnaObj,
group,
contrast,
design,
method = "edgeR",
logFC = 0,
pCutoff = 0.05,
pAdjustment = "fdr",
filterByExpr.args = list(),
calcNormFactors.args = list(),
estimateDisp.args = list(robust = TRUE),
glmQLFit.args = list(),
glmQLFTest.args = list(),
DESeq.args = list(),
useVoomWithQualityWeights = TRUE,
voom.args = list(),
lmFit.args = list(),
eBayes.args = list(),
useArrayWeights = TRUE,
useWsva = FALSE,
wsva.args = list(),
arrayWeights.args = list(),
useDuplicateCorrelation = FALSE,
correlationBlockVariable = NULL,
duplicateCorrelation.args = list()
)
Arguments
- mirnaObj
A
MirnaExperiment
object containing miRNA and gene data- group
The variable of interest for differential expression analysis. It must be the column name of a variable present in the metadata (colData) of a
MirnaExperiment
object. See the details section for additional information- contrast
A
character
object that specifies the groups to be compared during differential expression analysis, separated by a dash (e.g. 'disease-healthy'). Note that reference group must be the last one, for additional information see the details section- design
An R
formula
that indicates the model to fit. It must include the variable of interest (group
) together with eventual covariates (e.g. '~ 0 + disease + sex'). Please note thatgroup
variable must be the first one. See the details section for additional information- method
The statistical package used to compute differential expression. For NGS experiments, it must be one of
edgeR
(default),DESeq2
, andvoom
(for limma-voom). Instead, for microarray data, onlylimma
can be used- logFC
The minimum log2 fold change required to consider a gene as differentially expressed. Optional, default is 0
- pCutoff
The adjusted p-value cutoff to use for statistical significance. The default value is
0.05
- pAdjustment
The p-value correction method for multiple testing. It must be one of:
fdr
(default),BH
,none
,holm
,hochberg
,hommel
,bonferroni
,BY
- filterByExpr.args
A
list
object containing additional arguments passed toedgeR::filterByExpr()
function. It is used whenmethod
is set toedgeR
orvoom
- calcNormFactors.args
A
list
object containing additional arguments passed toedgeR::calcNormFactors()
function. It is used whenmethod
is set toedgeR
orvoom
- estimateDisp.args
A
list
object containing additional arguments passed toedgeR::estimateDisp()
function. It is used whenmethod
is set toedgeR
. Default islist(robust = TRUE)
to use the robust parameter- glmQLFit.args
A
list
object containing additional arguments passed toedgeR::glmQLFit()
function. It is used whenmethod
is set toedgeR
- glmQLFTest.args
A
list
object containing additional arguments passed toedgeR::glmQLFTest()
function. It is used whenmethod
is set toedgeR
- DESeq.args
A
list
object containing additional arguments passed toDESeq2::DESeq()
function. It is used whenmethod
is set toDESeq
- useVoomWithQualityWeights
Logical, whether to use the
limma::voomWithQualityWeights()
function or just thelimma::voom()
function. It is used whenmethod
is set tovoom
. Default is TRUE- voom.args
A
list
object containing additional arguments passed tolimma::voom()
function orlimma::voomWithQualityWeights()
function. It is used whenmethod
is set tovoom
- lmFit.args
A
list
object containing additional arguments passed tolimma::lmFit()
function. It is used whenmethod
is set tovoom
orlimma
- eBayes.args
A
list
object containing additional arguments passed tolimma::eBayes()
function. It is used whenmethod
is set tovoom
orlimma
- useArrayWeights
Logical, whether to use the
limma::arrayWeights()
function or not. It is used whenmethod
is set tolimma
. Default is TRUE- useWsva
Logical, whether to use the
limma::wsva()
function or not. It is used whenmethod
is set tolimma
. Default is FALSE- wsva.args
A
list
object containing additional arguments passed tolimma::wsva()
function. It is used whenmethod
is set tolimma
- arrayWeights.args
A
list
object containing additional arguments passed tolimma::arrayWeights()
function. It is used whenmethod
is set tolimma
- useDuplicateCorrelation
Logical, whether to use the
limma::duplicateCorrelation()
function or not. It is used whenmethod
is set tolimma
. Default is FALSE- correlationBlockVariable
It is the blocking variable to use for
limma::duplicateCorrelation()
. Default is NULL- duplicateCorrelation.args
A
list
object containing additional arguments passed tolimma::duplicateCorrelation()
function. It is used whenmethod
is set tolimma
Value
A MirnaExperiment
object containing differential
expression results. To access these results, the user may run the
mirnaDE()
and geneDE()
functions for miRNAs and genes, respectively.
Details
When performing differential expression for NGS experiments, count matrices
are detected and method
parameter must be one of edgeR
, DESeq2
,
and voom
. On the other hand, when dealing with microarray studies, only
limma
can be used.
To calculate differential expression, MIRit must be informed about the
variable of interest and the desired contrast. In particular, the group
parameter must be the name of a variable present in the metadata (colData)
of a MirnaExperiment
object, which specifies the
variable used to compute differential expression analysis, between the groups
indicated in contrast
. Specifically, contrast
must be a character vector
that defines the levels to compare separated by a dash. For example, if we
have a variable named 'condition', with two levels, namely 'disease' and
'healthy', we can identify differentially expressed genes in 'disease'
samples compared to 'healthy' subjects by specifying: group = 'condition'
and contrast = 'disease-healthy'
. Furthermore, the user needs to specify
the model to fit expression values. To do so, the user has to state the
model formula in the design
parameter. Please note that for a correct
inner working of these functions, the group
variable of interest must be
the first variable in model formula. Moreover, the user can include in the
design any other sources of variation by specifying covariates that will be
taken into account. For instance, if we want to compare 'disease' subjects
against 'healthy' individuals, without the influence of sex differences,
we may specify design = ~ condition + sex
, where 'sex' is also a
variable present in the metadata (colData) of mirnaObj
.
Notably, for all the methods available, the user can supply additional
arguments to the functions implemented in edgeR
, DESeq2
and limma
.
Therefore, the user has finer control over how the differential expression
analysis is performed. In this regard, for microarray studies, the user
may opt to include weighted surrogate variable analysis (WSVA) to correct
for unknown sources of variation (useWsva = TRUE
). Moreover, for
microarray data, the arrayWeights()
function in limma
can be used to
assess differential expression with respect to array qualities. Also, the
duplicateCorrelation()
function in limma
may be included in the pipeline
in order to block the effect of correlated samples. To do this, the user
must set useDuplicateCorrelation = TRUE
, and must specify the blocking
variable through the correlationBlockVariable
parameter. Additionally,
when using limma-voom
, the user may estimate voom transformation with or
without quality weights (by specifying useVoomWithQualityWeights = TRUE
).
Functions
performMirnaDE()
: Perform differential expression analysis for miRNAsperformGeneDE()
: Perform differential expression analysis for genes
References
Ritchie ME, Phipson B, Wu D, Hu Y, Law CW, Shi W, Smyth GK (2015). “limma powers differential expression analyses for RNA-sequencing and microarray studies.” Nucleic Acids Research, 43(7), e47. doi:10.1093/nar/gkv007.
Law, CW, Chen, Y, Shi, W, and Smyth, GK (2014). "Voom: precision weights unlock linear model analysis tools for RNA-seq read counts". Genome Biology 15, R29
Robinson MD, McCarthy DJ, Smyth GK (2010). “edgeR: a Bioconductor package for differential expression analysis of digital gene expression data.” Bioinformatics, 26(1), 139-140. doi:10.1093/bioinformatics/btp616.
Love MI, Huber W, Anders S (2014). “Moderated estimation of fold change and dispersion for RNA-seq data with DESeq2.” Genome Biology, 15, 550. doi:10.1186/s13059-014-0550-8.
Author
Jacopo Ronchi, jacopo.ronchi@unimib.it
Examples
# \donttest{
# load example MirnaExperiment object
obj <- loadExamples()
# perform miRNA DE with edgeR
obj <- performMirnaDE(obj,
group = "disease", contrast = "PTC-NTH",
design = ~ 0 + disease + patient, method = "edgeR"
)
#> Performing differential expression analysis with edgeR...
#> Differential expression analysis reported 61 significant miRNAs with p < 0.05 (correction: fdr). You can use the 'mirnaDE()' function to access results.
# perform miRNA DE with DESeq2
obj <- performMirnaDE(obj,
group = "disease", contrast = "PTC-NTH",
design = ~ 0 + disease + patient, method = "DESeq2"
)
#> Performing differential expression analysis with DESeq2...
#> Warning: some variables in design formula are characters, converting to factors
#> estimating size factors
#> estimating dispersions
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimates
#> fitting model and testing
#> Differential expression analysis reported 87 significant miRNAs with p < 0.05 (correction: fdr). You can use the 'mirnaDE()' function to access results.
# perform miRNA DE with limma-voom
obj <- performMirnaDE(obj,
group = "disease", contrast = "PTC-NTH",
design = ~ 0 + disease + patient, method = "voom"
)
#> Performing differential expression analysis with voom...
#> Differential expression analysis reported 55 significant miRNAs with p < 0.05 (correction: fdr). You can use the 'mirnaDE()' function to access results.
# }