mergeWindowsList: Consolidate window sizes

Description Usage Arguments Details Value Author(s) See Also Examples

View source: R/mergeWindowsList.R

Description

Consolidate windows of different sizes into a common set of regions.

Usage

1
2
3
mergeWindowsList(ranges.list, tol, ...)

findOverlapsList(ranges.list, regions, ...)

Arguments

ranges.list

A list of RangedSummarizedExperiment and/or GRanges objects.

tol

An integer scalar specifying the tolerance to use in mergeWindows.

...

For mergeWindowsList, further arguments to pass to mergeWindows.

For findOverlapsList, further arguments to pass to findOverlaps.

regions

A GRanges object specifying reference regions of interest for overlapping with windows.

Details

DB analyses with csaw can be repeated using windows of different sizes to provide comprehensive detection of DB at a range of spatial resolutions. This function merges together those sets of windows of differing sizes, typically generated by multiple calls to windowCounts using different width values. For mergeWindowsList, windows of all sizes are clustered together through mergeWindows. For findOverlapsList, windows of all sizes are overlapped with ref using findOverlaps.

The aim is to pass the output of this function to combineTests or combineOverlaps. This takes the statistical machinery used to combine p-values across windows for a cluster or reference region, and recycles it to provide a rigorous method for consolidating statistics across multiple analyses with different window sizes. However, it requires some care to balance the contribution of analyses with different window sizes to the combined p-value. Otherwise, analyses with many small windows will dominate the calculation for each cluster or reference region.

To compensate for this effect, both functions weight each window's contribution to a cluster or region. The weight is inversely proportional to the number of windows of the same size in the same cluster or region. For mergeWindowsList, this yields a weight for each window, while for findOverlapsList, this yields a weight for each overlap. Application of those weights in combineTests and related functions ensures that the combined p-value is not dominated by numerous small windows within a cluster.

Value

mergeWindowsList returns a named list containing:

ranges:

A GRanges object containing the concatenated intervals from all elements of x. The element-wise metadata of this GRanges contains origin, an integer field specifying the index of x from which each interval was obtained.

ids:

An integer vector of length equal to ranges. This specifies the cluster (i.e., entry of merged) to which each interval in ranges was assigned.

regions:

A GRanges object containing the genomic coordinates of the clusters of merged windows, as computed by mergeWindows.

weights:

An numeric vector of length equal to ranges, specifying the weight of each interval.

findOverlapsList returns a named list containing:

ranges:

A GRanges object containing the concatenated intervals from all elements of x. The element-wise metadata of this object contains the integer field origin, as described above.

olap:

A Hits object containing the overlaps between ref (query) and ranges (subject).

weights:

A numeric vector of length equal to olap, containing the weight of each overlap.

Author(s)

Aaron Lun

See Also

mergeWindows, for the functionality behind mergeWindowsList.

findOverlaps, for the functionality behind findOverlapsList.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
bamFiles <- system.file("exdata", c("rep1.bam", "rep2.bam"), package="csaw")
low <- windowCounts(bamFiles, width=1, filter=1)
med <- windowCounts(bamFiles, width=100, filter=1)
high <- windowCounts(bamFiles, width=500, filter=1)

cons <- mergeWindowsList(list(low, med, high), tol=100)
cons$ranges
cons$merged

of.interest <- GRanges(c('chrA', 'chrA', 'chrB', 'chrC'), 
    IRanges(c(1, 500, 100, 1000), c(200, 1000, 700, 1500)))
cons <- findOverlapsList(list(low, med, high), regions=of.interest)
cons$ranges
cons$olap

csaw documentation built on Nov. 12, 2020, 2:03 a.m.