R/consensusSeekeR.R
A549_NR3C1_CFQ_Peaks_partial.Rd
Sites representing the greatest evidence of enrichment for the NR3C1 transcription factor (DCC accession: ENCFF002CFQ) for regions chr2:40000000-50000000 and chr3:10000000-13000000 from the Encyclopedia of DNA Elements (ENCODE) data (Dunham I et al. 2012).
data(A549_NR3C1_CFQ_Peaks_partial)
A GRanges
containing one entry per site. The peaks are
surronded by ranges present in the dataset
A549_NR3C1_CFQ_NarrowPeaks_partial
.
The Encyclopedia of DNA Elements (ENCODE) (DCC accession: ENCFF002CFQ)
The peaks and ranges have been obtained using an optimal IDR analysis done on all replicates.
Dunham I, Kundaje A, Aldred SF, et al. An integrated encyclopedia of DNA elements in the human genome. Nature. 2012 Sep 6;489(7414):57-74.
A549_NR3C1_CFQ_NarrowPeaks_partial
the associate
genomic regions dataset.
findConsensusPeakRegions
for extracting regions
sharing the same features in more than one experiment.
## Loading datasets
data(A549_NR3C1_CFQ_NarrowPeaks_partial)
data(A549_NR3C1_CFQ_Peaks_partial)
data(A549_NR3C1_CFS_NarrowPeaks_partial)
data(A549_NR3C1_CFS_Peaks_partial)
data(A549_NR3C1_CFR_NarrowPeaks_partial)
data(A549_NR3C1_CFR_Peaks_partial)
## Assigning experiment name to each row of the dataset.
## NarrowPeak and Peak datasets from the same experiment must
## have identical names.
names(A549_NR3C1_CFQ_NarrowPeaks_partial) <- rep("NR3C1_CFQ",
length(A549_NR3C1_CFQ_NarrowPeaks_partial))
names(A549_NR3C1_CFQ_Peaks_partial) <- rep("NR3C1_CFQ",
length(A549_NR3C1_CFQ_Peaks_partial))
names(A549_NR3C1_CFS_NarrowPeaks_partial) <-rep("NR3C1_CFS",
length(A549_NR3C1_CFS_NarrowPeaks_partial))
names(A549_NR3C1_CFS_Peaks_partial) <- rep("NR3C1_CFS",
length(A549_NR3C1_CFS_Peaks_partial))
names(A549_NR3C1_CFR_NarrowPeaks_partial) <-rep("NR3C1_CFR",
length(A549_NR3C1_CFR_NarrowPeaks_partial))
names(A549_NR3C1_CFR_Peaks_partial) <- rep("NR3C1_CFR",
length(A549_NR3C1_CFR_Peaks_partial))
## Calculating consensus regions for chromosome 3
## with a default region size of 140 bp (2 * extendingSize)
## which is extended to include all genomic regions for the closest
## peak to the median position of all peaks included in the region (for
## each experiment).
## Peaks from at least 2 experiments must be present in a region to
## be retained as a consensus region.
chrList <- Seqinfo(c("chr3"), c(198022430), NA)
findConsensusPeakRegions(
narrowPeaks = c(A549_NR3C1_CFQ_NarrowPeaks_partial,
A549_NR3C1_CFS_NarrowPeaks_partial,
A549_NR3C1_CFR_NarrowPeaks_partial),
peaks = c(A549_NR3C1_CFQ_Peaks_partial,
A549_NR3C1_CFS_Peaks_partial,
A549_NR3C1_CFR_Peaks_partial),
chrInfo = chrList,
extendingSize = 70,
expandToFitPeakRegion = FALSE,
shrinkToFitPeakRegion = FALSE,
minNbrExp = 2,
nbrThreads = 1)
#> $call
#> findConsensusPeakRegions(narrowPeaks = c(A549_NR3C1_CFQ_NarrowPeaks_partial,
#> A549_NR3C1_CFS_NarrowPeaks_partial, A549_NR3C1_CFR_NarrowPeaks_partial),
#> peaks = c(A549_NR3C1_CFQ_Peaks_partial, A549_NR3C1_CFS_Peaks_partial,
#> A549_NR3C1_CFR_Peaks_partial), chrInfo = chrList, extendingSize = 70,
#> expandToFitPeakRegion = FALSE, shrinkToFitPeakRegion = FALSE,
#> minNbrExp = 2, nbrThreads = 1)
#>
#> $consensusRanges
#> GRanges object with 17 ranges and 0 metadata columns:
#> seqnames ranges strand
#> <Rle> <IRanges> <Rle>
#> [1] chr3 10018764-10018904 *
#> [2] chr3 10210879-10211019 *
#> [3] chr3 10233837-10233977 *
#> [4] chr3 10245710-10245850 *
#> [5] chr3 10247207-10247347 *
#> ... ... ... ...
#> [13] chr3 12447285-12447425 *
#> [14] chr3 12465942-12466082 *
#> [15] chr3 12549897-12550037 *
#> [16] chr3 12685950-12686090 *
#> [17] chr3 12721119-12721259 *
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths
#>
#> attr(,"class")
#> [1] "consensusRanges"