plotRLE: Plot relative log expression

View source: R/plotRLE.R

plotRLER Documentation

Plot relative log expression

Description

Produce a relative log expression (RLE) plot of one or more transformations of cell expression values.

Usage

plotRLE(
  object,
  exprs_values = "logcounts",
  exprs_logged = TRUE,
  style = "minimal",
  legend = TRUE,
  ordering = NULL,
  colour_by = color_by,
  by_exprs_values = exprs_values,
  BPPARAM = BiocParallel::bpparam(),
  color_by = NULL,
  assay.type = exprs_values,
  by.assay.type = by_exprs_values,
  assay_logged = exprs_logged,
  ...
)

Arguments

object

A SingleCellExperiment object.

exprs_values

Alias to assay.type.

exprs_logged

A logical scalar indicating whether the expression matrix is already log-transformed. If not, a log2-transformation (+1) will be performed prior to plotting.

style

String defining the boxplot style to use, either "minimal" (default) or "full"; see Details.

legend

Logical scalar specifying whether a legend should be shown.

ordering

A vector specifying the ordering of cells in the RLE plot. This can be useful for arranging cells by experimental conditions or batches.

colour_by

Specification of a column metadata field or a feature to colour by, see the by argument in ?retrieveCellInfo for possible values.

by_exprs_values

Alias to by.assay.type.

BPPARAM

A BiocParallelParam object to be used to parallelise operations using DelayedArray.

color_by

Alias to colour_by.

assay.type

A string or integer scalar specifying the expression matrix in object to use.

by.assay.type

A string or integer scalar specifying which assay to obtain expression values from, for use in point aesthetics - see the assay.type argument in ?retrieveCellInfo.

assay_logged

Alias to exprs_logged.

...

further arguments passed to geom_boxplot when style="full".

Details

Relative log expression (RLE) plots are a powerful tool for visualising unwanted variation in high dimensional data. These plots were originally devised for gene expression data from microarrays but can also be used on single-cell expression data. RLE plots are particularly useful for assessing whether a procedure aimed at removing unwanted variation (e.g., scaling normalisation) has been successful.

If style is “full”, the usual ggplot2 boxplot is created for each cell. Here, the box shows the inter-quartile range and whiskers extend no more than 1.5 times the IQR from the hinge (the 25th or 75th percentile). Data beyond the whiskers are called outliers and are plotted individually. The median (50th percentile) is shown with a white bar. This approach is detailed and flexible, but can take a long time to plot for large datasets.

If style is “minimal”, a Tufte-style boxplot is created for each cell. Here, the median is shown with a circle, the IQR in a grey line, and “whiskers” (as defined above) for the plots are shown with coloured lines. No outliers are shown for this plot style. This approach is more succinct and faster for large numbers of cells.

Value

A ggplot object

Author(s)

Davis McCarthy, with modifications by Aaron Lun

References

Gandolfo LC, Speed TP (2017). RLE plots: visualising unwanted variation in high dimensional data. arXiv.

Examples

example_sce <- mockSCE()
example_sce <- logNormCounts(example_sce)

plotRLE(example_sce, colour_by = "Mutation_Status", style = "minimal")

plotRLE(example_sce, colour_by = "Mutation_Status", style = "full",
       outlier.alpha = 0.1, outlier.shape = 3, outlier.size = 0)


davismcc/scater documentation built on July 19, 2024, 2:01 a.m.