bdiv_heatmap: Display beta diversities in an all vs all grid.

bdiv_heatmapR Documentation

Display beta diversities in an all vs all grid.

Description

Display beta diversities in an all vs all grid.

Usage

bdiv_heatmap(
  biom,
  bdiv = "Bray-Curtis",
  weighted = TRUE,
  tree = NULL,
  grid = "devon",
  color.by = NULL,
  order.by = NULL,
  limit.by = NULL,
  label = TRUE,
  label_size = NULL,
  rescale = "none",
  clust = "complete",
  trees = TRUE,
  asp = 1,
  tree_height = 10,
  track_height = 10,
  legend = "right",
  title = TRUE,
  xlab.angle = "auto",
  ...
)

Arguments

biom

An rbiom object, such as from as_rbiom(). Any value accepted by as_rbiom() can also be given here.

bdiv

Beta diversity distance algorithm(s) to use. Options are: "Bray-Curtis", "Manhattan", "Euclidean", "Jaccard", and "UniFrac". For "UniFrac", a phylogenetic tree must be present in biom or explicitly provided via ⁠tree=⁠. Default: "Bray-Curtis"

Multiple/abbreviated values allowed.

weighted

Take relative abundances into account. When weighted=FALSE, only presence/absence is considered. Default: TRUE

Multiple values allowed.

tree

A phylo object representing the phylogenetic relationships of the taxa in biom. Only required when computing UniFrac distances. Default: biom$tree

grid

Color palette name, or a list with entries for label, colors, range, bins, na.color, and/or guide. See the Track Definitions section for details. Default: "devon"

color.by

Add annotation tracks for these metadata column(s). See "Annotation Tracks" section below for details. Default: NULL

order.by

Which metadata column(s) to use for ordering the samples across the x and y axes. Overrides any clust argument. See "Ordering and Limiting" section below for details. Default: NULL

limit.by

Metadata definition(s) to use for sample subsetting prior to calculations. See "Ordering and Limiting" section below for details. Default: NULL

label

Label the matrix rows and columns. You can supply a list or logical vector of length two to control row labels and column labels separately, for example label = c(rows = TRUE, cols = FALSE), or simply label = c(T, F). Other valid options are "rows", "cols", "both", "bottom", "right", and "none". Default: TRUE

label_size

The font size to use for the row and column labels. You can supply a numeric vector of length two to control row label sizes and column label sizes separately, for example c(rows = 20, cols = 8), or simply c(20, 8). Default: NULL, which computes: pmax(8, pmin(20, 100 / dim(mtx)))

rescale

Rescale rows or columns to all have a common min/max. Options: "none", "rows", or "cols". Default: "none"

clust

Clustering algorithm for reordering the rows and columns by similarity. You can supply a list or character vector of length two to control the row and column clustering separately, for example clust = c(rows = "complete", cols = NA), or simply clust = c("complete", NA). Options are:

  • FALSE or NA - Disable reordering.

  • An hclust class object E.g. from stats::hclust().

  • A method name - "ward.D", "ward.D2", "single", "complete", "average", "mcquitty", "median", or "centroid".

Default: "complete"

trees

Draw a dendrogram for rows (left) and columns (top). You can supply a list or logical vector of length two to control the row tree and column tree separately, for example trees = c(rows = T, cols = F), or simply trees = c(T, F). Other valid options are "rows", "cols", "both", "left", "top", and "none". Default: TRUE

asp

Aspect ratio (height/width) for entire grid. Default: 1 (square)

tree_height, track_height

The height of the dendrogram or annotation tracks as a percentage of the overall grid size. Use a numeric vector of length two to assign c(top, left) independently. Default: 10 (10% of the grid's height)

legend

Where to place the legend. Options are: "right" or "bottom". Default: "right"

title

Plot title. Set to TRUE for a default title, NULL for no title, or any character string. Default: TRUE

xlab.angle

Angle of the labels at the bottom of the plot. Options are "auto", '0', '30', and '90'. Default: "auto".

...

Additional arguments to pass on to ggplot2::theme(). For example, labs.subtitle = "Plot subtitle".

Value

A ggplot2 plot.
The computed data points and ggplot command are available as ⁠$data⁠ and ⁠$code⁠, respectively.

Annotation Tracks

Metadata can be displayed as colored tracks above the heatmap. Common use cases are provided below, with more thorough documentation available at https://cmmr.github.io/rbiom .

## Categorical ----------------------------
color.by = "Body Site"
color.by = list('Body Site' = "bright")
color.by = list('Body Site' = c("Stool", "Saliva"), 'colors' = "bright")
color.by = list('Body Site' = c('Stool' = "blue", 'Saliva' = "green"))

## Numeric --------------------------------
color.by = "Age"
color.by = list('Age' = "reds")
color.by = list('Age' = c(20,NA), 'colors' = "reds") # at least 20 years old
color.by = list('Age' = c(20,40)) # between 20 and 40 years old (inclusive)

## Multiple Tracks ------------------------
color.by = c("Body Site", "Age")
color.by = list('Body Site' = "bright", 'Age' = "reds")
color.by = list(
  'Body Site' = c('Stool' = "blue", 'Saliva' = "green"),
  'Age'       = list(range = c(20,40), 'colors' = "reds") )

The following entries in the track definitions are understood:

  • colors - A pre-defined palette name or custom set of colors to map to.

  • range - The c(min,max) to use for scale values.

  • label - Label for this track. Defaults to the name of this list element.

  • side - Options are "top" (default) or "left".

  • na.color - The color to use for NA values.

  • bins - Bin a gradient into this many bins/steps.

  • guide - A list of arguments for guide_colorbar() or guide_legend().

All built-in color palettes are colorblind-friendly.

Categorical palette names: "okabe", "carto", "r4", "polychrome", "tol", "bright", "light", "muted", "vibrant", "tableau", "classic", "alphabet", "tableau20", "kelly", and "fishy".

Numeric palette names: "reds", "oranges", "greens", "purples", "grays", "acton", "bamako", "batlow", "bilbao", "buda", "davos", "devon", "grayC", "hawaii", "imola", "lajolla", "lapaz", "nuuk", "oslo", "tokyo", "turku", "bam", "berlin", "broc", "cork", "lisbon", "roma", "tofino", "vanimo", and "vik".

Ordering and Limiting

order.by controls which metadata column(s) are used to arrange samples on the plot. It also enables subsetting to a particular set or range of values. Prefix a column name with - to arrange values in descending order rather than ascending.

## Categorical ----------------------------
order.by = "Body Site"
order.by = list('Body Site' = c("Stool", "Saliva"))

## Numeric --------------------------------
order.by = "-Age"
order.by = list('Age'  = c(20,NA)) # at least 20 years old
order.by = list('-Age' = c(20,40)) # between 20 and 40 years old (inclusive)

## Multiple / Mixed -----------------------
order.by = c("-Body Site", "Age")
order.by = list("Body Site", '-Age' = c(20,40))

limit.by is used to specify a subset of samples without any side-effects on aesthetics. It is especially useful for limiting the data to a single categorical metadata value. Unlike the other *.by parameters, limit.by must always be a named list().

## Categorical ----------------------------
limit.by = list('Sex' = "Male")

## Numeric --------------------------------
limit.by = list('Age' = c(20,NA)) # at least 20 years old
limit.by = list('Age' = c(20,40)) # between 20 and 40 years old (inclusive)

## Multiple / Mixed -----------------------
limit.by = list(
  'Sex'       = "Male", 
  'Body Site' = c("Stool", "Saliva")
  'Age'       = c(20,40) )

See Also

Other beta_diversity: bdiv_boxplot(), bdiv_corrplot(), bdiv_ord_plot(), bdiv_ord_table(), bdiv_stats(), bdiv_table(), distmat_stats()

Other visualization: adiv_boxplot(), adiv_corrplot(), bdiv_boxplot(), bdiv_corrplot(), bdiv_ord_plot(), plot_heatmap(), rare_corrplot(), rare_multiplot(), rare_stacked(), stats_boxplot(), stats_corrplot(), taxa_boxplot(), taxa_corrplot(), taxa_heatmap(), taxa_stacked()

Examples

    library(rbiom)
    
    # Keep and rarefy the 10 most deeply sequenced samples.
    hmp10 <- rarefy(hmp50, n = 10)
    
    bdiv_heatmap(hmp10, color.by=c("Body Site", "Age"))
    
    bdiv_heatmap(hmp10, bdiv="uni", weighted=c(T,F), color.by="sex")
    

cmmr/rbiom documentation built on April 28, 2024, 6:38 a.m.