View source: R/makeConsensus.R
makeConsensus | R Documentation |
Make a set of consensus peaks based on the number of replicates
makeConsensus(
x,
p = 0,
var = NULL,
method = c("union", "coverage"),
ignore.strand = TRUE,
simplify = FALSE,
min_width = 0,
merge_within = 1L,
...
)
x |
A GRangesList |
p |
The minimum proportion of samples (i.e. elements of |
var |
Additional columns in the |
method |
Either return the union of all overlapping ranges, or the regions within the overlapping ranges which are covered by the specified proportion of replicates. When using p = 0, both methods will return identical results |
ignore.strand , simplify , ... |
Passed to reduceMC or intersectMC internally |
min_width |
Discard any regions below this width |
merge_within |
Passed to reduce as |
This takes a list of GRanges objects and forms a set of consensus peaks.
When using method = "union" the union ranges of all overlapping peaks will be returned, using the minimum proportion of replicates specified. When using method = "coverage", only the regions within each overlapping range which are 'covered' by the minimum proportion of replicates specified are returned. This will return narrower peaks in general, although some artefactual very small ranges may be included (e.g. 10bp). Careful setting of the min_width and merge_within parameters may be very helpful for these instances. It is also expected that setting method = "coverage" should return the region within each range which is more likely to contain the true binding site for the relevant ChIP targets
GRanges
object with mcols containing a logical vector for every element of
x, along with the column n
which adds all logical columns. These columns
denote which replicates contain an overlapping peak for each range
If any additional columns have been requested using var
, these will be
returned as CompressedList objects as produced by reduceMC()
or
intersectMC()
.
reduceMC intersectMC
data("peaks")
## The first three replicates are from the same treatment group
grl <- peaks[1:3]
names(grl) <- gsub("_peaks.+", "", names(grl))
makeConsensus(grl)
makeConsensus(grl, p = 2/3, var = "score")
## Using method = 'coverage' finds ranges based on the intersection
makeConsensus(grl, p = 2/3, var = "score", method = "coverage")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.