transformAssay: Transform assay

transformAssayR Documentation

Transform assay

Description

Variety of transformations for abundance data, stored in assay. See details for options.

Usage

transformAssay(
  x,
  assay.type = "counts",
  assay_name = NULL,
  method = c("alr", "chi.square", "clr", "frequency", "hellinger", "log", "log10",
    "log2", "max", "normalize", "pa", "range", "rank", "rclr", "relabundance", "rrank",
    "standardize", "total", "z"),
  MARGIN = "samples",
  name = method,
  pseudocount = FALSE,
  ...
)

## S4 method for signature 'SummarizedExperiment'
transformAssay(
  x,
  assay.type = "counts",
  assay_name = NULL,
  method = c("alr", "chi.square", "clr", "frequency", "hellinger", "log", "log10",
    "log2", "max", "normalize", "pa", "range", "rank", "rclr", "relabundance", "rrank",
    "standardize", "total", "z"),
  MARGIN = "samples",
  name = method,
  pseudocount = FALSE,
  ...
)

Arguments

x

TreeSummarizedExperiment.

assay.type

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

assay_name

Deprecated. Use assay.type instead.

method

Character scalar. Specifies the transformation method.

MARGIN

Character scalar. Determines whether the transformation is applied sample (column) or feature (row) wise. (Default: "samples")

name

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

pseudocount

Logical scalar or numeric scalar. When TRUE, automatically adds half of the minimum positive value of assay.type (missing values ignored by default: na.rm = TRUE). When FALSE, does not add any pseudocount (pseudocount = 0). Alternatively, a user-specified numeric value can be added as pseudocount. (Default: FALSE).

...

additional arguments passed on to vegan:decostand:

  • reference: Character scalar. use to to fill reference sample's column in returned assay when calculating alr. (Default: NA)

  • ref_vals Deprecated. Use reference instead.

Details

transformAssay function provides a variety of options for transforming abundance data. The transformed data is calculated and stored in a new assay.

The transformAssay provides sample-wise (column-wise) or feature-wise (row-wise) transformation to the abundance table (assay) based on specified MARGIN.

The available transformation methods include:

  • 'alr', 'chi.square', 'clr', 'frequency', 'hellinger', 'log', 'normalize', 'pa', 'rank', 'rclr' 'relabundance', 'rrank', 'standardize', 'total': please refer to decostand for details.

  • 'log10': log10 transformation can be used for reducing the skewness of the data.

    log10 = \log_{10} x

    where x is a single value of data.

  • 'log2': log2 transformation can be used for reducing the skewness of the data.

    log2 = \log_{2} x

    where x is a single value of data.

Value

transformAssay returns the input object x, with a new transformed abundance table named name added in the assay.

Author(s)

Leo Lahti and Tuomas Borman. Contact: microbiome.github.io

Examples

data(esophagus, package="mia")
tse <- esophagus

# By specifying 'method', it is possible to apply different transformations, 
# e.g. compositional transformation.
tse <- transformAssay(tse, method = "relabundance")

# The target of transformation can be specified with "assay.type"
# Pseudocount can be added by specifying 'pseudocount'.

# Perform CLR with smallest positive value as pseudocount
tse <- transformAssay(tse, assay.type = "relabundance", method = "clr", 
                     pseudocount = TRUE
                     )
                      
head(assay(tse, "clr"))

# With MARGIN, you can specify the if transformation is done for samples or
# for features. Here Z-transformation is done feature-wise.
tse <- transformAssay(tse, method = "standardize", MARGIN = "features")
head(assay(tse, "standardize"))

# Name of the stored table can be specified.
tse <- transformAssay(tse, method="hellinger", name="test")
head(assay(tse, "test"))

# pa returns presence absence table.
tse <- transformAssay(tse, method = "pa")
head(assay(tse, "pa"))

# rank returns ranks of taxa.
tse <- transformAssay(tse, method = "rank")
head(assay(tse, "rank"))

# In order to use other ranking variants, modify the chosen assay directly:
assay(tse, "rank_average", withDimnames = FALSE) <- colRanks(
    assay(tse, "counts"), ties.method = "average", preserveShape = TRUE)  


microbiome/mia documentation built on Aug. 14, 2024, 4:42 p.m.