| plotPrevalence {miaViz} | R Documentation |
plotPrevalence and plotTaxaPrevalence visualize prevalence
information.
plotPrevalence(x, ...) ## S4 method for signature 'SummarizedExperiment' plotPrevalence( x, detections = c(0.01, 0.1, 1, 2, 5, 10, 20)/100, prevalences = seq(0.1, 1, 0.1), abund_values = "counts", as_relative = TRUE, rank = NULL, BPPARAM = BiocParallel::SerialParam(), ... ) plotPrevalentAbundance(x, ...) ## S4 method for signature 'SummarizedExperiment' plotPrevalentAbundance( x, rank = taxonomyRanks(x)[1L], abund_values = "counts", as_relative = TRUE, colour_by = NULL, size_by = NULL, shape_by = NULL, label = NULL, facet_by = NULL, ... ) plotTaxaPrevalence(x, ...) ## S4 method for signature 'SummarizedExperiment' plotTaxaPrevalence( x, rank = taxonomyRanks(x)[1L], abund_values = "counts", detections = NULL, ndetections = 20, as_relative = TRUE, min_prevalence = 0, BPPARAM = BiocParallel::SerialParam(), ... )
x |
a
|
detections |
Detection thresholds for absence/presence. Either an
absolutes value compared directly to the values of |
prevalences |
Prevalence thresholds (in 0 to 1). The
required prevalence is strictly greater by default. To include the
limit, set |
abund_values |
a |
as_relative |
logical scalar: Should the detection threshold be applied
on compositional (relative) abundances? Passed onto
|
rank, ... |
additional arguments
|
BPPARAM |
A
|
colour_by |
Specification of a feature to colour points by, see the
|
size_by |
Specification of a feature to size points by, see the
|
shape_by |
Specification of a feature to shape points by, see the
|
label |
a |
facet_by |
Taxonomic rank to facet the plot by.
Value must be of |
ndetections |
If |
min_prevalence |
a single numeric value to apply as a threshold for
plotting. The threshold is applied per row and column.
(default: |
Whereas plotPrevalence procudes a line plot, plotTaxaPrevalence
returns a heatmap.
Agglomeration on different taxonomic levels is available through the
rank argument.
To exclude certain taxa, preprocess x to your liking, for example
with subsetting via getPrevalentTaxa or
agglomerateByPrevalence.
A ggplot2 object or plotly object, if more than one
prevalences was defined.
getPrevalence,
agglomerateByPrevalence,
agglomerateByRank
data(GlobalPatterns, package = "mia")
# plotting N of prevalence exceeding taxa on the Phylum level
plotPrevalence(GlobalPatterns, rank = "Phylum")
plotPrevalence(GlobalPatterns, rank = "Phylum") + scale_x_log10()
# plotting prevalence per taxa for different detectio thresholds as heatmap
plotTaxaPrevalence(GlobalPatterns, rank = "Phylum")
# by default a continous scale is used for different detections levels,
# but this can be adjusted
plotTaxaPrevalence(GlobalPatterns, rank = "Phylum",
detections = c(0, 0.001, 0.01, 0.1, 0.2))
# point layout for plotTaxaPrevalence can be used to visualize by additional
# information
plotPrevalentAbundance(GlobalPatterns, rank = "Family",
colour_by = "Phylum") +
scale_x_log10()
# When using function plotPrevalentAbundace, it is possible to create facets
# with 'facet_by'.
plotPrevalentAbundance(GlobalPatterns, rank = "Family",
colour_by = "Phylum", facet_by = "Kingdom") +
scale_x_log10()