R/readNarrowPeak.R
readNarrowPeakFile.Rd
Read a narrowPeak file and extract the narrow regions and/or the peaks, as specified by used. The narrowPeak file must fit the UCSC specifications. See https://genome.ucsc.edu/FAQ/FAQformat.html#format12 for more details. The file can have one or many header lines. However, the total number of header lines must be inferior to 250 lines.
readNarrowPeakFile(file_path, extractRegions = TRUE, extractPeaks = TRUE)
the name of the file.
a logical
indicating if the narrow regions must
be extracted. If TRUE
, a GRanges
containing the
narrow regions will be returned. Otherwise, NULL
is
returned. Default = TRUE
.
a logical
indicating if the peaks must
be extracted. If TRUE
, a GRanges
containing the peaks
will be returned. Otherwise, NULL
is
returned. Default = TRUE
.
a list
containing 2 entries:
narrowPeak a GRanges
containing
the narrow regions extracted from the file. NULL
when
not needed by user.
peak a GRanges
containing
the peaks extracted from the file. NULL
when not
## Set file information
test_narrowPeak <- system.file("extdata",
"A549_FOSL2_ENCSR000BQO_MZW_part_chr_1_and_12.narrowPeak",
package = "consensusSeekeR")
## Read file to extract peaks and regions
data <- readNarrowPeakFile(test_narrowPeak, extractRegions = TRUE,
extractPeaks = TRUE)
## To access peak data (GRanges format)
head(data$peak)
#> GRanges object with 6 ranges and 6 metadata columns:
#> seqnames ranges strand | name score signalValue
#> <Rle> <IRanges> <Rle> | <character> <numeric> <numeric>
#> [1] chr1 846687 * | peaks/Hosa_A549_FOSL.. 57 5.59984
#> [2] chr1 856112 * | peaks/Hosa_A549_FOSL.. 43 5.16770
#> [3] chr1 856643 * | peaks/Hosa_A549_FOSL.. 33 4.65093
#> [4] chr1 873964 * | peaks/Hosa_A549_FOSL.. 189 9.91576
#> [5] chr1 878672 * | peaks/Hosa_A549_FOSL.. 33 4.65093
#> [6] chr1 881440 * | peaks/Hosa_A549_FOSL.. 29 4.03975
#> pValue qValue peak
#> <numeric> <numeric> <integer>
#> [1] 8.75159 5.77648 98
#> [2] 7.21902 4.33609 108
#> [3] 6.18592 3.38310 96
#> [4] 22.37844 18.94888 314
#> [5] 6.18592 3.38310 217
#> [6] 5.68403 2.92251 143
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
## To access region data (GRanges format)
head(data$narrowPeak)
#> GRanges object with 6 ranges and 6 metadata columns:
#> seqnames ranges strand | name score
#> <Rle> <IRanges> <Rle> | <character> <numeric>
#> [1] chr1 846589-846847 * | peaks/Hosa_A549_FOSL.. 57
#> [2] chr1 856004-856159 * | peaks/Hosa_A549_FOSL.. 43
#> [3] chr1 856547-856768 * | peaks/Hosa_A549_FOSL.. 33
#> [4] chr1 873650-874067 * | peaks/Hosa_A549_FOSL.. 189
#> [5] chr1 878455-878768 * | peaks/Hosa_A549_FOSL.. 33
#> [6] chr1 881297-881473 * | peaks/Hosa_A549_FOSL.. 29
#> signalValue pValue qValue peak
#> <numeric> <numeric> <numeric> <integer>
#> [1] 5.59984 8.75159 5.77648 98
#> [2] 5.16770 7.21902 4.33609 108
#> [3] 4.65093 6.18592 3.38310 96
#> [4] 9.91576 22.37844 18.94888 314
#> [5] 4.65093 6.18592 3.38310 217
#> [6] 4.03975 5.68403 2.92251 143
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths