diff_expression | R Documentation |
The diff_expression()
function performs differential expression analysis
using a method of preference.
A corresponding autoplot()
method is visualizing the results as a volcano plot.
diff_expression(object, group, method = c("voom", "deseq2"), ...)
## S4 method for signature 'HermesDataDiffExpr'
autoplot(object, adj_p_val_thresh = 0.05, log2_fc_thresh = 2.5)
object |
( |
group |
( |
method |
( |
... |
additional arguments passed to the helper function associated with the selected method. |
adj_p_val_thresh |
( |
log2_fc_thresh |
( |
Possible method choices are:
voom: uses limma::voom()
, see h_diff_expr_voom()
for details.
deseq2
: uses DESeq2::DESeq()
, see h_diff_expr_deseq2()
for details.
A HermesDataDiffExpr
object which is a data frame with the following columns for each gene
in the HermesData
object:
log2_fc
(the estimate of the log2 fold change between the 2 levels of the
provided factor)
stat
(the test statistic, which one depends on the method used)
p_val
(the raw p-value)
adj_p_val
(the multiplicity adjusted p-value value)
autoplot(HermesDataDiffExpr)
: generates a volcano plot for a HermesDataDiffExpr
object.
We provide the df_cols_to_factor()
utility function that makes it easy to convert the
colData()
character and logical variables to factors, so that they can be subsequently
used as group
inputs. See the example.
In order to avoid a warning when using deseq2
, it can be necessary to specify
fitType = "local"
as additional argument. This could e.g. be the case when only few samples
are present in which case the default parametric dispersions estimation will not work.
object <- hermes_data %>%
add_quality_flags() %>%
filter()
# Convert character and logical to factor variables in `colData`,
# including the below used `group` variable.
colData(object) <- df_cols_to_factor(colData(object))
res1 <- diff_expression(object, group = "SEX", method = "voom")
head(res1)
res2 <- diff_expression(object, group = "SEX", method = "deseq2")
head(res2)
# Pass method arguments to the internally used helper functions.
res3 <- diff_expression(object, group = "SEX", method = "voom", robust = TRUE, trend = TRUE)
head(res3)
res4 <- diff_expression(object, group = "SEX", method = "deseq2", fitType = "local")
head(res4)
# Create the corresponding volcano plots.
autoplot(res1)
autoplot(res3)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.