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.
gi$border_color <- NA ## highlight some regions
gi$border_color[sample(seq_along(gi), 20)] <- sample(1:7, 20, replace=TRUE)
## The TADs will be drawn as lines at points start(first), center point, end(second).
tads <- GInteractions(
GRanges("chr6",
IRanges(c(51130001, 51130001, 51450001, 52210001), width = 20000)),
GRanges("chr6",
IRanges(c(51530001, 52170001, 52210001, 53210001), width = 20000)))
range <- GRanges("chr6", IRanges(51120000, 53200000))
heatmap <- gi2track(gi)
ctcf <- readRDS(system.file("extdata", "ctcf.sample.rds", package="trackViewer"))
viewTracks(trackList(ctcf, heatmap, heightDist = c(1, 3)),
gr=range, autoOptimizeStyle = TRUE)
## add TAD information
addInteractionAnnotation(tads, "heatmap", grid.lines, gp=gpar(col="#E69F00", lwd=3, lty=3))
## add highlight interested regions
gi_sub <- gi[order(gi$score, decreasing = TRUE)]
gi_sub <- head(gi_sub[distance(first(gi_sub), second(gi_sub))>200000], n=5)
start(regions(gi_sub)) <- start(regions(gi_sub))-40000
end(regions(gi_sub)) <- end(regions(gi_sub))+40000
addInteractionAnnotation(gi_sub, "heatmap", grid.polygon, gp=gpar(col="red", lwd=2, lty=2, fill=NA))
## add interesting anchor at giving coordinate.
addInteractionAnnotation(52900000, "heatmap", gp=gpar(col="blue", lwd=3))
addInteractionAnnotation(-52900000, "heatmap", gp=gpar(col="cyan", lwd=3, lty=4))