Map overlapping features in data onto intervals in scaffold and apply statistics and/or summary operations on those features.

map_bed(
  x,
  y,
  operation,
  max_gap = -1L,
  min_overlap = 0L,
  min_overlap_type = c("bp", "frac1", "frac2")
)

Arguments

x

A GRanges object.

y

A GRanges object containing intervals upon which you want to map data.

operation

List of functions for the statistics and summary operations. This is similar to merge_bed()

max_gap

The largest gap for two intervals to be considered as overlapping. Default is -1 (no gap allowed, adjacent intervals not allowd).

min_overlap

See intersect_bed().

min_overlap_type

intersect_bed().

Value

A mapped GRanges object.

References

Manual page of bedtools map: https://bedtools.readthedocs.io/en/latest/content/tools/map.html

See also

Examples

# Load BED tables
tbl_x <- read_bed(system.file("extdata", "example_merge.bed", package = "bedtorch"))
tbl_y <- read_bed(system.file("extdata", "example_intersect_y.bed", package = "bedtorch"))

# Basic usage
result <- map_bed(tbl_x, tbl_y, operation = list(score_mean = list(on = "score", func = mean)))
head(result)
#> GRanges object with 3 ranges and 1 metadata column:
#>       seqnames    ranges strand | score_mean
#>          <Rle> <IRanges>  <Rle> |  <numeric>
#>   [1]       21     23-36      * |    7.66667
#>   [2]       21     48-55      * |    2.66667
#>   [3]       22     41-60      * |    6.66667
#>   -------
#>   seqinfo: 2 sequences from an unspecified genome; no seqlengths

# Perform the mapping, requiring the minimum overlapping of 5bp
result <- map_bed(tbl_x, tbl_y, operation = list(score_mean = list(on = "score", func = mean)), 
                  min_overlap = 5, min_overlap_type = "bp")
head(result)
#> GRanges object with 2 ranges and 1 metadata column:
#>       seqnames    ranges strand | score_mean
#>          <Rle> <IRanges>  <Rle> |  <numeric>
#>   [1]       21     23-36      * |          9
#>   [2]       22     41-60      * |          9
#>   -------
#>   seqinfo: 2 sequences from an unspecified genome; no seqlengths