test_differential_cellularity | R Documentation |
test_differential_cellularity() takes as input A 'tbl' (with at least three columns for sample, feature and transcript abundance) or 'SummarizedExperiment' (more convenient if abstracted to tibble with library(tidySummarizedExperiment)) and returns a consistent object (to the input) with additional columns for the statistics from the hypothesis test.
test_differential_cellularity(
.data,
.formula,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "cibersort",
reference = X_cibersort,
significance_threshold = 0.05,
...
)
## S4 method for signature 'spec_tbl_df'
test_differential_cellularity(
.data,
.formula,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "cibersort",
reference = X_cibersort,
significance_threshold = 0.05,
...
)
## S4 method for signature 'tbl_df'
test_differential_cellularity(
.data,
.formula,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "cibersort",
reference = X_cibersort,
significance_threshold = 0.05,
...
)
## S4 method for signature 'tidybulk'
test_differential_cellularity(
.data,
.formula,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "cibersort",
reference = X_cibersort,
significance_threshold = 0.05,
...
)
## S4 method for signature 'SummarizedExperiment'
test_differential_cellularity(
.data,
.formula,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "cibersort",
reference = X_cibersort,
significance_threshold = 0.05,
...
)
## S4 method for signature 'RangedSummarizedExperiment'
test_differential_cellularity(
.data,
.formula,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "cibersort",
reference = X_cibersort,
significance_threshold = 0.05,
...
)
.data |
A 'tbl' (with at least three columns for sample, feature and transcript abundance) or 'SummarizedExperiment' (more convenient if abstracted to tibble with library(tidySummarizedExperiment)) |
.formula |
A formula representing the desired linear model. The formula can be of two forms: multivariable (recommended) or univariable Respectively: \"factor_of_interest ~ .\" or \". ~ factor_of_interest\". The dot represents cell-type proportions, and it is mandatory. If censored regression is desired (coxph) the formula should be of the form \"survival::Surv\(y, dead\) ~ .\" |
.sample |
The name of the sample column |
.transcript |
The name of the transcript/gene column |
.abundance |
The name of the transcript/gene abundance column |
method |
A string character. Either \"cibersort\", \"epic\" or \"llsr\". The regression method will be chosen based on being multivariable: lm or cox-regression (both on logit-transformed proportions); or univariable: beta or cox-regression (on logit-transformed proportions). See .formula for multi- or univariable choice. |
reference |
A data frame. The transcript/cell_type data frame of integer transcript abundance |
significance_threshold |
A real between 0 and 1 (usually 0.05). |
... |
Further parameters passed to the method deconvolve_cellularity |
'r lifecycle::badge("maturing")'
This routine applies a deconvolution method (e.g., Cibersort; DOI: 10.1038/nmeth.3337) and passes the proportions inferred into a generalised linear model (DOI:dx.doi.org/10.1007/s11749-010-0189-z) or a cox regression model (ISBN: 978-1-4757-3294-8)
Underlying method for the generalised linear model: data |> deconvolve_cellularity( !!.sample, !!.transcript, !!.abundance, method=method, reference = reference, action="get", ... ) [..] betareg::betareg(.my_formula, .)
Underlying method for the cox regression: data |> deconvolve_cellularity( !!.sample, !!.transcript, !!.abundance, method=method, reference = reference, action="get", ... ) [..] mutate(.proportion_0_corrected = .proportion_0_corrected |> boot::logit()) survival::coxph(.my_formula, .)
A consistent object (to the input) with additional columns for the statistics from the hypothesis test (e.g., log fold change, p-value and false discovery rate).
A 'SummarizedExperiment' object
A 'SummarizedExperiment' object
# Regular regression
test_differential_cellularity(
tidybulk::se_mini ,
. ~ condition,
cores = 1
)
# Cox regression - multiple
tidybulk::se_mini |>
# Test
test_differential_cellularity(
survival::Surv(days, dead) ~ .,
cores = 1
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.