Coverage of a GRanges or GRangesList object
coverage
methods for GRanges and
GRangesList objects.
NOTE: The coverage
generic function and methods
for IntegerRanges and IntegerRangesList
objects are defined and documented in the IRanges package.
Methods for GAlignments and
GAlignmentPairs objects are defined and
documented in the GenomicAlignments package.
## S4 method for signature 'GenomicRanges' coverage(x, shift=0L, width=NULL, weight=1L, method=c("auto", "sort", "hash", "naive")) ## S4 method for signature 'GRangesList' coverage(x, shift=0L, width=NULL, weight=1L, method=c("auto", "sort", "hash", "naive"))
x |
A GenomicRanges or GRangesList object. |
shift, weight |
A numeric vector or a list-like object. If numeric, it must be parallel
to Alternatively, each of these arguments can also be specified as a
single string naming a metadata column in See Note that when |
width |
Either See |
method |
See |
When x
is a GRangesList object, coverage(x, ...)
is equivalent to coverage(unlist(x), ...)
.
A named RleList object with one coverage vector per
seqlevel in x
.
H. Pagès and P. Aboyoun
coverage
in the IRanges package.
coverage-methods in the GenomicAlignments package.
RleList objects in the IRanges package.
GRanges, GPos, and GRangesList objects.
## Coverage of a GRanges object: gr <- GRanges( seqnames=Rle(c("chr1", "chr2", "chr1", "chr3"), c(1, 3, 2, 4)), ranges=IRanges(1:10, end=10), strand=Rle(strand(c("-", "+", "*", "+", "-")), c(1, 2, 2, 3, 2)), seqlengths=c(chr1=11, chr2=12, chr3=13)) cvg <- coverage(gr) pcvg <- coverage(gr[strand(gr) == "+"]) mcvg <- coverage(gr[strand(gr) == "-"]) scvg <- coverage(gr[strand(gr) == "*"]) stopifnot(identical(pcvg + mcvg + scvg, cvg)) ## Coverage of a GPos object: pos_runs <- GRanges(c("chr1", "chr1", "chr2"), IRanges(c(1, 5, 9), c(10, 8, 15))) gpos <- GPos(pos_runs) coverage(gpos) ## Coverage of a GRangesList object: gr1 <- GRanges(seqnames="chr2", ranges=IRanges(3, 6), strand = "+") gr2 <- GRanges(seqnames=c("chr1", "chr1"), ranges=IRanges(c(7,13), width=3), strand=c("+", "-")) gr3 <- GRanges(seqnames=c("chr1", "chr2"), ranges=IRanges(c(1, 4), c(3, 9)), strand=c("-", "-")) grl <- GRangesList(gr1=gr1, gr2=gr2, gr3=gr3) stopifnot(identical(coverage(grl), coverage(unlist(grl))))
Please choose more modern alternatives, such as Google Chrome or Mozilla Firefox.