geom_alignment | R Documentation |
Show interval data as alignment.
## S4 method for signature 'GRanges'
geom_alignment(data, ..., xlab, ylab, main, facets = NULL, stat =
c("stepping", "identity"), range.geom = c("rect",
"arrowrect"), gap.geom = c("chevron", "arrow",
"segment"), rect.height = NULL, group.selfish = TRUE,
label = TRUE)
## S4 method for signature 'TxDbOREnsDb'
geom_alignment(data, ..., which, columns = c("tx_id", "tx_name",
"gene_id"), names.expr = "tx_name", facets = NULL,
truncate.gaps = FALSE, truncate.fun = NULL, ratio =
0.0025)
## S4 method for signature 'GRangesList'
geom_alignment(data, ..., which = NULL,
cds.rect.h = 0.25,
exon.rect.h = cds.rect.h,
utr.rect.h = cds.rect.h/2,
xlab, ylab, main,
facets = NULL, geom = "alignment",
stat = c("identity", "reduce"),
range.geom = "rect",
gap.geom = "arrow",
utr.geom = "rect",
names.expr = NULL,
label = TRUE,
label.color = "gray40",
label.size = 3,
arrow.rate = 0.015,
length = unit(0.1, "cm"))
## S4 method for signature 'OrganismDb'
geom_alignment(data, ..., which,
columns = c("TXNAME", "SYMBOL", "TXID", "GENEID"),
names.expr = "SYMBOL",
facets = NULL,
truncate.gaps = FALSE,
truncate.fun = NULL, ratio = 0.0025
)
data |
A |
... |
Extra parameters such as aes() passed. |
which |
|
cds.rect.h |
cds heights. |
exon.rect.h |
exon heights. |
utr.rect.h |
utr heights. |
label.color |
label color. |
label.size |
label size. |
arrow.rate |
arrow rate. |
length |
arrow length. |
columns |
columns to get from object. |
xlab |
Label for x |
ylab |
Label for y |
main |
Title for plot. |
facets |
Faceting formula to use. |
stat |
For For |
gap.geom |
Geom for 'gap' computed from the data you passed based on the group information. |
rect.height |
Half height of the arrow body. |
group.selfish |
Passed to |
truncate.gaps |
logical value indicate to truncate gaps or not. |
truncate.fun |
shrinkage function. Please see |
ratio |
used in |
geom |
geometric object. only support "gene" now. |
range.geom |
geom for main intevals or exons. |
utr.geom |
geom for utr region. |
names.expr |
expression for showing y label. |
label |
logical value. Whether to label the intervals with names specified
by argument |
A 'Layer'.
Tengfei Yin
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_alignment()
## or
ggplot() + geom_alignment(gr)
## ======================================================================
## facetting and aesthetics
## ======================================================================
ggplot(gr) + geom_alignment(facets = sample ~ seqnames, aes(color = strand, fill = strand))
## ======================================================================
## stat:stepping
## ======================================================================
ggplot(gr) + geom_alignment(stat = "stepping", aes(group = pair))
## ======================================================================
## group.selfish controls when
## ======================================================================
ggplot(gr) + geom_alignment(stat = "stepping", aes(group = pair), group.selfish = FALSE)
## =======================================
## main/gap geom
## =======================================
ggplot(gr) + geom_alignment(range.geom = "arrowrect", gap.geom = "chevron")
## =======================================
## For TxDb
## =======================================
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
data(genesymbol, package = "biovizBase")
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
## made a track comparing full/reduce stat.
ggbio() + geom_alignment(data = txdb, which = genesymbol["RBM17"])
p1 <- ggplot(txdb) + geom_alignment(which = genesymbol["RBM17"])
p1
p2 <- ggplot(txdb) + geom_alignment(which = genesymbol["RBM17"], stat = "reduce")
tracks(full = p1, reduce = p2, heights = c(3, 1))
tracks(full = p1, reduce = p2, heights = c(3, 1)) + theme_tracks_sunset()
tracks(full = p1, reduce = p2, heights = c(3, 1)) +
theme_tracks_sunset(axis.line.color = NA)
## change y labels
ggplot(txdb) + geom_alignment(which = genesymbol["RBM17"], names.expr = "tx_id:::gene_id")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.