coverage.hist: Coverage histogram

Description Usage Arguments Value Author(s) See Also Examples

View source: R/coverage.hist.R

Description

Histogram and cumulative density of target base coverages

Usage

1
coverage.hist(coverageTarget, col.hist = "lightblue", col.line = "orange", covthreshold, outline = FALSE, breaks = "Sturges", xlab, ylab, main, lwd, ...)

Arguments

coverageTarget

RleList containing Rle vectors of per-target-base coverages for each chromosome, i.e. coverageTarget output from coverage.target

col.hist

histogram color

col.line

color of the cumulative density line

covthreshold

indicates with dashed vertical and horizontal lines, which fraction of bases has a coverage of at least covthreshold; if missing, no dashed lines are drawn

outline

if FALSE, outliers (according to boxplot.stats) are removed before plotting.

breaks

number of cells for the histogram, or string naming an algorithm to compute the number of cells, or function to compute the number of cells, or vector giving the breakpoints between histogram cells (see ?hist) but the latter option only with equidistant breakpoints

xlab, ylab

x- and y-axis labels

main

plot title

lwd

line width

...

further graphical parameters, passed to plot(histogram)

Value

Histogram of read coverages for bases within the target. Additionally, a line and the right axis indicate the cumulative fraction of target bases with coverage of at least x. If option covthreshold is specified, red dashed lines highlight the cumulative fraction of target bases with at least the specified coverage.

Author(s)

Manuela Hummel m.hummel@dkfz.de

See Also

coverage.target, coverage.uniformity, coverage.density, coverage.plot, coverage.targetlength.plot

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
## 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 histogram
coverage.hist(Coverage$coverageTarget, covthreshold=8)

Example output

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, sd, var, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colMeans, colSums, 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, rowMeans, rowSums, rownames, sapply,
    setdiff, sort, table, tapply, union, unique, unsplit, which,
    which.max, which.min

Loading required package: IRanges
Loading required package: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

The following object is masked from 'package:base':

    expand.grid

Loading required package: Rsamtools
Loading required package: GenomeInfoDb
Loading required package: GenomicRanges
Loading required package: Biostrings
Loading required package: XVector

Attaching package: 'Biostrings'

The following object is masked from 'package:base':

    strsplit

Loading required package: hwriter
[1] "read 19546 sequenced reads"
[1] "read 50 (non-overlapping) target regions"

TEQC documentation built on Nov. 8, 2020, 8:07 p.m.