plotPrevalence and plotRowPrevalence visualize prevalence information.

plotRowPrevalence(x, ...)

plotPrevalentAbundance(x, ...)

plotPrevalence(x, ...)

# S4 method for class 'SummarizedExperiment'
plotPrevalence(
  x,
  detection = detections,
  detections = c(0.01, 0.1, 1, 2, 5, 10, 20),
  prevalence = prevalences,
  prevalences = seq(0.1, 1, 0.1),
  assay.type = assay_name,
  assay_name = "counts",
  rank = NULL,
  BPPARAM = BiocParallel::SerialParam(),
  ...
)

# S4 method for class 'SummarizedExperiment'
plotPrevalentAbundance(
  x,
  rank = NULL,
  assay.type = assay_name,
  assay_name = "counts",
  colour.by = colour_by,
  colour_by = NULL,
  size.by = size_by,
  size_by = NULL,
  shape.by = shape_by,
  shape_by = NULL,
  show.label = label,
  label = NULL,
  facet.by = facet_by,
  facet_by = NULL,
  ...
)

# S4 method for class 'SummarizedExperiment'
plotRowPrevalence(
  x,
  rank = NULL,
  assay.type = assay_name,
  assay_name = "counts",
  detection = detections,
  detections = c(0.01, 0.1, 1, 2, 5, 10, 20),
  min.prevalence = min_prevalence,
  min_prevalence = 0,
  BPPARAM = BiocParallel::SerialParam(),
  ...
)

Arguments

x

a SummarizedExperiment object.

detection

Numeric scalar. Detection thresholds for absence/presence. Either an absolutes value compared directly to the values of x or a relative value between 0 and 1, if TRUE.

detections

Deprecated. Use detection instead.

prevalence

Numeric scalar. Prevalence thresholds (in 0 to 1). The required prevalence is strictly greater by default. To include the limit, set include.lowest to TRUE.

prevalences

Deprecated. Use prevalence instead.

assay.type

Character scalar. Defines which assay data to use. (Default: "relabundance")

assay_name

Deprecated. Use assay.type instead.

rank, ...

additional arguments

  • as.relative Logical scalar. Should the relative values be calculated? (Default: FALSE)

  • ndetection Integer scalar. Determines the number of breaks calculated detection thresholds when detection=NULL. When TRUE, as_relative is then also regarded as TRUE. (Default: 20)

  • If !is.null(rank) matching arguments are passed on to agglomerateByRank. See ?agglomerateByRank for more details.

  • additional arguments for plotting. See mia-plot-args for more details i.e. call help("mia-plot-args")

BPPARAM

A BiocParallelParam object specifying whether the UniFrac calculation should be parallelized.

colour.by

Character scalar. Specification of a feature to colour points by, see the by argument in ?retrieveFeatureInfo for possible values. Only used with layout = "point". (Default: NULL)

colour_by

Deprecated. Use colour.by instead.

size.by

Character scalar. Specification of a feature to size points by, see the by argument in ?retrieveFeatureInfo for possible values. Only used with layout = "point". (Default: NULL)

size_by

Deprecated. Use size.by instead.

shape.by

Character scalar. Specification of a feature to shape points by, see the by argument in ?retrieveFeatureInfo for possible values. Only used with layout = "point". (Default: NULL)

shape_by

Deprecated. Use shape.by instead.

show.label

Logical scalar, character scalar or integer vector for selecting labels from the rownames of x. If rank is not NULL the rownames might change. (Default: NULL)

label

Deprecated. Use show.label instead.

facet.by

Character scalar. Taxonomic rank to facet the plot by. Value must be of taxonomyRanks(x) Argument can only be used in function plotPrevalentAbundance.

facet_by

Deprecated. Use facet.by instead.

min.prevalence

Numeric scalar. Applied as a threshold for plotting. The threshold is applied per row and column. (Default: 0)

min_prevalence

Deprecated. Use min.prevalence instead.

Value

A ggplot2 object or plotly object, if more than one prevalence was defined.

Details

Whereas plotPrevalence produces a line plot, plotRowPrevalence 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 getPrevalent or agglomerateByPrevalence.

Examples

data(GlobalPatterns, package = "mia")

# Apply relative transformation
GlobalPatterns <- transformAssay(GlobalPatterns, method = "relabundance")

# 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 detection thresholds as heatmap
plotRowPrevalence(GlobalPatterns, rank = "Phylum")


# by default a continuous scale is used for different detection levels,
# but this can be adjusted
plotRowPrevalence(
    GlobalPatterns, rank = "Phylum", assay.type = "relabundance",
    detection = c(0, 0.001, 0.01, 0.1, 0.2))


# point layout for plotRowPrevalence 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()