estimateEvenness: Estimate Evenness measures

estimateEvennessR Documentation

Estimate Evenness measures

Description

This function calculates community evenness indices. These include the ‘Camargo’, ‘Pielou’, ‘Simpson’, ‘Evar’ and ‘Bulla’ evenness measures. See details for more information and references.

Usage

estimateEvenness(
  x,
  assay.type = assay_name,
  assay_name = "counts",
  index = c("pielou", "camargo", "simpson_evenness", "evar", "bulla"),
  name = index,
  ...
)

## S4 method for signature 'SummarizedExperiment'
estimateEvenness(
  x,
  assay.type = assay_name,
  assay_name = "counts",
  index = c("camargo", "pielou", "simpson_evenness", "evar", "bulla"),
  name = index,
  ...,
  BPPARAM = SerialParam()
)

Arguments

x

a SummarizedExperiment object

assay.type

A single character value for selecting the assay used for calculation of the sample-wise estimates.

assay_name

a single character value for specifying which assay to use for calculation. (Please use assay.type instead. At some point assay_name will be disabled.)

index

a character vector, specifying the evenness measures to be calculated.

name

a name for the column(s) of the colData the results should be stored in.

...

optional arguments:

  • threshold a numeric threshold. assay values below or equal to this threshold will be set to zero.

BPPARAM

A BiocParallelParam object specifying whether calculation of estimates should be parallelized.

Details

Evenness is a standard index in community ecology, and it quantifies how evenly the abundances of different species are distributed. The following evenness indices are provided:

By default, this function returns all indices.

The available evenness indices include the following (all in lowercase):

  • 'camargo' Camargo's evenness (Camargo 1992)

  • 'simpson_evenness' Simpson’s evenness is calculated as inverse Simpson diversity (1/lambda) divided by observed species richness S: (1/lambda)/S.

  • 'pielou' Pielou's evenness (Pielou, 1966), also known as Shannon or Shannon-Weaver/Wiener/Weiner evenness; H/ln(S). The Shannon-Weaver is the preferred term; see Spellerberg and Fedor (2003).

  • 'evar' Smith and Wilson’s Evar index (Smith & Wilson 1996).

  • 'bulla' Bulla’s index (O) (Bulla 1994).

Desirable statistical evenness metrics avoid strong bias towards very large or very small abundances; are independent of richness; and range within the unit interval with increasing evenness (Smith & Wilson 1996). Evenness metrics that fulfill these criteria include at least camargo, simpson, smith-wilson, and bulla. Also see Magurran & McGill (2011) and Beisel et al. (2003) for further details.

Value

x with additional colData named *name*

References

Beisel J-N. et al. (2003) A Comparative Analysis of Evenness Index Sensitivity. Internal Rev. Hydrobiol. 88(1):3-15. URL: https://portais.ufg.br/up/202/o/2003-comparative_evennes_index.pdf

Bulla L. (1994) An index of evenness and its associated diversity measure. Oikos 70:167–171.

Camargo, JA. (1992) New diversity index for assessing structural alterations in aquatic communities. Bull. Environ. Contam. Toxicol. 48:428–434.

Locey KJ and Lennon JT. (2016) Scaling laws predict global microbial diversity. PNAS 113(21):5970-5975; doi:10.1073/pnas.1521291113.

Magurran AE, McGill BJ, eds (2011) Biological Diversity: Frontiers in Measurement and Assessment (Oxford Univ Press, Oxford), Vol 12.

Pielou, EC. (1966) The measurement of diversity in different types of biological collections. J Theoretical Biology 13:131–144.

Smith B and Wilson JB. (1996) A Consumer's Guide to Evenness Indices. Oikos 76(1):70-82.

Spellerberg and Fedor (2003). A tribute to Claude Shannon (1916 –2001) and a plea for more rigorous use of species richness, species diversity and the ‘Shannon–Wiener’ Index. Alpha Ecology & Biogeography 12, 177–197.

See Also

plotColData

  • estimateRichness

  • estimateDominance

  • estimateDiversity

Examples

data(esophagus)
tse <- esophagus

# Specify index and their output names
index <- c("pielou", "camargo", "simpson_evenness", "evar", "bulla")
name  <- c("Pielou", "Camargo", "SimpsonEvenness",  "Evar", "Bulla")

# Estimate evenness and give polished names to be used in the output
tse <- estimateEvenness(tse, index = index, name = name)

# Check the output
head(colData(tse))


microbiome/mia documentation built on April 27, 2024, 4:04 a.m.