Abstract
Visualize chromatin interactions along with annotation as track layers. The interactions can be compared by back to back heatmaps. The interactions can be plot as heatmap and links.The chromatin interactions is involved in precise quantitative and spatiotemporal control of gene expression. The development of high-throughput experimental techniques, such as HiC-seq, HiCAR-seq, and InTAC-seq, for analyzing both the higher-order structure of chromatin and the interactions between protein and their nearby and remote regulatory elements has been developed to reveal how gene expression is controlled in genome-wide.
The interaction data will be saved in the format of paired genome coordinates with the interaction score. The popular format are .validPairs
, .hic
, and .cool
. The trackViewer
package can be used to handle those data to plot the heatmap or the interaction links.
Plot chromatin interactions tracks as heatmap.
library(trackViewer)
library(InteractionSet)
gi <- readRDS(system.file("extdata", "nij.chr6.51120000.53200000.gi.rds", package="trackViewer"))
head(gi)
## GInteractions object with 6 interactions and 1 metadata column:
## seqnames1 ranges1 seqnames2 ranges2 | score
## <Rle> <IRanges> <Rle> <IRanges> | <numeric>
## [1] chr6 51120000-51160000 --- chr6 51120000-51160000 | 45.1227
## [2] chr6 51120000-51160000 --- chr6 51160000-51200000 | 35.0006
## [3] chr6 51120000-51160000 --- chr6 51200000-51240000 | 44.7322
## [4] chr6 51120000-51160000 --- chr6 51240000-51280000 | 29.3507
## [5] chr6 51120000-51160000 --- chr6 51280000-51320000 | 38.8417
## [6] chr6 51120000-51160000 --- chr6 51320000-51360000 | 31.7063
## -------
## regions: 53 ranges and 0 metadata columns
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
## hicexplorer:hicConvertFormat tool can be used to convert other formats into GInteractions
## eg: hicConvertFormat -m mESC_rep.hic --inputFormat hic --outputFormat cool -o mESC_rep.mcool
## hicConvertFormat -m mESC_rep.mcool::resolutions/10000 --inputFormat cool --outputFormat ginteractions -o mESC_rep.ginteractions --resolutions 10000
## please note that metadata:score is used for plot.
range <- GRanges("chr6", IRanges(51120000, 53200000))
tr <- gi2track(gi)
ctcf <- readRDS(system.file("extdata", "ctcf.sample.rds", package="trackViewer"))
#viewTracks(trackList(ctcf, tr, heightDist = c(1, 3)),
# gr=range, autoOptimizeStyle = TRUE)
## view the interaction data back to back.
## Please make sure the data are normalized.
gi2 <- gi
set.seed(123)
gi2$score <- gi$score + rnorm(length(gi), sd = sd(gi$score))
back2back <- gi2track(gi, gi2)
## change the color
setTrackStyleParam(back2back, "breaks",
c(seq(from=0, to=50, by=10), 200))
setTrackStyleParam(back2back, "color",
c("lightblue", "yellow", "red"))
## chang the lim of y-axis (by default, [0, 1])
setTrackStyleParam(back2back, "ylim", c(0, .5))
viewTracks(trackList(ctcf, back2back, heightDist=c(1, 5)),
gr=range, autoOptimizeStyle = TRUE)
Plot chromatin interactions track as links.
setTrackStyleParam(tr, "tracktype", "link")
setTrackStyleParam(tr, "breaks",
c(seq(from=0, to=50, by=10), 200))
setTrackStyleParam(tr, "color",
c("lightblue", "yellow", "red"))
## filter the links to simulate the real data
keep <- distance(tr$dat, tr$dat2) > 5e5 & tr$dat$score>20
tr$dat <- tr$dat[keep]
tr$dat2 <- tr$dat2[keep]
viewTracks(trackList(tr), gr=range, autoOptimizeStyle = TRUE)
To import interactions data from “.hic” (reference to the script of hic-straw and the documentation). The function importGInteractions
(trackViewer version>=1.27.6) can be used to import data from .hic
format file.
hic <- system.file("extdata", "test_chr22.hic", package = "trackViewer",
mustWork=TRUE)
if(.Platform$OS.type!="windows"){
importGInteractions(file=hic, format="hic",
ranges=GRanges("22", IRanges(50000000, 100000000)),
out = "GInteractions")
}
## GInteractions object with 70 interactions and 1 metadata column:
## seqnames1 ranges1 seqnames2 ranges2 | score
## <Rle> <IRanges> <Rle> <IRanges> | <numeric>
## [1] 22 50000001-50100000 --- 22 50000001-50100000 | 26
## [2] 22 50000001-50100000 --- 22 50100001-50200000 | 2
## [3] 22 50100001-50200000 --- 22 50100001-50200000 | 22
## [4] 22 50100001-50200000 --- 22 50200001-50300000 | 7
## [5] 22 50200001-50300000 --- 22 50200001-50300000 | 31
## ... ... ... ... ... ... . ...
## [66] 22 50400001-50500000 --- 22 51200001-51300000 | 1
## [67] 22 50500001-50600000 --- 22 51200001-51300000 | 2
## [68] 22 50800001-50900000 --- 22 51200001-51300000 | 2
## [69] 22 51100001-51200000 --- 22 51200001-51300000 | 3
## [70] 22 51200001-51300000 --- 22 51200001-51300000 | 5
## -------
## regions: 13 ranges and 0 metadata columns
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
Another widely used genomic interaction data format is .cool
files and the cooler index contains analyzed HiC data for hg19 and mm9 from many different sources. Those files can be used as data resources for visualizations and annotations (see ChIPpeakAnno::findEnhancers). The importGInteractions
function can also be used to import data from .cool
format (trackViewer version>=1.27.6).
cool <- system.file("extdata", "test.mcool", package = "trackViewer",
mustWork=TRUE)
importGInteractions(file=cool, format="cool",
resolution = 2,
ranges=GRanges("chr1", IRanges(10, 28)),
out = "GInteractions")
Different from most of the available tools, plotGInteractions try to plot the data with the 2D structure. The nodes indicate the region with interactions and the edges indicates the interactions. The size of the nodes are relative to the width of the region. The features could be the enhancers, promoters or genes. The enhancer and promoter are shown as points with symbol 11 and 13.
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(InteractionSet)
gi <- readRDS(system.file("extdata", "gi.rds", package="trackViewer"))
range <- GRanges("chr2", IRanges(234500000, 235000000))
feature.gr <- suppressMessages(genes(TxDb.Hsapiens.UCSC.hg19.knownGene))
feature.gr <- subsetByOverlaps(feature.gr, regions(gi))
feature.gr$col <- sample(1:7, length(feature.gr), replace=TRUE)
feature.gr$type <- sample(c("promoter", "enhancer", "gene"),
length(feature.gr), replace=TRUE,
prob=c(0.1, 0.2, 0.7))
plotGInteractions(gi, range, feature.gr)
sessionInfo()
R version 4.2.0 (2022-04-22) Platform: x86_64-pc-linux-gnu (64-bit) Running under: Ubuntu 20.04.4 LTS
Matrix products: default BLAS: /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
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
attached base packages: [1] grid stats4 stats graphics grDevices utils datasets [8] methods base
other attached packages: [1] InteractionSet_1.24.0
[2] httr_1.4.3
[3] VariantAnnotation_1.42.1
[4] Rsamtools_2.12.0
[5] Biostrings_2.64.0
[6] XVector_0.36.0
[7] SummarizedExperiment_1.26.1
[8] MatrixGenerics_1.8.0
[9] matrixStats_0.62.0
[10] org.Hs.eg.db_3.15.0
[11] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2 [12] GenomicFeatures_1.48.1
[13] AnnotationDbi_1.58.0
[14] Biobase_2.56.0
[15] Gviz_1.40.1
[16] rtracklayer_1.56.0
[17] trackViewer_1.32.1
[18] Rcpp_1.0.8.3
[19] GenomicRanges_1.48.0
[20] GenomeInfoDb_1.32.2
[21] IRanges_2.30.0
[22] S4Vectors_0.34.0
[23] BiocGenerics_0.42.0
loaded via a namespace (and not attached): [1] colorspace_2.0-3 rjson_0.2.21 ellipsis_0.3.2
[4] biovizBase_1.44.0 htmlTable_2.4.0 base64enc_0.1-3
[7] dichromat_2.0-0.1 rstudioapi_0.13 bit64_4.0.5
[10] fansi_1.0.3 xml2_1.3.3 splines_4.2.0
[13] cachem_1.0.6 knitr_1.39 Formula_1.2-4
[16] jsonlite_1.8.0 cluster_2.1.3 dbplyr_2.1.1
[19] png_0.1-7 grImport_0.9-5 graph_1.74.0
[22] BiocManager_1.30.17 compiler_4.2.0 backports_1.4.1
[25] lazyeval_0.2.2 assertthat_0.2.1 Matrix_1.4-1
[28] fastmap_1.1.0 cli_3.3.0 htmltools_0.5.2
[31] prettyunits_1.1.1 tools_4.2.0 gtable_0.3.0
[34] glue_1.6.2 GenomeInfoDbData_1.2.8 dplyr_1.0.9
[37] rappdirs_0.3.3 jquerylib_0.1.4 rhdf5filters_1.8.0
[40] vctrs_0.4.1 xfun_0.31 stringr_1.4.0
[43] lifecycle_1.0.1 ensembldb_2.20.1 restfulr_0.0.13
[46] XML_3.99-0.9 zlibbioc_1.42.0 scales_1.2.0
[49] BiocStyle_2.24.0 BSgenome_1.64.0 ProtGenerics_1.28.0
[52] hms_1.1.1 parallel_4.2.0 rhdf5_2.40.0
[55] AnnotationFilter_1.20.0 RColorBrewer_1.1-3 yaml_2.3.5
[58] curl_4.3.2 memoise_2.0.1 gridExtra_2.3
[61] ggplot2_3.3.6 sass_0.4.1 biomaRt_2.52.0
[64] rpart_4.1.16 latticeExtra_0.6-29 stringi_1.7.6
[67] RSQLite_2.2.14 highr_0.9 BiocIO_1.6.0
[70] plotrix_3.8-2 checkmate_2.1.0 filelock_1.0.2
[73] BiocParallel_1.30.2 rlang_1.0.2 pkgconfig_2.0.3
[76] bitops_1.0-7 evaluate_0.15 lattice_0.20-45
[79] Rhdf5lib_1.18.2 purrr_0.3.4 GenomicAlignments_1.32.0 [82] htmlwidgets_1.5.4 bit_4.0.4 tidyselect_1.1.2
[85] magrittr_2.0.3 R6_2.5.1 generics_0.1.2
[88] Hmisc_4.7-0 DelayedArray_0.22.0 DBI_1.1.2
[91] pillar_1.7.0 foreign_0.8-82 survival_3.3-1
[94] KEGGREST_1.36.0 RCurl_1.98-1.6 nnet_7.3-17
[97] tibble_3.1.7 crayon_1.5.1 utf8_1.2.2
[100] BiocFileCache_2.4.0 rmarkdown_2.14 jpeg_0.1-9
[103] progress_1.2.2 data.table_1.14.2 Rgraphviz_2.40.0
[106] blob_1.2.3 digest_0.6.29 munsell_0.5.0
[109] bslib_0.3.1