Import functions for plus/minus pairs of bigWig
or bedGraph
files.
import_bigWig(
plus_file = NULL,
minus_file = NULL,
genome = NULL,
keep.X = TRUE,
keep.Y = TRUE,
keep.M = FALSE,
keep.nonstandard = FALSE,
makeBRG = TRUE,
ncores = getOption("mc.cores", 2L)
)
import_bedGraph(
plus_file = NULL,
minus_file = NULL,
genome = NULL,
keep.X = TRUE,
keep.Y = TRUE,
keep.M = FALSE,
keep.nonstandard = FALSE,
ncores = getOption("mc.cores", 2L)
)
Paths for strand-specific input files, or a vector of such paths. If vectors are given, the user should take care that the orders match!
Optional string for UCSC reference genome, e.g. "hg38". If given, non-standard chromosomes are trimmed, and options for sex and mitochondrial chromosomes are applied.
Logicals indicating which non-autosomes should be kept. By default, sex chromosomes are kept, but mitochondrial and non-standard chromosomes are removed.
If TRUE
(the default), the output ranges are made
single-width using makeGRangesBRG
Number of cores to use, if importing multiple objects simultaneously.
Imports a GRanges object containing base-pair resolution data, with
the score
metadata column indicating the number of reads represented
by each range.
For import_bigWig
, the output GRanges is formatted by
makeGRangesBRG
, such that all
ranges are disjoint and have width = 1, and the score
is single-base
coverage, i.e. the number of reads for each position.
import_bedGraph
is useful for when both 5'- and 3'-end information
is to be maintained for each sequenced molecule. For this purpose, one use
bedGraphs to store entire reads, with the score
representing the
number of reads sharing identical 5' and 3' ends. However,
import_bedGraph
doesn't modify the information in the bedGraph
files. If the bedGraph file represents basepair-resolution coverage data,
then users can coerce it to a basepair-resolution GRanges object by using
getStrandedCoverage
followed
by makeGRangesBRG
.
#--------------------------------------------------#
# Import PRO-seq bigWigs -> coverage of 3' bases
#--------------------------------------------------#
# get local address for included bigWig files
p.bw <- system.file("extdata", "PROseq_dm6_chr4_plus.bw",
package = "BRGenomics")
m.bw <- system.file("extdata", "PROseq_dm6_chr4_minus.bw",
package = "BRGenomics")
# import bigWigs (not supported on windows)
if (.Platform$OS.type == "unix") {
PROseq <- import_bigWig(p.bw, m.bw, genome = "dm6")
PROseq
}
#> GRanges object with 150 ranges and 1 metadata column:
#> seqnames ranges strand | score
#> <Rle> <IRanges> <Rle> | <integer>
#> [1] chr4 1295 + | 1
#> [2] chr4 41428 + | 1
#> [3] chr4 42588 + | 1
#> [4] chr4 42590 + | 2
#> [5] chr4 42593 + | 5
#> ... ... ... ... . ...
#> [146] chr4 1307742 - | 1
#> [147] chr4 1316537 - | 1
#> [148] chr4 1318960 - | 1
#> [149] chr4 1319004 - | 1
#> [150] chr4 1319369 - | 1
#> -------
#> seqinfo: 1 sequence from dm6 genome
#--------------------------------------------------#
# Import PRO-seq bedGraphs -> whole reads (matched 5' and 3' ends)
#--------------------------------------------------#
# get local address for included bedGraph files
p.bg <- system.file("extdata", "PROseq_dm6_chr4_plus.bedGraph",
package = "BRGenomics")
m.bg <- system.file("extdata", "PROseq_dm6_chr4_minus.bedGraph",
package = "BRGenomics")
# import bedGraphs
PROseq_paired <- import_bedGraph(p.bg, m.bg, genome = "dm6")
PROseq_paired
#> GRanges object with 164 ranges and 1 metadata column:
#> seqnames ranges strand | score
#> <Rle> <IRanges> <Rle> | <integer>
#> [1] chr4 1270-1295 + | 1
#> [2] chr4 41411-41428 + | 1
#> [3] chr4 42556-42590 + | 1
#> [4] chr4 42559-42588 + | 1
#> [5] chr4 42559-42593 + | 3
#> ... ... ... ... . ...
#> [160] chr4 1307742-1307776 - | 1
#> [161] chr4 1316537-1316563 - | 1
#> [162] chr4 1318960-1318994 - | 1
#> [163] chr4 1319004-1319038 - | 1
#> [164] chr4 1319369-1319403 - | 1
#> -------
#> seqinfo: 1 sequence from dm6 genome; no seqlengths