plotCCA | R Documentation |
plotRDA
and plotCCA
create an RDA/CCA plot starting from the
output of CCA and RDA
functions, two common methods
for supervised ordination of microbiome data.
plotCCA(x, ...)
## S4 method for signature 'SingleCellExperiment'
plotCCA(x, dimred, ...)
## S4 method for signature 'matrix'
plotCCA(x, ...)
plotRDA(x, ...)
## S4 method for signature 'SingleCellExperiment'
plotRDA(x, dimred, ...)
## S4 method for signature 'matrix'
plotRDA(x, ...)
x |
a
|
... |
additional parameters for plotting, inherited from
|
dimred |
|
plotRDA
and plotCCA
create an RDA/CCA plot starting from the
output of CCA and RDA
functions, two common methods
for supervised ordination of microbiome data. Either a
TreeSummarizedExperiment
or a matrix object is supported as input. When the input is a
TreeSummarizedExperiment, this should contain the output of addRDA
in the reducedDim slot and the argument dimred
needs to be defined.
When the input is a matrix, this should be returned as output from
getRDA. However, the first method is recommended because it provides
the option to adjust aesthetics to the colData variables through the
arguments inherited from plotReducedDim
.
A ggplot2
object
# Load dataset
library(miaViz)
data("enterotype", package = "mia")
tse <- enterotype
# Run RDA and store results into TreeSE
tse <- addRDA(
tse,
formula = assay ~ ClinicalStatus + Gender + Age,
FUN = getDissimilarity,
distance = "bray",
na.action = na.exclude
)
# Create RDA plot coloured by variable
plotRDA(tse, "RDA", colour.by = "ClinicalStatus")
# Create RDA plot with empty ellipses
plotRDA(tse, "RDA", colour.by = "ClinicalStatus", add.ellipse = "colour")
# Create RDA plot with text encased in labels
plotRDA(tse, "RDA", colour.by = "ClinicalStatus", vec.text = FALSE)
# Create RDA plot without repelling text
plotRDA(tse, "RDA", colour.by = "ClinicalStatus", repel.labels = FALSE)
# Create RDA plot without vectors
plotRDA(tse, "RDA", colour.by = "ClinicalStatus", add.vectors = FALSE)
# Calculate RDA as a separate object
rda_mat <- getRDA(
tse,
formula = assay ~ ClinicalStatus + Gender + Age,
FUN = getDissimilarity,
distance = "bray",
na.action = na.exclude
)
# Create RDA plot from RDA matrix
plotRDA(rda_mat)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.