trackViewer 1.16.1
There are two packages available in Bioconductor for visualizing genomic data: rtracklayer and Gviz. rtracklayer provides an interface to genome browsers and associated annotation tracks. Gviz plots data and annotation information along genomic coordinates. TrackViewer is a light-weighted visualization tool for generating neat figures for publication. It utilizes Gviz, is easy to use, and has a low memory and cpu consumption.
library(Gviz)
library(rtracklayer)
library(trackViewer)
extdata <- system.file("extdata", package="trackViewer",
mustWork=TRUE)
gr <- GRanges("chr11", IRanges(122929275, 122930122), strand="-")
fox2 <- importScore(file.path(extdata, "fox2.bed"), format="BED",
ranges=gr)
fox2$dat <- coverageGR(fox2$dat)
viewTracks(trackList(fox2), gr=gr, autoOptimizeStyle=TRUE, newpage=FALSE)
dt <- DataTrack(range=fox2$dat[strand(fox2$dat)=="-"] ,
genome="hg19", type="hist", name="fox2",
window=-1, chromosome="chr11",
fill.histogram="black", col.histogram="NA",
background.title="white",
col.frame="white", col.axis="black",
col="black", col.title="black")
plotTracks(dt, from=122929275, to=122930122, strand="-")
TrackViewer not only has the functionalities to plot the figures generated by Gviz, as shown in Figure above, but also provides additional plotting styles as shown in Figure below. The mimimalist design requires minimum input from users while retaining the flexibility to change output style easily.
gr <- GRanges("chr1", IRanges(c(1, 6, 10), c(3, 6, 12)), score=c(3, 4, 1))
dt <- DataTrack(range=gr, data="score", type="hist")
plotTracks(dt, from=2, to=11)
tr <- new("track", dat=gr, type="data", format="BED")
viewTracks(trackList(tr), chromosome="chr1", start=2, end=11)
It requires huge memory space to handle big wig files. To solve this problem, trackViewer rewrote the import function to import whole file first and parse it later when plot. trackViewer provides higher import speed (21 min vs. over 180 min) and acceptable memory cost (5.32G vs. over 10G) for a half giga wig file (GSM917672) comparing to Gviz.
Function importScore is used to import BED, WIG, bedGraph or BigWig files. Function importBam is employed to import bam file. Here is the example.
library(trackViewer)
extdata <- system.file("extdata", package="trackViewer",
mustWork=TRUE)
repA <- importScore(file.path(extdata, "cpsf160.repA_-.wig"),
file.path(extdata, "cpsf160.repA_+.wig"),
format="WIG")
## because the wig file does not contain strand info,
## we need to set it manually
strand(repA$dat) <- "-"
strand(repA$dat2) <- "+"
Function coverageGR could be used to calculate coverage after importing if needed.
fox2 <- importScore(file.path(extdata, "fox2.bed"), format="BED",
ranges=GRanges("chr11", IRanges(122929000, 122931000)))
dat <- coverageGR(fox2$dat)
## we can split the data by strand into two different track channels
## here we set the dat2 slot to save the negative strand info,
## reverse order as previous.
fox2$dat <- dat[strand(dat)=="+"]
fox2$dat2 <- dat[strand(dat)=="-"]
The gene model can be built for a given genomic range using geneModelFromTxdb function which uses TranscriptDb object as input.
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(org.Hs.eg.db)
gr <- GRanges("chr11", IRanges(122929275, 122930122), strand="-")
trs <- geneModelFromTxdb(TxDb.Hsapiens.UCSC.hg19.knownGene,
org.Hs.eg.db,
gr=gr)
Use viewTracks function to plot data and annotation information along genomic coordinates. addGuideLine or addArrowMark can be used to highlight the peaks.
viewerStyle <- trackViewerStyle()
setTrackViewerStyleParam(viewerStyle, "margin", c(.1, .05, .02, .02))
vp <- viewTracks(trackList(repA, fox2, trs),
gr=gr, viewerStyle=viewerStyle,
autoOptimizeStyle=TRUE)
addGuideLine(c(122929767, 122929969), vp=vp)
addArrowMark(list(x=122929650,
y=2), # 2 means track 2 from bottom.
label="label",
col="blue",
vp=vp)
In most cases, researchers are interested in the relative position of peaks in the gene. Sometimes, margin needs to be adjusted to be able to show the entire gene model. Figure below shows how to add an x-scale and remove x-axis using addGuideLine Function .
optSty <- optimizeStyle(trackList(repA, fox2, trs))
trackList <- optSty$tracks
viewerStyle <- optSty$style
setTrackViewerStyleParam(viewerStyle, "xaxis", FALSE)
setTrackViewerStyleParam(viewerStyle, "margin", c(.01, .05, .01, .01))
setTrackXscaleParam(trackList[[1]], "draw", TRUE)
setTrackXscaleParam(trackList[[1]], "gp", list(cex=.5))
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)