Abstract
OGRE calculates overlap between user defined annotated genomic region datasets. Any regions can be supplied such as public annotations (genes), genetic variation (SNPs, mutations), regulatory elements (TFBS, promoters, CpG islands) and basically all types of NGS output from sequencing experiments. After overlap calculation, key numbers help analyse the extend of overlaps which can also be visualized at a genomic level. To start OGRE’s GUI use function SHREC() in your R console. Find additional information and tutorials on github. OGRE package version: 1.9.0
Install OGRE using Bioconductor’s package installer.
if(!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("OGRE")
Load the OGRE package:
To start up OGRE you have to generate an OGREDataSet
that is used to store your datasets and additional information about the analysis that you are conducting. Query and subjects files can be conveniently stored in their own folders as GenomicRanges objects in form of stored .rds / .RDS files. We point OGRE to the correct location by supplying a path for each folder with the character vectors queryFolder
and subjectFolder
. In this vignette we are using lightweight query and subject example data sets to show OGRE’s functionality.
myQueryFolder <- file.path(system.file('extdata', package = 'OGRE'),"query")
mySubjectFolder <- file.path(system.file('extdata', package = 'OGRE'),"subject")
myOGRE <- OGREDataSetFromDir(queryFolder=myQueryFolder,
subjectFolder=mySubjectFolder)
## Initializing OGREDataSet...
By monitoring OGRE’s metadata information you can make sure the input paths you supplied are stored correctly.
## $queryFolder
## [1] "/tmp/RtmpxdWPsr/Rinstc9c49398afabd/OGRE/extdata/query"
##
## $subjectFolder
## [1] "/tmp/RtmpxdWPsr/Rinstc9c49398afabd/OGRE/extdata/subject"
##
## $outputFolder
## [1] "/tmp/RtmpxdWPsr/Rinstc9c49398afabd/OGRE/extdata/output"
##
## $gvizPlotsFolder
## [1] "/tmp/RtmpxdWPsr/Rinstc9c49398afabd/OGRE/extdata/gvizPlots"
##
## $summaryDT
## list()
##
## $itracks
## list()
Query and subject datasets are read by loadAnnotations()
and stored in the OGREDataSet
as GRanges
objects. We are going to read in the following example datasets:
## Reading query dataset...
## Reading subject datasets...
OGRE uses your dataset file names to label query and subjects internally, we can check these names by using the names()
function since every OGREDataSet
is a GRangesList
.
## [1] "genes" "CGI" "TFBS"
Let’s have a look at the stored datasets:
## GRangesList object of length 3:
## $genes
## GRanges object with 242 ranges and 3 metadata columns:
## seqnames ranges strand | ID name
## <Rle> <IRanges> <Rle> | <character> <character>
## [1] 21 10906201-11029719 - | ENSG00000166157 TPTE
## [2] 21 14741931-14745386 - | ENSG00000256715 AL050302.1
## [3] 21 14982498-15013906 + | ENSG00000166351 POTED
## [4] 21 15051621-15053459 - | ENSG00000269011 AL050303.1
## [5] 21 15481134-15583166 - | ENSG00000188992 LIPI
## ... ... ... ... . ... ...
## [238] 21 47720095-47743789 - | ENSG00000160298 C21orf58
## [239] 21 47744036-47865682 + | ENSG00000160299 PCNT
## [240] 21 47878812-47989926 + | ENSG00000160305 DIP2A
## [241] 21 48018875-48025121 - | ENSG00000160307 S100B
## [242] 21 48055079-48085036 + | ENSG00000160310 PRMT2
## score
## <numeric>
## [1] NA
## [2] NA
## [3] NA
## [4] NA
## [5] NA
## ... ...
## [238] NA
## [239] NA
## [240] NA
## [241] NA
## [242] NA
## -------
## seqinfo: 25 sequences (1 circular) from hg19 genome
##
## $CGI
## GRanges object with 365 ranges and 3 metadata columns:
## seqnames ranges strand | ID name score
## <Rle> <IRanges> <Rle> | <character> <character> <numeric>
## [1] 21 9437273-9439473 * | 26635 CpG:_285 NA
## [2] 21 9483486-9484663 * | 26636 CpG:_165 NA
## [3] 21 9647867-9648116 * | 26637 CpG:_18 NA
## [4] 21 9708936-9709231 * | 26638 CpG:_31 NA
## [5] 21 9825443-9826296 * | 26639 CpG:_120 NA
## ... ... ... ... . ... ... ...
## [361] 21 48018543-48018791 * | 26995 CpG:_21 NA
## [362] 21 48055200-48056060 * | 26996 CpG:_88 NA
## [363] 21 48068518-48068808 * | 26997 CpG:_24 NA
## [364] 21 48081242-48081849 * | 26998 CpG:_55 NA
## [365] 21 48087201-48088106 * | 26999 CpG:_93 NA
## -------
## seqinfo: 25 sequences (1 circular) from hg19 genome
##
## $TFBS
## GRanges object with 48761 ranges and 3 metadata columns:
## seqnames ranges strand | ID name
## <Rle> <IRanges> <Rle> | <character> <character>
## [1] 21 29884415-29884427 + | GATA1.85108 GATA1_04
## [2] 21 46923766-46923780 + | CDP.81529 CDP_02
## [3] 21 9491627-9491638 - | HFH1.46541 HFH1_01
## [4] 21 9491706-9491725 - | PPARA.24892 PPARA_01
## [5] 21 9491792-9491815 + | GFI1.35413 GFI1_01
## ... ... ... ... . ... ...
## [48757] 21 48083381-48083404 + | STAT5A.43326 STAT5A_02
## [48758] 21 48083400-48083419 + | ARNT.19751 ARNT_02
## [48759] 21 48084826-48084841 + | BRN2.40426 BRN2_01
## [48760] 21 48084830-48084847 + | FOXJ2.121681 FOXJ2_01
## [48761] 21 48084834-48084845 + | NKX3A.47953 NKX3A_01
## score
## <numeric>
## [1] 891
## [2] 831
## [3] 865
## [4] 757
## [5] 817
## ... ...
## [48757] 751
## [48758] 792
## [48759] 803
## [48760] 889
## [48761] 851
## -------
## seqinfo: 25 sequences (1 circular) from hg19 genome
To find overlaps between your query and subject datasets we call fOverlaps()
. Internally OGRE makes use of the GenomicRanges
package to calculate full and partial overlap as schematically shown.
Any existing subject - query hits are then listed in detailDT
and stored as a data.table
.
## queryID queryType subjID subjType queryChr queryStart queryEnd
## <char> <char> <char> <char> <char> <int> <int>
## 1: ENSG00000166157 genes 26649 CGI 21 10906201 11029719
## 2: ENSG00000269011 genes 26654 CGI 21 15051621 15053459
## queryStrand subjChr subjStart subjEnd subjStrand overlapWidth overlapRatio
## <char> <char> <int> <int> <char> <int> <num>
## 1: - 21 10989914 10991413 * 1500 0.01214388
## 2: - 21 15052411 15052644 * 234 0.12724307
The summary plot provides us with useful information about the number of overlaps between your datasets.
Using the Gviz
visualization each query can be displayed with all overlapping subject elements. Choose labels for all region tracks by supplying a trackRegionLabels
vector. Plots are stored in the same location as your dataset files.
myOGRE <- gvizPlot(myOGRE,"ENSG00000142168",showPlot = TRUE,
trackRegionLabels = setNames(c("name","name"),c("genes","CGI")))
## Plotting query: ENSG00000142168
The overlap distribution can be generated with summarizeOverlap(myOGRE)
and outputs a table with informative statistics such as minimum, lower quantile, mean, median, upper quantile, and maximum number of overlaps per region and per dataset. Overlap distribution can also be displayed as histograms using plotHist(myOGRE)
and accessed by metadata(myOGRE)$hist
and metadata(myOGRE)$summaryDT
. Two tables / plots are generated. The first one showing numbers for regions with and without overlap and the second one showing numbers only for regions with overlap by excluding all others. Next, we generate an histogram with the number of TFBS per gene (x-axis, log scale) and the TFBS frequency (y-axis). When focusing only on regions with overlap, we see that genes have on average (median) 54 TFBS overlaps (black dashed line).
## $includes0
## CGI TFBS
## Min. 0.000000 0.0000
## 1st Qu. 0.000000 8.0000
## Median 1.000000 36.0000
## Mean 1.210744 119.6116
## 3rd Qu. 1.750000 129.7500
## Max. 14.000000 3136.0000
##
## $excludes0
## CGI TFBS
## Min. 1.00000 1.0000
## 1st Qu. 1.00000 15.0000
## Median 1.00000 54.0000
## Mean 2.02069 139.8357
## 3rd Qu. 2.00000 159.5000
## Max. 14.00000 3136.0000
## NA's 97.00000 35.0000
## `stat_bin()` using `bins = 30`. Pick better value with `binwidth`.
It is possible to create an average coverage profile of all gene-TFBS overlaps, split in 100 bins, which represent gene bodies of all 242 genes. Both, forward and reverse coding genes are arranged on the x-Axis and peaks indicate an TFBS overlap enrichment. Overlap coverage is calculated as the sum of all gene TFBS overlaps in 5’-3’direction. Generated plots can be accessed by metadata(myOGRE)$covPlot$TFBS
and the resulting profile shows an accumulation of TFBS around gene start and end positions.
## Generating coverage plot(s), this might take a while...
## Excluding regions with nucleotides<nbin
## `geom_smooth()` using method = 'loess' and formula = 'y ~ x'
AnnotationHub offers a wide range of annotated datasets which can be manually aquired but need some parsing to work with OGRE as detailed in vignette section Frequently Asked Questions(FAQ). For convenience addDataSetFromHub()
adds one of the predefined human datasets of listPredefinedDataSets()
to an OGREDataSet. Those are taken from AnnotationHub and are ready to use for OGRE. We start by creating an empty OGREDataSet and attaching one dataset after another, whereby one query and two subjects are added. The datasets are now ready for further analysis.
myOGRE <- OGREDataSet()
listPredefinedDataSets()
myOGRE <- addDataSetFromHub(myOGRE,"protCodingGenes","query")
myOGRE <- addDataSetFromHub(myOGRE,"CGI","subject")
myOGRE <- addDataSetFromHub(myOGRE,"TFBS","subject")
names(myOGRE)
As you can see, the three datasets proteinCodingGenes, CGI and TFBS are stored within OGRE. You can then continue with overlap analysis using fOverlaps()
.
To offer more flexibility addGRanges()
enables the user to attach additional datasets to OGRE in form of GenomicRanges objects. Again we start by creating an empty OGREDataSet and generate an example GenomicRanges object which is then added to your OGREDataSet either as “query” or “subject”.
Use AnnotationHub()
to connect to AnnotationHub. Each dataset is stored under a unique ID and can be accessed in a list like fashion i.e. aH[["AH5086"]]
. Queries like c("GRanges","Homo sapiens", "CpG")
enable browsing through datasets. In this case we are searching for human CpG islands ranges stored as GenomicRanges objects. For more information refer to ?AnnotationHub
To make those datasets compatible with OGRE additional parsing is needed as stated in How to add custom GenomicRanges datasets?
Any GenomicRanges datasets can be added that fulfill basic compatibility requirements. GenomicRanges objects must:
Use GenomeInfoDb::genome()
on any GenomicRanges object to get/set genome information
Use GenomeInfoDb::seqinfo()
on any GenomicRanges object to get/set chromosome information
Use S4Vectors::mcols()
on any GenomicRanges object to get/set metadata information
Datasets from external sources are often stored as .gff (v2&v3) files. Once those files exist in the query or subject folder and their attribute columns contain “ID” and “name” information, OGRE tries to load them. Working example .gff files can be found on OGRE’s github page in folder: inst/extdata/gffTest.
Datasets stored as tabular files like .csv or .bed may need some preprocessing for them work with OGRE. We recommend reading them in with read.table()
or data.table::fread()
to obtain a data frame object. After making sure the dataset complies with the requirements in section How to add custom GenomicRanges datasets?, GenomicRanges::makeGRangesFromDataFrame()
offers a convenient way to generate GenomicRanges object from data frames.
Both, partial overlap, where only a part of two (or more) regions are overlapping and complete overlap, where one region is completely overlapped by another, are reported.
OGRE automatically infers dataset names based on your file names. You can either change your file names before you start OGRE or you can use names(myOGRE) <- c("NewName1", "NewName2","...")
after you read in your datasets.
## R version 4.4.0 RC (2024-04-16 r86468)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## time zone: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] OGRE_1.9.0 S4Vectors_0.43.0 BiocGenerics_0.51.0
##
## loaded via a namespace (and not attached):
## [1] RColorBrewer_1.1-3 rstudioapi_0.16.0
## [3] jsonlite_1.8.8 magrittr_2.0.3
## [5] GenomicFeatures_1.57.0 farver_2.1.1
## [7] rmarkdown_2.26 ragg_1.3.0
## [9] fs_1.6.4 BiocIO_1.15.0
## [11] zlibbioc_1.51.0 vctrs_0.6.5
## [13] memoise_2.0.1 Rsamtools_2.21.0
## [15] RCurl_1.98-1.14 base64enc_0.1-3
## [17] htmltools_0.5.8.1 S4Arrays_1.5.0
## [19] progress_1.2.3 AnnotationHub_3.13.0
## [21] curl_5.2.1 SparseArray_1.5.0
## [23] Formula_1.2-5 shinyFiles_0.9.3
## [25] sass_0.4.9 bslib_0.7.0
## [27] htmlwidgets_1.6.4 Gviz_1.49.0
## [29] httr2_1.0.1 cachem_1.0.8
## [31] GenomicAlignments_1.41.0 mime_0.12
## [33] lifecycle_1.0.4 pkgconfig_2.0.3
## [35] Matrix_1.7-0 R6_2.5.1
## [37] fastmap_1.1.1 GenomeInfoDbData_1.2.12
## [39] MatrixGenerics_1.17.0 shiny_1.8.1.1
## [41] digest_0.6.35 colorspace_2.1-0
## [43] AnnotationDbi_1.67.0 textshaping_0.3.7
## [45] Hmisc_5.1-2 GenomicRanges_1.57.0
## [47] RSQLite_2.3.6 labeling_0.4.3
## [49] filelock_1.0.3 fansi_1.0.6
## [51] mgcv_1.9-1 httr_1.4.7
## [53] abind_1.4-5 compiler_4.4.0
## [55] withr_3.0.0 bit64_4.0.5
## [57] htmlTable_2.4.2 backports_1.4.1
## [59] BiocParallel_1.39.0 DBI_1.2.2
## [61] highr_0.10 biomaRt_2.61.0
## [63] rappdirs_0.3.3 DelayedArray_0.31.0
## [65] rjson_0.2.21 tools_4.4.0
## [67] foreign_0.8-86 httpuv_1.6.15
## [69] nnet_7.3-19 glue_1.7.0
## [71] restfulr_0.0.15 nlme_3.1-164
## [73] promises_1.3.0 grid_4.4.0
## [75] checkmate_2.3.1 cluster_2.1.6
## [77] generics_0.1.3 gtable_0.3.5
## [79] BSgenome_1.73.0 shinyBS_0.61.1
## [81] tidyr_1.3.1 ensembldb_2.29.0
## [83] data.table_1.15.4 hms_1.1.3
## [85] xml2_1.3.6 utf8_1.2.4
## [87] XVector_0.45.0 BiocVersion_3.20.0
## [89] pillar_1.9.0 stringr_1.5.1
## [91] later_1.3.2 splines_4.4.0
## [93] dplyr_1.1.4 BiocFileCache_2.13.0
## [95] lattice_0.22-6 rtracklayer_1.65.0
## [97] bit_4.0.5 deldir_2.0-4
## [99] biovizBase_1.53.0 tidyselect_1.2.1
## [101] Biostrings_2.73.0 knitr_1.46
## [103] gridExtra_2.3 IRanges_2.39.0
## [105] ProtGenerics_1.37.0 SummarizedExperiment_1.35.0
## [107] shinydashboard_0.7.2 xfun_0.43
## [109] Biobase_2.65.0 matrixStats_1.3.0
## [111] DT_0.33 stringi_1.8.3
## [113] UCSC.utils_1.1.0 lazyeval_0.2.2
## [115] yaml_2.3.8 evaluate_0.23
## [117] codetools_0.2-20 interp_1.1-6
## [119] tibble_3.2.1 BiocManager_1.30.22
## [121] cli_3.6.2 rpart_4.1.23
## [123] systemfonts_1.0.6 xtable_1.8-4
## [125] munsell_0.5.1 jquerylib_0.1.4
## [127] dichromat_2.0-0.1 Rcpp_1.0.12
## [129] GenomeInfoDb_1.41.0 dbplyr_2.5.0
## [131] png_0.1-8 XML_3.99-0.16.1
## [133] parallel_4.4.0 assertthat_0.2.1
## [135] ggplot2_3.5.1 blob_1.2.4
## [137] prettyunits_1.2.0 latticeExtra_0.6-30
## [139] jpeg_0.1-10 AnnotationFilter_1.29.0
## [141] bitops_1.0-7 VariantAnnotation_1.51.0
## [143] scales_1.3.0 purrr_1.0.2
## [145] crayon_1.5.2 rlang_1.1.3
## [147] KEGGREST_1.45.0