Abstract
Visualize mapped reads along with annotation as track layers for NGS dataset such as ChIP-seq, RNA-seq, miRNA-seq, DNA-seq, SNPs and methylation data.Plot high dense methylation/mutation/variant data by dandelion plot
Plot genomic interaction data (HiC, cooler, etc) by viewTracks
First, we import the sample data.
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 any strand info,
## we need to set it manually.
strand(repA$dat) <- "-"
strand(repA$dat2) <- "+"
fox2 <- importScore(file.path(extdata, "fox2.bed"), format="BED",
ranges=GRanges("chr11", IRanges(122830799, 123116707)))
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.
fox2$dat <- dat[strand(dat)=="+"]
fox2$dat2 <- dat[strand(dat)=="-"]
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)
In most cases, researchers are interested in the relative position of the peaks in the gene. Sometimes, margin needs to be adjusted to be able to show the entire gene model. The Figure below shows how to add an X scale (x-scale) and remove the x-axis using the setTrackXscaleParam and setTrackViewerStyleParam functions.
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=0.8))
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
setTrackXscaleParam(trackList[[1]], attr="position",
value=list(x=122929700, y=3, label=200))
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
The y-axis can be put to the right side of the track by setting the main slot to FALSE in the y-axis slot of each track. In addition, the limit of y-axis (ylim) can be set by setTrackStyleParam.
setTrackViewerStyleParam(viewerStyle, "margin", c(.01, .05, .01, .05))
for(i in 1:2){
setTrackYaxisParam(trackList[[i]], "main", FALSE)
}
## Adjust the limit of y-axis
setTrackStyleParam(trackList[[1]], "ylim", c(0, 25))
setTrackStyleParam(trackList[[2]], "ylim", c(-25, 0))
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
The function setTrackYaxisParam
can be sued to adjust the track color/size/etc of y-axis and the y-axis labels.
## change the
setTrackYaxisParam(trackList[[1]], "gp", list(cex=.8, col="green"))
setTrackYaxisParam(trackList[[2]], "gp", list(cex=.8, col="blue"))
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
The style of y-axis can be changed by setting the ylabgp slot in the style of each track.
setTrackStyleParam(trackList[[1]], "ylabgp", list(cex=.8, col="green"))
## set cex to avoid automatic adjust
setTrackStyleParam(trackList[[2]], "ylabgp", list(cex=.8, col="blue"))
setTrackStyleParam(trackList[[2]], "marginBottom", .2)
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
The y-axis label can be put at the top or the bottom of each track.
setTrackStyleParam(trackList[[1]], "ylabpos", "bottomleft")
setTrackStyleParam(trackList[[2]], "ylabpos", "topright")
setTrackStyleParam(trackList[[2]], "marginTop", .2)
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
For each transcript, the transcript name can be put on the upstream or downstream of the transcript.
trackListN <- trackList
setTrackStyleParam(trackListN[[3]], "ylabpos", "upstream")
setTrackStyleParam(trackListN[[4]], "ylabpos", "downstream")
## set cex to avoid automatic adjust
setTrackStyleParam(trackListN[[3]], "ylabgp", list(cex=.6))
setTrackStyleParam(trackListN[[4]], "ylabgp", list(cex=.6))
gr1 <- range(unname(unlist(GRangesList(sapply(trs, function(.ele) .ele$dat)))))
start(gr1) <- start(gr1) - 2000
end(gr1) <- end(gr1) + 2000
viewTracks(trackListN, gr=gr1, viewerStyle=viewerStyle)
The track color can be changed by setting the color slot in the style of each track. The first color is for the dat slot of track and the second color is for the dat2 slot.
setTrackStyleParam(trackList[[1]], "color", c("green", "black"))
setTrackStyleParam(trackList[[2]], "color", c("black", "blue"))
for(i in 3:length(trackList))
setTrackStyleParam(trackList[[i]], "color", "black")
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
The track height can be changed by setting the height slot in the style of each track. However, the total height for all the tracks should be 1.
trackListH <- trackList
setTrackStyleParam(trackListH[[1]], "height", .1)
setTrackStyleParam(trackListH[[2]], "height", .44)
for(i in 3:length(trackListH)){
setTrackStyleParam(trackListH[[i]], "height",
(1-(0.1+0.44))/(length(trackListH)-2))
}
viewTracks(trackListH, gr=gr, viewerStyle=viewerStyle)
The track names such as gene model names can be edited easily by changing the names of trackList.
names(trackList) <- c("cpsf160", "fox2", rep("HSPA8", 5))
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
trackViewer can be used to show to-be-compared data in the same track side by side.
cpsf160 <- importScore(file.path(extdata, "cpsf160.repA_-.wig"),
file.path(extdata, "cpsf160.repB_-.wig"),
format="WIG")
strand(cpsf160$dat) <- strand(cpsf160$dat2) <- "-"
setTrackStyleParam(cpsf160, "color", c("black", "red"))
viewTracks(trackList(trs, cpsf160), gr=gr, viewerStyle=viewerStyle)
The x-axis can be horizotally flipped for the genes in the negative strand.
viewerStyleF <- viewerStyle
setTrackViewerStyleParam(viewerStyleF, "flip", TRUE)
setTrackViewerStyleParam(viewerStyleF, "xaxis", TRUE)
setTrackViewerStyleParam(viewerStyleF, "margin", c(.1, .05, .01, .01))
vp <- viewTracks(trackList, gr=gr, viewerStyle=viewerStyleF)
addGuideLine(c(122929767, 122929969), vp=vp)
addArrowMark(list(x=122929650,
y=2),
label="label",
col="blue",
vp=vp)
Currently, we support two themes: bw (black and white) and col (colored).
optSty <- optimizeStyle(trackList(repA, fox2, trs), theme="bw")
trackList <- optSty$tracks
viewerStyle <- optSty$style
vp <- viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
optSty <- optimizeStyle(trackList(repA, fox2, trs), theme="col")
trackList <- optSty$tracks
viewerStyle <- optSty$style
vp <- viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
optSty <- optimizeStyle(trackList(repA, fox2, trs), theme="safe")
trackList <- optSty$tracks
viewerStyle <- optSty$style
vp <- viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
We could plot the tracks with breaks by setting multiple genomic ranges.
gr.breaks <- GRanges("chr11",
IRanges(c(122929275, 122929575, 122929775),
c(122929555, 122929725, 122930122)),
strand="-", percentage=c(.4, .2, .4))
vp <- viewTracks(trackList, gr=gr.breaks, viewerStyle=viewerStyle)
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] httr_1.4.3
[2] VariantAnnotation_1.42.1
[3] Rsamtools_2.12.0
[4] Biostrings_2.64.0
[5] XVector_0.36.0
[6] SummarizedExperiment_1.26.1
[7] MatrixGenerics_1.8.0
[8] matrixStats_0.62.0
[9] org.Hs.eg.db_3.15.0
[10] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2 [11] GenomicFeatures_1.48.1
[12] AnnotationDbi_1.58.0
[13] Biobase_2.56.0
[14] Gviz_1.40.1
[15] rtracklayer_1.56.0
[16] trackViewer_1.32.1
[17] Rcpp_1.0.8.3
[18] GenomicRanges_1.48.0
[19] GenomeInfoDb_1.32.2
[20] IRanges_2.30.0
[21] S4Vectors_0.34.0
[22] 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] InteractionSet_1.24.0 XML_3.99-0.9 zlibbioc_1.42.0
[49] scales_1.2.0 BiocStyle_2.24.0 BSgenome_1.64.0
[52] ProtGenerics_1.28.0 hms_1.1.1 parallel_4.2.0
[55] rhdf5_2.40.0 AnnotationFilter_1.20.0 RColorBrewer_1.1-3
[58] yaml_2.3.5 curl_4.3.2 memoise_2.0.1
[61] gridExtra_2.3 ggplot2_3.3.6 sass_0.4.1
[64] biomaRt_2.52.0 rpart_4.1.16 latticeExtra_0.6-29
[67] stringi_1.7.6 RSQLite_2.2.14 highr_0.9
[70] BiocIO_1.6.0 plotrix_3.8-2 checkmate_2.1.0
[73] filelock_1.0.2 BiocParallel_1.30.2 rlang_1.0.2
[76] pkgconfig_2.0.3 bitops_1.0-7 evaluate_0.15
[79] lattice_0.20-45 Rhdf5lib_1.18.2 purrr_0.3.4
[82] GenomicAlignments_1.32.0 htmlwidgets_1.5.4 bit_4.0.4
[85] tidyselect_1.1.2 magrittr_2.0.3 R6_2.5.1
[88] generics_0.1.2 Hmisc_4.7-0 DelayedArray_0.22.0
[91] DBI_1.1.2 pillar_1.7.0 foreign_0.8-82
[94] survival_3.3-1 KEGGREST_1.36.0 RCurl_1.98-1.6
[97] nnet_7.3-17 tibble_3.1.7 crayon_1.5.1
[100] utf8_1.2.2 BiocFileCache_2.4.0 rmarkdown_2.14
[103] jpeg_0.1-9 progress_1.2.2 data.table_1.14.2
[106] Rgraphviz_2.40.0 blob_1.2.3 digest_0.6.29
[109] munsell_0.5.0 bslib_0.3.1