Partition a set of Genomic Ranges by another

partitionRanges(x, y, ...)

# S4 method for class 'GRanges,GRanges'
partitionRanges(
  x,
  y,
  y_as_both = TRUE,
  ignore.strand = FALSE,
  simplify = TRUE,
  suffix = c(".x", ".y"),
  ...
)

Arguments

x, y

GenomicRanges objects

...

Not used

y_as_both

logical(1) If there are any unstranded regions in y, should these be assigned to both strands. If TRUE unstranded regions can be used to partition stranded regions

ignore.strand

If set to TRUE, then the strand of x and y is set to "*" prior to any computation.

simplify

Pass to chopMC and simplify mcols in the output

suffix

Added to any shared column names in the provided objects

Value

A GRanges object

Details

The query set of ranges can be broken in regions which strictly overlap a second set of ranges. The complete set of mcols from both initial objects will included in the set of partitioned ranges

Examples

x <- GRanges(c("chr1:1-10", "chr1:6-15"))
x$id <- paste0("range", seq_along(x))
x
#> GRanges object with 2 ranges and 1 metadata column:
#>       seqnames    ranges strand |          id
#>          <Rle> <IRanges>  <Rle> | <character>
#>   [1]     chr1      1-10      * |      range1
#>   [2]     chr1      6-15      * |      range2
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
y <- GRanges(c("chr1:2-5", "chr1:6-12"))
y$id <- paste0("range", seq_along(y))
y
#> GRanges object with 2 ranges and 1 metadata column:
#>       seqnames    ranges strand |          id
#>          <Rle> <IRanges>  <Rle> | <character>
#>   [1]     chr1       2-5      * |      range1
#>   [2]     chr1      6-12      * |      range2
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
partitionRanges(x, y)
#> GRanges object with 5 ranges and 2 metadata columns:
#>       seqnames    ranges strand |        id.y        id.x
#>          <Rle> <IRanges>  <Rle> | <character> <character>
#>   [1]     chr1         1      * |        <NA>      range1
#>   [2]     chr1       2-5      * |      range1      range1
#>   [3]     chr1      6-10      * |      range2      range1
#>   [4]     chr1      6-12      * |      range2      range2
#>   [5]     chr1     13-15      * |        <NA>      range2
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths