View source: R/ranges_helpers.R
extendTrailers | R Documentation |
Will extend the trailers or transcripts downstream (3' end) by extension.
The extension is general not relative, that means splicing
will not be taken into account.
Requires the grl
to be sorted beforehand,
use sortPerGroup
to get sorted grl.
extendTrailers(
grl,
extension = 1000L,
is.circular = all(isCircular(grl) %in% TRUE)
)
grl |
usually a |
extension |
an integer, how much to extend downstream (3' end). Eiter single value that will apply for all, or same as length of grl which will give 1 update value per grl object. Or a GRangesList where start / stops sites by strand are the positions to use as new starts. |
is.circular |
logical, default FALSE if not any is: all(isCircular(grl) Where grl is the ranges checked. If TRUE, allow ranges to extend below position 1 on chromosome. Since circular genomes can have negative coordinates. |
an extended GRangeslist
Other ExtendGenomicRanges:
asTX()
,
coveragePerTiling()
,
extendLeaders()
,
reduceKeepAttr()
,
tile1()
,
txSeqsFromFa()
,
windowPerGroup()
library(GenomicFeatures)
samplefile <- system.file("extdata", "hg19_knownGene_sample.sqlite",
package = "GenomicFeatures")
txdb <- loadDb(samplefile)
threeUTRs <- threeUTRsByTranscript(txdb) # <- extract only 5' leaders
tx <- exonsBy(txdb, by = "tx", use.names = TRUE)
## now try(extend downstream 1000):
extendTrailers(threeUTRs, extension = 1000)
## Or on transcripts
extendTrailers(tx, extension = 1000)
## Circular genome (allow negative coordinates)
circular_three <- threeUTRs
isCircular(circular_three) <- rep(TRUE, length(isCircular(circular_three)))
extendTrailers(circular_three, extension = 126200008L)[41] # <- negative stop coordinate
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.