getMediation
and addMediation
provide a wrapper of
mediate
for
SummarizedExperiment
.
addMediation(x, ...)
# S4 method for class 'SummarizedExperiment'
addMediation(
x,
outcome,
treatment,
name = "mediation",
mediator = NULL,
assay.type = NULL,
dimred = NULL,
family = gaussian(),
covariates = NULL,
p.adj.method = "holm",
add.metadata = FALSE,
verbose = TRUE,
...
)
getMediation(x, ...)
# S4 method for class 'SummarizedExperiment'
getMediation(
x,
outcome,
treatment,
mediator = NULL,
assay.type = "counts",
dimred = NULL,
family = gaussian(),
covariates = NULL,
p.adj.method = "holm",
add.metadata = FALSE,
verbose = TRUE,
...
)
additional parameters that can be passed to
mediate
.
Character scalar
. Indicates the colData variable used
as outcome in the model.
Character scalar
. Indicates the colData variable
used as treatment in the model.
Character scalar
. A name for the column of the
colData
where results will be stored. (Default: "mediation"
)
Character scalar
. Indicates the colData variable used
as mediator in the model. (Default: NULL
)
Character scalar
. Specifies the assay used for
feature-wise mediation analysis. (Default: "counts"
)
Character scalar
. Indicates the reduced dimension
result in reducedDims(object)
for component-wise mediation analysis.
(Default: NULL
)
Character scalar
. A specification for the outcome model link function.
(Default: gaussian("identity")
)
Character scalar
or character vector
. Indicates the colData
variables used as covariates in the model.
(Default: NULL
)
Character scalar
. Selects adjustment method
of p-values. Passed to p.adjust
function.
(Default: "holm"
)
Logical scalar
. Should the model metadata be returned.
(Default: FALSE
)
Logical scalar
. Should execution messages be printed.
(Default: TRUE
)
getMediation
returns a data.frame
of adjusted p-values and
effect sizes for the ACMEs and ADEs of every mediator given as input, whereas
addMediation
returns an updated
SummarizedExperiment
instance with the same data.frame
stored in the metadata as
"mediation". Its columns include:
the mediator variable
the Average Causal Mediation Effect (ACME) estimate
the Average Direct Effect (ADE) estimate
the adjusted p-value for the ACME estimate
the adjusted p-value for the ADE estimate
This wrapper of mediate
for
SummarizedExperiment
provides a simple method to analyse the effect of a treatment variable on an
outcome variable found in colData(se)
through the mediation of either
another variable in colData (argument mediator
) or an assay
(argument assay.type
) or a reducedDim (argument dimred
).
Importantly, those three arguments are mutually exclusive.
if (FALSE) { # \dontrun{
# Import libraries
library(mia)
library(scater)
# Load dataset
data(hitchip1006, package = "miaTime")
tse <- hitchip1006
# Agglomerate features by family (merely to speed up execution)
tse <- agglomerateByRank(tse, rank = "Phylum")
# Convert BMI variable to numeric
tse$bmi_group <- as.numeric(tse$bmi_group)
# Analyse mediated effect of nationality on BMI via alpha diversity
# 100 permutations were done to speed up execution, but ~1000 are recommended
med_df <- getMediation(tse,
outcome = "bmi_group",
treatment = "nationality",
mediator = "diversity",
covariates = c("sex", "age"),
treat.value = "Scandinavia",
control.value = "CentralEurope",
boot = TRUE, sims = 100,
add.metadata = TRUE)
# Visualise model statistics for 1st mediator
plot(attr(med_df, "metadata")[[1]])
# Apply clr transformation to counts assay
tse <- transformAssay(tse,
method = "clr",
pseudocount = 1)
# Analyse mediated effect of nationality on BMI via clr-transformed features
# 100 permutations were done to speed up execution, but ~1000 are recommended
tse <- addMediation(tse, name = "assay_mediation",
outcome = "bmi_group",
treatment = "nationality",
assay.type = "clr",
covariates = c("sex", "age"),
treat.value = "Scandinavia",
control.value = "CentralEurope",
boot = TRUE, sims = 100,
p.adj.method = "fdr")
# Show results for first 5 mediators
head(metadata(tse)$assay_mediation, 5)
# Perform ordination
tse <- runMDS(tse, name = "MDS",
method = "euclidean",
assay.type = "clr",
ncomponents = 3)
# Analyse mediated effect of nationality on BMI via NMDS components
# 100 permutations were done to speed up execution, but ~1000 are recommended
tse <- addMediation(tse, name = "reddim_mediation",
outcome = "bmi_group",
treatment = "nationality",
dimred = "MDS",
covariates = c("sex", "age"),
treat.value = "Scandinavia",
control.value = "CentralEurope",
boot = TRUE, sims = 100,
p.adj.method = "fdr")
# Show results for first 5 mediators
head(metadata(tse)$reddim_mediation, 5)
} # }