getAbundant: Determine abundant and rare taxa. Rare taxa can be further...

getAbundantR Documentation

Determine abundant and rare taxa. Rare taxa can be further classified to conditionally rare and permanently rare.

Description

These functions determine abundant and rare taxa based on the abundances of taxa. Compared to getPrevalent and getRare, these functions determine abundant and rare taxa based on abundance while the first mentioned are based on prevalence.

Usage

getAbundant(x, ...)

getLowAbundant(x, ...)

getConditionallyLowAbundant(x, ...)

getPermanentlyLowAbundant(x, ...)

getAbundanceClass(x, ...)

addAbundanceClass(x, ...)

## S4 method for signature 'SingleCellExperiment'
getAbundant(x, ...)

## S4 method for signature 'SummarizedExperiment'
getAbundant(x, assay.type = "relabundance", ...)

## S4 method for signature 'ANY'
getAbundant(x, abundant.th = 1/100, ...)

## S4 method for signature 'SingleCellExperiment'
getLowAbundant(x, ...)

## S4 method for signature 'SummarizedExperiment'
getLowAbundant(x, assay.type = "relabundance", abundant.th = 1/100, ...)

## S4 method for signature 'ANY'
getLowAbundant(x, abundant.th = 1/100, ...)

## S4 method for signature 'SingleCellExperiment'
getConditionallyLowAbundant(x, ...)

## S4 method for signature 'SummarizedExperiment'
getConditionallyLowAbundant(x, assay.type = "relabundance", ...)

## S4 method for signature 'ANY'
getConditionallyLowAbundant(x, abundant.th = 1/100, crt.th = 100, ...)

## S4 method for signature 'SingleCellExperiment'
getPermanentlyLowAbundant(x, ...)

## S4 method for signature 'SummarizedExperiment'
getPermanentlyLowAbundant(x, assay.type = "relabundance", ...)

## S4 method for signature 'ANY'
getPermanentlyLowAbundant(x, abundant.th = 1/100, prt.th = 5, ...)

## S4 method for signature 'SingleCellExperiment'
getAbundanceClass(x, ...)

## S4 method for signature 'SummarizedExperiment'
getAbundanceClass(x, assay.type = "relabundance", ...)

## S4 method for signature 'ANY'
getAbundanceClass(x, abundant.th = 1/100, crt.th = 100, prt.th = 5, ...)

## S4 method for signature 'SingleCellExperiment'
addAbundanceClass(x, ...)

## S4 method for signature 'SummarizedExperiment'
addAbundanceClass(x, name = "abundance_class", ...)

Arguments

x

a SummarizedExperiment object.

...

additional arguments.

assay.type

Character scalar. Specifies the name of assay used in calculation. (Default: "relabundance")

abundant.th

Numeric scalar. Specifies threshold that is used to separate abundant features from rare. (Default: 1/100)

crt.th

Numeric scalar. Specifies threshold that is used to separate conditionally rare features from other rare features. (Default: 100)

prt.th

Numeric scalar. Specifies threshold that is used to separate permanently rare features from other rare features. (Default: 5)

name

Character scalar. Specifies name of column in rowData where the results will be stored. (Default: "abundance_class")

Details

These functions identify abundant and rare taxa in a dataset. Abundant taxa are characterized by high average abundance across the dataset, while rare taxa are characterized by consistently low abundance.

Conditionally rare taxa exhibit variable abundance, being abundant in some samples and rare in others. In contrast, permanently rare taxa consistently maintain low abundance across all samples.

  • Abundant taxa: Taxa with an average abundance exceeding abundant.th.

  • Low abundant / rare taxa: Taxa with an average abundance not exceeding abundant.th. Optionally, if specified, they must also satisfy the condition crt.th >= \frac{abundance_{max}}{abundance_{min}} > prt.th.

  • Conditionally rare or low abundant taxa (CRT): Taxa with an average abundance not exceeding abundant.th and with a maximum-to-minimum abundance ratio (\frac{abundance_{max}}{abundance_{min}}) greater than crt.th.

  • Permanently rare or low abundant taxa (PRT): Taxa with an average abundance not exceeding abundant.th and with a maximum-to-minimum abundance ratio (\frac{abundance_{max}}{abundance_{min}}) less than or equal to prt.th.

Value

For getAbundant, getLowAbundant, getConditionallyLowAbundant, and getPermanentlyLowAbundant a vector of taxa. For getAbudanceClass a vector of abundance classes for each feature. For addAbudanceClass, a SummarizedExperiment object.

References

Sizhong Y. et al. (2017) Community structure of rare methanogenic archaea: insight from a single functional group- FEMS Microbiol. Ecol. 93(11). https://doi.org/10.1093/femsec/fix126

See Also

getPrevalent and getRare

Examples


data(GlobalPatterns)
tse <- GlobalPatterns

# Agglomerate to family level
tse <- agglomerateByRank(tse, rank = "Family")
# Transform to relative abundances. Note that we add pseudocount. This is
# because otherwise we cannot calculate CRT and PRT due to zeroes and
# zero division in calculating abundance ratio.
tse <- transformAssay(tse, method = "relabundance", pseudocount = TRUE)

# Get abundant taxa
abundant <- getAbundant(tse, assay.type = "relabundance")
abundant |> head()

# Get all rare taxa that have average relative abundance below 10%
rare <- getLowAbundant(
    tse, assay.type = "relabundance", abundant.th = 10/100)
rare |> head()

# Get rare taxa that are not permanently or conditionally rare
rare <- getLowAbundant(
    tse, assay.type = "relabundance", prt.th = 5, crt.th = 100)
rare |> head()

# Get permanently rare taxa
prt <- getPermanentlyLowAbundant(
    tse, assay.type = "relabundance", prt.th = 5)
prt |> head()

# Get conditionally rare taxa
prt <- getConditionallyLowAbundant(
    tse, assay.type = "relabundance", crt.th = 100)
prt |> head()

# To classify all features, one can use *AbundantClass function
tse <- addAbundanceClass(tse)
# When one uses add* function, the results are stored to rowData
rowData(tse)


FelixErnst/mia documentation built on Jan. 19, 2025, 2:30 a.m.