partitionRegions-methods | R Documentation |
Partitions genomic regions into windows of fixed size
## S4 method for signature 'GRanges'
partitionRegions(x, chrs=character(), width=5000, overlap=0.5)
## S4 method for signature 'GRangesList'
partitionRegions(x, chrs=character(), width=5000, overlap=0.5)
## S4 method for signature 'MaskedBSgenome'
partitionRegions(x, chrs=character(), width=5000, overlap=0.5, ...)
x |
an object of class |
chrs |
a character vector (possibly empty) with names of chromosomes to limit to |
width |
window size |
overlap |
amount of overlap; a zero value corresponds to non-overlapping windows and the default 0.5 corresponds to 50% overlap. The largest possible value is 0.8 which corresponds to an overlap of 80%. |
... |
further arguments are passed on to
|
For a GRanges
object x
, this method
partitions each genomic region into possibly overlapping, equally large
windows of size width
. The amount of overlap is controlled with
the overlap
parameter. The windows are placed such that
possible overhangs are balanced at the beginning and end of the
region. As an example, suppose we have a region from bases 1 to 14,000 and
that we want to cover it with windows of 10,000bp length and
50% overlap. The straightforward approach would be to have two
windows 1-10,000 and 5,001-15,000, and to crop the latter to
5,001-14,000. As said, the partitionRegions
balances the
overhangs, so it will return two windows 1-9,500 and 4,501-14,000
instead.
If chrs
is not empty, partitionRegions
will only
consider regions from those chromosomes (i.e. regions in the
GRanges
object whose seqnames
occur
in chrs
).
If called for a GRangesList
object, all
componentes of the GRangesList
object are
partitioned separately as described above.
For convenience, this function can also be called
for a MaskedBSgenome
object. In this case,
unmaskedRegions
is called before partitioning.
If x
is a GRanges
object, the function
also returns a GRanges
object. In the other two
cases, a GRangesList
object is returned.
Ulrich Bodenhofer
https://github.com/UBod/podkat
assocTest
, unmaskedRegions
,
unmasked-datasets
, GRangesList
,
GRanges
## create a toy example
gr <- GRanges(seqnames="chr1", ranges=IRanges(start=1, end=14000))
partitionRegions(gr, width=10000, overlap=0.5)
partitionRegions(gr, width=10000, overlap=0.8)
partitionRegions(gr, width=10000, overlap=0)
## a toy example of a 'GRangesList'
grL <- GRangesList(gr, GRanges(seqnames="chr2",
ranges=IRanges(start=1, end=22000)))
partitionRegions(grL, width=10000, overlap=0.5)
partitionRegions(grL, width=10000, overlap=0.8)
## real-world example
data(hg38Unmasked)
partitionRegions(hg38Unmasked, chrs="chr22", width=20000)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.