Skip to contents

For each interval in the `regions` GRanges object, this function finds all overlapping fragments in `binding_data` and computes a weighted mean of their signal values. The mean is weighted by the length of the overlap between each fragment and the region.

Usage

gr_occupancy(
  binding_data,
  regions,
  buffer = 0,
  BPPARAM = BiocParallel::bpparam()
)

Arguments

binding_data

A data.frame as returned by `build_dataframes()`. It must contain columns 'chr', 'start', 'end', followed by numeric sample columns.

regions

A GRanges object of genomic intervals (e.g., reduced peaks) over which to calculate occupancy.

buffer

Optional integer. Number of base pairs to expand each interval in `regions` on both sides before calculating occupancy. Default is 0.

BPPARAM

A BiocParallel parameter object for parallel computation. Default is `BiocParallel::bpparam()`.

Value

A data.frame with one row per region from the input `regions` object. Columns include 'name', 'nfrags' (number of overlapping fragments), and one column for the calculated weighted mean occupancy for each sample.

Examples

# 1. Create a set of regions (e.g., reduced peaks)
regions_gr <- GenomicRanges::GRanges("chrX", IRanges::IRanges(start = c(100, 500), width = 100))
S4Vectors::mcols(regions_gr)$name <- paste0(
    GenomicRanges::seqnames(regions_gr), ":",
    GenomicRanges::start(regions_gr), "-", GenomicRanges::end(regions_gr)
)

# 2. Create a mock binding data data.frame
binding_df <- data.frame(
    chr = "chrX",
    start = c(90, 150, 480, 550),
    end = c(110, 170, 520, 580),
    sampleA = c(1.2, 0.8, 2.5, 3.0),
    sampleB = c(1.0, 0.9, 2.8, 2.9)
)

# 3. Calculate occupancy over the regions
# Use BiocParallel::SerialParam() for deterministic execution in examples
if (requireNamespace("BiocParallel", quietly = TRUE)) {
    occupancy_data <- gr_occupancy(binding_df, regions_gr,
        BPPARAM = BiocParallel::SerialParam()
    )
    print(occupancy_data)
}
#> Calculating average occupancy per region ...
#>                      name nfrags  sampleA  sampleB
#> chrX:100-199 chrX:100-199      2 0.937500 0.934375
#> chrX:500-599 chrX:500-599      2 2.798077 2.859615