View source: R/coverage_helpers.R
scaledWindowPositions | R Documentation |
For example scale a coverage table of a all human CDS to width 100
scaledWindowPositions(
grl,
reads,
scaleTo = 100,
scoring = "meanPos",
weight = "score",
is.sorted = FALSE,
drop.zero.dt = FALSE
)
grl |
a |
reads |
a |
scaleTo |
an integer (100), if windows have different size, a meta window can not directly be created, since a meta window must have equal size for all windows. Rescale all windows to scaleTo. i.e c(1,2,3) -> size 2 -> c(1, mean(2,3)) etc. Can also be a vector, 1 number per grl group. |
scoring |
a character, one of (zscore, transcriptNormalized, mean, median, sum, log2sum, log10sum, sumLength, meanPos and frameSum, periodic, NULL). More info in details |
weight |
(default: 'score'), if defined a character name of valid meta column in subject. GRanges("chr1", 1, "+", score = 5), would mean score column tells that this alignment region was found 5 times. Formats which loads a score column like this: Bigwig, wig, ORFik ofst, collapsed bam, bedoc and .bedo. As do CAGEr CAGE files and many other package formats. You can also assign a score column manually. |
is.sorted |
logical (FALSE), is grl sorted. That is + strand groups in increasing ranges (1,2,3), and - strand groups in decreasing ranges (3,2,1) |
drop.zero.dt |
logical FALSE, if TRUE and as.data.table is TRUE, remove all 0 count positions. This greatly speeds up and most importantly, greatly reduces memory usage. Will not change any plots, unless 0 positions are used in some sense. (mean, median, zscore coverage will only scale differently) |
Nice for making metaplots, the score will be mean of merged positions.
A data.table with scored counts (counts) of reads mapped to positions (position) specified in windows along with frame (frame).
Other coverage:
coverageScorings()
,
metaWindow()
,
regionPerReadLength()
,
windowPerReadLength()
library(GenomicRanges)
windows <- GRangesList(GRanges("chr1", IRanges(1, 200), "-"))
x <- GenomicRanges::GRanges(
seqnames = "chr1",
ranges = IRanges::IRanges(c(1, 100, 199), c(2, 101, 200)),
strand = "-")
scaledWindowPositions(windows, x, scaleTo = 100)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.