tile-methods | R Documentation |
tile
and slidingWindows
methods for GenomicRanges. tile
partitions each range
into a set of tiles, which are defined in terms of their number or
width. slidingWindows
generates sliding windows of a specified
width and frequency.
## S4 method for signature 'GenomicRanges'
tile(x, n, width)
## S4 method for signature 'GenomicRanges'
slidingWindows(x, width, step=1L)
x |
A GenomicRanges object, like a |
n |
The number of tiles to generate.
See |
width |
The (maximum) width of each tile.
See |
step |
The distance between the start positions of the sliding windows. |
The tile
function splits x
into a GRangesList
,
each element of which corresponds
to a tile, or partition, of x
. Specify the tile geometry with either
n
or width
(not both). Passing n
creates n
tiles
of approximately equal width, truncated by sequence end, while passing
width
tiles the region with ranges of the given width, again truncated
by sequence end.
The slidingWindows
function generates sliding windows within
each range of x
, according to width
and step
,
returning a GRangesList
. If the sliding windows do not exactly
cover a range in x
, the last window is partial.
A GRangesList
object, each element of which corresponds to a window.
M. Lawrence
tile
in the IRanges package.
gr <- GRanges(
seqnames=Rle(c("chr1", "chr2", "chr1", "chr3"), c(1, 3, 2, 4)),
ranges=IRanges(1:10, end=11),
strand=Rle(strand(c("-", "+", "*", "+", "-")), c(1, 2, 2, 3, 2)),
seqlengths=c(chr1=11, chr2=12, chr3=13))
# split every range in half
tiles <- tile(gr, n = 2L)
stopifnot(all(elementNROWS(tiles) == 2L))
# split ranges into subranges of width 2
# odd width ranges must contain one subrange of width 1
tiles <- tile(gr, width = 2L)
stopifnot(all(all(width(tiles) %in% c(1L, 2L))))
windows <- slidingWindows(gr, width=3L, step=2L)
width(windows[[1L]]) # last range is truncated
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.