These functions perform Canonical Correspondence Analysis on data stored
in a SummarizedExperiment
.
getCCA(x, ...)
addCCA(x, ...)
getRDA(x, ...)
addRDA(x, ...)
calculateCCA(x, ...)
runCCA(x, ...)
# S4 method for class 'ANY'
getCCA(x, formula, data, ...)
# S4 method for class 'SummarizedExperiment'
getCCA(
x,
formula = NULL,
col.var = variables,
variables = NULL,
test.signif = TRUE,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
...
)
# S4 method for class 'SingleCellExperiment'
addCCA(x, altexp = NULL, name = "CCA", ...)
calculateRDA(x, ...)
runRDA(x, ...)
# S4 method for class 'ANY'
getRDA(x, formula, data, ...)
# S4 method for class 'SummarizedExperiment'
getRDA(
x,
formula = NULL,
col.var = variables,
variables = NULL,
test.signif = TRUE,
assay.type = assay_name,
assay_name = exprs_values,
exprs_values = "counts",
...
)
# S4 method for class 'SingleCellExperiment'
addRDA(x, altexp = NULL, name = "RDA", ...)
additional arguments passed to vegan::cca or vegan::dbrda and other internal functions.
method
a dissimilarity measure to be applied in dbRDA and
possible following homogeneity test. (By default:
method="euclidean"
)
scale
: Logical scalar
. Should the expression values be
standardized? scale
is disabled when using *RDA
functions.
Please scale before performing RDA. (Default: TRUE
)
na.action
: function
. Action to take when missing
values for any of the variables in formula
are encountered.
(Default: na.fail
)
full
Logical scalar
. should all the results from the
significance calculations be returned. When full=FALSE
, only
summary tables are returned. (Default: FALSE
)
homogeneity.test
: Character scalar
. Specifies
the significance test used to analyse vegan::betadisper
results.
Options include 'permanova' (vegan::permutest
), 'anova'
(stats::anova
) and 'tukeyhsd' (stats::TukeyHSD
).
(By default: homogeneity.test="permanova"
)
permutations
a numeric value specifying the number of
permutations for significance testing in vegan::anova.cca
.
(By default: permutations=999
)
formula
. If x
is a
SummarizedExperiment
a formula can be supplied. Based on the right-hand side of the given formula
colData
is subset to col.var
.
col.var
and formula
can be missing, which turns the CCA
analysis into a CA analysis and dbRDA into PCoA/MDS.
data.frame
or coarcible to one. The covariance table
including covariates defined by formula
.
Character scalar
. When x
is a
SummarizedExperiment
,col.var
can be used to specify variables
from colData
.
Deprecated. Use "col.var"
instead.
Logical scalar
. Should the PERMANOVA and analysis
of multivariate homogeneity of group dispersions be performed.
(Default: TRUE
)
Character scalar
. Specifies which assay to use for
calculation. (Default: "counts"
)
Deprecated. Use assay.type
instead.
Deprecated. Use assay.type
instead.
Character scalar
or integer scalar
. Specifies an
alternative experiment containing the input data.
Character scalar
. A name for the column of the
colData
where results will be stored. (Default: "CCA"
)
For getCCA
a matrix with samples as rows and CCA dimensions as
columns. Attributes include output from scores
,
eigenvalues, the cca
/rda
object and significance analysis
results.
For addCCA
a modified x
with the results stored in
reducedDim
as the given name
.
*CCA functions utilize vegan:cca
and *RDA functions
vegan:dbRDA
. By default, dbRDA is done with euclidean distances, which
is equivalent to RDA. col.var
and formula
can be missing,
which turns the CCA analysis into a CA analysis and dbRDA into PCoA/MDS.
Significance tests are done with vegan:anova.cca
(PERMANOVA). Group
dispersion, i.e., homogeneity within groups is analyzed with
vegan:betadisper
(multivariate homogeneity of groups dispersions
(variances)) and statistical significance of homogeneity is tested with a
test specified by homogeneity.test
parameter.
library(miaViz)
#> Loading required package: ggraph
#>
#> Attaching package: ‘miaViz’
#> The following object is masked from ‘package:mia’:
#>
#> plotNMDS
data("enterotype", package = "mia")
tse <- enterotype
# Perform CCA and exclude any sample with missing ClinicalStatus
tse <- addCCA(
tse,
formula = data ~ ClinicalStatus,
na.action = na.exclude
)
# Plot CCA
plotCCA(tse, "CCA", colour_by = "ClinicalStatus")
#> Too few points to calculate an ellipse
#> Too few points to calculate an ellipse
# Fetch significance results
attr(reducedDim(tse, "CCA"), "significance")
#> $permanova
#> Df ChiSquare F Pr(>F) Total variance Explained variance
#> Model 4 0.150386 0.7379888 0.676 1.8825 0.07988631
#> ClinicalStatus 4 0.150386 0.7379888 0.651 1.8825 0.07988631
#> Residual 34 1.732114 NA NA 1.8825 0.92011369
#>
#> $homogeneity
#> Df Sum Sq Mean Sq F N.Perm Pr(>F) Total variance
#> ClinicalStatus 4 0.2184213 0.05460533 1.757098 999 0.152 1.275039
#> Explained variance
#> ClinicalStatus 0.1713056
#>
tse <- transformAssay(tse, method = "relabundance")
# Specify dissimilarity measure
tse <- addRDA(
tse,
formula = data ~ ClinicalStatus,
assay.type = "relabundance",
method = "bray",
name = "RDA_bray",
na.action = na.exclude
)
# To scale values when using *RDA functions, use
# transformAssay(MARGIN = "features", ...)
tse <- transformAssay(tse, method = "standardize", MARGIN = "features")
# Data might include taxa that do not vary. Remove those because after
# z-transform their value is NA
tse <- tse[rowSums(is.na(assay(tse, "standardize"))) == 0, ]
# Calculate RDA
tse <- addRDA(
tse,
formula = data ~ ClinicalStatus,
assay.type = "standardize",
name = "rda_scaled",
na.action = na.omit
)
#> Warning: Certain samples are removed from the result because they did not include sufficient metadata.
# Plot RDA
plotRDA(tse, "rda_scaled", colour_by = "ClinicalStatus")
#> Too few points to calculate an ellipse
#> Too few points to calculate an ellipse
# A common choice along with PERMANOVA is ANOVA when statistical significance
# of homogeneity of groups is analysed. Moreover, full significance test
# results can be returned.
tse <- addRDA(
tse,
formula = data ~ ClinicalStatus,
homogeneity.test = "anova",
full = TRUE
)
# Example showing how to pass extra parameters, such as 'permutations',
# to anova.cca
tse <- addRDA(
tse,
formula = data ~ ClinicalStatus,
permutations = 500
)