geom_chevron {ggbio}R Documentation

Chevron geoms for GRanges object

Description

Break normal intervals stroed in GRanges object and show them as chevron, useful for showing model or splice summary.

Usage

## S4 method for signature 'GRanges'
geom_chevron(data, ..., xlab, ylab, main,
             offset = 0.1,
             facets = NULL,
             stat = c("stepping", "identity"),
             chevron.height.rescale = c(0.1, 0.8),
             group.selfish = TRUE)

Arguments

data

A GRanges object.

...

Extra parameters passed to autoplot function.

xlab

Label for x

ylab

Label for y

main

Title for plot.

offset

A nunmeric value or characters. If it's numeric value, indicate how much you want the chevron to wiggle, usually the rectangle for drawing GRanges is of height unit 1, so it's better between -0.5 and 0.5 to make it nice looking. Unless you specify offset as one of those columns, this will use height of the chevron to indicate the columns. Of course you could use size of the chevron to indicate the column variable easily, please see the examples.

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.

chevron.height.rescale

A numeric vector of length 2. When the offset parameters is a character which is one of the data columns, this parameter rescale the offset.

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.

Details

To draw a normal GRanges as Chevron, we need to provide a special geom for this purpose. Chevron is popular in gene viewer or genomoe browser, when they try to show isoforms or gene model.geom_chevron, just like any other geom_* function in ggplot2, you can pass aes() to it to use height of chevron or width of chevron to show statistics summary.

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_chevron()

plot of chunk unnamed-chunk-1

## or
ggplot() + geom_chevron(gr)

plot of chunk unnamed-chunk-1

## ======================================================================
##  facetting and aesthetics
## ======================================================================
ggplot(gr) + geom_chevron(facets = sample ~ seqnames, aes(color = strand))

plot of chunk unnamed-chunk-1

## ======================================================================
##  stat:identity
## ======================================================================
ggplot(gr) + geom_chevron(stat = "identity", aes(y = value))

plot of chunk unnamed-chunk-1

## ======================================================================
##  stat:stepping
## ======================================================================
ggplot(gr) + geom_chevron(stat = "stepping", aes(group = pair))

plot of chunk unnamed-chunk-1

## ======================================================================
##  group.selfish controls when 
## ======================================================================
ggplot(gr) + geom_chevron(stat = "stepping", aes(group = pair), group.selfish = FALSE,
                        xlab = "xlab", ylab = "ylab", main = "main")

plot of chunk unnamed-chunk-1

p <- qplot(x = mpg, y = cyl, data = mtcars)
## ======================================================================
##  offset
## ======================================================================
gr2 <- GRanges("chr1", IRanges(c(1, 10, 20), width = 5))
gr2.p <- gaps(gr2)
## resize to connect them
gr2.p <- resize(gr2.p, fix = "center", width = width(gr2.p)+2)
ggplot(gr2) + geom_rect() + geom_chevron(gr2.p)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

## notice the rectangle height is 0.8
## offset = 0 just like a line
ggplot(gr2) + geom_rect() + geom_chevron(gr2.p, offset = 0)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

## equal height
ggplot(gr2) + geom_rect() + geom_chevron(gr2.p, offset = 0.4)
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

## ======================================================================
##  chevron.height
## ======================================================================
values(gr2.p)$score <- c(100, 200)
ggplot(gr2) + geom_rect() + geom_chevron(gr2.p, offset = "score")
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

## chevron.height
ggplot(gr2) + geom_rect() + geom_chevron(gr2.p, offset = "score",
                                         chevron.height.rescale = c(0.4, 10))
## Scale for 'y' is already present. Adding another scale for 'y', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1



[Package ggbio version 1.5.20 ]