View source: R/clusterWindowsList.R
clusterWindowsList | R Documentation |
Consolidate DB results from multiple analyses with cluster-level FDR control.
clusterWindowsList(ranges.list, tab.list, equiweight=TRUE, ...)
ranges.list |
A list of GRanges or RangedSummarizedExperiment objects,
usually containing windows of varying sizes from separate calls to |
tab.list |
A list of data.frames of differential binding results,
usually from separate analyses at differing window sizes.
Each should contain one row per interval for the corresponding entry of |
equiweight |
a logical scalar indicating whether equal weighting from each analysis should be enforced |
... |
arguments to be passed to |
This function consolidates DB results from multiple analyses, typically involving different window sizes.
The aim is to provide comprehensive detection of DB at a range of spatial resolutions.
Significant windows from each analysis are identified and used for clustering with clusterWindows
.
This represents the post-hoc counterpart to mergeResultsList
.
Some effort is required to equalize the contribution of the results from each analysis.
This is done by setting equiweight=TRUE
,
where the weight of each window is inversely proportional to the number of windows from that analysis.
These weights are used as frequency weights for window-level FDR control (to identify DB windows prior to clustering).
Otherwise, the final results would be dominated by large number of small windows.
Users can cluster by the sign of log-fold changes, to obtain clusters of DB windows of the same sign.
However, note that nested windows with opposite signs may give unintuitive results - see mergeWindows
for details.
A named list is returned 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
,
an integer field specifying the index of x
from which each interval was obtained.
ids
:An integer vectors indicating the cluster ID for each window in ranges
.
Non-significant windows that were not assigned to a cluster have IDs of NA
, as described in ?clusterWindows
.
regions
:A GRanges object containing the coordinates for each cluster.
FDR
:A numeric scalar containing the cluster-level FDR estimate.
Aaron Lun
clusterWindows
, the equivalent function for a single GRanges input.
mergeResultsList
, for a more rigorous approach to clustering windows.
# Making up some GRanges.
low <- GRanges("chrA", IRanges(runif(100, 1, 1000), width=5))
med <- GRanges("chrA", IRanges(runif(40, 1, 1000), width=10))
high <- GRanges("chrA", IRanges(runif(10, 1, 1000), width=20))
# Making up some DB results.
dbl <- data.frame(logFC=rnorm(length(low)), PValue=rbeta(length(low), 1, 20))
dbm <- data.frame(logFC=rnorm(length(med)), PValue=rbeta(length(med), 1, 20))
dbh <- data.frame(logFC=rnorm(length(high)), PValue=rbeta(length(high), 1, 20))
result.list <- list(dbl, dbm, dbh)
# Consolidating.
cons <- clusterWindowsList(list(low, med, high), result.list, tol=20)
cons$region
cons$id
cons$FDR
# Without weights.
cons <- clusterWindowsList(list(low, med, high), result.list, tol=20,
equiweight=FALSE)
cons$FDR
# Using the signs.
cons <- clusterWindowsList(list(low, med, high), result.list, tol=20,
fc.col="logFC")
cons$FDR
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.