Last data update: 2014.03.03
R: Arrowrect geoms for GRanges object
geom_arrowrect R Documentation
Arrowrect geoms for GRanges object
Description
Show interval data as rectangle with a arrow head.
Usage
## S4 method for signature 'GRanges'
geom_arrowrect(data, ..., xlab, ylab, main,
facets = NULL, stat = c("stepping", "identity"),
rect.height = NULL, arrow.head = 0.06,
arrow.head.rate = arrow.head, arrow.head.fix = NULL,
group.selfish = TRUE)
Arguments
data
A GRanges
object.
...
Extra parameters such as aes() passed.
xlab
Label for x
ylab
Label for y
main
Title for plot.
facets
Faceting formula to use.
stat
Character vector specifying statistics to use. "stepping" with
randomly assigned stepping levels as y varialbe. "identity" allow
users to specify y
value in aes
.
rect.height
Half height of the arrow body.
arrow.head
Arrow head to body ratio.
arrow.head.rate
Arrow head to body ratio. same with arrow.head.
arrow.head.fix
fixed length of arrow head.
group.selfish
Passed to addStepping
, control whether to show each group as
unique level or not. If set to FALSE
, if two groups are not
overlapped with each other, they will probably be layout in the same
level to save space.
Value
A 'Layer'.
Author(s)
Tengfei Yin
Examples
set.seed(1)
N <- 100
require(GenomicRanges)
## ======================================================================
## 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(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))
## ======================================================================
## default
## ======================================================================
ggplot(gr) + geom_arrowrect()
## or
ggplot() + geom_arrowrect(gr)
## ======================================================================
## facetting and aesthetics
## ======================================================================
ggplot(gr) + geom_arrowrect(facets = sample ~ seqnames, aes(color = strand, fill = strand))
## ======================================================================
## stat:identity
## ======================================================================
ggplot(gr) + geom_arrowrect(stat = "identity", aes(y = value))
## ======================================================================
## stat:stepping
## ======================================================================
ggplot(gr) + geom_arrowrect(stat = "stepping", aes(y = value, group = pair))
## ======================================================================
## group.selfish controls when
## ======================================================================
ggplot(gr) + geom_arrowrect(gr, stat = "stepping", aes(y = value, group = pair), group.selfish = FALSE)
Results
R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)
R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.
R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.
Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.
> library(ggbio)
Loading required package: BiocGenerics
Loading required package: parallel
Attaching package: 'BiocGenerics'
The following objects are masked from 'package:parallel':
clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
clusterExport, clusterMap, parApply, parCapply, parLapply,
parLapplyLB, parRapply, parSapply, parSapplyLB
The following objects are masked from 'package:stats':
IQR, mad, xtabs
The following objects are masked from 'package:base':
Filter, Find, Map, Position, Reduce, anyDuplicated, append,
as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
rbind, rownames, sapply, setdiff, sort, table, tapply, union,
unique, unsplit
Loading required package: ggplot2
Need specific help about ggbio? try mailing
the maintainer or visit http://tengfei.github.com/ggbio/
Attaching package: 'ggbio'
The following objects are masked from 'package:ggplot2':
geom_bar, geom_rect, geom_segment, ggsave, stat_bin, stat_identity,
xlim
Warning message:
replacing previous import 'ggplot2::Position' by 'BiocGenerics::Position' when loading 'ggbio'
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/ggbio/geom_arrowrect-method.Rd_%03d_medium.png", width=480, height=480)
> ### Name: geom_arrowrect
> ### Title: Arrowrect geoms for GRanges object
> ### Aliases: geom_arrowrect geom_arrowrect,GRanges-method
> ### geom_arrowrect,missing-method geom_arrowrect,uneval-method
>
> ### ** Examples
>
> set.seed(1)
> N <- 100
> require(GenomicRanges)
Loading required package: GenomicRanges
Loading required package: S4Vectors
Loading required package: stats4
Attaching package: 'S4Vectors'
The following objects are masked from 'package:base':
colMeans, colSums, expand.grid, rowMeans, rowSums
Loading required package: IRanges
Loading required package: GenomeInfoDb
> ## ======================================================================
> ## 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(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))
>
>
> ## ======================================================================
> ## default
> ## ======================================================================
> ggplot(gr) + geom_arrowrect()
> ## or
> ggplot() + geom_arrowrect(gr)
>
> ## ======================================================================
> ## facetting and aesthetics
> ## ======================================================================
> ggplot(gr) + geom_arrowrect(facets = sample ~ seqnames, aes(color = strand, fill = strand))
>
>
> ## ======================================================================
> ## stat:identity
> ## ======================================================================
> ggplot(gr) + geom_arrowrect(stat = "identity", aes(y = value))
>
>
> ## ======================================================================
> ## stat:stepping
> ## ======================================================================
> ggplot(gr) + geom_arrowrect(stat = "stepping", aes(y = value, group = pair))
>
>
> ## ======================================================================
> ## group.selfish controls when
> ## ======================================================================
> ggplot(gr) + geom_arrowrect(gr, stat = "stepping", aes(y = value, group = pair), group.selfish = FALSE)
>
>
>
>
>
> dev.off()
null device
1
>