runCCA: Canonical Correspondence Analysis and Redundancy Analysis

View source: R/runCCA.R

runCCAR Documentation

Canonical Correspondence Analysis and Redundancy Analysis

Description

These functions perform Canonical Correspondence Analysis on data stored in a SummarizedExperiment.

Usage

getCCA(x, ...)

addCCA(x, ...)

getRDA(x, ...)

addRDA(x, ...)

## S4 method for signature 'ANY'
getCCA(x, ...)

## S4 method for signature 'SummarizedExperiment'
getCCA(
  x,
  formula,
  variables,
  test.signif = TRUE,
  assay.type = assay_name,
  assay_name = exprs_values,
  exprs_values = "counts",
  scores = "wa",
  ...
)

calculateCCA(x, ...)

## S4 method for signature 'SingleCellExperiment'
addCCA(x, formula, variables, altexp = NULL, name = "CCA", ...)

runCCA(x, ...)

## S4 method for signature 'ANY'
getRDA(x, ...)

## S4 method for signature 'SummarizedExperiment'
getRDA(
  x,
  formula,
  variables,
  test.signif = TRUE,
  assay.type = assay_name,
  assay_name = exprs_values,
  exprs_values = "counts",
  scores = "wa",
  ...
)

calculateRDA(x, ...)

## S4 method for signature 'SingleCellExperiment'
addRDA(x, formula, variables, altexp = NULL, name = "RDA", ...)

runRDA(x, ...)

Arguments

x

TreeSummarizedExperiment.

...

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 variables.

variables and formula can be missing, which turns the CCA analysis into a CA analysis and dbRDA into PCoA/MDS.

variables

Character scalar. When x is a SummarizedExperiment, variables can be used to specify variables from colData.

When x is a matrix, variables is a data.frame or an object coercible to one containing the variables to use.

All variables are used. Please subset, if you want to consider only some of them. variables and formula can be missing, which turns the CCA analysis into a CA analysis and dbRDA into PCoA/MDS.

test.signif

Logical scalar. Should the PERMANOVA and analysis of multivariate homogeneity of group dispersions be performed. (Default: TRUE)

assay.type

Character scalar. Specifies which assay to use for calculation. (Default: "counts")

assay_name

Deprecated. Use assay.type instead.

exprs_values

Deprecated. Use assay.type instead.

scores

Character scalar. Specifies scores to be returned. Must be 'wa' (site scores found as weighted averages (cca) or weighted sums (rda) of v with weights Xbar, but the multiplying effect of eigenvalues removed) or 'u' ((weighted) orthonormal site scores). (Default: 'wa')

altexp

Character scalar or integer scalar. Specifies an alternative experiment containing the input data.

name

Character scalar. A name for the column of the colData where results will be stored. (Default: "CCA")

Details

For run* a SingleCellExperiment or a derived object.

*CCA functions utilize vegan:cca and *RDA functions vegan:dbRDA. By default, dbRDA is done with euclidean distances, which is equivalent to RDA.

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.

Value

For getCCA a matrix with samples as rows and CCA dimensions as columns. Attributes include calculated cca/rda object and significance analysis results.

For addCCA a modified x with the results stored in reducedDim as the given name.

See Also

For more details on the actual implementation see cca and dbrda

Examples

library(miaViz)
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")

# Fetch significance results
attr(reducedDim(tse, "CCA"), "significance")

tse <- transformAssay(tse, method = "relabundance")

# Specify dissimilarity measure
tse <- addRDA(tse,
              formula = data ~ ClinicalStatus,
              assay.type = "relabundance",
              method = "bray",
              name = "RDA_bray")

# 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)

# Plot RDA
plotRDA(tse, "rda_scaled",
        colour_by = "ClinicalStatus")

# 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)


microbiome/mia documentation built on Sept. 19, 2024, 11:17 p.m.