getMediation | R Documentation |
getMediation
and addMediation
provide a wrapper of
mediate
for
SummarizedExperiment
.
addMediation(x, ...)
## S4 method for signature '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 signature '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,
...
)
x |
a |
... |
additional parameters that can be passed to
|
outcome |
|
treatment |
|
name |
|
mediator |
|
assay.type |
|
dimred |
|
family |
|
covariates |
|
p.adj.method |
|
add.metadata |
|
verbose |
|
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.
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
## Not run:
# 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)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.