getPrevalence: Calculation prevalence information for features across...

getPrevalenceR Documentation

Calculation prevalence information for features across samples

Description

These functions calculate the population prevalence for taxonomic ranks in a SummarizedExperiment-class object.

Usage

getPrevalence(x, ...)

## S4 method for signature 'ANY'
getPrevalence(
  x,
  detection = 0,
  include_lowest = FALSE,
  sort = FALSE,
  na.rm = TRUE,
  ...
)

## S4 method for signature 'SummarizedExperiment'
getPrevalence(
  x,
  assay.type = assay_name,
  assay_name = "counts",
  as_relative = FALSE,
  rank = NULL,
  ...
)

getPrevalentFeatures(x, ...)

## S4 method for signature 'ANY'
getPrevalentFeatures(x, prevalence = 50/100, include_lowest = FALSE, ...)

## S4 method for signature 'SummarizedExperiment'
getPrevalentFeatures(
  x,
  rank = NULL,
  prevalence = 50/100,
  include_lowest = FALSE,
  ...
)

getPrevalentTaxa(x, ...)

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

getRareFeatures(x, ...)

## S4 method for signature 'ANY'
getRareFeatures(x, prevalence = 50/100, include_lowest = FALSE, ...)

## S4 method for signature 'SummarizedExperiment'
getRareFeatures(
  x,
  rank = NULL,
  prevalence = 50/100,
  include_lowest = FALSE,
  ...
)

getRareTaxa(x, ...)

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

subsetByPrevalentFeatures(x, ...)

## S4 method for signature 'SummarizedExperiment'
subsetByPrevalentFeatures(x, rank = NULL, ...)

subsetByPrevalentTaxa(x, ...)

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

subsetByRareFeatures(x, ...)

## S4 method for signature 'SummarizedExperiment'
subsetByRareFeatures(x, rank = NULL, ...)

subsetByRareTaxa(x, ...)

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

getPrevalentAbundance(
  x,
  assay.type = assay_name,
  assay_name = "relabundance",
  ...
)

## S4 method for signature 'ANY'
getPrevalentAbundance(
  x,
  assay.type = assay_name,
  assay_name = "relabundance",
  ...
)

## S4 method for signature 'SummarizedExperiment'
getPrevalentAbundance(x, assay.type = assay_name, assay_name = "counts", ...)

Arguments

x

a SummarizedExperiment object

...

additional arguments

  • If !is.null(rank) arguments are passed on to agglomerateByRank. See ?agglomerateByRank for more details. Note that you can specify whether to remove empty ranks with agg.na.rm instead of na.rm. (default: FALSE)

  • for getPrevalentFeatures, getRareFeatures, subsetByPrevalentFeatures and subsetByRareFeatures additional parameters passed to getPrevalence

  • for getPrevalentAbundance additional parameters passed to getPrevalentFeatures

detection

Detection threshold for absence/presence. Either an absolute value compared directly to the values of x or a relative value between 0 and 1, if as_relative = FALSE.

include_lowest

logical scalar: Should the lower boundary of the detection and prevalence cutoffs be included? (default: FALSE)

sort

logical scalar: Should the result be sorted by prevalence? (default: FALSE)

na.rm

logical scalar: Should NA values be omitted when calculating prevalence? (default: na.rm = TRUE)

assay.type

A single character value for selecting the assay to use for prevalence calculation.

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

as_relative

logical scalar: Should the detection threshold be applied on compositional (relative) abundances? (default: FALSE)

rank

a single character defining a taxonomic rank. Must be a value of taxonomyRanks() function.

prevalence

Prevalence threshold (in 0 to 1). The required prevalence is strictly greater by default. To include the limit, set include_lowest to TRUE.

Details

getPrevalence calculates the relative frequency of samples that exceed the detection threshold. For SummarizedExperiment objects, the prevalence is calculated for the selected taxonomic rank, otherwise for the rows. The absolute population prevalence can be obtained by multiplying the prevalence by the number of samples (ncol(x)). If as_relative = FALSE the relative frequency (between 0 and 1) is used to check against the detection threshold.

The core abundance index from getPrevalentAbundance gives the relative proportion of the core species (in between 0 and 1). The core taxa are defined as those that exceed the given population prevalence threshold at the given detection level as set for getPrevalentFeatures.

subsetPrevalentFeatures and subsetRareFeatures return a subset of x. The subset includes the most prevalent or rare taxa that are calculated with getPrevalentFeatures or getRareFeatures respectively.

getPrevalentFeatures returns taxa that are more prevalent with the given detection threshold for the selected taxonomic rank.

getRareFeatures returns complement of getPrevalentTaxa.

Value

subsetPrevalentFeatures and subsetRareFeatures return subset of x.

All other functions return a named vectors:

  • getPrevalence returns a numeric vector with the names being set to either the row names of x or the names after agglomeration.

  • getPrevalentAbundance returns a numeric vector with the names corresponding to the column name of x and include the joint abundance of prevalent taxa.

  • getPrevalentTaxa and getRareFeatures return a character vector with only the names exceeding the threshold set by prevalence, if the rownames of x is set. Otherwise an integer vector is returned matching the rows in x.

Author(s)

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

References

A Salonen et al. The adult intestinal core microbiota is determined by analysis depth and health status. Clinical Microbiology and Infection 18(S4):16 20, 2012. To cite the R package, see citation('mia')

See Also

agglomerateByRank, getTopTaxa

Examples

data(GlobalPatterns)
tse <- GlobalPatterns
# Get prevalence estimates for individual ASV/OTU
prevalence.frequency <- getPrevalence(tse,
                                      detection = 0,
                                      sort = TRUE,
                                      as_relative = TRUE)
head(prevalence.frequency)

# Get prevalence estimates for phylums
# - the getPrevalence function itself always returns population frequencies
prevalence.frequency <- getPrevalence(tse,
                                      rank = "Phylum",
                                      detection = 0,
                                      sort = TRUE,
                                      as_relative = TRUE)
head(prevalence.frequency)

# - to obtain population counts, multiply frequencies with the sample size,
# which answers the question "In how many samples is this phylum detectable"
prevalence.count <- prevalence.frequency * ncol(tse)
head(prevalence.count)

# Detection threshold 1 (strictly greater by default);
# Note that the data (GlobalPatterns) is here in absolute counts
# (and not compositional, relative abundances)
# Prevalence threshold 50 percent (strictly greater by default)
prevalent <- getPrevalentFeatures(tse,
                              rank = "Phylum",
                              detection = 10,
                              prevalence = 50/100,
                              as_relative = FALSE)
head(prevalent)

# Gets a subset of object that includes prevalent taxa
altExp(tse, "prevalent") <- subsetByPrevalentFeatures(tse,
                                       rank = "Family",
                                       detection = 0.001,
                                       prevalence = 0.55,
                                       as_relative = TRUE)
altExp(tse, "prevalent")                                 

# getRareFeatures returns the inverse
rare <- getRareFeatures(tse,
                    rank = "Phylum",
                    detection = 1/100,
                    prevalence = 50/100,
                    as_relative = TRUE)
head(rare)

# Gets a subset of object that includes rare taxa
altExp(tse, "rare") <- subsetByRareFeatures(tse,
                             rank = "Class",
                             detection = 0.001,
                             prevalence = 0.001,
                             as_relative = TRUE)
altExp(tse, "rare")      

# Names of both experiments, prevalent and rare, can be found from slot altExpNames
tse
                         
data(esophagus)
getPrevalentAbundance(esophagus, assay.type = "counts")


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