Analyze with GREAT

Zuguang Gu ( z.gu@dkfz.de )

2019-09-13

Note: On Aug 19 2019 GREAT released version 4 where it supports hg38 genome and removes some ontologies such pathways. submitGreatJob() still takes hg19 as default. hg38 can be specified by the species = "hg38" argument. To use the older versions such as 3.0.0, specify as submitGreatJob(..., version = "3.0.0").

GREAT (Genomic Regions Enrichment of Annotations Tool) is a popular web-based tool to associate biological functions to genomic regions. The rGREAT package makes GREAT anlaysis automatic by first constructing a HTTP POST request according to user’s input and retrieving results from GREAT web server afterwards.

Load the package:

library(rGREAT)

The input data is either a GRanges object or a BED-format data frame, no matter it is sorted or not. In following example, we use a data frame which is randomly generated.

set.seed(123)
bed = circlize::generateRandomBed(nr = 1000, nc = 0)
bed[1:2, ]
##    chr   start      end
## 1 chr1 7634457  9204434
## 2 chr1 9853594 10435028

Submit genomic regions by submitGreatJob(). Before submitting, genomic regions will be sorted and overlapping regions will be merged.

The returned variable job is a GreatJob class instance which can be used to retrieve results from GREAT server and stored results which are already downloaded.

job = submitGreatJob(bed)

You can get the summary of your job by directly calling job variable.

job
## Submit time: 2019-09-13 20:59:47 
## Version: 4.0.4 
## Species: hg19 
## Inputs: 1005 regions
## Background: wholeGenome 
## Model: Basal plus extension 
##   Proximal: 5 kb upstream, 1 kb downstream,
##   plus Distal: up to 1000 kb
## Include curated regulatory domains
## 
## Enrichment tables for following ontologies have been downloaded:
##   None

More parameters can be set for the job:

job = submitGreatJob(bed, species = "mm9")
job = submitGreatJob(bed, bg, species = "mm9")
job = submitGreatJob(bed, adv_upstream = 10, adv_downstream = 2, adv_span = 2000)
job = submitGreatJob(bed, rule = "twoClosest", adv_twoDistance = 2000)
job = submitGreatJob(bed, rule = "oneClosest", adv_oneDistance = 2000)

Also you can choose different versions of GREAT for the analysis.

job = submitGreatJob(bed, version = "3.0")
job = submitGreatJob(bed, version = "2.0")

Available parameters are (following content is copied from GREAT website):

With job, we can now retrieve results from GREAT. The first and the primary results are the tables which contain enrichment statistics for the analysis. By default it will retrieve results from three GO Ontologies and all pathway ontologies. All tables contains statistics for all terms no matter they are significant or not. Users can then make filtering yb self-defined cutoff.

There is a column for adjusted p-values by “BH” method. Other p-value adjustment methods can be applied by p.adjust().

The returned value of getEnrichmentTables() is a list of data frames in which each one corresponds to tables for single ontology. The structure of data frames are same as the tables on GREAT website.

tb = getEnrichmentTables(job)
names(tb)
## [1] "GO Molecular Function" "GO Biological Process" "GO Cellular Component"
tb[[1]][1:2, ]
##           ID                                         name Binom_Genome_Fraction Binom_Expected
## 1 GO:0016798 hydrolase activity, acting on glycosyl bonds          0.0100418100     10.0920200
## 2 GO:0051787                    misfolded protein binding          0.0007646973      0.7685208
##   Binom_Observed_Region_Hits Binom_Fold_Enrichment Binom_Region_Set_Coverage Binom_Raw_PValue
## 1                         20              1.981763                 0.0199005      0.003633144
## 2                          4              5.204804                 0.0039801      0.007894598
##   Binom_Adjp_BH Hyper_Total_Genes Hyper_Expected Hyper_Observed_Gene_Hits Hyper_Fold_Enrichment
## 1             1               127      10.756210                       17              1.580482
## 2             1                12       1.016335                        4              3.935710
##   Hyper_Gene_Set_Coverage Hyper_Term_Gene_Coverage Hyper_Raw_PValue Hyper_Adjp_BH
## 1             0.010821130                0.1338583       0.03955741     1.0000000
## 2             0.002546149                0.3333333       0.01458949     0.6280924

Information stored in job will be updated after retrieving enrichment tables.

job
## Submit time: 2019-09-13 20:59:47 
## Version: 4.0.4 
## Species: hg19 
## Inputs: 1005 regions
## Background: wholeGenome 
## Model: Basal plus extension 
##   Proximal: 5 kb upstream, 1 kb downstream,
##   plus Distal: up to 1000 kb
## Include curated regulatory domains
## 
## Enrichment tables for following ontologies have been downloaded:
##   GO Biological Process
##   GO Cellular Component
##   GO Molecular Function

You can get results by either specifying the ontologies or by the pre-defined categories (categories already contains pre-defined sets of ontologies):

tb = getEnrichmentTables(job, ontology = c("GO Molecular Function", "BioCyc Pathway"))
tb = getEnrichmentTables(job, category = c("GO"))

All available ontology names for given species can be get by availableOntologies() and all available ontology categories can be get by availableCategories(). Here you do not need to provide species information because job already contains it.

availableOntologies(job)
## [1] "GO Molecular Function"     "GO Biological Process"     "GO Cellular Component"    
## [4] "Mouse Phenotype"           "Mouse Phenotype Single KO" "Human Phenotype"          
## [7] "Ensembl Genes"
availableCategories(job)
## [1] "GO"        "Phenotype" "Genes"
availableOntologies(job, category = "GO")
## [1] "GO Molecular Function" "GO Biological Process" "GO Cellular Component"

Association between genomic regions and genes can be get by plotRegionGeneAssociationGraphs(). The function will make the three plots which are same as on GREAT website and returns a GRanges object which contains the gene-region associations.

res = plotRegionGeneAssociationGraphs(job)

res[1:2, ]
## GRanges object with 2 ranges and 2 metadata columns:
##       seqnames          ranges strand |        gene   distTSS
##          <Rle>       <IRanges>  <Rle> | <character> <numeric>
##   [1]     chr1 7634457-9204434      * |     SLC45A1     41560
##   [2]     chr1 7634457-9204434      * |        RERE    458078
##   -------
##   seqinfo: 24 sequences from an unspecified genome; no seqlengths

For those regions that are not associated with any genes under current settings, the corresponding gene and distTSS columns will be NA.

You can also choose only plotting one of the three figures.

plotRegionGeneAssociationGraphs(job, type = 1)

By specifying ontology and term ID, you can get the association in a certain term. Here the term ID is from the first column of the data frame which is returned by getEnrichmentTables().

res = plotRegionGeneAssociationGraphs(job, ontology = "GO Molecular Function",
    termID = "GO:0004984")
## The webpage for 'GOMolecularFunction:GO:0004984' is available at:
##   http://great.stanford.edu/public-4.0.4/cgi-bin/showTermDetails.php?termId=GO:0004984&ontoName=GOMolecularFunction&ontoUiName=GO Molecular Function&sessionName=20190913-public-4.0.4-tW0RI7&species=hg19&foreName=file75b5298cd139.gz&backName=&table=region

res[1:2, ]
## GRanges object with 2 ranges and 2 metadata columns:
##       seqnames              ranges strand |        gene   distTSS
##          <Rle>           <IRanges>  <Rle> | <character> <numeric>
##   [1]     chr1 247871555-248931674      * |       OR2M4      -616
##   [2]     chr3   96368632-99050767      * |      OR5AC2    -96317
##   -------
##   seqinfo: 7 sequences from an unspecified genome; no seqlengths

Session info

sessionInfo()
## R version 3.6.1 (2019-07-05)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.9-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.9-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C               LC_TIME=en_US.UTF-8       
##  [4] LC_COLLATE=C               LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                  LC_ADDRESS=C              
## [10] LC_TELEPHONE=C             LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] parallel  stats4    stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] rGREAT_1.16.1        GenomicRanges_1.36.1 GenomeInfoDb_1.20.0  IRanges_2.18.2      
## [5] S4Vectors_0.22.1     BiocGenerics_0.30.0  knitr_1.24          
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_1.0.2             XVector_0.24.0         magrittr_1.5           zlibbioc_1.30.0       
##  [5] colorspace_1.4-1       rjson_0.2.20           stringr_1.4.0          tools_3.6.1           
##  [9] grid_3.6.1             circlize_0.4.8         xfun_0.9               htmltools_0.3.6       
## [13] yaml_2.2.0             digest_0.6.20          GenomeInfoDbData_1.2.1 GlobalOptions_0.1.0   
## [17] bitops_1.0-6           RCurl_1.95-4.12        shape_1.4.4            evaluate_0.14         
## [21] rmarkdown_1.15         stringi_1.4.3          compiler_3.6.1         GetoptLong_0.1.7