Last data update: 2014.03.03

R: Heatmap representing signals in given ranges
featureAlignedHeatmapR Documentation

Heatmap representing signals in given ranges

Description

plot heatmap in the given feature ranges

Usage

featureAlignedHeatmap(cvglists, feature.gr, upstream, downstream, 
             zeroAt, n.tile=100, 
             annoMcols=c(), sortBy=names(cvglists)[1],
             color=colorRampPalette(c("yellow", "red"))(50),
             lower.extreme, upper.extreme,
             margin=c(0.1, 0.01, 0.15, 0.1), gap=0.01, 
             newpage=TRUE, gp=gpar(fontsize=10),
             ...)

Arguments

cvglists

Output of featureAlignedSignal or a list of SimpleRleList or RleList

feature.gr

An object of GRanges with identical width. If the width equal to 1, you can use upstream and downstream to set the range for plot. If the width not equal to 1, you can use zeroAt to set the zero point of the heatmap.

upstream, downstream

upstream or dwonstream from the feature.gr. It must keep same as featureAlignedSignal. It is used for x-axis label.

zeroAt

zero point position of feature.gr

n.tile

The number of tiles to generate for each element of feature.gr, default is 100

annoMcols

The columns of metadata of feature.gr that specifies the annotations shown of the right side of the heatmap.

sortBy

Sort the feature.gr by columns by annoMcols and then the signals of the given samples. Default is the first sample. Set to NULL to disable sort.

color

vector of colors used in heatmap

lower.extreme, upper.extreme

The lower and upper boundary value of each samples

margin

Margin for of the plot region.

gap

Gap between each heatmap columns.

newpage

Call grid.newpage or not. Default, TRUE

gp

A gpar object can be used for text.

...

Not used.

Value

invisible gList object.

Author(s)

Jianhong Ou

See Also

See Also as featureAlignedSignal, featureAlignedDistribution

Examples

  cvglists <- list(A=RleList(chr1=Rle(sample.int(5000, 100), 
                                      sample.int(300, 100))), 
                   B=RleList(chr1=Rle(sample.int(5000, 100), 
                                      sample.int(300, 100))))
  feature.gr <- GRanges("chr1", IRanges(seq(1, 4900, 100), width=100))
  feature.gr$anno <- rep(c("type1", "type2"), c(25, 24))
  featureAlignedHeatmap(cvglists, feature.gr, zeroAt=50, annoMcols="anno")

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(ChIPpeakAnno)
Loading required package: grid
Loading required package: IRanges
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: 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: Biostrings
Loading required package: XVector
Loading required package: GenomicRanges
Loading required package: GenomeInfoDb
Loading required package: VennDiagram
Loading required package: futile.logger

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/ChIPpeakAnno/featureAlignedHeatmap.Rd_%03d_medium.png", width=480, height=480)
> ### Name: featureAlignedHeatmap
> ### Title: Heatmap representing signals in given ranges
> ### Aliases: featureAlignedHeatmap
> ### Keywords: misc
> 
> ### ** Examples
> 
>   cvglists <- list(A=RleList(chr1=Rle(sample.int(5000, 100), 
+                                       sample.int(300, 100))), 
+                    B=RleList(chr1=Rle(sample.int(5000, 100), 
+                                       sample.int(300, 100))))
>   feature.gr <- GRanges("chr1", IRanges(seq(1, 4900, 100), width=100))
>   feature.gr$anno <- rep(c("type1", "type2"), c(25, 24))
>   featureAlignedHeatmap(cvglists, feature.gr, zeroAt=50, annoMcols="anno")
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>