Title: | ChIPexoQual |
---|---|
Description: | Package with a quality control pipeline for ChIP-exo/nexus data. |
Authors: | Rene Welch, Dongjun Chung, Sunduz Keles |
Maintainer: | Rene Welch <[email protected]> |
License: | GPL (>=2) |
Version: | 1.31.0 |
Built: | 2024-11-18 03:12:42 UTC |
Source: | https://github.com/bioc/ChIPexoQual |
ARCvURCplot
returns a ggplot
object with the ARC vs
URC plot to analyze enrichment and library complexity in ChIP-exo data.
ARCvURCplot(..., names.input = NULL, both.strand = FALSE)
ARCvURCplot(..., names.input = NULL, both.strand = FALSE)
... |
a |
names.input |
a character vector with the names to use in the
plot. If it is empty |
both.strand |
A logical value indicating if the |
A ggplot2
object with the ARC vs URC plot.
data(exoExample) ARCvURCplot(exoExample)
data(exoExample) ARCvURCplot(exoExample)
beta1
returns a vector with all the estimated values of the
models fitted by
ChIPexoQual
beta1(object) ## S4 method for signature 'ExoData' beta1(object)
beta1(object) ## S4 method for signature 'ExoData' beta1(object)
object |
a |
A numeric vector with estimated values for .
data(exoExample) beta1(exoExample)
data(exoExample) beta1(exoExample)
beta2
returns a vector with all the estimated values of the
models fitted by
ChIPexoQual
beta2(object) ## S4 method for signature 'ExoData' beta2(object)
beta2(object) ## S4 method for signature 'ExoData' beta2(object)
object |
a |
A numeric vector with estimated values for .
data(exoExample) beta2(exoExample)
data(exoExample) beta2(exoExample)
list
of GRanges
objects with the blacklists generated by the ENCODE and
modENCODE projects.list
of GRanges
objects with the blacklists generated by the ENCODE and
modENCODE projects.
data(blacklists)
data(blacklists)
list
of GRanges
objects.
A list
with the blacklists listed in https://sites.google.com/site/anshulkundaje/projects/blacklists.
calculateParamDist
calculates the quality parameters of one iteration.
This function samples nregions
rows from the stat matrix and fits
the linear model lm(d ~ 0 + u + w)
calculateParamDist
calculateParamDist
calculates the quality parameters of one iteration.
This function samples nregions
rows from the stat matrix and fits
the linear model lm(d ~ 0 + u + w)
calculateParamDist(i, stats, nregions)
calculateParamDist(i, stats, nregions)
i |
a numeric value indicating the current iteration. |
stats |
a |
nregions |
a numeric value indicating the number of regions sampled. |
a data.table
with both parameters and some extra info
data("exoExample") DT <- formatRegions(exoExample) calculateParamDist(1,DT,100)
data("exoExample") DT <- formatRegions(exoExample) calculateParamDist(1,DT,100)
ExoData
is a subclass of GenomicRanges
, used to asses the
quality of ChIP-exo/nexus sample.
ExoData(file = NULL, reads = NULL, height = 1, mc.cores = getOption("mc.cores", 2L), save.reads = FALSE, nregions = 1000, ntimes = 100, verbose = TRUE)
ExoData(file = NULL, reads = NULL, height = 1, mc.cores = getOption("mc.cores", 2L), save.reads = FALSE, nregions = 1000, ntimes = 100, verbose = TRUE)
file |
a character value with location of the bam file with the aligned reads. |
reads |
a |
height |
a numeric value indicating the value used to slice the coverage of the experiment into a set of regions. |
mc.cores |
a numeric value with the number of cores to use, i.e. at most how many child processes will be run simultaneously. |
save.reads |
a logical value to indicate if the reads are stored in the
|
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. |
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.
files <- list.files(system.file("extdata",package = "ChIPexoQualExample"), full.names = TRUE) ExoData(files[5],mc.cores = 2L)
files <- list.files(system.file("extdata",package = "ChIPexoQualExample"), full.names = TRUE) ExoData(files[5],mc.cores = 2L)
ExoDataBlacklist
separates the regions in an ExoData
object by overlapping them
with a set of blacklisted regions and calculates the quality parameters in both collections of
islands.
ExoDataBlacklist(exo, blacklist, which.param = "beta1", nregions = NULL, ntimes = NULL)
ExoDataBlacklist(exo, blacklist, which.param = "beta1", nregions = NULL, ntimes = NULL)
exo |
a |
blacklist |
a |
which.param |
a character value with either |
nregions |
a numeric value indicating the number of regions sampled to
estimate the quality parameter distributions. The default value is extracted from |
ntimes |
a numeric value indicating the number of times that regions are
sampled to estimate the quality parameter distributions. The default value
is extracted from |
A ggplot
object with a boxplot that compares the quality scores distribution when the regions
overlap a pre-defined collection of blacklists.
data(exoExample) data(blacklists) ExoDataBlacklist(exoExample,blacklists[["mm9"]],ntimes = 10,nregions = 500)
data(exoExample) data(blacklists) ExoDataBlacklist(exoExample,blacklists[["mm9"]],ntimes = 10,nregions = 500)
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))
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))
file |
a character value with location of the bam file with the aligned reads. |
reads |
a |
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
|
mc.cores |
a numeric value with the number of cores to use, i.e. at most how many child processes will be run simultaneously. |
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.
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)
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)
ExoData
results for FoxA1 ChIP-exo experimentExoData
object, generated with ChIPexoQual
and the
file:
data(exoExample)
data(exoExample)
ExoData object, which are GRanges
with additional columns.
ChIPexo_carroll_FoxA1_mouse_rep3_chr1.bam
An ExoData
object with the 3rd replicate of the FoxA1
experiment from ChIPExoQualExample
.
formatRegions
separates the width, depth and uniquePos summary statistics from the ExoData
object to calculate the quality parameters/formatRegions
formatRegions
separates the width, depth and uniquePos summary statistics from the ExoData
object to calculate the quality parameters/
formatRegions(exo)
formatRegions(exo)
exo |
a |
a data.table
with the width, depth and uniquePos of the regions in exo
.
data("exoExample") formatRegions(exoExample)
data("exoExample") formatRegions(exoExample)
FSRDistplot
returns a ggplot
object with the Forward
Strand Ratio distribution plot to analyze strand imbalance in
ChIP-exo data.
FSRDistplot(..., names.input = NULL, quantiles = c(0, 0.25, 0.5, 0.75, 1), depth.values = seq_len(30), both.strand = FALSE)
FSRDistplot(..., names.input = NULL, quantiles = c(0, 0.25, 0.5, 0.75, 1), depth.values = seq_len(30), both.strand = FALSE)
... |
a |
names.input |
a character vector with the names to use in the
plot. If it is empty |
quantiles |
a numeric vector with the quantiles used to estimate the
FSR distribution at a given depth. The default value is
|
depth.values |
a numeric vector indicating the regions with depth
less or equal to, that are going to be filtered out. The defaulta values
are |
both.strand |
a logical value indicating if the |
A ggplot2
object with the FSR distribution plot.
data(exoExample) FSRDistplot(exoExample)
data(exoExample) FSRDistplot(exoExample)
MAplot
returns a ggplot
object with the MA plot to
analyze the strand imbalance in ChIP-exo data.
MAplot(..., names.input = NULL)
MAplot(..., names.input = NULL)
... |
a |
names.input |
a character vector with the names to use in the
plot. If it is empty |
A ggplot2
object with the MA plot.
data(exoExample) MAplot(exoExample)
data(exoExample) MAplot(exoExample)
nreads
returns the number of reads in the object.
nreads(object) ## S4 method for signature 'ExoData' nreads(object)
nreads(object) ## S4 method for signature 'ExoData' nreads(object)
object |
A |
The number of reads in the ExoData
object.
data(exoExample) nreads(exoExample)
data(exoExample) nreads(exoExample)
paramDist
returns a DataFrame
with all the estimated
coefficients in the models
fitted by
ChIPexoQual
paramDist(object) ## S4 method for signature 'ExoData' paramDist(object = "ExoData")
paramDist(object) ## S4 method for signature 'ExoData' paramDist(object = "ExoData")
object |
a |
A DataFrame
with the fitted values of and
.
data(exoExample) paramDist(exoExample)
data(exoExample) paramDist(exoExample)
paramDistBoxplot
returns a ggplot
object with a
boxplot comparing the ntimes
estimations of the chosen
parameter.
paramDistBoxplot(..., names.input = NULL, which.param = "beta1", sort.as.numeric = FALSE)
paramDistBoxplot(..., names.input = NULL, which.param = "beta1", sort.as.numeric = FALSE)
... |
a |
names.input |
a character vector with the names to use in the
plot. If it is empty |
which.param |
a character value with either |
sort.as.numeric |
a logical value indicating if the values of
|
A ggplot2
object with the boxplot of the chosen
parameter
data(exoExample) paramDistBoxplot(exoExample)
data(exoExample) paramDistBoxplot(exoExample)
regionCompplot
returns a ggplot
object with the
Region Composition plot to analyze strand imbalance in ChIP-exo data.
regionCompplot(..., names.input = NULL, depth.values = seq_len(15))
regionCompplot(..., names.input = NULL, depth.values = seq_len(15))
... |
a |
names.input |
a character vector with the names to use in the
plot. If it is empty |
depth.values |
a numeric vector indicating the regions with depth
less or equal to, that are going to be filtered out. The defaulta values
are |
A ggplot2
object with the Region Composition plot.
data(exoExample) regionCompplot(exoExample)
data(exoExample) regionCompplot(exoExample)