getBaselineDivergence: Beta diversity between the baseline and later time steps

View source: R/getBaselineDivergence.R

getBaselineDivergenceR Documentation

Beta diversity between the baseline and later time steps

Description

Calculates sample dissimilarity between the given baseline and other time points, optionally 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

getBaselineDivergence(
  x,
  group = NULL,
  time_field,
  name_divergence = "divergence_from_baseline",
  name_timedifference = "time_from_baseline",
  assay.type = "counts",
  FUN = vegan::vegdist,
  method = "bray",
  baseline_sample = NULL,
  altexp = NULL,
  dimred = NULL,
  n_dimred = NULL,
  ...
)

Arguments

x

A SummarizedExperiment object.

group

optional; a single character value for specifying the grouping factor (name of a colData field). If given, the divergence is calculated per group. e.g. subject, chamber, group etc.).

time_field

a single character value, specifying the name of the time series field in colData.

name_divergence

a column vector showing beta diversity between samples (default: name_divergence = "time_divergence")

name_timedifference

field name for adding the time difference between samples used to calculate beta diversity (default: name_timedifference = "time_difference")

assay.type

character indicating which assay values are used in the dissimilarity estimation (default: assay.type = "counts").

FUN

a function for dissimilarity calculation. The function must expect the input matrix as its first argument. With rows as samples and columns as features. By default, FUN is vegan::vegdist.

method

a method that is used to calculate the distance. Method is passed to the function that is specified by FUN. By default, method is "bray".

baseline_sample

Optional. A character vector specifying the baseline sample(s) to be used. If the "group" argument is given, this must be a named vector; one element per group.

altexp

String or integer scalar specifying the alternative experiment containing the input data.

dimred

A string or integer scalar indicating the reduced dimension result in reducedDims to use in the estimation.

n_dimred

Integer scalar or vector specifying the dimensions to use if dimred is specified.

...

Arguments to be passed

Details

The group argument allows calculating divergence per group. Otherwise, this is done across all samples at once.

The baseline sample/s always need to belong to the data object i.e. they can be merged into it before applying this function. The reason is that they need to have comparable sample data, at least some time point information for calculating time differences w.r.t. baseline sample.

The baseline time point is by default defined as the smallest time point (per group). Alternatively, the user can provide the baseline vector, or a list of baseline vectors per group (named list per group).

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)
library(dplyr)

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

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

tse2 <- getBaselineDivergence(tse,
                              group = "subject",
                              time_field = "time",
                              name_divergence = "divergence_from_baseline",
                              name_timedifference = "time_from_baseline",
                              assay.type="relabundance",
                              FUN = vegan::vegdist,
                              method="bray")

tse2 <- getBaselineDivergence(tse,
                              baseline_sample = "Sample-875",
                              group = "subject",
                              time_field = "time",
                              name_divergence = "divergence_from_baseline",
                              name_timedifference = "time_from_baseline",
                              assay.type="relabundance",
                              FUN = vegan::vegdist,
                              method="bray")


microbiome/miaTime documentation built on May 7, 2023, 2:06 p.m.