stdPeakLocation: Peak density with respect to closest gene.

Description Usage Arguments Value Methods Examples

Description

stdPeakLocation plots the density of peaks with respect to the genomic feature (e.g. gene) in standardized gene coordinates so that genes with different lengths are comparable.

PeakLocation produces the same plot in non-standardized coordinates (i.e. distances are measured in base pairs).

plotMeanCoverage plots the mean coverage in a series of selected genomic regions.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
stdPeakLocation(x, peakDistance=1000, startpos='start_position', endpos='end_position',
strand='strand', distance, main='', xlab='Distance relative to feature length', xaxt='n',
  xlim=c(-1,2), densityType="kernel", nbreaks=10, ...)

PeakLocation(x, peakDistance=1000, startpos='start_position', endpos='end_position',
strand='strand', distance, main='', xlab='Distance (bp)',
  densityType="kernel", breaks, ...)

plotMeanCoverage(cover, x, upstreambp=1000, downstreambp=5000,
startpos='start_position', endpos='end_position', normalize=FALSE,
smooth=FALSE, span=0.05, main='', xlab='(bp)', ylab='Average coverage', ...)

Arguments

x

A RangedData or data.frame indicating peak start and end in start and end, and start and end of the closest genomic feature (e.g. gene) in startpos and endpos.

cover

An RleList object containing the coverage, as returned by the function coverage.

peakDistance

Peaks more than peakDistance bases upstream or more than 3*peakDistance downstream of the closest feature are discarded.

startpos

Name of the variable storing the start position of the closest genomic feature.

endpos

Name of the variable storing the end position of the closest genomic feature.

strand

Name of the variable storing the strand for the closest genomic feature.

distance

Name of the variable indicating the distance between the peak and the closest genomic feature. If left missing the distance between the feature start and the mid-point of the peak is computed.

main

Graphical parameter passed on to plot.

xlab

Graphical parameter passed on to plot.

ylab

Graphical parameter passed on to plot.

xaxt

Graphical parameter passed on to plot.

xlim

In stdPeakLocation the x-axis limit is set to xlim*peakDistance.

densityType

If we eant a density plot or a histogram. Has to be one of "kernel" (for the density plot) or "hist" for the histogram.

nbreaks

Number of breaks to be used. It will not be used if densityType is different from "hist".

breaks

This parameter will be passed to the hist plotting function. It will not be used if densityType is different from "hist".

upstreambp

Number of bp upstream of the TSS where the coverage should be computed

downstreambp

Number of bp downstream of the TSS where the coverage should be computed

normalize

When set to TRUE the average coverage in each position is divided by the average across all positions. This is useful when trying to super-impose data from several experiments that had different read coverage.

smooth

If set to TRUE, the average coverage is smooth by calling loess.

span

Parameter controlling smoothing, passed on to loess. Larger values indicate more smoothing.

...

Further parameters passed on to plot.

Value

This function produces a density plot.

Methods

Methods for stdPeakLocation, PeakLocation

signature(x = "data.frame")

The data frame should contain columns named start and end indicating the peak location, txStart, txEnd indicating transcription start/end of the closest gene and strand indicating the strand.

signature(x = "RangedData")

start(x) and end(x) indicate the peak location. x should contain variables x[['txStart']], x[['txEnd']] indicating the transcription start/end of the closest gene and x[['strand']] indicating the strand.

Methods for plotMeanCoverage

signature(cover="RleList", x="RangedData")

cover contains the coverage and x the genomic regions of interest.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
#Generate synthetic peaks
set.seed(1)
st <- runif(100,1,1000)
en <- st+runif(length(st),25,100)
peaks <- RangedData(IRanges(st,en),space='chr1')

#Assign distance to closest gene
#(typically one would call annotatePeakInBatch
#from package ChIPpeakAnno to do this)
peaks[['start_position']] <- start(peaks) + runif(nrow(peaks),-500,1000)
peaks[['end_position']] <- peaks[['start_position']] + 500
peaks[['distance']] <- peaks[['start_position']] - start(peaks)
peaks[['strand']] <- sample(c('+','-'),nrow(peaks),replace=TRUE)
PeakLocation(peaks,peakDistance=1000)

htSeqTools documentation built on May 6, 2019, 3:39 a.m.