getStability: Estimate stability

getStabilityR Documentation

Estimate stability

Description

Quantify intermediate stability with respect to a given reference point.

Usage

getStability(x, ...)

addStability(x, ...)

## S4 method for signature 'SummarizedExperiment'
addStability(x, time.col, name = "stability", ...)

## S4 method for signature 'SummarizedExperiment'
getStability(
  x,
  time.col,
  assay.type = "counts",
  reference = NULL,
  group = NULL,
  ...
)

Arguments

x

A SummarizedExperiment object.

...

additional arguments.

  • time.interval: Integer scalar. Indicates the increment between time steps. By default, the function compares each sample to the previous one. If you need to take every second, every third, or so, time step, then increase this accordingly. (Default: 1L)

  • calc.separately: Logical scalar. Specifies whether to calculate stability separately for data points with abundance below or at/above the reference point. (Defaul: FALSE)

  • mode: Character scalar. Specifies whether to calculate stability using correlation ("correlation") or a linear model ("lm"). (Default: "correlation")

time.col

Character scalar. Specifies a name of the column from colData that identifies the sampling time points for the samples.

name

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

assay.type

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

reference

Character scalar. Specifies a column from rowData(x) that determines the reference points. If NULL, the default reference is used. (Default: NULL)

group

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

Details

These methods estimate intermediate stability described in Lahti et al. 2014. The method is heuristic and makes many simplifying assumptions. However, the stability estimation can be useful for exploration, and can provide a baseline for more advanced measures.

The stability is calculated by first defining reference point, R_{f}, for each feature. User can define the reference points with reference. If reference points are not defined, they are calculated by taking median

R_{f} = median(x_{f})

where f denotes a single feature and x abundance.

Then difference between consecutive time points, \Delta_{f, x},

\Delta_{f, x} = ∣x_{f, t} - x_{f, t-1}∣

and difference between the previous time point and reference, \Delta_{f. R}, are calculated:

\Delta_{f. R} = ∣x_{f, t-1} - R_{f}∣

where t denotes time point.

Optionally, time difference \Delta_{f, t} is calculated

\Delta_{f, t} = t_{f, t-1} - t_{f, t-1}

The stability coefficient s_{f} is calculated with correlation

s_{f} = corr(\Delta_{f, x}, \Delta_{f, R})

or with linear model as follows

s_{f} = \frac{ \Delta_{f, x} - \beta_{0} - \beta_{2}*\Delta_{f, t} - \epsilon_{f} }{ \Delta_{f, R} }

Value

getStability returns DataFrame while addStability returns results added to its rowData(x).

References

Lahti L, et al. (2014) Tipping elements in the human intestinal ecosystem. Nat Commun. doi: 10.1038/ncomms5344

See Also

getBimodality()

Examples

library(miaTime)

# Load time series data
data(minimalgut)
tse <- minimalgut

# Apply clr transformation
tse <- transformAssay(tse, method = "rclr")

# Calculate stability single system
tse_sub <- tse[, tse[["StudyIdentifier"]] == "Bioreactor A"]
tse_sub <- addStability(tse_sub, assay.type = "rclr", time.col = "Time.hr")
rowData(tse_sub)

# Add custom reference values
rowData(tse)[["ref_col"]] <- 1
# Calculate stability for each system simultaneously by taking time
# difference into account
tse <- addStability(
    tse, assay.type = "rclr", time.col = "Time.hr",
    group = "StudyIdentifier", ref_col = "ref_col", mode = "lm")
rowData(tse)


microbiome/miaTime documentation built on Feb. 8, 2025, 7:33 p.m.