ExoDataSubsampling samples sample.reads from the ChIP-exo experiment and creates a list of ExoData objects

ExoDataSubsampling(file = NULL, reads = NULL, sample.depth = NULL,
  height = 1, nregions = 1000, ntimes = 1000, verbose = TRUE,
  save.reads = FALSE, mc.cores = getOption("mc.cores", 2L))

Arguments

file
a character value with location of the bam file with the aligned reads.
reads
a GAlignments object with the aligned reads of a ChIP-exo sample. It is meant to be used instead of file.
sample.depth
a numeric vector with the number of reads to be sampled.
height
a numeric value indicating the value used to slice the coverage of the experiment into a set of regions.
nregions
a numeric value indicating the number of regions sampled to estimate the quality parameter distributions. The default value is 1e3.
ntimes
a numeric value indicating the number of times that regions are sampled to estimate the quality parameter distributions. The default value is 1e2.
verbose
a logical value indicating if the user want to receive progress details. The default value is FALSE.
save.reads
a logical value to indicate if the reads are stored in the ExoData object. The default value is FALSE.
mc.cores
a numeric value with the number of cores to use, i.e. at most how many child processes will be run simultaneously.

Value

It returns an ExoData object with the regions obtained after partitioning the genome and the summary statistics for each region. If the save.reads parameter is TRUE then it contains a GRanges object with the reads of the ChIP-exo experiment.

Examples

files <- list.files(system.file("extdata",package = "ChIPexoQualExample"), full.names = TRUE) sample.depth <- seq(1e5,2e5,5e4) ExoDataSubsampling(file = files[5],sample.depth = sample.depth)
#> Creating 'ExoData' for 100,000 reads
#> Sampling reads.
#> Using 'reads' argument
#> Keeping reads in object: No
#> Calculating summary statistics
#> Calculating quality scores distribution
#> Done!
#> Creating 'ExoData' for 150,000 reads
#> Sampling reads.
#> Using 'reads' argument
#> Keeping reads in object: No
#> Calculating summary statistics
#> Calculating quality scores distribution
#> Done!
#> Creating 'ExoData' for 200,000 reads
#> Sampling reads.
#> Using 'reads' argument
#> Keeping reads in object: No
#> Calculating summary statistics
#> Calculating quality scores distribution
#> Done!
#> $`100,000 reads` #> ExoData object with 72562 ranges and 11 metadata columns: #> seqnames ranges strand | fwdReads revReads #> <Rle> <IRanges> <Rle> | <integer> <integer> #> [1] chr1 [3001926, 3001961] * | 1 0 #> [2] chr1 [3002412, 3002447] * | 0 2 #> [3] chr1 [3003987, 3004022] * | 0 1 #> [4] chr1 [3006189, 3006224] * | 1 0 #> [5] chr1 [3014235, 3014270] * | 1 0 #> ... ... ... ... . ... ... #> [72558] chr1 [197171309, 197171344] * | 1 0 #> [72559] chr1 [197175350, 197175385] * | 1 0 #> [72560] chr1 [197184957, 197184992] * | 0 1 #> [72561] chr1 [197188843, 197188878] * | 1 0 #> [72562] chr1 [197189679, 197189714] * | 1 0 #> fwdPos revPos depth uniquePos ARC URC #> <integer> <integer> <integer> <integer> <numeric> <numeric> #> [1] 1 0 1 1 0.0277777777777778 1 #> [2] 0 1 2 1 0.0555555555555556 0.5 #> [3] 0 1 1 1 0.0277777777777778 1 #> [4] 1 0 1 1 0.0277777777777778 1 #> [5] 1 0 1 1 0.0277777777777778 1 #> ... ... ... ... ... ... ... #> [72558] 1 0 1 1 0.0277777777777778 1 #> [72559] 1 0 1 1 0.0277777777777778 1 #> [72560] 0 1 1 1 0.0277777777777778 1 #> [72561] 1 0 1 1 0.0277777777777778 1 #> [72562] 1 0 1 1 0.0277777777777778 1 #> FSR M A #> <numeric> <numeric> <numeric> #> [1] 1 -Inf Inf #> [2] 0 -Inf -Inf #> [3] 0 -Inf -Inf #> [4] 1 -Inf Inf #> [5] 1 -Inf Inf #> ... ... ... ... #> [72558] 1 -Inf Inf #> [72559] 1 -Inf Inf #> [72560] 0 -Inf -Inf #> [72561] 1 -Inf Inf #> [72562] 1 -Inf Inf #> ------- #> seqinfo: 1 sequence from an unspecified genome; no seqlengths #> #> $`150,000 reads` #> ExoData object with 96484 ranges and 11 metadata columns: #> seqnames ranges strand | fwdReads revReads #> <Rle> <IRanges> <Rle> | <integer> <integer> #> [1] chr1 [3002412, 3002447] * | 0 1 #> [2] chr1 [3014235, 3014270] * | 1 0 #> [3] chr1 [3014692, 3014727] * | 1 0 #> [4] chr1 [3016288, 3016323] * | 1 0 #> [5] chr1 [3016377, 3016412] * | 1 0 #> ... ... ... ... . ... ... #> [96480] chr1 [197184957, 197184992] * | 0 1 #> [96481] chr1 [197188843, 197188878] * | 1 0 #> [96482] chr1 [197189055, 197189090] * | 2 0 #> [96483] chr1 [197192288, 197192323] * | 0 2 #> [96484] chr1 [197192797, 197192832] * | 1 0 #> fwdPos revPos depth uniquePos ARC URC #> <integer> <integer> <integer> <integer> <numeric> <numeric> #> [1] 0 1 1 1 0.0277777777777778 1 #> [2] 1 0 1 1 0.0277777777777778 1 #> [3] 1 0 1 1 0.0277777777777778 1 #> [4] 1 0 1 1 0.0277777777777778 1 #> [5] 1 0 1 1 0.0277777777777778 1 #> ... ... ... ... ... ... ... #> [96480] 0 1 1 1 0.0277777777777778 1 #> [96481] 1 0 1 1 0.0277777777777778 1 #> [96482] 1 0 2 1 0.0555555555555556 0.5 #> [96483] 0 1 2 1 0.0555555555555556 0.5 #> [96484] 1 0 1 1 0.0277777777777778 1 #> FSR M A #> <numeric> <numeric> <numeric> #> [1] 0 -Inf -Inf #> [2] 1 -Inf Inf #> [3] 1 -Inf Inf #> [4] 1 -Inf Inf #> [5] 1 -Inf Inf #> ... ... ... ... #> [96480] 0 -Inf -Inf #> [96481] 1 -Inf Inf #> [96482] 1 -Inf Inf #> [96483] 0 -Inf -Inf #> [96484] 1 -Inf Inf #> ------- #> seqinfo: 1 sequence from an unspecified genome; no seqlengths #> #> $`200,000 reads` #> ExoData object with 115669 ranges and 11 metadata columns: #> seqnames ranges strand | fwdReads revReads #> <Rle> <IRanges> <Rle> | <integer> <integer> #> [1] chr1 [3002412, 3002447] * | 0 1 #> [2] chr1 [3003987, 3004022] * | 0 1 #> [3] chr1 [3014692, 3014727] * | 2 0 #> [4] chr1 [3015843, 3015878] * | 0 1 #> [5] chr1 [3016288, 3016340] * | 2 0 #> ... ... ... ... . ... ... #> [115665] chr1 [197183268, 197183303] * | 1 0 #> [115666] chr1 [197188637, 197188672] * | 0 1 #> [115667] chr1 [197189055, 197189090] * | 1 0 #> [115668] chr1 [197189679, 197189714] * | 2 0 #> [115669] chr1 [197192288, 197192323] * | 0 1 #> fwdPos revPos depth uniquePos ARC URC #> <integer> <integer> <integer> <integer> <numeric> <numeric> #> [1] 0 1 1 1 0.0277777777777778 1 #> [2] 0 1 1 1 0.0277777777777778 1 #> [3] 1 0 2 1 0.0555555555555556 0.5 #> [4] 0 1 1 1 0.0277777777777778 1 #> [5] 2 0 2 2 0.0377358490566038 1 #> ... ... ... ... ... ... ... #> [115665] 1 0 1 1 0.0277777777777778 1 #> [115666] 0 1 1 1 0.0277777777777778 1 #> [115667] 1 0 1 1 0.0277777777777778 1 #> [115668] 1 0 2 1 0.0555555555555556 0.5 #> [115669] 0 1 1 1 0.0277777777777778 1 #> FSR M A #> <numeric> <numeric> <numeric> #> [1] 0 -Inf -Inf #> [2] 0 -Inf -Inf #> [3] 1 -Inf Inf #> [4] 0 -Inf -Inf #> [5] 1 -Inf Inf #> ... ... ... ... #> [115665] 1 -Inf Inf #> [115666] 0 -Inf -Inf #> [115667] 1 -Inf Inf #> [115668] 1 -Inf Inf #> [115669] 0 -Inf -Inf #> ------- #> seqinfo: 1 sequence from an unspecified genome; no seqlengths #>