plotLoadings: Plot feature loadings for TreeSummarizedExperiment objects or...

plotLoadingsR Documentation

Plot feature loadings for TreeSummarizedExperiment objects or feature loadings numeric matrix.

Description

This function is used after performing a reduction method. If TreeSE is given it retrieves the feature loadings matrix to plot values. A tree from rowTree can be added to heatmap layout.

Usage

plotLoadings(x, ...)

## S4 method for signature 'TreeSummarizedExperiment'
plotLoadings(
  x,
  dimred,
  layout = "barplot",
  ncomponents = 5,
  tree.name = "phylo",
  row.var = NULL,
  add.tree = FALSE,
  ...
)

## S4 method for signature 'SingleCellExperiment'
plotLoadings(x, dimred, layout = "barplot", ncomponents = 5, ...)

## S4 method for signature 'matrix'
plotLoadings(x, layout = "barplot", ncomponents = 5, ...)

Arguments

x

a TreeSummarizedExperiment x.

...

additional parameters for plotting.

  • n: Integer scalar. Number of features to be plotted. Applicable when layout="barplot". (Default: 10))

dimred

Character scalar. Determines the reduced dimension to plot.

layout

Character scalar. Determines the layout of plot. Must be either "barplot" or "heatmap". (Default: "barplot")

ncomponents

Numeric scalar. Number of components must be lower or equal to the number of components chosen in the reduction method. (Default: 5)

tree.name

Character scalar. Specifies a rowTree/colTree from x. (Default: tree.name = "phylo")

row.var

NULL or Character scalar. Specifies a variable from rowData to plot with tree heatmap layout. (Default: NULL)

add.tree

Logical scalar. Whether to add tree to heatmap layout. (Default: FALSE)

Details

These method visualize feature loadings of dimension reduction results. Inspired by the plotASVcircular method using phyloseq. TreeSummarizedExperiment object is expected to have content in reducedDim slot calculated with standardized methods from mia or scater package.

Value

A ggplot2 object.

Examples


library(mia)
library(scater)
data("GlobalPatterns", package = "mia")
tse <- GlobalPatterns

# Calculate PCA
tse <- agglomerateByPrevalence(tse, rank="Phylum", update.tree = TRUE)
tse <- transformAssay(tse, method = "clr", pseudocount = 1)
tse <- runPCA(tse, ncomponents = 5, assay.type = "clr")

#' # Plotting feature loadings with tree
plotLoadings(tse, dimred = "PCA", layout = "heatmap", add.tree = TRUE)

# Plotting matrix as a barplot
loadings_matrix <- attr(reducedDim(tse, "PCA"), "rotation")
plotLoadings(loadings_matrix)

# Plotting more features but less components
plotLoadings(tse, dimred = "PCA", ncomponents = 2, n = 12)

# Plotting matrix as heatmap without tree
plotLoadings(loadings_matrix, layout = "heatmap")

# Plot with less components
plotLoadings(tse, "PCA", layout = "heatmap", ncomponents = 2)


microbiome/miaViz documentation built on Aug. 13, 2024, 8:57 p.m.