R/AllGenerics.R
, R/getBaselineDivergence.R
addBaselineDivergence.Rd
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.
getBaselineDivergence(x, ...)
addBaselineDivergence(x, ...)
# S4 method for class 'SummarizedExperiment'
getBaselineDivergence(
x,
time.col,
assay.type = "counts",
reference = NULL,
group = NULL,
method = "bray",
...
)
# S4 method for class 'SummarizedExperiment'
addBaselineDivergence(
x,
name = c("divergence", "time_diff", "ref_samples"),
...
)
A
SummarizedExperiment
object.
Optional arguments passed into
mia::addDivergence()
.
Character scalar
. Specifies a name of the column from
colData
that identifies the sampling time points for the samples.
Character scalar
. Specifies which assay values are
used in the dissimilarity estimation. (Default: "counts"
)
Character scalar
. Specifies a name of the column from
colData
that identifies the baseline samples to be used.
(Default: NULL
)
Character scalar
. Specifies a name of the column from
colData
that identifies the grouping of the samples.
(Default: NULL
)
Character scalar
. Used to calculate the dissimilarity
Method is passed to the function that is specified by dis.fun
.
(Default: "bray"
)
Character vector
. Specifies a column name for storing
divergence results.
(Default: c("divergence", "time_diff", "ref_samples")
)
getBaselineDivergence
returns DataFrame
object
containing the sample dissimilarity and corresponding time difference between
samples. addBaselineDivergence
, on the other hand, returns a
SummarizedExperiment
object with these results in its colData
.
The group argument allows calculating divergence per group. If given, the divergence is calculated per group. e.g. subject, chamber, group etc. 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).
library(miaTime)
library(mia)
data(hitchip1006)
tse <- transformAssay(hitchip1006, method = "relabundance")
# By default, reference samples are the samples from the first timepoint
tse <- addBaselineDivergence(
tse,
group = "subject",
time.col = "time",
assay.type = "relabundance",
method = "bray")
# Add reference samples to colData, if you want to specify reference
# samples manually
colData(tse)[["reference"]] <- "Sample-875"
tse <- addBaselineDivergence(
tse,
reference = "reference",
group = "subject",
time.col = "time",
name = c("divergence_from_baseline",
"time_from_baseline", "reference_samples"),
assay.type = "relabundance",
method = "bray")