View source: R/getBaselineDivergence.R
getBaselineDivergence | R Documentation |
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
.
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,
...
)
x |
A
|
group |
|
time_field |
|
name_divergence |
|
name_timedifference |
|
assay.type |
|
FUN |
|
method |
|
baseline_sample |
|
altexp |
|
dimred |
|
n_dimred |
|
... |
Arguments to be passed |
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).
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.
#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")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.