Last data update: 2014.03.03

R: Coverage density plot
coverage.densityR Documentation

Coverage density plot

Description

Visualization of target coverage density for one or more samples.

Usage

coverage.density(coveragelist, normalized = TRUE, legend, main, xlab, col, lwd, lty, xlim, ylim, ...)

Arguments

coveragelist

Output of function coverage.target, where option perBase had to be set to TRUE, i.e. a list with elements coverageTarget and avgTargetCoverage. Or, when density of several samples shall be visualized, a list with respective outputs of coverage.target.

normalized

if TRUE, densities of normalized coverages will be shown; original coverages otherwise

legend

legend text. If missing, names of coveragelist will be taken. If NULL, no legend will be drawn.

main

main title

xlab

x-axis label

col

line color(s)

lwd

line width(s)

lty

line style(s)

xlim, ylim

x- and y-axis coordinate ranges

...

further graphical parameters passed to plot

Details

If normalized = TRUE, the function calculates normalized coverages: per-base coverages divided by average coverage over all targeted bases. Normalized coverages are not dependent on the absolute quantity of reads and are hence better comparable between different samples or even different experiments.

Value

Line plot(s) showing densities.

Author(s)

Manuela Hummel m.hummel@dkfz.de

See Also

coverage.target, covered.k, coverage.hist, coverage.uniformity, coverage.correlation, coverage.plot

Examples

## get reads and targets
exptPath <- system.file("extdata", package="TEQC")
readsfile <- file.path(exptPath, "ExampleSet_Reads.bed")
reads <- get.reads(readsfile, idcol=4, skip=0)
targetsfile <- file.path(exptPath, "ExampleSet_Targets.bed")
targets <- get.targets(targetsfile, skip=0)

## calculate per-base coverages
Coverage <- coverage.target(reads, targets, perBase=TRUE)

## coverage density
coverage.density(Coverage)

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(TEQC)
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: IRanges
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: Rsamtools
Loading required package: GenomeInfoDb
Loading required package: GenomicRanges
Loading required package: Biostrings
Loading required package: XVector
Loading required package: hwriter
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/TEQC/coverage.density.Rd_%03d_medium.png", width=480, height=480)
> ### Name: coverage.density
> ### Title: Coverage density plot
> ### Aliases: coverage.density
> ### Keywords: hplot
> 
> ### ** Examples
> 
> ## get reads and targets
> exptPath <- system.file("extdata", package="TEQC")
> readsfile <- file.path(exptPath, "ExampleSet_Reads.bed")
> reads <- get.reads(readsfile, idcol=4, skip=0)
[1] "read 19546 sequenced reads"
> targetsfile <- file.path(exptPath, "ExampleSet_Targets.bed")
> targets <- get.targets(targetsfile, skip=0)
[1] "read 50 (non-overlapping) target regions"
Warning message:
the "reduce" method for RangedData object is deprecated 
> 
> ## calculate per-base coverages
> Coverage <- coverage.target(reads, targets, perBase=TRUE)
> 
> ## coverage density
> coverage.density(Coverage)
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>