translationalEff | R Documentation |
Uses RnaSeq and RiboSeq to get translational efficiency of every element in 'grl'. Translational efficiency is defined as:
(density of RPF within ORF) / (RNA expression of ORFs transcript)
translationalEff(
grl,
RNA,
RFP,
tx,
with.fpkm = FALSE,
pseudoCount = 0,
librarySize = "full",
weight.RFP = 1L,
weight.RNA = 1L
)
grl |
a |
RNA |
RnaSeq reads as |
RFP |
RiboSeq reads as |
tx |
a GRangesList of the transcripts. If you used cage data, then the tss for the the leaders have changed, therefor the tx lengths have changed. To account for that call: ' translationalEff(grl, RNA, RFP, tx = extendLeaders(tx, cageFiveUTRs)) ' where cageFiveUTRs are the reannotated by CageSeq data leaders. |
with.fpkm |
logical, default: FALSE, if true return the fpkm values together with translational efficiency as a data.table |
pseudoCount |
a numeric, default 0, set it to 1 if you want to avoid NA and inf values. |
librarySize |
either numeric value or character vector. Default ("full"), number of alignments in library (reads). If you just have a subset, you can give the value by librarySize = length(wholeLib) or sum(wholeLib$score), if you want lib size to be only number of reads overlapping grl, do: librarySize = "overlapping" sum(countOverlaps(reads, grl) > 0), if reads[1] has 3 hits in grl, and reads[2] has 2 hits, librarySize will be 2, not 5. You can also get the inverse overlap, if you want lib size to be total number of overlaps, do: librarySize = "DESeq" This is standard fpkm way of DESeq2::fpkm(robust = FALSE) sum(countOverlaps(grl, reads)) if grl[1] has 3 reads and grl[2] has 2 reads, librarySize is 5, not 2. |
weight.RFP |
a vector (default: 1L). Can also be character name of column in RFP. As in translationalEff(weight = "score") for: GRanges("chr1", 1, "+", score = 5), would mean score column tells that this alignment region was found 5 times. |
weight.RNA |
Same as weightRFP but for RNA weights. (default: 1L) |
a numeric vector of fpkm ratios, if with.fpkm is TRUE, return a data.table with te and fpkm values (total 3 columns then)
doi: 10.1126/science.1168978
Other features:
computeFeatures()
,
computeFeaturesCage()
,
countOverlapsW()
,
disengagementScore()
,
distToCds()
,
distToTSS()
,
entropy()
,
floss()
,
fpkm()
,
fpkm_calc()
,
fractionLength()
,
initiationScore()
,
insideOutsideORF()
,
isInFrame()
,
isOverlapping()
,
kozakSequenceScore()
,
orfScore()
,
rankOrder()
,
ribosomeReleaseScore()
,
ribosomeStallingScore()
,
startRegion()
,
startRegionCoverage()
,
stopRegion()
,
subsetCoverage()
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), "+")
RNA <- GRanges("1", IRanges(1, 50), "+")
tx <- GRangesList(tx1 = GRanges("1", IRanges(1, 50), "+"))
# grl must have same names as cds + _1 etc, so that they can be matched.
te <- translationalEff(grl, RNA, RFP, tx, with.fpkm = TRUE, pseudoCount = 1)
te$fpkmRFP
te$te
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.