merge_enrich_terms {ViSEAGO}R Documentation

Merge enriched GO terms.

Description

combine results from GO enrichment tests (obtained with topGO package) or from fgsea (obtained with runfgsea method), for a given ontology (MF, BP, or CC).

Usage

merge_enrich_terms(Input, cutoff = 0.01, envir = .GlobalEnv)

## S4 method for signature 'list'
merge_enrich_terms(Input, cutoff = 0.01, envir = .GlobalEnv)

Arguments

Input

a list containing named elements. Each element must contain the name of:

cutoff

default pvalue cutoff (default to 0.01). Several cutoff can be use in the same order as list elements.

envir

objects environment (default to .GlobalEnv).

Details

This method extracts for each result of GO enrichment test: informations about GO term (identifiant, name, and description), gene frequency (number of significant genes / Annotated genes), pvalue, -log10(pvalue), significant genes identifiants (GeneID, or Ensembl ID, or uniprot accession), and gene symbols. At the last, this method builds a merged data.table of enriched GO terms at least once and provides all mentionned columns.

Value

an enrich_GO_terms-class object.

References

Matt Dowle and Arun Srinivasan (2017). data.table: Extension of data.frame. R package version 1.10.4. https://CRAN.R-project.org/package=data.table

Herve Pages, Marc Carlson, Seth Falcon and Nianhua Li (2017). AnnotationDbi: Annotation Database Interface. R package version 1.38.0.

See Also

Other GO_terms: GOcount(), GOterms_heatmap(), annotate(), create_topGOdata(), gene2GO-class, runfgsea()

Examples

## topGO terms enrichment

# load genes identifiants (GeneID,ENS...) universe/background (Expressed genes)
background_L<-scan(
    system.file(
        "extdata/data/input",
        "background_L.txt",
        package = "ViSEAGO"
    ),
    quiet=TRUE,
    what=""
)

# load Differentialy Expressed (DE) gene identifiants from files
PregnantvslactateDE<-scan(
    system.file(
        "extdata/data/input",
        "pregnantvslactateDE.txt",
        package = "ViSEAGO"
    ),
    quiet=TRUE,
    what=""
)

VirginvslactateDE<-scan(
    system.file(
        "extdata/data/input",
        "virginvslactateDE.txt",
        package = "ViSEAGO"
    ),
    quiet=TRUE,
    what=""
)

VirginvspregnantDE<-scan(
    system.file(
        "extdata/data/input",
        "virginvspregnantDE.txt",
        package="ViSEAGO"
    ),
    quiet=TRUE,
    what=""
)

## Not run: 
# connect to Bioconductor
Bioconductor<-ViSEAGO::Bioconductor2GO()

# load GO annotations from Bioconductor
myGENE2GO<-ViSEAGO::annotate(
    "org.Mm.eg.db",
    Bioconductor
)

# create topGOdata for BP for each list of DE genes
BP_Pregnantvslactate<-ViSEAGO::create_topGOdata(
    geneSel=PregnantvslactateDE,
    allGenes=background_L,
    gene2GO=myGENE2GO,
    ont="BP",
    nodeSize=5
)

BP_Virginvslactate<-ViSEAGO::create_topGOdata(
    geneSel=VirginvslactateDE,
    allGenes=background_L,
    gene2GO=myGENE2GO,
    ont="BP",
    nodeSize=5
)

BP_Virginvspregnant<-ViSEAGO::create_topGOdata(
    geneSel=VirginvspregnantDE,
    allGenes=background_L,
    gene2GO=myGENE2GO,
    ont="BP",
    nodeSize=5
)

# perform TopGO tests
elim_BP_Pregnantvslactate<-topGO::runTest(
    BP_L_pregnantvslactate,
    algorithm ="elim",
    statistic = "fisher"
)

elim_BP_Virginvslactate<-topGO::runTest(
    BP_L_virginvslactate,
    algorithm ="elim",
    statistic = "fisher"
)

elim_BP_Virginvspregnant<-topGO::runTest(
    BP_L_virginvspregnant,
    algorithm ="elim",
    statistic = "fisher"
)

# merge topGO results
BP_sResults<-ViSEAGO::merge_enrich_terms(
    Input=list(
        Pregnantvslactate=c("BP_Pregnantvslactate","elim_BP_Pregnantvslactate"),
        Virginvslactate=c("BP_Virginvslactate","elim_BP_Virginvslactate"),
        Virginvspregnant=c("BP_Virginvspregnant","elim_BP_Virginvspregnant")
    )
)

## End(Not run)

## fgsea analysis

# load gene identifiants and padj test results from Differential Analysis complete tables
PregnantvsLactate<-data.table::fread(
    system.file(
        "extdata/data/input",
        "pregnantvslactate.complete.txt",
        package = "ViSEAGO"
    ),
    select = c("Id","padj")
)

VirginvsLactate<-data.table::fread(
    system.file(
        "extdata/data/input",
        "virginvslactate.complete.txt",
        package = "ViSEAGO"
   ),
   select = c("Id","padj")
)

VirginvsPregnant<-data.table::fread(
    system.file(
       "extdata/data/input",
       "virginvspregnant.complete.txt",
        package = "ViSEAGO"
    ),
    select = c("Id","padj")
)

# rank Id based on statistical value (padj)
PregnantvsLactate<-data.table::setorder(PregnantvsLactate,padj)

VirginvsLactate<-data.table::setorder(VirginvsLactate,padj)

VirginvsPregnant<-data.table::setorder(VirginvsPregnant,padj)

## Not run: 
# connect to Bioconductor
Bioconductor<-ViSEAGO::Bioconductor2GO()

# load GO annotations from Bioconductor
myGENE2GO<-ViSEAGO::annotate(
    "org.Mm.eg.db",
    Bioconductor
)

# perform fgseaMultilevel tests
BP_PregnantvsLactate<-runfgsea(
    geneSel=PregnantvsLactate,
    gene2GO=myGENE2GO, 
    ont="BP",
    params = list(
        scoreType = "pos",
        minSize=5
    )
)

BP_VirginvsLactate<-runfgsea(
    geneSel=VirginvsLactate,
    gene2GO=myGENE2GO, 
    ont="BP",
    params = list(
        scoreType = "pos",
        minSize=5
    )
)

BP_VirginvsPregnant<-runfgsea(
    geneSel=VirginvsPregnant,
    gene2GO=myGENE2GO, 
    ont="BP",
    params = list(
        scoreType = "pos",
        minSize=5
    )
)

# merge fgsea results
BP_sResults<-merge_enrich_terms(
    cutoff=0.01,
    Input=list(
        PregnantvsLactate="BP_PregnantvsLactate",
        VirginvsLactate="BP_VirginvsLactate",
        VirginvsPregnant="BP_VirginvsPregnant"
    )
)

## End(Not run)

[Package ViSEAGO version 1.8.0 Index]