Display beta diversities in an all vs all grid.
bdiv_heatmap(
biom,
bdiv = "Bray-Curtis",
weighted = TRUE,
tree = NULL,
tracks = NULL,
grid = "devon",
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",
underscores = FALSE,
...
)
A ggplot2
plot. The computed data points and ggplot
command are available as $data
and $code
,
respectively.
An rbiom object, such as from as_rbiom()
.
Any value accepted by as_rbiom()
can also be given here.
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=
. Multiple/abbreviated values allowed.
Default: "Bray-Curtis"
Take relative abundances into account. When
weighted=FALSE
, only presence/absence is considered.
Multiple values allowed. Default: TRUE
A phylo
object representing the phylogenetic
relationships of the taxa in biom
. Only required when
computing UniFrac distances. Default: biom$tree
A character vector of metadata fields to display as tracks
at the top of the plot. Or, a list as expected by the tracks
argument of plot_heatmap()
. Default: NULL
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"
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(TRUE, FALSE)
. Other valid options are "rows"
,
"cols"
, "both"
, "bottom"
, "right"
,
and "none"
.
Default: TRUE
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 rows or columns to all have a common min/max.
Options: "none"
, "rows"
, or "cols"
.
Default: "none"
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.
hclust
class objectE.g. from stats::hclust()
.
"ward.D"
,
"ward.D2"
, "single"
, "complete"
,
"average"
, "mcquitty"
, "median"
, or
"centroid"
.
Default: "complete"
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 = TRUE, cols = FALSE)
,
or simply trees = c(TRUE, FALSE)
.
Other valid options are "rows"
, "cols"
, "both"
,
"left"
, "top"
, and "none"
.
Default: TRUE
Aspect ratio (height/width) for entire grid.
Default: 1
(square)
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)
Where to place the legend. Options are: "right"
or
"bottom"
. Default: "right"
Plot title. Set to TRUE
for a default title, NULL
for
no title, or any character string. Default: TRUE
Angle of the labels at the bottom of the plot.
Options are "auto"
, '0'
, '30'
, and '90'
.
Default: "auto"
.
When parsing the tree, should underscores be kept as
is? By default they will be converted to spaces (unless the entire ID
is quoted). Default FALSE
Additional arguments to pass on to ggplot2::theme().
For example, labs.subtitle = "Plot subtitle"
.
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 ----------------------------
tracks = "Body Site"
tracks = list('Body Site' = "bright")
tracks = list('Body Site' = c('Stool' = "blue", 'Saliva' = "green"))## Numeric --------------------------------
tracks = "Age"
tracks = list('Age' = "reds")
## Multiple Tracks ------------------------
tracks = c("Body Site", "Age")
tracks = list('Body Site' = "bright", 'Age' = "reds")
tracks = list(
'Body Site' = c('Stool' = "blue", 'Saliva' = "green"),
'Age' = list('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"
.
Other beta_diversity:
bdiv_boxplot()
,
bdiv_clusters()
,
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()
library(rbiom)
# Keep and rarefy the 10 most deeply sequenced samples.
hmp10 <- rarefy(hmp50, n = 10)
bdiv_heatmap(hmp10, tracks=c("Body Site", "Age"))
bdiv_heatmap(hmp10, bdiv="uni", weighted=c(TRUE,FALSE), tracks="sex")
Run the code above in your browser using DataLab