Skip to contents

Given a GRanges object, pack produces a GRangesList of the same ranges grouped and re-ordered.

Usage

# S4 method for class 'GRanges'
pack(x, ..., range_len = 1e9, inter_range_len = 1e7)

Arguments

x

A GRanges object.

range_len

A numeric specifying the max length allowed for ranges in x.

inter_range_len

A numeric specifying the max length allowed between ranges in x.

...

Arguments passed to other methods.

Details

Packing ranges

The pack method attempts to re-package ranges in optimal form for extracting data from files. Ranges are not modified (made shorter or longer) but re-ordered and / or re-grouped according to the following criteria.

  • order: Ranges are ordered by genomic position within chromosomes.

  • distance: Ranges separted by a distance greater than the inter_range_len are packed in groups around the gap separating the distant ranges.

  • length: Ranges longer than range_len are packed `individually' (i.e., retrived from the file as a single range vs grouped with other ranges).

Utilities

isPacked(x, ...): Returns a logical indicating if the ranges in x are packed. x must be a GRangesList object.

Value

A GRanges object.

See also

  • unpack for unpacking the result obtained with `packed' ranges.

Examples

  ## Ranges are ordered by position within chromosome.
  gr1 <- GRanges("chr1", IRanges(5:1*5, width = 3)) 
  pack(gr1)
#> GRangesList object of length 1:
#> [[1]]
#> GRanges object with 5 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr1       5-7      *
#>   [2]     chr1     10-12      *
#>   [3]     chr1     15-17      *
#>   [4]     chr1     20-22      *
#>   [5]     chr1     25-27      *
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#> 
  
  ## Ranges separated by > inter_range_len are partitioned
  ## into groups defined by the endpoints of the gap.
  gr2 <- GRanges("chr2", IRanges(c(1:3, 30000:30003), width = 1000))
  pack(gr2, inter_range_len = 20000)
#> GRangesList object of length 2:
#> [[1]]
#> GRanges object with 3 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr2    1-1000      *
#>   [2]     chr2    2-1001      *
#>   [3]     chr2    3-1002      *
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#> 
#> [[2]]
#> GRanges object with 4 ranges and 0 metadata columns:
#>       seqnames      ranges strand
#>          <Rle>   <IRanges>  <Rle>
#>   [1]     chr2 30000-30999      *
#>   [2]     chr2 30001-31000      *
#>   [3]     chr2 30002-31001      *
#>   [4]     chr2 30003-31002      *
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#> 
  
  ## Ranges exceeding 'range_len' are isolated in a single element
  ## of the GRangesList.
  gr3 <- GRanges("chr3", IRanges(c(1:4), width=c(45, 1e8, 45, 45)))
  width(gr3)
#> [1]        45 100000000        45        45
  pack(gr3, range_len = 1e7)
#> GRangesList object of length 3:
#> [[1]]
#> GRanges object with 1 range and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr3      1-45      *
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#> 
#> [[2]]
#> GRanges object with 1 range and 0 metadata columns:
#>       seqnames      ranges strand
#>          <Rle>   <IRanges>  <Rle>
#>   [1]     chr3 2-100000001      *
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#> 
#> [[3]]
#> GRanges object with 2 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr3      3-47      *
#>   [2]     chr3      4-48      *
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
#>