coveragePerTiling | R Documentation |
It tiles each GRangesList group to width 1, and finds hits per position.
A range from 1:5 will split into c(1,2,3,4,5) and count hits on each.
This is a safer speedup of coverageByTranscript from GenomicFeatures.
It also gives the possibility to return as data.table, for faster
computations.
coveragePerTiling(
grl,
reads,
is.sorted = FALSE,
keep.names = TRUE,
as.data.table = FALSE,
withFrames = FALSE,
weight = "score",
drop.zero.dt = FALSE,
fraction = NULL
)
grl |
a |
reads |
a |
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) |
keep.names |
logical (TRUE), keep names or not. If as.data.table is TRUE, names (genes column) will be a factor column, if FALSE it will be an integer column (index of gene), so first input grl element is 1. Dropping names gives ~ 20 % speedup. If drop.zero.dt is FALSE, data.table will not return names, will use index (to avoid memory explosion). |
as.data.table |
a logical (FALSE), return as data.table with 2 columns, position and count. |
withFrames |
a logical (FALSE), only available if as.data.table is TRUE, return the ORF frame, 1,2,3, where position 1 is 1, 2 is 2 and 4 is 1 etc. |
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. ORFik ofst, bedoc and .bedo files contains a score column like this. As do CAGEr CAGE files and many other package formats. You can also assign a score column manually. |
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) |
fraction |
integer or character, a description column. Useful for grouping multiple outputs together. If returned as Rle, this is added as: metadata(coverage) <- list(fraction = fraction). If as.data.table it will be added as an additional column. |
NOTE: If reads contains a $score column, it will presume that this is the number of replicates per reads, weights for the coverage() function. So delete the score column or set weight to something else if this is not wanted.
a numeric RleList, one numeric-Rle per group with # of hits per position. Or data.table if as.data.table is TRUE, with column names c("count" [numeric or integer], "genes" [integer], "position" [integer])
Other ExtendGenomicRanges:
asTX()
,
extendLeaders()
,
extendTrailers()
,
reduceKeepAttr()
,
tile1()
,
txSeqsFromFa()
,
windowPerGroup()
ORF <- GRanges(seqnames = "1",
ranges = IRanges(start = c(1, 10, 20),
end = c(5, 15, 25)),
strand = "+")
grl <- GRangesList(tx1_1 = ORF)
RFP <- GRanges("1", IRanges(25, 25), "+")
coveragePerTiling(grl, RFP, is.sorted = TRUE)
# now as data.table with frames
coveragePerTiling(grl, RFP, is.sorted = TRUE, as.data.table = TRUE,
withFrames = TRUE)
# With score column (usually replicated reads on that position)
RFP <- GRanges("1", IRanges(25, 25), "+", score = 5)
dt <- coveragePerTiling(grl, RFP, is.sorted = TRUE,
as.data.table = TRUE, withFrames = TRUE)
class(dt$count) # numeric
# With integer score column (faster and less space usage)
RFP <- GRanges("1", IRanges(25, 25), "+", score = 5L)
dt <- coveragePerTiling(grl, RFP, is.sorted = TRUE,
as.data.table = TRUE, withFrames = TRUE)
class(dt$count) # integer
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.