getStepwiseDivergence: Beta diversity between consecutive time steps

View source: R/getStepwiseDivergence.R

getStepwiseDivergenceR Documentation

Beta diversity between consecutive time steps

Description

Calculates sample dissimilarity between consecutive time points (t, t+i), within a group (subject, reaction chamber, or similar). The corresponding time difference is returned as well. The method operates on SummarizedExperiment objects, and the results are stored in colData.

Usage

getStepwiseDivergence(
  x,
  group = NULL,
  time_field,
  time_interval = 1,
  name_divergence = "time_divergence",
  name_timedifference = "time_difference",
  assay.type = "counts",
  FUN = vegan::vegdist,
  method = "bray",
  altexp = NULL,
  dimred = NULL,
  n_dimred = NULL,
  ...
)

getTimeDivergence(x, ...)

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

Arguments

x

A SummarizedExperiment object.

group

Character scalar. Specifies the grouping factor (name of a colData field). If given, the divergence is calculated per group. e.g. subject, chamber, group etc.). (Default: NULL)

time_field

Character scalar. Specifies the name of the time series field in colData.

time_interval

Integer scalar. Indicates the increment between time steps. If you need to take every second, every third, or so, time step only, then increase this accordingly. (Default: 1)

name_divergence

Character scalar. Shows beta diversity between samples. (Default: "time_divergence")

name_timedifference

Character scalar. Field name for adding the time difference between samples used to calculate beta diversity. (Default: "time_difference")

assay.type

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

FUN

Function for dissimilarity calculation. The function must expect the input matrix as its first argument. With rows as samples and columns as features. (Default: vegan::vegdist)

method

Character scalar. Used to calculate the distance. Method is passed to the function that is specified by FUN. (Default: "bray")

altexp

Character scalar or integer scalar. Specifies the alternative experiment containing the input data. (Default: NULL)

dimred

Character scalar or integer scalar. indicates the reduced dimension result in reducedDims to use in the estimation. (Default: NULL)

n_dimred

Integer vector. Specifies the dimensions to use if dimred is specified. (Default: NULL)

...

Arguments to be passed

Value

a SummarizedExperiment or TreeSummarizedExperiment containing the sample dissimilarity and corresponding time difference between samples (across n time steps), within each level of the grouping factor.

Examples

#library(miaTime)
library(TreeSummarizedExperiment)

data(hitchip1006)
tse <- mia::transformCounts(hitchip1006, method = "relabundance")

# Subset to speed up example
tse <- tse[, colData(tse)$subject %in% c("900", "934", "843", "875")]

# Using vegdist for divergence calculation, one can pass
# the dissimilarity method from the vegan::vegdist options
# via the "method" argument
tse2 <- getStepwiseDivergence(tse, group = "subject",
                              time_interval = 1,
                              time_field = "time",
                              assay.type="relabundance",
                              FUN = vegan::vegdist,
                              method="bray")


microbiome/miaTime documentation built on Oct. 3, 2024, 5:04 a.m.