plotAbundance.Rd
plotAbundance
plots the abundance on a selected taxonomic rank.
Since this probably makes sense only for relative abundance data, the
assay used by default is expected to be in the slot ‘relabundance’.
If only ‘counts’ is present, the relative abundance is computed.
plotAbundance(x, ...)
# S4 method for SummarizedExperiment
plotAbundance(
x,
rank = taxonomyRanks(x)[1],
features = NULL,
order_rank_by = c("name", "abund", "revabund"),
order_sample_by = NULL,
decreasing = TRUE,
use_relative = TRUE,
layout = c("bar", "point"),
one_facet = TRUE,
ncol = 2,
scales = "fixed",
assay.type = assay_name,
assay_name = "counts",
...
)
a
SummarizedExperiment
object.
additional parameters for plotting. See
mia-plot-args
for more details i.e. call help("mia-plot-args")
a single character
value defining the taxonomic rank to
use. Must be a value of taxonomyRanks(x)
.
a single character
value defining a column from
colData
to be plotted below the abundance plot.
Continuous numeric values will be plotted as point, whereas factors and
character will be plotted as colour-code bar. (default: features =
NULL
)
How to order abundance value: By name (“name”) for sorting the taxonomic labels alphabetically, by abundance (“abund”) to sort by abundance values or by a reverse order of abundance values (“revabund”).
A single character value from the chosen rank of abundance
data or from colData
to select values to order the abundance
plot by. (default: order_sample_by = NULL
)
TRUE or FALSE: If the order_sample_by
is defined and the
values are numeric, should the values used to order in decreasing or
increasing fashion? (default: decreasing = FALSE
)
TRUE
or FALSE
: Should the relative values
be calculated? (default: use_relative = TRUE
)
Either “bar” or “point”.
Should the plot be returned in on facet or split into
different facet, one facet per different value detect in rank
. If
features
or order_sample_by
is not NULL
, this setting will
be disregarded.
if one_facet = FALSE
, ncol
defines many
columns should be for plotting the different facets and scales
is
used to define the behavior of the scales of each facet. Both values are
passed onto facet_wrap
.
a character
value defining which assay data to
use. (default: assay.type = "relabundance"
)
a single character
value for specifying which
assay to use for calculation.
(Please use assay.type
instead. At some point assay_name
will be disabled.)
Subsetting to rows of interested and ordering of those is expected to be done
outside of this functions, e.g. x[1:2,]
. This will plot data of all
features present.
data(GlobalPatterns, package="mia")
se <- GlobalPatterns
## Plotting abundance using the first taxonomic rank as default
plotAbundance(se, assay.type="counts")
## Using "Phylum" as rank
plotAbundance(se, assay.type="counts", rank = "Phylum", add_legend = FALSE)
## If rank is set to NULL plotAbundance behaves like plotExpression
plotAbundance(se, assay.type="counts", rank = NULL,
features = head(rownames(se)))
## A feature from colData or taxon from chosen rank can be used for ordering samples.
plotAbundance(se, assay.type="counts", rank = "Phylum",
order_sample_by = "Bacteroidetes")
## Features from colData can be plotted together with abundance plot.
# Returned object is a list that includes two plot; other visualizes abundance
# other features.
plot <- plotAbundance(se, assay.type = "counts", rank = "Phylum",
features = "SampleType")
# \donttest{
# These two plots can be combined with wrap_plots function from patchwork package
library(patchwork)
wrap_plots(plot, ncol = 1)
# }
## Same plot as above but showing sample IDs as labels for the x axis on the top plot
plot[[1]] <- plotAbundance(se, assay.type = "counts", rank = "Phylum",
features = "SampleType", add_legend = FALSE,
add_x_text = TRUE)[[1]] +
theme(axis.text.x = element_text(angle = 90))
# \donttest{
wrap_plots(plot, ncol = 1, heights = c(0.8,0.2))
# }
## Compositional barplot with top 5 taxa and samples sorted by "Bacteroidetes"
# Getting top taxa on a Phylum level
se <- transformAssay(se, method="relabundance")
se_phylum <- agglomerateByRank(se, rank ="Phylum", onRankOnly=TRUE)
top_taxa <- getTopFeatures(se_phylum,top = 5, assay.type = "relabundance")
# Renaming the "Phylum" rank to keep only top taxa and the rest to "Other"
phylum_renamed <- lapply(rowData(se)$Phylum,
function(x){if (x %in% top_taxa) {x} else {"Other"}})
rowData(se)$Phylum <- as.character(phylum_renamed)
# Compositional barplot
plotAbundance(se, assay.type="relabundance", rank = "Phylum",
order_rank_by="abund", order_sample_by = "Bacteroidetes")