View source: R/getFeatureRanges.R
getFeatureRanges | R Documentation |
Generate a GRangesList object with genomic ranges for (any combination of) spliced transcripts, unspliced transcripts and introns.
getFeatureRanges(
gtf,
featureType = c("spliced", "intron"),
intronType = "separate",
flankLength = 90L,
joinOverlappingIntrons = FALSE,
collapseIntronsByGene = FALSE,
keepIntronInFeature = FALSE,
verbose = TRUE
)
gtf |
Path to gtf file. |
featureType |
Character vector indicating the type(s) of features to
extract, any subset of |
intronType |
Character vector indicating how to define the introns
(only used if "intron" is part of |
flankLength |
Integer scalar indicating the length of the flanking
sequence added to each side of each extracted intron (only used if
"intron" is included in |
joinOverlappingIntrons |
Logical scalar indicating whether two introns that overlap (after adding the flanking sequences) should be joined into one feature. |
collapseIntronsByGene |
Logical scalar indicating whether to collapse overlapping intron ranges by gene after extraction. |
keepIntronInFeature |
Logical scalar indicating whether introns
(after adding the flank length) should be restricted to stay within the
boundaries of the feature they were generated from. For backwards
compatibility, the default is |
verbose |
Logical scalar, whether to print out progress messages. |
Returns a GRangesList
object where each element represents
one extracted feature. The metadata of this object contains two
data.frame
s mapping corresponding identifiers between the
different feature types, as well as a list of all features for each type.
Charlotte Soneson
## Get feature ranges
grl <- getFeatureRanges(
gtf = system.file("extdata/small_example.gtf", package = "eisaR"),
featureType = c("spliced", "intron"),
intronType = "separate",
flankLength = 5L,
joinOverlappingIntrons = FALSE,
collapseIntronsByGene = FALSE,
verbose = TRUE
)
## GRangesList
grl
## Corresponding transcript/gene IDs
S4Vectors::metadata(grl)$corrtx
S4Vectors::metadata(grl)$corrgene
## List of features of different types
S4Vectors::metadata(grl)$featurelist
## Get feature sequences
if (requireNamespace("BSgenome", quietly = TRUE) &&
requireNamespace("GenomicFeatures", quietly = TRUE)) {
library(BSgenome)
genome <- Biostrings::readDNAStringSet(
system.file("extdata/small_example_genome.fa", package = "eisaR"))
seqs <- GenomicFeatures::extractTranscriptSeqs(x = genome,
transcripts = grl)
seqs
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.