Skip to contents

Merge GRanges or a list of GRanges

Usage

MergeGRanges(..., sortRanges = FALSE, reduceRanges = FALSE)

Arguments

...

<GRanges or GRangesList or listGRanges>: Some GRanges or a list of GRanges or a GRangesList.

sortRanges

: Whether the result should be sorted. (Default FALSE)

reduceRanges

: Whether the result should be reduced. See GenomicRanges::reduce for more details. (Default FALSE)

Value

a GRange object.

Details

MergeGRanges

Examples

GRange_1.grn <- GenomicRanges::GRanges(
    seqnames = S4Vectors::Rle(c("chr1", "chr2", "chr1"), c(1, 3, 1)),
    ranges = IRanges::IRanges(101:105, end = 111:115,
        names = letters[seq_len(5)]),
    strand = S4Vectors::Rle(BiocGenerics::strand(c("-", "+", "*", "+")),
        c(1, 1, 2, 1)),
    score = seq_len(5)
)
GRange_2.grn <- GenomicRanges::GRanges(
    seqnames = S4Vectors::Rle(c("chr1", "chr3"), c(1, 4)),
    ranges = IRanges::IRanges(106:110, end = 116:120, names = letters[6:10]),
    strand = S4Vectors::Rle(BiocGenerics::strand(c("*", "+", "-")),
        c(2, 1, 2)),
    score = 6:10
)
GRange_1.grn
#> GRanges object with 5 ranges and 1 metadata column:
#>     seqnames    ranges strand |     score
#>        <Rle> <IRanges>  <Rle> | <integer>
#>   a     chr1   101-111      - |         1
#>   b     chr2   102-112      + |         2
#>   c     chr2   103-113      * |         3
#>   d     chr2   104-114      * |         4
#>   e     chr1   105-115      + |         5
#>   -------
#>   seqinfo: 2 sequences from an unspecified genome; no seqlengths
GRange_2.grn
#> GRanges object with 5 ranges and 1 metadata column:
#>     seqnames    ranges strand |     score
#>        <Rle> <IRanges>  <Rle> | <integer>
#>   f     chr1   106-116      * |         6
#>   g     chr3   107-117      * |         7
#>   h     chr3   108-118      + |         8
#>   i     chr3   109-119      - |         9
#>   j     chr3   110-120      - |        10
#>   -------
#>   seqinfo: 2 sequences from an unspecified genome; no seqlengths
MergeGRanges(GRange_1.grn, GRange_2.grn)
#> GRanges object with 10 ranges and 1 metadata column:
#>     seqnames    ranges strand |     score
#>        <Rle> <IRanges>  <Rle> | <integer>
#>   a     chr1   101-111      - |         1
#>   b     chr2   102-112      + |         2
#>   c     chr2   103-113      * |         3
#>   d     chr2   104-114      * |         4
#>   e     chr1   105-115      + |         5
#>   f     chr1   106-116      * |         6
#>   g     chr3   107-117      * |         7
#>   h     chr3   108-118      + |         8
#>   i     chr3   109-119      - |         9
#>   j     chr3   110-120      - |        10
#>   -------
#>   seqinfo: 3 sequences from an unspecified genome; no seqlengths
GRange.lst <- list(GRange_1.grn, GRange_2.grn)
MergeGRanges(GRange.lst)
#> GRanges object with 10 ranges and 1 metadata column:
#>     seqnames    ranges strand |     score
#>        <Rle> <IRanges>  <Rle> | <integer>
#>   a     chr1   101-111      - |         1
#>   b     chr2   102-112      + |         2
#>   c     chr2   103-113      * |         3
#>   d     chr2   104-114      * |         4
#>   e     chr1   105-115      + |         5
#>   f     chr1   106-116      * |         6
#>   g     chr3   107-117      * |         7
#>   h     chr3   108-118      + |         8
#>   i     chr3   109-119      - |         9
#>   j     chr3   110-120      - |        10
#>   -------
#>   seqinfo: 3 sequences from an unspecified genome; no seqlengths
MergeGRanges(GRange.lst, reduceRanges = TRUE)
#> GRanges object with 8 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr1   105-115      +
#>   [2]     chr1   101-111      -
#>   [3]     chr1   106-116      *
#>   [4]     chr2   102-112      +
#>   [5]     chr2   103-114      *
#>   [6]     chr3   108-118      +
#>   [7]     chr3   109-120      -
#>   [8]     chr3   107-117      *
#>   -------
#>   seqinfo: 3 sequences from an unspecified genome; no seqlengths