R/methods.R
, R/methods_SE.R
test_gene_enrichment-methods.Rd
test_gene_enrichment() 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 `tbl` of gene set information
test_gene_enrichment(
.data,
.formula,
.sample = NULL,
.entrez,
.abundance = NULL,
contrasts = NULL,
methods = c("camera", "roast", "safe", "gage", "padog", "globaltest", "ora"),
gene_sets = c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease",
"kegg_metabolism", "kegg_signaling"),
species,
cores = 10,
method = NULL,
.contrasts = NULL
)
# S4 method for class 'spec_tbl_df'
test_gene_enrichment(
.data,
.formula,
.sample = NULL,
.entrez,
.abundance = NULL,
contrasts = NULL,
methods = c("camera", "roast", "safe", "gage", "padog", "globaltest", "ora"),
gene_sets = c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease",
"kegg_metabolism", "kegg_signaling"),
species,
cores = 10,
method = NULL,
.contrasts = NULL
)
# S4 method for class 'tbl_df'
test_gene_enrichment(
.data,
.formula,
.sample = NULL,
.entrez,
.abundance = NULL,
contrasts = NULL,
methods = c("camera", "roast", "safe", "gage", "padog", "globaltest", "ora"),
gene_sets = c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease",
"kegg_metabolism", "kegg_signaling"),
species,
cores = 10,
method = NULL,
.contrasts = NULL
)
# S4 method for class 'tidybulk'
test_gene_enrichment(
.data,
.formula,
.sample = NULL,
.entrez,
.abundance = NULL,
contrasts = NULL,
methods = c("camera", "roast", "safe", "gage", "padog", "globaltest", "ora"),
gene_sets = c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease",
"kegg_metabolism", "kegg_signaling"),
species,
cores = 10,
method = NULL,
.contrasts = NULL
)
# S4 method for class 'SummarizedExperiment'
test_gene_enrichment(
.data,
.formula,
.sample = NULL,
.entrez,
.abundance = NULL,
contrasts = NULL,
methods = c("camera", "roast", "safe", "gage", "padog", "globaltest", "ora"),
gene_sets = c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease",
"kegg_metabolism", "kegg_signaling"),
species,
cores = 10,
method = NULL,
.contrasts = NULL
)
# S4 method for class 'RangedSummarizedExperiment'
test_gene_enrichment(
.data,
.formula,
.sample = NULL,
.entrez,
.abundance = NULL,
contrasts = NULL,
methods = c("camera", "roast", "safe", "gage", "padog", "globaltest", "ora"),
gene_sets = c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease",
"kegg_metabolism", "kegg_signaling"),
species,
cores = 10,
method = NULL,
.contrasts = NULL
)
A `tbl` (with at least three columns for sample, feature and transcript abundance) or `SummarizedExperiment` (more convenient if abstracted to tibble with library(tidySummarizedExperiment))
A formula with no response variable, representing the desired linear model
The name of the sample column
The ENTREZ ID of the transcripts/genes
The name of the transcript/gene abundance column
This parameter takes the format of the contrast parameter of the method of choice. For edgeR and limma-voom is a character vector. For DESeq2 is a list including a character vector of length three. The first covariate is the one the model is tested against (e.g., ~ factor_of_interest)
A character vector. One or 3 or more methods to use in the testing (currently EGSEA errors if 2 are used). Type EGSEA::egsea.base() to see the supported GSE methods.
A character vector or a list. It can take one or more of the following built-in collections as a character vector: c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease", "kegg_metabolism", "kegg_signaling"), to be used with EGSEA buildIdx. c1 is human specific. Alternatively, a list of user-supplied gene sets can be provided, to be used with EGSEA buildCustomIdx. In that case, each gene set is a character vector of Entrez IDs and the names of the list are the gene set names.
A character. It can be human, mouse or rat.
An integer. The number of cores available
DEPRECATED. Please use methods.
DEPRECATED - This parameter takes the format of the contrast parameter of the method of choice. For edgeR and limma-voom is a character vector. For DESeq2 is a list including a character vector of length three. The first covariate is the one the model is tested against (e.g., ~ factor_of_interest)
A consistent object (to the input)
A consistent object (to the input)
A consistent object (to the input)
A consistent object (to the input)
A consistent object (to the input)
A consistent object (to the input)
`r lifecycle::badge("maturing")`
This wrapper executes ensemble gene enrichment analyses of the dataset using EGSEA (DOI:0.12688/f1000research.12544.1)
dge = data |> keep_abundant( factor_of_interest = !!as.symbol(parse_formula(.formula)[[1]]), !!.sample, !!.entrez, !!.abundance )
# Make sure transcript names are adjacent [...] as_matrix(rownames = !!.entrez) edgeR::DGEList(counts = .)
idx = buildIdx(entrezIDs = rownames(dge), species = species, msigdb.gsets = msigdb.gsets, kegg.exclude = kegg.exclude)
dge |>
# Calculate weights limma::voom(design, plot = FALSE) |>
# Execute EGSEA egsea( contrasts = my_contrasts, baseGSEAs = methods, gs.annots = idx, sort.by = "med.rank", num.threads = cores, report = FALSE )
if (FALSE) { # \dontrun{
library(SummarizedExperiment)
se = tidybulk::se_mini
rowData( se)$entrez = rownames(se )
df_entrez = aggregate_duplicates(se,.transcript = entrez )
library("EGSEA")
test_gene_enrichment(
df_entrez,
~ condition,
.sample = sample,
.entrez = entrez,
.abundance = count,
methods = c("roast" , "safe", "gage" , "padog" , "globaltest", "ora" ),
gene_sets = c("h", "c1", "c2", "c3", "c4", "c5", "c6", "c7", "kegg_disease", "kegg_metabolism", "kegg_signaling"),
species="human",
cores = 2
)
} # }