This function calculates coefficient of bimodality for each taxa.

getBimodality(x, ...)

addBimodality(x, ...)

# S4 method for class 'SummarizedExperiment'
addBimodality(x, name = "bimodality", ...)

# S4 method for class 'SummarizedExperiment'
getBimodality(x, assay.type = "counts", ...)

Arguments

x

A SummarizedExperiment object.

...

additional arguments.

  • group: Character scalar. Specifies a name of the column from colData that identifies the grouping of the samples. (Default: NULL)

name

Character scalar. Specifies a column name for storing bimodality results. (Default: "bimodality")

assay.type

Character scalar. Specifies which assay values are used in the dissimilarity estimation. (Default: "counts")

Value

DataFrame or x with results added to its rowData.

Details

This function calculates coefficient of bimodality for each taxa. If the dataset includes grouping, for instance, individual systems or patients, the coefficient is calculated for each group separately.

The coefficient of bimodality measures whether a taxon has bimodal abundance. For instance, certain taxon can be high-abundant in some timepoints while in others it might be rare. The coefficient can help to determine these taxa.

The coefficient of bimodality (b) is defined as follows:

$$b = \frac{1+skewness^{2}}{kurtosis+3}$$

where skewness is calculated as follows

$$skewness = \frac{\sum_{i=1}^{n}(x_{i}-\overline{x})^{3}/n}{ \sum_{i=1}^{n}(x_{i}-\overline{x})^{2}/n]^{3/2}}$$

and kurtosis as follows

$$kurtosis = \frac{\sum_{i=1}^{n}(x_{i}-\overline{x})^{4}/n}{ (\sum_{i=1}^{n}(x_{i}-\overline{x})^{2}/n)^{2}}$$

The coefficient ranges from 0-1, where 1 means high bimodality. The coefficient was introduced in the paper Shade A et al. 2014, where they used bimodality and abundance to determine conditionally rare taxa (CRT).

References

Shade A, et al. (2014) Conditionally Rare Taxa Disproportionately Contribute to Temporal Changes in Microbial Diversity. doi: 10.1128/mbio.01371-14

Examples

library(miaTime)

data(SilvermanAGutData)
tse <- SilvermanAGutData

# In this example, we are only interested on vessel 1.
tse <- tse[, tse[["Vessel"]] == 1]
tse <- transformAssay(tse, method = "relabundance")

# Calculate bimodality
b <- getBimodality(tse)
# Determine taxa with high bimodality
bimodal_taxa <- names(b)[ which(b[[1]] > 0.95) ]

if (FALSE) { # \dontrun{
# Determine taxa with abundance > 0.5%
abundant <- getAbundant(
    tse, assay.type = "relabundance", abundant.th = 0.5/100)

# The detected CRT
crt <- intersect(bimodal_taxa, abundant)
head(crt)
} # }