railMatrix | R Documentation |
Rail (available at http://rail.bio) generates coverage BigWig files. These files are faster to load in R and to process. Rail creates an un-adjusted coverage BigWig file per sample and an adjusted summary coverage BigWig file by chromosome (median or mean). railMatrix reads in the mean (or median) coverage BigWig file and applies a threshold cutoff to identify expressed regions (ERs). Then it goes back to the sample coverage BigWig files and extracts the base level coverage for each sample. Finally it summarizes this information in a matrix with one row per ERs and one column per sample. This function is similar to regionMatrix but is faster thanks to the advantages presented by BigWig files.
railMatrix(
chrs,
summaryFiles,
sampleFiles,
L,
cutoff = NULL,
maxClusterGap = 300L,
totalMapped = NULL,
targetSize = 4e+07,
file.cores = 1L,
chrlens = NULL,
...
)
chrs |
A character vector with the names of the chromosomes. |
summaryFiles |
A character vector (or BigWigFileList) with the paths to
the summary BigWig files created by Rail. Either mean or median files. These
are library size adjusted by default in Rail. The order of the files in this
vector must match the order in |
sampleFiles |
A character vector with the paths to the BigWig files by sample. These files are unadjusted for library size. |
L |
The width of the reads used. Either a vector of length 1 or length equal to the number of samples. |
cutoff |
The base-pair level cutoff to use. It's behavior is controlled
by |
maxClusterGap |
This determines the maximum gap between candidate ERs. |
totalMapped |
A vector with the total number of reads mapped for each
sample. The vector should be in the same order as the samples in |
targetSize |
The target library size to adjust the coverage to. Used
only when |
file.cores |
Number of cores used for loading the BigWig files per chr. |
chrlens |
The chromosome lengths in base pairs. If it's |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
Passed to filterData, findRegions and define_cluster. |
Given a set of un-filtered coverage data (see fullCoverage), create candidate regions by applying a cutoff on the coverage values, and obtain a count matrix where the number of rows corresponds to the number of candidate regions and the number of columns corresponds to the number of samples. The values are the mean coverage for a given sample for a given region.
This function uses several other derfinder-package functions. Inspect the code if interested.
You should use at most one core per chromosome.
A list with one entry per chromosome. Then per chromosome, a list with two components.
A set of regions based on the coverage filter cutoff as returned by findRegions.
A matrix with the mean coverage by sample for each candidate region.
Leonardo Collado-Torres
## BigWig files are not supported in Windows
if (.Platform$OS.type != "windows") {
## Get data
library("derfinderData")
## Identify sample files
sampleFiles <- rawFiles(system.file("extdata", "AMY",
package =
"derfinderData"
), samplepatt = "bw", fileterm = NULL)
names(sampleFiles) <- gsub(".bw", "", names(sampleFiles))
## Create the mean bigwig file. This file is normally created by Rail
## but in this example we'll create it manually.
library("GenomicRanges")
fullCov <- fullCoverage(files = sampleFiles, chrs = "chr21")
meanCov <- Reduce("+", fullCov$chr21) / ncol(fullCov$chr21)
createBw(list("chr21" = DataFrame("meanChr21" = meanCov)),
keepGR =
FALSE
)
summaryFile <- "meanChr21.bw"
## Get the regions
regionMat <- railMatrix(
chrs = "chr21", summaryFiles = summaryFile,
sampleFiles = sampleFiles, L = 76, cutoff = 5.1,
maxClusterGap = 3000L
)
## Explore results
names(regionMat$chr21)
regionMat$chr21$regions
dim(regionMat$chr21$coverageMatrix)
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.