enrichedIn: This function builds a cross-tabulation of enriched (TRUE)...

View source: R/enrichedIn.R

enrichedInR Documentation

This function builds a cross-tabulation of enriched (TRUE) and non-enriched (FALSE) GO terms vs. gene lists

Description

This function builds a cross-tabulation of enriched (TRUE) and non-enriched (FALSE) GO terms vs. gene lists

Usage

enrichedIn(x, ...)

## Default S3 method:
enrichedIn(
  x,
  geneUniverse,
  orgPackg,
  onto,
  GOLevel,
  pAdjustMeth = "BH",
  pvalCutoff = 0.01,
  qvalCutoff = 0.05,
  parallel = FALSE,
  nOfCores = 1,
  ...
)

## S3 method for class 'character'
enrichedIn(
  x,
  geneUniverse,
  orgPackg,
  onto,
  GOLevel,
  pAdjustMeth = "BH",
  pvalCutoff = 0.01,
  qvalCutoff = 0.05,
  parallel = FALSE,
  nOfCores = 1,
  ...
)

## S3 method for class 'list'
enrichedIn(
  x,
  geneUniverse,
  orgPackg,
  onto,
  GOLevel,
  pAdjustMeth = "BH",
  pvalCutoff = 0.01,
  qvalCutoff = 0.05,
  parallel = FALSE,
  nOfCores = min(detectCores() - 1, length(x)),
  ...
)

Arguments

x

either an object of class "character" (or coerzable to "character") or "list". In the "character" interface, these values should represent Entrez gene (or, in general, feature) identifiers. In the "list" interface, each element of the list must be a "character" vector of Entrez identifiers

...

Additional parameters

geneUniverse

character vector containing the universe of genes from where gene lists have been extracted. This vector must be obtained from the annotation package declared in orgPackg. For more details see README File.

orgPackg

A string with the name of the genomic annotation package corresponding to a specific species to be analyzed, which must be previously installed and activated. For more details see README File.

onto

string describing the ontology. Belongs to c('BP', 'MF', 'CC')

GOLevel

GO level, an integer

pAdjustMeth

string describing the adjust method. Belongs to c('BH', 'BY', 'Bonf')

pvalCutoff

adjusted pvalue cutoff on enrichment tests to report

qvalCutoff

qvalue cutoff on enrichment tests to report as significant. Tests must pass i) pvalueCutoff on unadjusted pvalues, ii) pvalueCutoff on adjusted pvalues and iii) qvalueCutoff on qvalues to be reported

parallel

Logical. Only in "list" interface. Defaults to FALSE but put it at TRUE for parallel computation

nOfCores

Number of cores for parallel computations. Only in "list" interface

Details

The arguments 'parallel' and 'nOfCores' are ignored in the 'default' and "character" interfaces because (in the present implementation) parallelisation is only applied to repeated calls to function 'clusterProfiler::enrichGO' which, in turn, does not provide for the possibility of parallelisation. They only apply to the "list" interface.

Value

In the "character" interface, a length k vector of TRUE/FALSE values corresponding to enrichment or not, where k stands for the total number of GO terms at level 'GOLev' in ontology 'onto'. In the "list" interface, a boolean matrix of TRUE/FALSE values indicating enrichment or not, with k rows and s columns, where k corresponds to the total number of GO terms at level 'GOLev' in ontology 'onto' and s corresponds to the length of "list" 'x'.

Methods (by class)

  • enrichedIn(default): S3 default method

  • enrichedIn(character): S3 method for class "character"

  • enrichedIn(list): S3 method for class "list"

Examples

# Obtaining ENTREZ identifiers for the gene universe of humans:
library(org.Hs.eg.db)
humanEntrezIDs <- keys(org.Hs.eg.db, keytype = "ENTREZID")

# Gene lists to be explored for enrichment:
data(allOncoGeneLists)
?allOncoGeneLists

# Computing the cross table:
enrichd <- enrichedIn(allOncoGeneLists[["Vogelstein"]],
                      geneUniverse = humanEntrezIDs, orgPackg = "org.Hs.eg.db",
                      onto = "MF", GOLevel = 6)
enrichd

# Cross table of enriched GO terms (GO ontology MF, level 6) for all gene 
# lists in 'allOncoGeneLists':
enrichedAllOncoMF.6 <- enrichedIn(allOncoGeneLists,
                          geneUniverse = humanEntrezIDs, orgPackg = "org.Hs.eg.db",
                          onto = "MF", GOLevel = 6)
enrichedAllOncoMF.6

pablof1988/goSorensen documentation built on July 31, 2024, 10:26 p.m.