View source: R/getRegionCoverage.R
getRegionCoverage | R Documentation |
This function extracts the raw coverage information calculated by fullCoverage at each base for a set of regions found with calculatePvalues. It can further calculate the mean coverage per sample for each region.
getRegionCoverage(
fullCov = NULL,
regions,
totalMapped = NULL,
targetSize = 8e+07,
files = NULL,
...
)
fullCov |
A list where each element is the result from
loadCoverage used with |
regions |
The |
totalMapped |
The total number of reads mapped for each sample.
Providing this data adjusts the coverage to reads in |
targetSize |
The target library size to adjust the coverage to. Used
only when |
files |
A character vector with the full path to the sample BAM files
(or BigWig files).
The names are used for the column names of the DataFrame. Check
rawFiles for constructing |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
Passed to extendedMapSeqlevels and define_cluster. When |
When fullCov
is the output of loadCoverage with
cutoff
non-NULL, getRegionCoverage assumes that the regions
come from the same data. Meaning that filterData was not used again.
This ensures that the regions are a subset of the data available in
fullCov
.
If fullCov
is NULL
and files
is specified, this function
will attempt to read the coverage from the files. Note that if you used
'totalMapped' and 'targetSize' before, you will have to specify them again
to get the same results.
You should use at most one core per chromosome.
a list of data.frame where each data.frame has the coverage
information (nrow = width of region, ncol = number of samples) for a given
region. The names of the list correspond to the region indexes in
regions
Andrew Jaffe, Leonardo Collado-Torres
fullCoverage, calculatePvalues
## Obtain fullCov object
fullCov <- list("21" = genomeDataRaw$coverage)
## Assign chr lengths using hg19 information, use only first two regions
library("GenomicRanges")
regions <- genomeRegions$regions[1:2]
seqlengths(regions) <- seqlengths(getChromInfoFromUCSC("hg19",
as.Seqinfo = TRUE
))[
mapSeqlevels(names(seqlengths(regions)), "UCSC")
]
## Finally, get the region coverage
regionCov <- getRegionCoverage(fullCov = fullCov, regions = regions)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.