plot_deseq2_pq | R Documentation |
Graphical representation of DESeq2 analysis.
plot_deseq2_pq(
data,
contrast = NULL,
tax_table = NULL,
pval = 0.05,
taxolev = "Genus",
select_taxa = NULL,
color_tax = "Phylum",
tax_depth = NULL,
verbose = TRUE,
jitter_width = 0.1,
...
)
data |
(required) a |
contrast |
(required) contrast specifies what comparison to extract
from the object to build a results table. See |
tax_table |
Required if data is a
|
pval |
(default: 0.05) the significance cutoff used for optimizing the independent filtering. If the adjusted p-value cutoff (FDR) will be a value other than 0.05, pval should be set to that value. |
taxolev |
taxonomic level of interest |
select_taxa |
Either the name of the taxa (in the form of |
color_tax |
taxonomic level used for color or a color vector. |
tax_depth |
Taxonomic depth to test for differential
distribution among contrast. If Null the analysis is done at the OTU
(i.e. Species) level. If not Null, data need to be a column name in
the |
verbose |
whether the function print some information during the computation |
jitter_width |
width for the jitter positioning |
... |
Additional arguments passed on to |
Please cite DESeq2
package if you use chis function.
A ggplot
2 plot representing DESeq2 results
Adrien Taudière
DESeq
results
plot_edgeR_pq
data("GlobalPatterns", package = "phyloseq")
GP <- subset_taxa(GlobalPatterns, GlobalPatterns@tax_table[, 1] == "Archaea")
GP <- subset_samples(GP, SampleType %in% c("Soil", "Skin"))
if (requireNamespace("DESeq2")) {
res <- DESeq2::DESeq(phyloseq_to_deseq2(GP, ~SampleType),
test = "Wald", fitType = "local"
)
plot_deseq2_pq(res, c("SampleType", "Soil", "Skin"),
tax_table = GP@tax_table, color_tax = "Kingdom"
)
plot_deseq2_pq(res, c("SampleType", "Soil", "Skin"),
tax_table = GP@tax_table, color_tax = "Kingdom",
pval = 0.7
)
plot_deseq2_pq(res, c("SampleType", "Soil", "Skin"),
tax_table = GP@tax_table, color_tax = "Class",
select_taxa = c("522457", "271582")
)
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.