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.
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)
y-axis can be put to right side of the track by setting main slot to FALSE in y-axis slot of each track. And 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 y scale
setTrackStyleParam(trackList[[1]], "ylim", c(0, 25))
setTrackStyleParam(trackList[[2]], "ylim", c(-25, 0))
viewTracks(trackList, gr=gr, viewerStyle=viewerStyle)
Y label style can be changed by setting the ylabgp slot in 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)
Y label can be also put to top or 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 to 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(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 style of each track. The first color is for dat slot of track and seconde color is for 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 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 also 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 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)
We support two themes now: bw and col.
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)
If there are two tracks and we want to draw the two track by adding or substract one from another, we can try operators.
newtrack <- repA
## must keep same format for dat and dat2
newtrack <- parseWIG(newtrack, "chr11", 122929275, 122930122)
newtrack$dat2 <- newtrack$dat
newtrack$dat <- fox2$dat2
setTrackStyleParam(newtrack, "color", c("blue", "red"))
viewTracks(trackList(newtrack, trs),
gr=gr, viewerStyle=viewerStyle, operator="+")
viewTracks(trackList(newtrack, trs), gr=gr, viewerStyle=viewerStyle, operator="-")
Or try GRoperator before view tracks.
newtrack$dat <- GRoperator(newtrack$dat, newtrack$dat2, col="score", operator="-")
newtrack$dat2 <- GRanges()
viewTracks(trackList(newtrack, trs), gr=gr, viewerStyle=viewerStyle)
Lolliplot is a mutation distribution graphics tool.
SNP <- c(10, 12, 1400, 1402)
sample.gr <- GRanges("chr1", IRanges(SNP, width=1, names=paste0("snp", SNP)))
features <- GRanges("chr1", IRanges(c(1, 501, 1001),
width=c(120, 400, 405),
names=paste0("block", 1:3)))
lolliplot(sample.gr, features)
## more SNPs
SNP <- c(10, 100, 105, 108, 400, 410, 420, 600, 700, 805, 840, 1400, 1402)
sample.gr <- GRanges("chr1", IRanges(SNP, width=1, names=paste0("snp", SNP)))
lolliplot(sample.gr, features)
## define the range
lolliplot(sample.gr, features, ranges = GRanges("chr1", IRanges(104, 109)))
features$fill <- c("#FF8833", "#51C6E6", "#DFA32D")
lolliplot(sample.gr, features)
sample.gr$color <- sample.int(6, length(SNP), replace=TRUE)
sample.gr$border <- sample(c("gray80", "gray30"), length(SNP), replace=TRUE)
lolliplot(sample.gr, features)
sample.gr$label <- as.character(1:length(sample.gr))
sample.gr$label.col <- "white"
lolliplot(sample.gr, features)
features$height <- c(0.02, 0.05, 0.08)
lolliplot(sample.gr, features)
## keep the height by giving the unit
features$height <- list(unit(1/16, "inches"),
unit(3, "mm"),
unit(12, "points"))
lolliplot(sample.gr, features)
The metadata ‘featureLayerID’ are used for drawing features in different layers.
features.mul <- rep(features, 2)
features.mul$height[4:6] <- list(unit(1/8, "inches"),
unit(0.5, "lines"),
unit(.2, "char"))
features.mul$fill <- c("#FF8833", "#F9712A", "#DFA32D",
"#51C6E6", "#009DDA", "#4B9CDF")
end(features.mul)[5] <- end(features.mul[5])+50
features.mul$featureLayerID <-
paste("tx", rep(1:2, each=length(features)), sep="_")
names(features.mul) <-
paste(features.mul$featureLayerID,
rep(1:length(features), 2), sep="_")
lolliplot(sample.gr, features.mul)
## one name per transcripts
names(features.mul) <- features.mul$featureLayerID
lolliplot(sample.gr, features.mul)
#Note: the score value is integer less than 10
sample.gr$score <- sample.int(5, length(sample.gr), replace = TRUE)
lolliplot(sample.gr, features)
##remove yaxis
lolliplot(sample.gr, features, yaxis=FALSE)
#Try score value greater than 10
sample.gr$score <- sample.int(20, length(sample.gr), replace=TRUE)
lolliplot(sample.gr, features)
#Try float numeric score
sample.gr$score <- runif(length(sample.gr))*10
lolliplot(sample.gr, features)
# score should not be smaller than 1
xaxis <- c(1, 200, 400, 701, 1000, 1200, 1402) ## define the position
lolliplot(sample.gr, features, xaxis=xaxis)
names(xaxis) <- xaxis # define the labels
names(xaxis)[4] <- "center"
lolliplot(sample.gr, features, xaxis=xaxis)
yaxis <- c(0, 5) ## define the position
lolliplot(sample.gr, features, yaxis=yaxis)
yaxis <- c(0, 5, 10, 15)
names(yaxis) <- yaxis # define the labels
names(yaxis)[3] <- "yaxis"
lolliplot(sample.gr, features, yaxis=yaxis)
legend <- 1:6 ## legend fill color
names(legend) <- paste0("legend", letters[1:6]) ## legend labels
lolliplot(sample.gr, features, legend=legend)
## use list to define more attributes. see ?grid::gpar to get more details.
legend <- list(labels=paste0("legend", LETTERS[1:6]),
col=palette()[6:1],
fill=palette()[legend])
lolliplot(sample.gr, features, legend=legend)
## if you have multiple tracks, try to set the legend by list.
## see more in section [Plot multiple samples](#plot-multiple-samples)
legend <- list(legend)
lolliplot(sample.gr, features, legend=legend)
Use can control the paramters of labels by name the metadata start as label.parameter.
sample.gr.rot <- sample.gr
sample.gr.rot$label.parameter.rot <- 45
lolliplot(sample.gr.rot, features, legend=legend)
sample.gr.rot$label.parameter.rot <- 60
sample.gr.rot$label.parameter.gp <- gpar(col="brown")
lolliplot(sample.gr.rot, features, legend=legend)
If you want to change the text in ylab, please try to set the labels in ylab. lolliplot does not support the parameters for title and xlab. If you want to add title and xlab, please try to add them by .
lolliplot(sample.gr.rot, features, legend=legend, ylab="y label here")
grid.text("x label here", x=.5, y=.01, just="bottom")
grid.text("title here", x=.5, y=.98, just="top",
gp=gpar(cex=1.5, fontface="bold"))
lolliplot(sample.gr, features, type="pin")
sample.gr$color <- lapply(sample.gr$color, function(.ele) c(.ele, sample.int(6, 1)))
sample.gr$border <- sample.int(6, length(SNP), replace=TRUE)
lolliplot(sample.gr, features, type="pin")
sample.gr$score <- NULL ## must be removed, because pie will consider all the numeric columns except column "color", "fill", "lwd", "id" and "id.col".
sample.gr$label <- NULL
sample.gr$label.col <- NULL
x <- sample.int(100, length(SNP))
sample.gr$value1 <- x
sample.gr$value2 <- 100 - x
## the length of color should be no less than the values number
sample.gr$color <- rep(list(c("#87CEFA", '#98CE31')), length(SNP))
sample.gr$border <- "gray30"
lolliplot(sample.gr, features, type="pie")
SNP2 <- sample(4000:8000, 30)
x2 <- sample.int(100, length(SNP2), replace=TRUE)
sample2.gr <- GRanges("chr3", IRanges(SNP2, width=1, names=paste0("snp", SNP2)),
value1=x2, value2=100-x2)
sample2.gr$color <- rep(list(c('#DB7575', '#FFD700')), length(SNP2))
sample2.gr$border <- "gray30"
features2 <- GRanges("chr3", IRanges(c(5001, 5801, 7001),
width=c(500, 500, 405),
names=paste0("block", 4:6)),
fill=c("orange", "gray30", "lightblue"),
height=unit(c(0.5, 0.3, 0.8), "cm"))
legends <- list(list(labels=c("WT", "MUT"), fill=c("#87CEFA", '#98CE31')),
list(labels=c("WT", "MUT"), fill=c('#DB7575', '#FFD700')))
lolliplot(list(A=sample.gr, B=sample2.gr),
list(x=features, y=features2),
type="pie", legend=legends)
Different layouts are also be possible.
sample2.gr$score <- sample2.gr$value1 ## circle layout need score column
lolliplot(list(A=sample.gr, B=sample2.gr),
list(x=features, y=features2),
type=c("pie", "circle"), legend=legends)
rand.id <- sample.int(length(sample.gr), 3*length(sample.gr), replace=TRUE)
rand.id <- sort(rand.id)
sample.gr.mul.patient <- sample.gr[rand.id]
## pie.stack require metadata "stack.factor", and the metadata can not be
## stack.factor.order or stack.factor.first
len.max <- max(table(rand.id))
stack.factors <- paste0("patient", formatC(1:len.max,
width=nchar(as.character(len.max)),
flag="0"))
sample.gr.mul.patient$stack.factor <-
unlist(lapply(table(rand.id), sample, x=stack.factors))
sample.gr.mul.patient$value1 <-
sample.int(100, length(sample.gr.mul.patient), replace=TRUE)
sample.gr.mul.patient$value2 <- 100 - sample.gr.mul.patient$value1
patient.color.set <- as.list(as.data.frame(rbind(rainbow(length(stack.factors)),
"#FFFFFFFF"),
stringsAsFactors=FALSE))
names(patient.color.set) <- stack.factors
sample.gr.mul.patient$color <-
patient.color.set[sample.gr.mul.patient$stack.factor]
legend <- list(labels=stack.factors, col="gray80",
fill=sapply(patient.color.set, `[`, 1))
lolliplot(sample.gr.mul.patient, features, type="pie.stack",
legend=legend, dashline.col="gray")
Metadata SNPsideID is used to trigger caterpillar layout. SNPsideID must be ‘top’ or ‘bottom’.
sample.gr$SNPsideID <- sample(c("top", "bottom"),
length(sample.gr),
replace=TRUE)
lolliplot(sample.gr, features, type="pie",
legend=legends[[1]])
## two layers
sample2.gr$SNPsideID <- "top"
idx <- sample.int(length(sample2.gr), 15)
sample2.gr$SNPsideID[idx] <- "bottom"
sample2.gr$color[idx] <- '#FFD700'
lolliplot(list(A=sample.gr, B=sample2.gr),
list(x=features.mul, y=features2),
type=c("pie", "circle"), legend=legends)
library(VariantAnnotation)
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(org.Hs.eg.db)
fl <- system.file("extdata", "chr22.vcf.gz", package="VariantAnnotation")
gr <- GRanges("22", IRanges(50968014, 50970514, names="TYMP"))
tab <- TabixFile(fl)
vcf <- readVcf(fl, "hg19", param=gr)
mutation.frequency <- rowRanges(vcf)
mcols(mutation.frequency) <- cbind(mcols(mutation.frequency),
VariantAnnotation::info(vcf))
mutation.frequency$border <- "gray30"
mutation.frequency$color <-
ifelse(grepl("^rs", names(mutation.frequency)), "lightcyan", "lavender")
## plot Global Allele Frequency based on AC/AN
mutation.frequency$score <- mutation.frequency$AF*100
seqlevelsStyle(gr) <- seqlevelsStyle(mutation.frequency) <- "UCSC"
trs <- geneModelFromTxdb(TxDb.Hsapiens.UCSC.hg19.knownGene,
org.Hs.eg.db,
gr=gr)
features <- c(range(trs[[1]]$dat), range(trs[[5]]$dat))
names(features) <- c(trs[[1]]$name, trs[[5]]$name)
features$fill <- c("lightblue", "mistyrose")
features$height <- c(.02, .04)
lolliplot(mutation.frequency, features, ranges=gr)
library(rtracklayer)
session <- browserSession()
query <- ucscTableQuery(session, track="HAIB Methyl RRBS",
range=GRangesForUCSCGenome("hg19", seqnames(gr), ranges(gr)))
tableName(query) <- tableNames(query)[1]
methy <- track(query)
methy <- GRanges(methy)
lolliplot(methy, features, ranges=gr, type="pin")
In above sample, some of the nodes overlap each other. To change the node size, cex prameter could be used.
methy$lwd <- .5
lolliplot(methy, features, ranges=gr, type="pin", cex=.5)
lolliplot(methy, features, ranges=gr, type="circle", cex=.5)
methy$score2 <- max(methy$score) - methy$score
lolliplot(methy, features, ranges=gr, type="pie", cex=.5)
Sometimes, there are too many SNPs invoved in one gene. If you want to plot such as more than handred SNPs, Dandelion plot may be the good chioce. Please note, the height of the dandelion means the desity of the SNPs.
dandelion.plot(methy, features, ranges=gr, type="pin")
There are one more type for dandelion plot, type “fan”. The area of the fan indicate the percentage of methylation or rate of mutation.
methy$color <- 3
methy$border <- "gray"
## score info is required, the score must be a number in [0, 1]
m <- max(methy$score)
methy$score <- methy$score/m
dandelion.plot(methy, features, ranges=gr, type="fan")
methy$color <- rep(list(c(3, 5)), length(methy))
methy$score2 <- methy$score2/m
legends <- list(list(labels=c("s1", "s2"), fill=c(3, 5)))
dandelion.plot(methy, features, ranges=gr, type="pie", legend=legends)
## want less
dandelion.plot(methy, features, ranges=gr, type="circle", maxgaps=1/10)
## want more
dandelion.plot(methy, features, ranges=gr, type="circle", maxgaps=1/100)
sessionInfo()
R version 3.3.2 (2016-10-31) Platform: x86_64-pc-linux-gnu (64-bit) Running under: Ubuntu 16.04.1 LTS
locale: [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_US.UTF-8 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 parallel stats4 stats graphics grDevices utils
[8] datasets methods base
other attached packages: [1] VariantAnnotation_1.20.2
[2] Rsamtools_1.26.1
[3] Biostrings_2.42.1
[4] XVector_0.14.0
[5] SummarizedExperiment_1.4.0
[6] org.Hs.eg.db_3.4.0
[7] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2 [8] GenomicFeatures_1.26.0
[9] AnnotationDbi_1.36.0
[10] Biobase_2.34.0
[11] Gviz_1.18.1
[12] rtracklayer_1.34.1
[13] trackViewer_1.10.2
[14] GenomicRanges_1.26.1
[15] GenomeInfoDb_1.10.1
[16] IRanges_2.8.1
[17] S4Vectors_0.12.1
[18] BiocGenerics_0.20.0
[19] BiocStyle_2.2.1
loaded via a namespace (and not attached): [1] httr_1.2.1 AnnotationHub_2.6.4
[3] splines_3.3.2 Formula_1.2-1
[5] shiny_0.14.2 assertthat_0.1
[7] interactiveDisplayBase_1.12.0 highr_0.6
[9] latticeExtra_0.6-28 BSgenome_1.42.0
[11] yaml_2.1.14 RSQLite_1.1
[13] backports_1.0.4 lattice_0.20-34
[15] biovizBase_1.22.0 digest_0.6.10
[17] RColorBrewer_1.1-2 colorspace_1.3-1
[19] htmltools_0.3.5 httpuv_1.3.3
[21] Matrix_1.2-7.1 plyr_1.8.4
[23] XML_3.98-1.5 biomaRt_2.30.0
[25] zlibbioc_1.20.0 xtable_1.8-2
[27] scales_0.4.1 BiocParallel_1.8.1
[29] htmlTable_1.7 tibble_1.2
[31] ggplot2_2.2.0 pbapply_1.3-1
[33] nnet_7.3-12 lazyeval_0.2.0
[35] survival_2.40-1 magrittr_1.5
[37] mime_0.5 memoise_1.0.0
[39] evaluate_0.10 foreign_0.8-67
[41] grImport_0.9-0 BiocInstaller_1.24.0
[43] tools_3.3.2 data.table_1.9.8
[45] matrixStats_0.51.0 stringr_1.1.0
[47] munsell_0.4.3 cluster_2.0.5
[49] ensembldb_1.6.2 RCurl_1.95-4.8
[51] dichromat_2.0-0 bitops_1.0-6
[53] rmarkdown_1.2 gtable_0.2.0
[55] DBI_0.5-1 R6_2.2.0
[57] GenomicAlignments_1.10.0 gridExtra_2.2.1
[59] knitr_1.15.1 Hmisc_4.0-0
[61] rprojroot_1.1 stringi_1.1.2
[63] Rcpp_0.12.8 rpart_4.1-10
[65] acepack_1.4.1