autoplot {ggbio} | R Documentation |
autoplot
is a generic function to visualize various data
object, it tries to give better default graphics and customized
choices for each data type, quick and convenient to explore your
genomic data compare to low level ggplot
method, it is much
simpler and easy to produce fairly complicate graphics, though you may
lose some flexibility for each layer.
## S4 method for signature 'GRanges' autoplot(object, ..., xlab, ylab, main, truncate.gaps = FALSE, truncate.fun = NULL, ratio = 0.0025, space.skip = 0.1, legend = TRUE, geom = NULL, stat = NULL, coord = c("default", "genome", "truncate_gaps"), layout = c("linear", "karyogram", "circle")) ## S4 method for signature 'GRangesList' autoplot(object, ..., xlab, ylab, main, indName = "grl_name", geom = NULL, stat = NULL, type = c("none", "sashimi"), coverage.col = "gray50", coverage.fill = coverage.col, group.selfish = FALSE, arch.offset = 1.3) ## S4 method for signature 'IRanges' autoplot(object, ..., xlab, ylab, main) ## S4 method for signature 'Seqinfo' autoplot(object, single.ideo = TRUE, ... ) ## S4 method for signature 'GappedAlignments' autoplot(object, ..., xlab, ylab, main, which, geom = NULL, stat = NULL) ## S4 method for signature 'BamFile' autoplot(object, ..., which, xlab, ylab, main, bsgenome, geom = "line", stat = "coverage", method = c("estimate", "raw"), coord = c("linear", "genome"), resize.extra = 10, space.skip = 0.1, show.coverage = TRUE) ## S4 method for signature 'character' autoplot(object, ..., xlab, ylab, main, which, asRangedData = FALSE) ## S4 method for signature 'TranscriptDb' autoplot(object, which, ..., xlab, ylab, main, truncate.gaps = FALSE, truncate.fun = NULL, ratio = 0.0025, geom = c("alignment"), stat = c("identity", "reduce"), names.expr = "tx_name(gene_id)") ## S4 method for signature 'BSgenome' autoplot(object, which, ..., xlab, ylab, main, geom = c("text", "segment", "point", "rect")) ## S4 method for signature 'Rle' autoplot(object, ..., xlab, ylab, main, binwidth, nbin = 30, geom = NULL, stat = c("bin", "identity", "slice"), type = c("viewSums", "viewMins", "viewMaxs", "viewMeans")) ## S4 method for signature 'RleList' autoplot(object, ..., xlab, ylab, main, nbin = 30, binwidth, facetByRow = TRUE, stat = c("bin", "identity", "slice"), geom = NULL, type = c("viewSums", "viewMins", "viewMaxs", "viewMeans")) ## S4 method for signature 'matrix' autoplot(object, ..., xlab, ylab, main, geom = c("tile", "raster"), axis.text.angle = NULL, hjust = 0.5, na.value = NULL, rownames.label = TRUE, colnames.label = TRUE, axis.text.x = TRUE, axis.text.y = TRUE) ## S4 method for signature 'ExpressionSet' autoplot(object, ..., type = c("heatmap", "none", "scatterplot.matrix", "pcp", "MA", "boxplot", "mean-sd", "volcano", "NUSE", "RLE"), test.method = "t", rotate = FALSE, pheno.plot = FALSE, main_to_pheno = 4.5, padding = 0.2) ## S4 method for signature 'SummarizedExperiment' autoplot(object, ..., type = c("heatmap", "link", "pcp", "boxplot", "scatterplot.matrix"), pheno.plot = FALSE, main_to_pheno = 4.5, padding = 0.2, assay.id = 1) ## S4 method for signature 'VCF' autoplot(object, ..., xlab, ylab, main, assay.id, type = c("geno", "info", "fixed"), full.string = FALSE, ref.show = TRUE, genome.axis = TRUE, transpose = TRUE)
object |
object to be plot. |
transpose |
logical value, defaut TRUE, always make features from VCF as x, so we can use it to map to genomic position. |
axis.text.angle |
axis text angle. |
axis.text.x |
logical value indicates whether to show x axis and labels or not. |
axis.text.y |
logical value indicates whether to show y axis and labels or not. |
hjust |
horizontal just for axis text. |
rownames.label |
logical value indicates whether to show rownames of matrix as y label or not. |
colnames.label |
logical value indicates whether to show colnames of matrix as y label or not. |
na.value |
color for NA value. |
rotate |
|
pheno.plot |
show pheno plot or not. |
main_to_pheno |
main matrix plot width to pheno plot width ratio. |
padding |
padding between plots. |
assay.id |
index for assay you are going to use. |
geom |
Geom to use (Single character for now). Please see section Geometry for details. |
truncate.gaps |
logical value indicate to truncate gaps or not. |
truncate.fun |
shrinkage function. Please see |
ratio |
used in |
space.skip |
space ratio between chromosome spaces in coordate genome. |
coord |
Coodinate system. |
legend |
A logical value indicates whether to show legend or not. Default is
|
which |
A |
show.coverage |
A logical value indicates whether to show coverage or not. This is used for geom "mismatch.summary". |
resize.extra |
A numeric value used to add buffer to intervals to compute stepping levels on. |
bsgenome |
A BSgenome object. Only need for geom "mismatch.summary". |
xlab |
x label. |
ylab |
y label. |
facetByRow |
A logical value, default is TRUE ,facet RleList by row. If FALSE, facet by column. |
type |
For Rle/RleList, "raw" plot everything, so be careful, that would be
pretty slow if you have too much data. For "viewMins", "viewMaxs",
"viewMeans", "viewSums", require extra arguments to slice the
object. so users need to at least provide For ExpreesionSet, ploting types. |
layout |
Layout including linear, circular and karyogram. for |
method |
method used for parsing coverage from bam files. 'estimate' use fast esitmated method and 'raw' use relatively slow parsing method. |
asRangedData |
when object is character, it may be imported by |
test.method |
test method |
... |
Extra parameters. Usually are those parameters used in autoplot to control aesthetics or geometries. |
main |
title. |
stat |
statistical transformation. |
indName |
When coerce |
coverage.col |
coverage stroke color. |
coverage.fill |
coverage fill color. |
group.selfish |
Passed to |
arch.offset |
arch.offset. |
names.expr |
names expression used for creating labels. |
binwidth |
width of the bins. |
nbin |
number of bins. |
genome.axis |
logical value, if TRUE, whenever possible, try to parse genomic postition for each column(e.g. SummarizedExperiment), show column as exatcly the genomic position instead of showing them side by side and indexed from 1. |
full.string |
logical value. If TRUE, show full string of indels in plot for VCF. |
ref.show |
logical value. If TRUE, show REF in VCF at bottom track. |
single.ideo |
for object 'Seqinfo', if length of object is 1, then this control
whether to use |
A ggplot
object, so you can use common features from ggplot2
package to manipulate the plot.
autoplot
is redefined as generic s4 method inside this package,
user could use autoplot
in the way they are familiar with, and
we are also setting limitation inside this package, like
scales X scales is always genomic coordinates in most cases, x could be specified as start/end/midpoint when it's special geoms for interval data like point/line
colors Try to use default color scheme defined in biovizBase package as possible as it can
We have developed new geom
for different objects, some of
them may require extra parameters you need to provide. Some of the
geom are more like geom + stat in ggplot2
package. e.g. "coverage.line" and "coverage.polygon".We simply combine
them together, but in the future, we plan to make it more general.
This package is designed for only biological data, especially genomic
data if users want to explore the data in a more flexible way, you
could simply coerce the GRanges
to a data.frame, then
just use formal autoplot
function in ggplot2, or autoplot
generic for data.frame
.
Some objects share the same geom so we introduce all the geom together in this section
Showing all the intervals as stepped rectangle, colored by strand automatically.
For TranscripDb
object, showing full model.
Showing all the intervals as stepped segments, colored by strand automatically.
For object BSgenome
, show nucleotides as colored segment.
For Rle/RleList, show histogram-like segments.
Showing interval as line, the interval data could also be just single position when start = end, x is one of start/end/midpoint, y value is unquoted name in elementMetadata column names. y value is required.
Showing interval as point, the interval data could also be just single position when start = end, x is one of start/end/midpoint, y value is unquoted name in elementMetadata column names. y value is required.
For object BSgenome
, show nucleotides as colored point.
Coverage showing as lines for interval data.
Coverage showing as polygon for interval data.
Splicing summary. The size and width of the line and rectangle should represent the counts in each model. Need to provide model.
For TranscripDb
object, showing single(reduced) model only.
For TranscripDb
object, showing transcirpts isoforms.
Showing color coded mismatched stacked bar to indicate the proportion of mismatching at each position, the reference is set to gray.
For object BSgenome
, show nucleotides as colored text.
For object BSgenome
, show nucleotides as colored rectangle.
Faceting in ggbio package is a little differnt from ggplot2 in several ways
The faceted column could only be seqnames or regions on the genome. So we limited the formula passing to facet argument, e.g something \~ seqnames, is accepted formula, you can change "something" to variable name in the elementMetadata. But you can not change the second part.
Sometime, we need to view different regions, so we also have a
facet_gr argument which accept a GRanges
. If
this is provided, it will override the default seqnames and use
provided region to facet the graphics, this might be useful for
different gene centric views.
Tengfei Yin
### R code from vignette source 'autoplot.Rnw'
###################################################
### code chunk number 1: load
###################################################
library(ggbio)
###################################################
### code chunk number 2: simul
###################################################
set.seed(1)
N <- 1000
library(GenomicRanges)
gr <- GRanges(seqnames =
sample(c("chr1", "chr2", "chr3"),
size = N, replace = TRUE),
IRanges(
start = sample(1:300, size = N, replace = TRUE),
width = sample(70:75, size = N,replace = TRUE)),
strand = sample(c("+", "-", "*"), size = N,
replace = TRUE),
value = rnorm(N, 10, 3), score = rnorm(N, 100, 30),
sample = sample(c("Normal", "Tumor"),
size = N, replace = TRUE),
pair = sample(letters, size = N,
replace = TRUE))
idx <- sample(1:length(gr), size = 50)
###################################################
### code chunk number 3: default
###################################################
autoplot(gr[idx])
###################################################
### code chunk number 4: bar-default-pre
###################################################
set.seed(123)
gr.b <- GRanges(seqnames = "chr1", IRanges(start = seq(1, 100, by = 10),
width = sample(4:9, size = 10, replace = TRUE)),
score = rnorm(10, 10, 3), value = runif(10, 1, 100))
gr.b2 <- GRanges(seqnames = "chr2", IRanges(start = seq(1, 100, by = 10),
width = sample(4:9, size = 10, replace = TRUE)),
score = rnorm(10, 10, 3), value = runif(10, 1, 100))
gr.b <- c(gr.b, gr.b2)
## Warning: Each of the 2 combined objects has sequence levels not in the
## other: - in 'x': chr1 - in 'y': chr2 Make sure to always combine/compare
## objects based on the same reference genome (use suppressWarnings() to
## suppress this warning).
head(gr.b)
## GRanges with 6 ranges and 2 metadata columns:
## seqnames ranges strand | score value
## <Rle> <IRanges> <Rle> | <numeric> <numeric>
## [1] chr1 [ 1, 5] * | 15.1451949606498 96.3393990211189
## [2] chr1 [11, 18] * | 11.3827486179676 90.327605466824
## [3] chr1 [21, 26] * | 6.2048162961804 69.3798225638457
## [4] chr1 [31, 39] * | 7.93944144431942 79.7512743510306
## [5] chr1 [41, 49] * | 8.66301408970013 3.436754766386
## [6] chr1 [51, 54] * | 13.6722453923184 48.3018011380918
## ---
## seqlengths:
## chr1 chr2
## NA NA
###################################################
### code chunk number 5: bar-default
###################################################
p1 <- autoplot(gr.b, geom = "bar")
## use score as y by default
## use value to fill the bar
p2 <- autoplot(gr.b, geom = "bar", aes(fill = value))
## use score as y by default
tracks(default = p1, fill = p2)
###################################################
### code chunk number 6: autoplot.Rnw:236-237
###################################################
autoplot(gr[idx], geom = "arch", aes(color = value), facets = sample ~ seqnames)
###################################################
### code chunk number 7: gr-group
###################################################
gra <- GRanges("chr1", IRanges(c(1,7,20), end = c(4,9,30)), group = c("a", "a", "b"))
## if you desn't specify group, then group based on stepping levels, and gaps are computed without
## considering extra group method
p1 <- autoplot(gra, aes(fill = group), geom = "alignment")
## when use group method, gaps only computed for grouped intervals.
## default is group.selfish = TRUE, each group keep one row.
## in this way, group labels could be shown as y axis.
p2 <- autoplot(gra, aes(fill = group, group = group), geom = "alignment")
## group.selfish = FALSE, save space
p3 <- autoplot(gra, aes(fill = group, group = group), geom = "alignment", group.selfish = FALSE)
tracks('non-group' = p1,'group.selfish = TRUE' = p2 , 'group.selfish = FALSE' = p3)
###################################################
### code chunk number 8: gr-facet-strand
###################################################
autoplot(gr, stat = "coverage", geom = "area",
facets = strand ~ seqnames, aes(fill = strand))
###################################################
### code chunk number 9: gr-autoplot-circle
###################################################
autoplot(gr[idx], layout = 'circle')
###################################################
### code chunk number 10: gr-circle
###################################################
seqlengths(gr) <- c(400, 500, 700)
values(gr)$to.gr <- gr[sample(1:length(gr), size = length(gr))]
idx <- sample(1:length(gr), size = 50)
gr <- gr[idx]
ggplot() + layout_circle(gr, geom = "ideo", fill = "gray70", radius = 7, trackWidth = 3) +
layout_circle(gr, geom = "bar", radius = 10, trackWidth = 4,
aes(fill = score, y = score)) +
layout_circle(gr, geom = "point", color = "red", radius = 14,
trackWidth = 3, grid = TRUE, aes(y = score)) +
layout_circle(gr, geom = "link", linked.to = "to.gr", radius = 6, trackWidth = 1)
###################################################
### code chunk number 11: seqinfo-src
###################################################
data(hg19Ideogram, package = "biovizBase")
sq <- seqinfo(hg19Ideogram)
sq
## Seqinfo of length 93
## seqnames seqlengths isCircular genome
## chr1 249250621 <NA> hg19
## chr1_gl000191_random 106433 <NA> hg19
## chr1_gl000192_random 547496 <NA> hg19
## chr2 243199373 <NA> hg19
## chr3 198022430 <NA> hg19
## chr4 191154276 <NA> hg19
## chr4_ctg9_hap1 590426 <NA> hg19
## chr4_gl000193_random 189789 <NA> hg19
## chr4_gl000194_random 191469 <NA> hg19
## ... ... ... ...
## chrUn_gl000242 43523 <NA> hg19
## chrUn_gl000243 43341 <NA> hg19
## chrUn_gl000244 39929 <NA> hg19
## chrUn_gl000245 36651 <NA> hg19
## chrUn_gl000246 38154 <NA> hg19
## chrUn_gl000247 36422 <NA> hg19
## chrUn_gl000248 39786 <NA> hg19
## chrUn_gl000249 38502 <NA> hg19
## chrM 16571 <NA> hg19
###################################################
### code chunk number 12: seqinfo
###################################################
autoplot(sq[paste0("chr", c(1:22, "X"))])
###################################################
### code chunk number 13: ir-load
###################################################
set.seed(1)
N <- 100
ir <- IRanges(start = sample(1:300, size = N, replace = TRUE),
width = sample(70:75, size = N,replace = TRUE))
## add meta data
df <- DataFrame(value = rnorm(N, 10, 3), score = rnorm(N, 100, 30),
sample = sample(c("Normal", "Tumor"),
size = N, replace = TRUE),
pair = sample(letters, size = N,
replace = TRUE))
values(ir) <- df
ir
## IRanges of length 100
## start end width
## [1] 80 152 73
## [2] 112 183 72
## [3] 172 242 71
## [4] 273 347 75
## [5] 61 133 73
## [6] 270 340 71
## [7] 284 353 70
## [8] 199 270 72
## [9] 189 263 75
## ... ... ... ...
## [92] 18 89 72
## [93] 193 262 70
## [94] 263 337 75
## [95] 234 304 71
## [96] 240 312 73
## [97] 137 206 70
## [98] 124 198 75
## [99] 244 314 71
## [100] 182 255 74
###################################################
### code chunk number 14: ir-exp
###################################################
p1 <- autoplot(ir)
p2 <- autoplot(ir, aes(fill = pair)) + theme(legend.position = "none")
p3 <- autoplot(ir, stat = "coverage", geom = "line", facets = sample ~. )
p4 <- autoplot(ir, stat = "reduce")
tracks(p1, p2, p3, p4)
###################################################
### code chunk number 15: grl-simul
###################################################
set.seed(1)
N <- 100
## ======================================================================
## simmulated GRanges
## ======================================================================
gr <- GRanges(seqnames =
sample(c("chr1", "chr2", "chr3"),
size = N, replace = TRUE),
IRanges(
start = sample(1:300, size = N, replace = TRUE),
width = sample(30:40, size = N,replace = TRUE)),
strand = sample(c("+", "-", "*"), size = N,
replace = TRUE),
value = rnorm(N, 10, 3), score = rnorm(N, 100, 30),
sample = sample(c("Normal", "Tumor"),
size = N, replace = TRUE),
pair = sample(letters, size = N,
replace = TRUE))
grl <- split(gr, values(gr)$pair)
###################################################
### code chunk number 16: grl-exp
###################################################
## default gap.geom is 'chevron'
p1 <- autoplot(grl, group.selfish = TRUE)
p2 <- autoplot(grl, group.selfish = TRUE, main.geom = "arrowrect", gap.geom = "segment")
tracks(p1, p2)
###################################################
### code chunk number 17: grl-name
###################################################
autoplot(grl, aes(fill = ..grl_name..))
## equal to
## autoplot(grl, aes(fill = grl_name))
###################################################
### code chunk number 18: rle-simul
###################################################
library(IRanges)
library(ggbio)
set.seed(1)
lambda <- c(rep(0.001, 4500), seq(0.001, 10, length = 500),
seq(10, 0.001, length = 500))
## @knitr create
xVector <- rpois(1e4, lambda)
xRle <- Rle(xVector)
xRle
## numeric-Rle of length 10000 with 823 runs
## Lengths: 779 1 208 1 1599 1 ... 1 5 2 9 1 4507
## Values : 0 1 0 1 0 1 ... 1 0 1 0 1 0
###################################################
### code chunk number 19: rle-bin
###################################################
p1 <- autoplot(xRle)
## Default use binwidth: range/30
p2 <- autoplot(xRle, nbin = 80)
## Default use binwidth: range/80
p3 <- autoplot(xRle, geom = "heatmap", nbin = 200)
## Default use binwidth: range/200
tracks('nbin = 30' = p1, "nbin = 80" = p2, "nbin = 200(heatmap)" = p3)
###################################################
### code chunk number 20: rle-id
###################################################
p1 <- autoplot(xRle, stat = "identity")
p2 <- autoplot(xRle, stat = "identity", geom = "point", color = "red")
tracks('line' = p1, "point" = p2)
###################################################
### code chunk number 21: rle-slice
###################################################
p1 <- autoplot(xRle, type = "viewMaxs", stat = "slice", lower = 5)
p2 <- autoplot(xRle, type = "viewMaxs", stat = "slice", lower = 5, geom = "heatmap")
tracks('bar' = p1, "heatmap" = p2)
###################################################
### code chunk number 22: rlel-simul
###################################################
xRleList <- RleList(xRle, 2L * xRle)
xRleList
## SimpleRleList of length 2
## [[1]]
## numeric-Rle of length 10000 with 823 runs
## Lengths: 779 1 208 1 1599 1 ... 1 5 2 9 1 4507
## Values : 0 1 0 1 0 1 ... 1 0 1 0 1 0
##
## [[2]]
## numeric-Rle of length 10000 with 823 runs
## Lengths: 779 1 208 1 1599 1 ... 1 5 2 9 1 4507
## Values : 0 2 0 2 0 2 ... 2 0 2 0 2 0
###################################################
### code chunk number 23: rlel-bin
###################################################
p1 <- autoplot(xRleList)
## Default use binwidth: range/30
p2 <- autoplot(xRleList, nbin = 80)
## Default use binwidth: range/80
p3 <- autoplot(xRleList, geom = "heatmap", nbin = 200)
## Default use binwidth: range/200
tracks('nbin = 30' = p1, "nbin = 80" = p2, "nbin = 200(heatmap)" = p3)
###################################################
### code chunk number 24: rlel-id
###################################################
p1 <- autoplot(xRleList, stat = "identity")
p2 <- autoplot(xRleList, stat = "identity", geom = "point", color = "red")
tracks('line' = p1, "point" = p2)
###################################################
### code chunk number 25: rlel-slice
###################################################
p1 <- autoplot(xRleList, type = "viewMaxs", stat = "slice", lower = 5)
p2 <- autoplot(xRleList, type = "viewMaxs", stat = "slice", lower = 5, geom = "heatmap")
tracks('bar' = p1, "heatmap" = p2)
###################################################
### code chunk number 26: txdb
###################################################
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
## Loading required package: GenomicFeatures
## Loading required package: AnnotationDbi
## Loading required package: Biobase
## Welcome to Bioconductor
##
## Vignettes contain introductory material; view with 'browseVignettes()'. To
## cite Bioconductor, see 'citation("Biobase")', and for packages
## 'citation("pkgname")'.
data(genesymbol, package = "biovizBase")
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
###################################################
### code chunk number 27: txdb-visual
###################################################
p1 <- autoplot(txdb, which = genesymbol["ALDOA"], names.expr = "tx_name:::gene_id")
## Warning: 'function (data, ...) standardGeneric("stat_gene")' is
## deprecated. Use 'geom_alignment' instead. See help("Deprecated")
## Aggregating TranscriptDb...
## Parsing exons...
## Parsing cds...
## Parsing transcripts...
## Aggregating...
## Done
## Constructing graphics...
p2 <- autoplot(txdb, which = genesymbol["ALDOA"], stat = "reduce", color = "brown",
fill = "brown")
## Warning: 'function (data, ...) standardGeneric("stat_gene")' is
## deprecated. Use 'geom_alignment' instead. See help("Deprecated")
## Aggregating TranscriptDb...
## Parsing exons...
## Parsing cds...
## Parsing transcripts...
## Aggregating...
## Done
## Constructing graphics...
tracks(full = p1, reduce = p2, heights = c(5, 1)) + ylab("")
###################################################
### code chunk number 28: ga-load
###################################################
library(Rsamtools)
## Loading required package: Biostrings
data("genesymbol", package = "biovizBase")
bamfile <- system.file("extdata", "SRR027894subRBM17.bam", package="biovizBase")
## need to set use.names = TRUE
ga <- readBamGappedAlignments(bamfile,
param = ScanBamParam(which = genesymbol["RBM17"]),
use.names = TRUE)
###################################################
### code chunk number 29: ga-exp
###################################################
p1 <- autoplot(ga)
p2 <- autoplot(ga, geom = "rect")
## extracting information...
p3 <- autoplot(ga, geom = "line", stat = "coverage")
## extracting information...
tracks(default = p1, rect = p2, coverage = p3)
###################################################
### code chunk number 30: bf-load (eval = FALSE)
###################################################
## library(Rsamtools)
## bamfile <- "./wgEncodeCaltechRnaSeqK562R1x75dAlignsRep1V2.bam"
## bf <- BamFile(bamfile)
###################################################
### code chunk number 31: bf-est-cov (eval = FALSE)
###################################################
## autoplot(bamfile)
## autoplot(bamfile, which = c("chr1", "chr2"))
## autoplot(bf)
## autoplot(bf, which = c("chr1", "chr2"))
##
## data(genesymbol, package = "biovizBase")
## autoplot(bamfile, method = "raw", which = genesymbol["ALDOA"])
##
## library(BSgenome.Hsapiens.UCSC.hg19)
## autoplot(bf, stat = "mismatch", which = genesymbol["ALDOA"], bsgenome = Hsapiens)
###################################################
### code chunk number 32: char-bam (eval = FALSE)
###################################################
## bamfile <- "./wgEncodeCaltechRnaSeqK562R1x75dAlignsRep1V2.bam"
## autoplot(bamfile)
###################################################
### code chunk number 33: char-gr
###################################################
library(rtracklayer)
test_path <- system.file("tests", package = "rtracklayer")
test_bed <- file.path(test_path, "test.bed")
autoplot(test_bed, aes(fill = name))
## reading in
## use score as y by default
###################################################
### matrix
###################################################
volcano <- volcano[20:70, 20:60] - 150
autoplot(volcano)
autoplot(volcano, xlab = "xlab", main = "main", ylab = "ylab")
## special scale theme for 0-centered values
autoplot(volcano, geom = "raster")+scale_fill_fold_change()
## when a matrix has colnames and rownames label them by default
colnames(volcano) <- sort(sample(1:300, size = ncol(volcano), replace = FALSE))
autoplot(volcano)+scale_fill_fold_change()
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
rownames(volcano) <- letters[sample(1:24, size = nrow(volcano), replace = TRUE)]
autoplot(volcano)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## even with row/col names, you could also disable it and just use numeric index
autoplot(volcano, colnames.label = FALSE)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
autoplot(volcano, rownames.label = FALSE, colnames.label = FALSE)
## don't want the axis has label??
autoplot(volcano, axis.text.x = FALSE)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(volcano, axis.text.y = FALSE)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
# or totally remove axis
colnames(volcano) <- lapply(letters[sample(1:24, size = ncol(volcano),
replace = TRUE)],
function(x){
paste(rep(x, 7), collapse = "")
})
## Oops, overlapped
autoplot(volcano)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## tweak with it.
autoplot(volcano, axis.text.angle = -45, hjust = 0)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## when character is the value
x <- sample(c(letters[1:3], NA), size = 100, replace = TRUE)
mx <- matrix(x, nrow = 5)
autoplot(mx)
## tile gives you a white margin
rownames(mx) <- LETTERS[1:5]
autoplot(mx, color = "white")
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
colnames(mx) <- LETTERS[1:20]
autoplot(mx, color = "white")
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(mx, color = "white", size = 2)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## weird in aes(), though works
## default tile is flexible
autoplot(mx, aes(width = 0.6, height = 0.6))
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(mx, aes(width = 0.6, height = 0.6), na.value = "white")
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(mx, aes(width = 0.6, height = 0.6)) + theme_clear()
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
###################################################
### Views
###################################################
lambda <- c(rep(0.001, 4500), seq(0.001, 10, length = 500),
seq(10, 0.001, length = 500))
xVector <- dnorm(1:5e3, mean = 1e3, sd = 200)
xRle <- Rle(xVector)
v1 <- Views(xRle, start = sample(.4e3:.6e3, size = 50, replace = FALSE), width =1000)
autoplot(v1)
names(v1) <- letters[sample(1:24, size = length(v1), replace = TRUE)]
autoplot(v1)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
autoplot(v1, geom = "tile", aes(width = 0.5, height = 0.5))
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
autoplot(v1, geom = "line")
autoplot(v1, geom = "line", aes(color = row)) + theme(legend.position = "none")
autoplot(v1, geom = "line", facets = NULL)
autoplot(v1, geom = "line", facets = NULL, alpha = 0.1)
###################################################
### ExpressionSet
###################################################
library(Biobase)
data(sample.ExpressionSet)
sample.ExpressionSet
## ExpressionSet (storageMode: lockedEnvironment)
## assayData: 500 features, 26 samples
## element names: exprs, se.exprs
## protocolData: none
## phenoData
## sampleNames: A B ... Z (26 total)
## varLabels: sex type score
## varMetadata: labelDescription
## featureData: none
## experimentData: use 'experimentData(object)'
## Annotation: hgu95av2
set.seed(1)
## select 50 features
idx <- sample(seq_len(dim(sample.ExpressionSet)[1]), size = 50)
eset <- sample.ExpressionSet[idx,]
eset
## ExpressionSet (storageMode: lockedEnvironment)
## assayData: 50 features, 26 samples
## element names: exprs, se.exprs
## protocolData: none
## phenoData
## sampleNames: A B ... Z (26 total)
## varLabels: sex type score
## varMetadata: labelDescription
## featureData: none
## experimentData: use 'experimentData(object)'
## Annotation: hgu95av2
autoplot(as.matrix(pData(eset)))
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## default heatmap
p1 <- autoplot(eset)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
p2 <- p1 + scale_fill_fold_change()
p2
autoplot(eset)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(eset, geom = "tile", color = "white", size = 2)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(eset, geom = "tile", aes(width = 0.6, height = 0.6))
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(eset, pheno.plot = TRUE)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
idx <- order(pData(eset)[,1])
eset2 <- eset[,idx]
autoplot(eset2, pheno.plot = TRUE)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## parallel coordainte plot
autoplot(eset, type = "pcp")
## boxplot
autoplot(eset, type = "boxplot")
## scatterplot.matrix
## slow, be carefull
autoplot(eset[, 1:7], type = "scatterplot.matrix")
## mean-sd
autoplot(eset, type = "mean-sd")
## Loading required package: vsn
## volcano
autoplot(eset, type = "volcano", fac = pData(sample.ExpressionSet)$type)
## Loading required package: genefilter
## genefilter::rowttests used
###################################################
### SummarizedExperiment
###################################################
library(GenomicRanges)
nrows <- 200; ncols <- 6
counts <- matrix(runif(nrows * ncols, 1, 1e4), nrows)
counts2 <- matrix(runif(nrows * ncols, 1, 1e4), nrows)
rowData <- GRanges(rep(c("chr1", "chr2"), c(50, 150)),
IRanges(floor(runif(200, 1e5, 1e6)), width=100),
strand=sample(c("+", "-"), 200, TRUE))
colData <- DataFrame(Treatment=rep(c("ChIP", "Input"), 3),
row.names=LETTERS[1:6])
sset <- SummarizedExperiment(assays=SimpleList(counts=counts,
counts2 = counts2),
rowData=rowData, colData=colData)
autoplot(sset) + scale_fill_fold_change()
## Assay index: 1 used
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(sset, pheno.plot = TRUE)
## Assay index: 1 used
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
###################################################
### pcp
###################################################
autoplot(sset, type = "pcp")
## Assay index: 1 used
###################################################
### boxplot
###################################################
autoplot(sset, type = "boxplot")
## Assay index: 1 used
###################################################
### scatterplot matrix
###################################################
autoplot(sset, type = "scatterplot.matrix")
## Assay index: 1 used
###################################################
### vcf
###################################################
library(VariantAnnotation)
## Attaching package: 'VariantAnnotation'
## The following object(s) are masked from 'package:Biobase':
##
## samples
vcffile <- system.file("extdata", "chr22.vcf.gz", package="VariantAnnotation")
vcf <- readVcf(vcffile, "hg19")
## default use type 'geno'
## default use genome position
autoplot(vcf)
## GT,DS,GL could be used for 'geno' type
## use GT for type geno as default
## Index: 744,3604,4738,7096,7287,10037,10198 snp with duplicated start
## position may be masked by each other
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## or disable it
autoplot(vcf, genome.axis = FALSE)
## GT,DS,GL could be used for 'geno' type
## use GT for type geno as default
## Index: 744,3604,4738,7096,7287,10037,10198 snp with duplicated start
## position may be masked by each other
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## not transpose
autoplot(vcf, genome.axis = FALSE, transpose = FALSE, rownames.label = FALSE)
## GT,DS,GL could be used for 'geno' type
## use GT for type geno as default
## Index: 744,3604,4738,7096,7287,10037,10198 snp with duplicated start
## position may be masked by each other
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
autoplot(vcf)
## GT,DS,GL could be used for 'geno' type
## use GT for type geno as default
## Index: 744,3604,4738,7096,7287,10037,10198 snp with duplicated start
## position may be masked by each other
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## use
autoplot(vcf, assay.id = "DS")
## GT,DS,GL could be used for 'geno' type
## Index: 744,3604,4738,7096,7287,10037,10198 snp with duplicated start
## position may be masked by each other
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## equivalent to
autoplot(vcf, assay.id = 2)
## GT,DS,GL could be used for 'geno' type
## Index: 744,3604,4738,7096,7287,10037,10198 snp with duplicated start
## position may be masked by each other
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
## Not run:
##D ## doesn't work when assay.id cannot find
##D autoplot(vcf, assay.id = "NO")
## End(Not run)
## use AF or first
autoplot(vcf, type = "info")
## use AF for type info as default
## Other options for potential mapping(only keep
## numeric/integer/character/factor variable):
## geom bar
autoplot(vcf, type = "info", aes(y = THETA))
## Other options for potential mapping(only keep
## numeric/integer/character/factor variable):
autoplot(vcf, type = "info", aes(y = THETA, fill = VT, color = VT))
## Other options for potential mapping(only keep
## numeric/integer/character/factor variable):
autoplot(vcf, type = "fixed")
autoplot(vcf, type = "fixed", size = 10) + xlim(c(50310860, 50310890)) + ylim(0.75, 1.25)
## Warning: Removed 7 rows containing missing values (geom_text).
## Warning: Removed 7 rows containing missing values (geom_text).
p1 <- autoplot(vcf, type = "fixed") + xlim(50310860, 50310890)
p2 <- autoplot(vcf, type = "fixed", full.string = TRUE) + xlim(50310860, 50310890)
tracks("full.string = FALSE" = p1, "full.string = TRUE" = p2)+
scale_y_continuous(breaks = NULL, limits = c(0, 3))
## Warning: Removed 1 rows containing missing values (geom_text).
p3 <- autoplot(vcf, type = "fixed", ref.show = FALSE) + xlim(50310860, 50310890) +
scale_y_continuous(breaks = NULL, limits = c(0, 2))
p3
###################################################
### code chunk number 56: bs-v
###################################################
library(BSgenome.Hsapiens.UCSC.hg19)
## Loading required package: BSgenome
## Attaching package: 'BSgenome'
## The following object(s) are masked from 'package:AnnotationDbi':
##
## species
data(genesymbol, package = "biovizBase")
p1 <- autoplot(Hsapiens, which = resize(genesymbol["ALDOA"], width = 50))
p2 <- autoplot(Hsapiens, which = resize(genesymbol["ALDOA"], width = 50), geom = "rect")
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.
tracks(text = p1, rect = p2)
###################################################
### code chunk number 57: sessionInfo
###################################################
sessionInfo()
## R Under development (unstable) (2012-07-17 r59871)
## Platform: x86_64-unknown-linux-gnu (64-bit)
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=en_US.UTF-8
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=C LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] methods stats graphics grDevices utils datasets base
##
## other attached packages:
## [1] BSgenome.Hsapiens.UCSC.hg19_1.3.19
## [2] BSgenome_1.25.9
## [3] VariantAnnotation_1.3.32
## [4] genefilter_1.39.0
## [5] vsn_3.25.0
## [6] rtracklayer_1.17.21
## [7] Rsamtools_1.9.31
## [8] Biostrings_2.25.12
## [9] TxDb.Hsapiens.UCSC.hg19.knownGene_2.8.0
## [10] GenomicFeatures_1.9.44
## [11] AnnotationDbi_1.19.46
## [12] Biobase_2.17.8
## [13] GenomicRanges_1.9.66
## [14] IRanges_1.15.48
## [15] BiocGenerics_0.3.2
## [16] ggbio_1.5.20
## [17] ggplot2_0.9.2.1
## [18] knitr_0.8
##
## loaded via a namespace (and not attached):
## [1] affy_1.35.1 affyio_1.25.0 annotate_1.35.5
## [4] BiocInstaller_1.5.12 biomaRt_2.13.2 biovizBase_1.5.9
## [7] bitops_1.0-4.1 cluster_1.14.2 colorspace_1.1-1
## [10] DBI_0.2-5 dichromat_1.2-4 digest_0.5.1
## [13] evaluate_0.4.2 formatR_0.6 grid_2.16.0
## [16] gridExtra_0.9.1 gtable_0.1.1 Hmisc_3.9-3
## [19] labeling_0.1 lattice_0.20-6 limma_3.13.24
## [22] markdown_0.5.2 MASS_7.3-19 memoise_0.1
## [25] munsell_0.4 parallel_2.16.0 plyr_1.7.1
## [28] preprocessCore_1.19.0 proto_0.3-9.2 RColorBrewer_1.0-5
## [31] RCurl_1.95-0 reshape2_1.2.1 RSQLite_0.11.2
## [34] scales_0.2.2 splines_2.16.0 stats4_2.16.0
## [37] stringr_0.6.1 survival_2.36-14 tools_2.16.0
## [40] XML_3.95-0 xtable_1.7-0 zlibbioc_1.3.0