View source: R/mod_understand_fct_enrichment.R
runMultiEnrichment | R Documentation |
Launches functional enrichment of features provided for every databases furnished. Run by utils function runMultiOmicEnrichment
runMultiEnrichment( omicSignature, databasesChosen, organismDb, pvAdjust = "BH", minGSSize = 5, maxGSSize = 800, pvStouffer = 0.1 )
omicSignature |
Feature lists from each omic data block. |
databasesChosen |
Which biological database c(reactome, kegg, wikiPathways, MF, CC, BP) |
organismDb |
Organism provided by user in Home tab. |
pvAdjust |
pv adjust method (e.g "BH" for Benjamini-Hochberg) |
minGSSize, maxGSSize, pvStouffer |
Numeric values chosen by user in ui. |
Wraps in obj enrichment results for all databases and all omics.
data("omic2", package = "multiSight") splitData <- splitDatatoTrainTest(omic2, 0.8) data.train <- splitData$data.train data.test <- splitData$data.test diabloRes <- runSPLSDA(data.train) diabloModels <- diabloRes$model #sPLS-DA model using all omics. diabloFeats <- diabloRes$biosignature #selected features for each omic. id_db <- list(omic1 = "ENSEMBL", omic2 = "ENSEMBL") if (requireNamespace("org.Mm.eg.db", quietly = TRUE)) { library(org.Mm.eg.db, warn.conflicts = FALSE) #Organism's database featList <- list(Omic1 = c("ENSMUSG00000039621", "ENSMUSG00000038733", "ENSMUSG00000062031"), Omic2 = c("ENSMUSG00000031170", "ENSMUSG00000077495", "ENSMUSG00000042992")) dbList <- list(Omic1 = "ENSEMBL", Omic2 = "ENSEMBL") convFeat <- convertToEntrezid(featList, dbList, "org.Mm.eg.db") ## To enrich features database <- c("reactome", "MF") #runMultiEnrichment_result <- runMultiEnrichment(databasesChosen = database, # omicSignature = convFeat, # organismDb = "org.Mm.eg.db") }
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.