Title: | Tools for spliced gene structure manipulation and analysis |
---|---|
Description: | GeneStructureTools can be used to create in silico alternative splicing events, and analyse potential effects this has on functional gene products. |
Authors: | Beth Signal |
Maintainer: | Beth Signal <[email protected]> |
License: | BSD_3_clause + file LICENSE |
Version: | 1.27.0 |
Built: | 2024-11-14 05:50:26 UTC |
Source: | https://github.com/bioc/GeneStructureTools |
Change transcript biotypes to a broader set in a GRanges GTF object
addBroadTypes(gtf)
addBroadTypes(gtf)
gtf |
GRanges object of the GTF |
GRanges object of the GTF with new transcript types
Beth Signal
Other gtf manipulation: UTR2UTR53
,
exonsToTranscripts
,
filterGtfOverlap
,
removeDuplicateTranscripts
,
removeSameExon
,
reorderExonNumbers
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- addBroadTypes(gtf)
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- addBroadTypes(gtf)
Add a retained intron to the transcripts it is skipped by
addIntronInTranscript(flankingExons, exons, whippetDataSet = NULL, match = "exact", glueExons = TRUE)
addIntronInTranscript(flankingExons, exons, whippetDataSet = NULL, match = "exact", glueExons = TRUE)
flankingExons |
data.frame generataed by findIntronContainingTranscripts() |
exons |
GRanges object made from a GTF with ONLY exon annotations (no gene, transcript, CDS etc.) |
whippetDataSet |
whippetDataSet generated from |
match |
what type of match replacement should be done? exact: exact matches to the intron only retain: keep non-exact intron match coordinates in spliced sets, and retain them in retained sets replace: replace non-exact intron match coordinates with event coordinates in spliced sets, and retain in retained sets |
glueExons |
Join together exons that are not seperated by introns? |
GRanges with transcripts containing retained introns
Beth Signal
Other whippet splicing isoform creation: findExonContainingTranscripts
,
findIntronContainingTranscripts
,
findJunctionPairs
,
replaceJunction
,
skipExonInTranscript
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.intronRetention <- filterWhippetEvents(wds, eventTypes="RI") exons.intronRetention <- findIntronContainingTranscripts(wds.intronRetention, exons) IntronRetentionTranscripts <- addIntronInTranscript(exons.intronRetention, exons, whippetDataSet=wds.intronRetention) exonsFromGRanges <- exons[exons$transcript_id=="ENSMUST00000139129.8" & exons$exon_number %in% c(3,4)] intronFromGRanges <- exonsFromGRanges[1] GenomicRanges::start(intronFromGRanges) <- GenomicRanges::end(exonsFromGRanges[exonsFromGRanges$exon_number==3]) GenomicRanges::end(intronFromGRanges) <- GenomicRanges::start(exonsFromGRanges[exonsFromGRanges$exon_number==4]) exons.intronRetention <- findIntronContainingTranscripts(intronFromGRanges, exons) IntronRetentionTranscripts <- addIntronInTranscript(exons.intronRetention, exons, match="retain")
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.intronRetention <- filterWhippetEvents(wds, eventTypes="RI") exons.intronRetention <- findIntronContainingTranscripts(wds.intronRetention, exons) IntronRetentionTranscripts <- addIntronInTranscript(exons.intronRetention, exons, whippetDataSet=wds.intronRetention) exonsFromGRanges <- exons[exons$transcript_id=="ENSMUST00000139129.8" & exons$exon_number %in% c(3,4)] intronFromGRanges <- exonsFromGRanges[1] GenomicRanges::start(intronFromGRanges) <- GenomicRanges::end(exonsFromGRanges[exonsFromGRanges$exon_number==3]) GenomicRanges::end(intronFromGRanges) <- GenomicRanges::start(exonsFromGRanges[exonsFromGRanges$exon_number==4]) exons.intronRetention <- findIntronContainingTranscripts(intronFromGRanges, exons) IntronRetentionTranscripts <- addIntronInTranscript(exons.intronRetention, exons, match="retain")
Creates transcript isoforms from alternative intron usage tested by leafcutter
alternativeIntronUsage(altIntronLocs, exons)
alternativeIntronUsage(altIntronLocs, exons)
altIntronLocs |
data.frame containing information from the per_intron_results.tab file output from leafcutter. Note that only one cluster of alternative introns can be processed at a time. |
exons |
GRanges object made from a GTF with ONLY exon annotations (no gene, transcript, CDS etc.) |
GRanges object with all potential alternative isoforms skipping the introns specified in either the upregulated or downregulated locations
Beth Signal
leafcutterFiles <- list.files(system.file("extdata","leafcutter/", package = "GeneStructureTools"), full.names = TRUE) leafcutterIntrons <- read.delim(leafcutterFiles[grep("intron_results", leafcutterFiles)],stringsAsFactors=FALSE) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] # single cluster processing cluster <- leafcutterIntrons[leafcutterIntrons$cluster=="chr16:clu_1396",] altIsoforms1396 <- alternativeIntronUsage(cluster, exons) unique(altIsoforms1396$transcript_id) cluster <- leafcutterIntrons[leafcutterIntrons$cluster=="chr16:clu_1395",] altIsoforms1395 <- alternativeIntronUsage(cluster, exons) unique(altIsoforms1395$transcript_id) # multiple cluster processing altIsoforms1396plus1395 <- alternativeIntronUsage(cluster, c(exons, altIsoforms1396)) unique(altIsoforms1396plus1395$transcript_id)
leafcutterFiles <- list.files(system.file("extdata","leafcutter/", package = "GeneStructureTools"), full.names = TRUE) leafcutterIntrons <- read.delim(leafcutterFiles[grep("intron_results", leafcutterFiles)],stringsAsFactors=FALSE) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] # single cluster processing cluster <- leafcutterIntrons[leafcutterIntrons$cluster=="chr16:clu_1396",] altIsoforms1396 <- alternativeIntronUsage(cluster, exons) unique(altIsoforms1396$transcript_id) cluster <- leafcutterIntrons[leafcutterIntrons$cluster=="chr16:clu_1395",] altIsoforms1395 <- alternativeIntronUsage(cluster, exons) unique(altIsoforms1395$transcript_id) # multiple cluster processing altIsoforms1396plus1395 <- alternativeIntronUsage(cluster, c(exons, altIsoforms1396)) unique(altIsoforms1396plus1395$transcript_id)
Annotate a GRanges gene model with ORF boundries for visualisation with Gviz
annotateGeneModel(transcripts, orfs)
annotateGeneModel(transcripts, orfs)
transcripts |
GRanges of gene model to be visualised |
orfs |
ORF predictions. Created by getORFs() |
data.frame of a gene model for visualisation
Beth Signal
Other Gviz gene structure visualisation: makeGeneModel
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 # longest ORF for each transcripts orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = TRUE) geneModelAnnotated <- annotateGeneModel(transcript, orfs)
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 # longest ORF for each transcripts orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = TRUE) geneModelAnnotated <- annotateGeneModel(transcript, orfs)
Evaluate the change in an attribute between a set of 'normal' transcripts and 'alternative' transcripts
attrChangeAltSpliced(orfsX, orfsY, attribute = "orf_length", compareBy = "gene", useMax = TRUE, compareUTR = FALSE)
attrChangeAltSpliced(orfsX, orfsY, attribute = "orf_length", compareBy = "gene", useMax = TRUE, compareUTR = FALSE)
orfsX |
orf information for 'normal' transcripts. Generated by getOrfs() |
orfsY |
orf information for 'alternative' transcripts. Generated by getOrfs() |
attribute |
attribute to compare |
compareBy |
compare by 'transcript' isoforms or by 'gene' groups |
useMax |
use max as the summary function when multiple isoforms are aggregated? If FALSE, will use min instead. |
compareUTR |
compare the UTR lengths between transcripts? Only runs if attribute = "orf_length" |
data.frame with attribute changes
Beth Signal
Other transcript isoform comparisons: orfDiff
,
transcriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) orfsSkipped <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="skipped_exon"], BSgenome = g) orfsIncluded <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="included_exon"], BSgenome = g) attrChangeAltSpliced(orfsSkipped, orfsIncluded,attribute = "orf_length")
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) orfsSkipped <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="skipped_exon"], BSgenome = g) orfsIncluded <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="included_exon"], BSgenome = g) attrChangeAltSpliced(orfsSkipped, orfsIncluded,attribute = "orf_length")
Method coordinates
coordinates(whippetDataSet) ## S4 method for signature 'whippetDataSet' coordinates(whippetDataSet)
coordinates(whippetDataSet) ## S4 method for signature 'whippetDataSet' coordinates(whippetDataSet)
whippetDataSet |
whippetDataSet generated from |
whippet splicing event coordinates as a GRanges object
Other whippet data processing: diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) coordinates <- coordinates(wds)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) coordinates <- coordinates(wds)
Convert DEXSeq ids to gene ids
DEXSeqIdsToGeneIds(DEXSeqIds, removeVersion = FALSE, containsE = TRUE)
DEXSeqIdsToGeneIds(DEXSeqIds, removeVersion = FALSE, containsE = TRUE)
DEXSeqIds |
vector of DEXSeq group or exon ids |
removeVersion |
remove the version (.xx) of the gene? |
containsE |
do the DEXSeq exons ids contain :E00X? |
vector of unique gene ids
Beth Signal
Other DEXSeq processing methods: findDEXexonType
,
summariseExonTypes
# multiple genes in name DEXSeqId <- "ENSMUSG00000027618.17+ENSMUSG00000098950.7+ENSMUSG00000089824.10+ENSMUSG00000074643.12" DEXSeqIdsToGeneIds(DEXSeqId) # exonic part number in id DEXSeqIdsToGeneIds("ENSMUSG00000001017.15:E013", removeVersion=TRUE)
# multiple genes in name DEXSeqId <- "ENSMUSG00000027618.17+ENSMUSG00000098950.7+ENSMUSG00000089824.10+ENSMUSG00000074643.12" DEXSeqIdsToGeneIds(DEXSeqId) # exonic part number in id DEXSeqIdsToGeneIds("ENSMUSG00000001017.15:E013", removeVersion=TRUE)
Method diffSplicingResults
diffSplicingResults(whippetDataSet) ## S4 method for signature 'whippetDataSet' diffSplicingResults(whippetDataSet)
diffSplicingResults(whippetDataSet) ## S4 method for signature 'whippetDataSet' diffSplicingResults(whippetDataSet)
whippetDataSet |
whippetDataSet generated from |
differential splicing results data.frame (originally from a whippet .diff file)
Other whippet data processing: coordinates
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) diffSplicingResults <- diffSplicingResults(wds)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) diffSplicingResults <- diffSplicingResults(wds)
Convert an exon-level gtf annotation to a transcript-level gtf annotation
exonsToTranscripts(exons)
exonsToTranscripts(exons)
exons |
GRanges object with exons |
GRanges object with transcripts
Beth Signal
Other gtf manipulation: UTR2UTR53
,
addBroadTypes
,
filterGtfOverlap
,
removeDuplicateTranscripts
,
removeSameExon
,
reorderExonNumbers
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon" & gtf$transcript_id=="ENSMUST00000126412.1"] exons transcripts <- exonsToTranscripts(exons) transcripts
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon" & gtf$transcript_id=="ENSMUST00000126412.1"] exons transcripts <- exonsToTranscripts(exons) transcripts
Filter a GTF overlap to remove exons when exon is annotated as a CDS/UTR
filterGtfOverlap(gtf.from)
filterGtfOverlap(gtf.from)
gtf.from |
GRanges object of the GTF produced from an overlap |
GRanges object of the GTF with redundant exons removed
Beth Signal
Other gtf manipulation: UTR2UTR53
,
addBroadTypes
,
exonsToTranscripts
,
removeDuplicateTranscripts
,
removeSameExon
,
reorderExonNumbers
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) overlap <- as.data.frame(GenomicRanges::findOverlaps(gtf[which(gtf$type=="CDS")[1]], gtf)) table(gtf$type[overlap$subjectHits]) overlapFiltered <- filterGtfOverlap(gtf[overlap$subjectHits]) table(overlapFiltered$type[overlap$subjectHits]) overlap <- as.data.frame(GenomicRanges::findOverlaps(gtf[which( gtf$transcript_type=="retained_intron")[1]],gtf)) table(gtf$type[overlap$subjectHits]) overlapFiltered <- filterGtfOverlap(gtf[overlap$subjectHits]) table(overlapFiltered$type[overlap$subjectHits])
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) overlap <- as.data.frame(GenomicRanges::findOverlaps(gtf[which(gtf$type=="CDS")[1]], gtf)) table(gtf$type[overlap$subjectHits]) overlapFiltered <- filterGtfOverlap(gtf[overlap$subjectHits]) table(overlapFiltered$type[overlap$subjectHits]) overlap <- as.data.frame(GenomicRanges::findOverlaps(gtf[which( gtf$transcript_type=="retained_intron")[1]],gtf)) table(gtf$type[overlap$subjectHits]) overlapFiltered <- filterGtfOverlap(gtf[overlap$subjectHits]) table(overlapFiltered$type[overlap$subjectHits])
Filter out significant events from a whippet diff comparison
filterWhippetEvents(whippetDataSet, probability = 0.95, psiDelta = 0.1, eventTypes = "all", idList = NA, minCounts = NA, medianCounts = NA, sampleTable)
filterWhippetEvents(whippetDataSet, probability = 0.95, psiDelta = 0.1, eventTypes = "all", idList = NA, minCounts = NA, medianCounts = NA, sampleTable)
whippetDataSet |
whippetDataSet generated from |
probability |
minimum probability required to call event as significant |
psiDelta |
minimum change in psi required to call an event as significant |
eventTypes |
which event type to filter for? default = |
idList |
(optional) list of gene ids to filter for |
minCounts |
minumum number of counts for all replicates in at least one condition to call an event as significant |
medianCounts |
median count for all replicates in at least one condition to call an event as significant |
sampleTable |
data.frame with sample names and conditions. Only needed if filtering with counts. |
filtered whippet differential comparison data.frame
Beth Signal
Other whippet data processing: coordinates
,
diffSplicingResults
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds)
Find a DEXSeq exons' biotype
findDEXexonType(DEXSeqExonId, DEXSeqGtf, gtf, set = "overlap")
findDEXexonType(DEXSeqExonId, DEXSeqGtf, gtf, set = "overlap")
DEXSeqExonId |
vector of DEXSeq exon ids |
DEXSeqGtf |
GRanges object of the DEXSeq formatted gtf |
gtf |
GRanges object of the GTF annotated with exon biotypes - i.e. exon, CDS, UTR |
set |
which overlapping set of exon biotypes to return - to, from, and/or overlap |
overlaping types
Beth Signal
Other DEXSeq processing methods: DEXSeqIdsToGeneIds
,
summariseExonTypes
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") DEXSeqGtfFile <- system.file("extdata","gencode.vM14.dexseq.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) DEXSeqGtf <- rtracklayer::import(DEXSeqGtfFile) findDEXexonType("ENSMUSG00000032366.15:E028", DEXSeqGtf, gtf) DEXSeqResultsFile <- system.file("extdata","dexseq_results_significant.txt", package = "GeneStructureTools") DEXSeqResults <- read.table(DEXSeqResultsFile, sep="\t") findDEXexonType(rownames(DEXSeqResults), DEXSeqGtf, gtf)
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") DEXSeqGtfFile <- system.file("extdata","gencode.vM14.dexseq.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) DEXSeqGtf <- rtracklayer::import(DEXSeqGtfFile) findDEXexonType("ENSMUSG00000032366.15:E028", DEXSeqGtf, gtf) DEXSeqResultsFile <- system.file("extdata","dexseq_results_significant.txt", package = "GeneStructureTools") DEXSeqResults <- read.table(DEXSeqResultsFile, sep="\t") findDEXexonType(rownames(DEXSeqResults), DEXSeqGtf, gtf)
Given the location of a whole spliced in exon, find transcripts which can splice out this exon
findExonContainingTranscripts(input, exons, variableWidth = 0, findIntrons = FALSE, transcripts)
findExonContainingTranscripts(input, exons, variableWidth = 0, findIntrons = FALSE, transcripts)
input |
whippetDataSet generated from |
exons |
GRanges object made from a GTF containing exon coordinates |
variableWidth |
How many nts overhang is allowed for finding matching exons (default = 0, i.e. complete match) |
findIntrons |
Find transcripts where the event occurs within the intron? |
transcripts |
GRanges object made from a GTF containing transcript coordinates (only required if findIntrons=TRUE) |
data.frame with all overlapping exons
Beth Signal
Other whippet splicing isoform creation: addIntronInTranscript
,
findIntronContainingTranscripts
,
findJunctionPairs
,
replaceJunction
,
skipExonInTranscript
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) exonFromGRanges <- exons[exons$exon_id == "ENSMUSE00001271768.1"] exons.exonSkip <- findExonContainingTranscripts(exonFromGRanges, exons, variableWidth=0, findIntrons=FALSE, transcripts)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) exonFromGRanges <- exons[exons$exon_id == "ENSMUSE00001271768.1"] exons.exonSkip <- findExonContainingTranscripts(exonFromGRanges, exons, variableWidth=0, findIntrons=FALSE, transcripts)
Given the location of a whole retained intron, find transcripts which splice out this intron
findIntronContainingTranscripts(input, exons, match = "exact")
findIntronContainingTranscripts(input, exons, match = "exact")
input |
whippetDataSet generated from |
exons |
GRanges object made from a GTF with ONLY exon annotations (no gene, transcript, CDS etc.) |
match |
what type of matching to perform? exact = only exons which bound the intron exactly, introns = any exon pairs which overlap the intron, all = any exon pairs AND single exons which overlap the intron |
data.frame with all flanking exon pairs
Beth Signal
Other whippet splicing isoform creation: addIntronInTranscript
,
findExonContainingTranscripts
,
findJunctionPairs
,
replaceJunction
,
skipExonInTranscript
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.intronRetention <- filterWhippetEvents(wds, eventTypes="RI") exons.intronRetention <- findIntronContainingTranscripts(input=wds.intronRetention, exons) exonsFromGRanges <- exons[exons$transcript_id=="ENSMUST00000139129.8" & exons$exon_number %in% c(3,4)] intronFromGRanges <- exonsFromGRanges[1] GenomicRanges::start(intronFromGRanges) <- GenomicRanges::end(exonsFromGRanges[exonsFromGRanges$exon_number==3]) GenomicRanges::end(intronFromGRanges) <- GenomicRanges::start(exonsFromGRanges[exonsFromGRanges$exon_number==4]) exons.intronRetention <- findIntronContainingTranscripts(intronFromGRanges, exons)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.intronRetention <- filterWhippetEvents(wds, eventTypes="RI") exons.intronRetention <- findIntronContainingTranscripts(input=wds.intronRetention, exons) exonsFromGRanges <- exons[exons$transcript_id=="ENSMUST00000139129.8" & exons$exon_number %in% c(3,4)] intronFromGRanges <- exonsFromGRanges[1] GenomicRanges::start(intronFromGRanges) <- GenomicRanges::end(exonsFromGRanges[exonsFromGRanges$exon_number==3]) GenomicRanges::end(intronFromGRanges) <- GenomicRanges::start(exonsFromGRanges[exonsFromGRanges$exon_number==4]) exons.intronRetention <- findIntronContainingTranscripts(intronFromGRanges, exons)
Find junctions that pair with each end of an AA (alt. acceptor) or AD (alt. donor) whippet range Find junctions that pair with the upsteam/downstream exon of an AF (alt. first exon) or an AL (alt. last exon)
findJunctionPairs(whippetDataSet, jncCoords, type = NA)
findJunctionPairs(whippetDataSet, jncCoords, type = NA)
whippetDataSet |
whippetDataSet generated from |
jncCoords |
GRanges object with Whippet junctions. Generated by readWhippetJNCfiles() |
type |
type of Whippet event (AA/AD/AF/AL). Note only one event type should be processed at a time. |
GRanges object with alternative junctions. Each event should have a set of X (for which the psi measurement is reported) junctions, and alternative Y junctions.
Beth Signal
Other whippet splicing isoform creation: addIntronInTranscript
,
findExonContainingTranscripts
,
findIntronContainingTranscripts
,
replaceJunction
,
skipExonInTranscript
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.altAce <- filterWhippetEvents(wds, eventTypes="AA") jncPairs.altAce <- findJunctionPairs(wds.altAce, type="AA") wds.altDon <- filterWhippetEvents(wds, eventTypes="AD") jncPairs.altDon <- findJunctionPairs(wds.altDon, type="AD") wds.altFirst <- filterWhippetEvents(wds, eventTypes="AF", psiDelta=0.2) jncPairs.altFirst <- findJunctionPairs(wds.altFirst, type="AF") wds.altLast <- filterWhippetEvents(wds, eventTypes="AL", psiDelta=0.2) jncPairs.altLast <- findJunctionPairs(wds.altLast, type="AL")
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.altAce <- filterWhippetEvents(wds, eventTypes="AA") jncPairs.altAce <- findJunctionPairs(wds.altAce, type="AA") wds.altDon <- filterWhippetEvents(wds, eventTypes="AD") jncPairs.altDon <- findJunctionPairs(wds.altDon, type="AD") wds.altFirst <- filterWhippetEvents(wds, eventTypes="AF", psiDelta=0.2) jncPairs.altFirst <- findJunctionPairs(wds.altFirst, type="AF") wds.altLast <- filterWhippetEvents(wds, eventTypes="AL", psiDelta=0.2) jncPairs.altLast <- findJunctionPairs(wds.altLast, type="AL")
Format Whippet co-ordinates as a GRanges object
formatWhippetEvents(whippet)
formatWhippetEvents(whippet)
whippet |
data.frame containing event location information. May be generated by readWhippetDIFFfiles() |
GRanges object with events
Beth Signal
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) diffFiles <- whippetFiles[grep(".diff", whippetFiles)] whippetDiffSplice <- readWhippetDIFFfiles(diffFiles) whippetCoords <- formatWhippetEvents(whippetDiffSplice)
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) diffFiles <- whippetFiles[grep(".diff", whippetFiles)] whippetDiffSplice <- readWhippetDIFFfiles(diffFiles) whippetCoords <- formatWhippetEvents(whippetDiffSplice)
Get open reading frames for transcripts
getOrfs(transcripts, BSgenome = NULL, returnLongestOnly = TRUE, allFrames = FALSE, longest = 1, exportFasta = FALSE, fastaFile = NULL, uORFs = FALSE)
getOrfs(transcripts, BSgenome = NULL, returnLongestOnly = TRUE, allFrames = FALSE, longest = 1, exportFasta = FALSE, fastaFile = NULL, uORFs = FALSE)
transcripts |
GRanges object with ONLY exon annotations (no gene, transcript, CDS etc.) with all transcripts for orf retrevial |
BSgenome |
BSgenome object |
returnLongestOnly |
only return longest ORF? |
allFrames |
return longest ORF for all 3 frames? |
longest |
return x longest ORFs (regardless of frames) |
exportFasta |
export a .fa.gz file with nucleotide sequences for each transcript? |
fastaFile |
file name for .fa.gz export |
uORFs |
get uORF summaries? |
data.frame with longest orf details
Beth Signal
Other ORF annotation: getUOrfs
,
maxLocation
, orfSimilarity
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 # longest ORF for each transcripts orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = TRUE) # longest ORF in all 3 frames for each transcript orfs <- getOrfs(transcript, BSgenome = g, allFrames = TRUE) # longest 3 ORFS in eacht transcript orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = FALSE, longest=3)
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 # longest ORF for each transcripts orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = TRUE) # longest ORF in all 3 frames for each transcript orfs <- getOrfs(transcript, BSgenome = g, allFrames = TRUE) # longest 3 ORFS in eacht transcript orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = FALSE, longest=3)
Get upstream open reading frames for transcripts with annotated main ORFs
getUOrfs(transcripts, BSgenome = NULL, orfs, findExonB = FALSE)
getUOrfs(transcripts, BSgenome = NULL, orfs, findExonB = FALSE)
transcripts |
GRanges object with ONLY exon annotations (no gene, transcript, CDS etc.) with all transcripts for orf retrevial |
BSgenome |
BSgenome object |
orfs |
orf annotation for the transcripts object. Generated by getOrfs(transcripts, ...) |
findExonB |
find the distance to and exon number of the downstream (B) junction? |
data.frame with all upstream ORF details.
Beth Signal
Other ORF annotation: getOrfs
,
maxLocation
, orfSimilarity
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 # longest ORF for each transcripts orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = FALSE) uORFS <- getUOrfs(transcript, BSgenome = g, orfs = orfs, findExonB = TRUE)
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 # longest ORF for each transcripts orfs <- getOrfs(transcript, BSgenome = g, returnLongestOnly = FALSE) uORFS <- getUOrfs(transcript, BSgenome = g, orfs = orfs, findExonB = TRUE)
Method junctions
junctions(whippetDataSet) ## S4 method for signature 'whippetDataSet' junctions(whippetDataSet)
junctions(whippetDataSet) ## S4 method for signature 'whippetDataSet' junctions(whippetDataSet)
whippetDataSet |
whippetDataSet generated from |
junctions GRanges object (originally from a whippet .jnc file)
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) junctions <- junctions(wds)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) junctions <- junctions(wds)
Compare open reading frames for whippet differentially spliced events
leafcutterTranscriptChangeSummary(significantEvents, combineGeneEvents = FALSE, exons, BSgenome, NMD = FALSE, showProgressBar = TRUE, exportGTF = NULL)
leafcutterTranscriptChangeSummary(significantEvents, combineGeneEvents = FALSE, exons, BSgenome, NMD = FALSE, showProgressBar = TRUE, exportGTF = NULL)
significantEvents |
data.frame containing information from the per_intron_results.tab file output from leafcutter. |
combineGeneEvents |
combine clusters occuring in the same gene? Currently not reccomended. |
exons |
GRanges gtf annotation of exons |
BSgenome |
BSGenome object containing the genome for the species analysed |
NMD |
Use NMD predictions? (Note: notNMD must be installed to use this feature) |
showProgressBar |
show a progress bar of alternative isoform generation? |
exportGTF |
file name to export alternative isoform GTFs (default=NULL) |
data.frame containing signficant whippet diff data and ORF change summaries
Beth Signal
leafcutterFiles <- list.files(system.file("extdata","leafcutter/", package = "GeneStructureTools"), full.names = TRUE) leafcutterIntrons <- read.delim(leafcutterFiles[ grep("intron_results", leafcutterFiles)],stringsAsFactors=FALSE) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 leafcutterTranscriptChangeSummary(significantEvents = leafcutterIntrons, exons=exons,BSgenome = g,NMD=FALSE)
leafcutterFiles <- list.files(system.file("extdata","leafcutter/", package = "GeneStructureTools"), full.names = TRUE) leafcutterIntrons <- read.delim(leafcutterFiles[ grep("intron_results", leafcutterFiles)],stringsAsFactors=FALSE) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 leafcutterTranscriptChangeSummary(significantEvents = leafcutterIntrons, exons=exons,BSgenome = g,NMD=FALSE)
Convert GRanges gene model to data.frame for visualisation with Gviz
makeGeneModel(transcript)
makeGeneModel(transcript)
transcript |
GRanges of gene model to be visualised |
data.frame of a gene model for visualisation
Beth Signal
Other Gviz gene structure visualisation: annotateGeneModel
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] geneModel <- makeGeneModel(transcript)
gtf <- rtracklayer::import(system.file("extdata", "example_gtf.gtf", package="GeneStructureTools")) transcript <- gtf[gtf$type=="exon" & gtf$gene_name=="Neurl1a"] geneModel <- makeGeneModel(transcript)
Find the largest distance between two vectors of numbers Helper function for get_orfs
maxLocation(startSite, stopSite, longest = 1)
maxLocation(startSite, stopSite, longest = 1)
startSite |
vector of start sites - i.e Met amino acid positions |
stopSite |
vector of stop sites - i.e Stop (*) amino acid positions |
longest |
which pair to return (1 = longest pair, 2= 2nd longest pair etc.) |
sequential start site and end site with the greatest difference
Beth Signal
Other ORF annotation: getOrfs
,
getUOrfs
, orfSimilarity
starts <- c(1,10,15,25) stops <- c(4,16,50,55) # longest start site = 25, longest stop site = 50 maxLocation(starts, stops, longest = 1) starts <- c(1,10,15,25) stops <- c(4,14,50,55) # longest start site = 15, longest stop site = 50 maxLocation(starts, stops, longest = 1) # 2nd longest start site = 10, 2nd longest stop site = 14 maxLocation(starts, stops, longest = 2)
starts <- c(1,10,15,25) stops <- c(4,16,50,55) # longest start site = 25, longest stop site = 50 maxLocation(starts, stops, longest = 1) starts <- c(1,10,15,25) stops <- c(4,14,50,55) # longest start site = 15, longest stop site = 50 maxLocation(starts, stops, longest = 1) # 2nd longest start site = 10, 2nd longest stop site = 14 maxLocation(starts, stops, longest = 2)
Evaluate changes to ORFs caused by alternative splicing
orfDiff(orfsX, orfsY, filterNMD = TRUE, geneSimilarity = TRUE, compareUTR = TRUE, compareBy = "gene", allORFs = NULL)
orfDiff(orfsX, orfsY, filterNMD = TRUE, geneSimilarity = TRUE, compareUTR = TRUE, compareBy = "gene", allORFs = NULL)
orfsX |
orf information for 'normal' transcripts. Generated by getOrfs() |
orfsY |
orf information for 'alternative' transcripts. Generated by getOrfs() |
filterNMD |
filter orf information for transcripts not targeted by nmd first? |
geneSimilarity |
compare orf to all orfs in gene? |
compareUTR |
compare UTRs? |
compareBy |
compare by 'transcript' isoforms or by 'gene' groups |
allORFs |
orf information for all transcripts for novel sequence comparisons. Generated by getOrfs() |
data.frame with orf changes
Beth Signal
Other transcript isoform comparisons: attrChangeAltSpliced
,
transcriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 orfsProteinCoding <- getOrfs(exons[exons$gene_name=="Prex2" & exons$transcript_type=="protein_coding"], BSgenome = g) orfsNMD <- getOrfs(exons[exons$gene_name=="Prex2" & exons$transcript_type=="nonsense_mediated_decay"], BSgenome = g) orfDiff(orfsProteinCoding, orfsNMD, filterNMD=FALSE) wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) orfsSkipped <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="skipped_exon"], BSgenome = g) orfsIncluded <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="included_exon"], BSgenome = g) orfDiff(orfsSkipped, orfsIncluded, filterNMD=FALSE)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 orfsProteinCoding <- getOrfs(exons[exons$gene_name=="Prex2" & exons$transcript_type=="protein_coding"], BSgenome = g) orfsNMD <- getOrfs(exons[exons$gene_name=="Prex2" & exons$transcript_type=="nonsense_mediated_decay"], BSgenome = g) orfDiff(orfsProteinCoding, orfsNMD, filterNMD=FALSE) wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) orfsSkipped <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="skipped_exon"], BSgenome = g) orfsIncluded <- getOrfs(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="included_exon"], BSgenome = g) orfDiff(orfsSkipped, orfsIncluded, filterNMD=FALSE)
calculate percentage of orfB contained in orfA
orfSimilarity(orfA, orfB, substitutionCost = 100)
orfSimilarity(orfA, orfB, substitutionCost = 100)
orfA |
character string of ORF amino acid sequence |
orfB |
character string of ORF amino acid sequence |
substitutionCost |
cost for substitutions in ORF sequences. Set to 1 if substitutions should be weighted equally to insertions and deletions. |
percentage of orfB contained in orfA
Beth Signal
Other ORF annotation: getOrfs
,
getUOrfs
, maxLocation
orfSimilarity("MFGLDIYAGTRSSFRQFSLT","MFGLDIYAGTRSSFRQFSLT") orfSimilarity("MFGLDIYAGTRSSFRQFSLT","MFGLDIYAFRQFSLT") orfSimilarity("MFGLDIYAFRQFSLT","MFGLDIYAGTRSSFRQFSLT") orfSimilarity("MFGLDIYAGTRXXFRQFSLT","MFGLDIYAGTRSSFRQFSLT") orfSimilarity("MFGLDIYAGTRXXFSLT","MFGLDIYAGTRSSFRQFSLT", 1)
orfSimilarity("MFGLDIYAGTRSSFRQFSLT","MFGLDIYAGTRSSFRQFSLT") orfSimilarity("MFGLDIYAGTRSSFRQFSLT","MFGLDIYAFRQFSLT") orfSimilarity("MFGLDIYAFRQFSLT","MFGLDIYAGTRSSFRQFSLT") orfSimilarity("MFGLDIYAGTRXXFRQFSLT","MFGLDIYAGTRSSFRQFSLT") orfSimilarity("MFGLDIYAGTRXXFSLT","MFGLDIYAGTRSSFRQFSLT", 1)
Annotate introns and exonic parts by overlaping exon biotype
overlapTypes(queryCoords, gtf, set = c("from", "to", "overlap"))
overlapTypes(queryCoords, gtf, set = c("from", "to", "overlap"))
queryCoords |
GRanges object of the query regions |
gtf |
GRanges object of the GTF annotated with exon biotypes - i.e. exon, CDS, UTR |
set |
which overlapping set of exon biotypes to return - to, from, and/or overlap |
overlaping types in a data.frame
Beth Signal
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") DEXSeqGtfFile <- system.file("extdata","gencode.vM14.dexseq.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) DEXSeqGtf <- rtracklayer::import(DEXSeqGtfFile) overlapTypes(DEXSeqGtf[1:10], gtf)
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") DEXSeqGtfFile <- system.file("extdata","gencode.vM14.dexseq.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) DEXSeqGtf <- rtracklayer::import(DEXSeqGtfFile) overlapTypes(DEXSeqGtf[1:10], gtf)
Method readCounts
readCounts(whippetDataSet) ## S4 method for signature 'whippetDataSet' readCounts(whippetDataSet)
readCounts(whippetDataSet) ## S4 method for signature 'whippetDataSet' readCounts(whippetDataSet)
whippetDataSet |
whippetDataSet generated from |
whippet read count data.frame (originally from a whippet .psi file)
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) readCounts <- readCounts(wds)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) readCounts <- readCounts(wds)
Import whippet results files as a whippetDataSet
readWhippetDataSet(filePath = ".")
readWhippetDataSet(filePath = ".")
filePath |
path to whippet output files |
whippetDataSet
Beth Signal
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles)
Read in a list of whippet .diff.gz files and format as a data.frame
readWhippetDIFFfiles(files)
readWhippetDIFFfiles(files)
files |
vector of *.diff.gz file names |
data.frame with junction counts for all files
Beth Signal
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) diffFiles <- whippetFiles[grep(".diff", whippetFiles)] whippetDiffSplice <- readWhippetDIFFfiles(diffFiles)
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) diffFiles <- whippetFiles[grep(".diff", whippetFiles)] whippetDiffSplice <- readWhippetDIFFfiles(diffFiles)
Read in a list of whippet .jnc.gz files and format as a GRanges object
readWhippetJNCfiles(files)
readWhippetJNCfiles(files)
files |
vector of *.jnc.gz file names |
GRanges object with junctions
Beth Signal
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetPSIfiles
,
whippetTranscriptChangeSummary
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) jncFiles <- whippetFiles[grep(".jnc", whippetFiles)] whippetJNC <- readWhippetJNCfiles(jncFiles)
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) jncFiles <- whippetFiles[grep(".jnc", whippetFiles)] whippetJNC <- readWhippetJNCfiles(jncFiles)
Read in a list of whippet .psi.gz files and format as a data.frame
readWhippetPSIfiles(files, attribute = "Total_Reads", maxNA = NA)
readWhippetPSIfiles(files, attribute = "Total_Reads", maxNA = NA)
files |
vector of *.psi.gz file names |
attribute |
which attribute from the PSI files to use (Total_Reads, Psi, CI_width) |
maxNA |
maximum number of NA values allowed before a site is removed |
data.frame with junction counts for all files
Beth Signal
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
whippetTranscriptChangeSummary
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) psiFiles <- whippetFiles[grep(".psi", whippetFiles)] whippetPSI <- readWhippetPSIfiles(psiFiles)
whippetFiles <- list.files(system.file("extdata","whippet/", package = "GeneStructureTools"), full.names = TRUE) psiFiles <- whippetFiles[grep(".psi", whippetFiles)] whippetPSI <- readWhippetPSIfiles(psiFiles)
Removes Structural duplicates of transcripts in a GRanges object Note that duplicates must have different transcript ids.
removeDuplicateTranscripts(transcripts)
removeDuplicateTranscripts(transcripts)
transcripts |
GRanges object with transcript structures in exon form |
GRanges object with unique transcript structures in exon form
Beth Signal
Other gtf manipulation: UTR2UTR53
,
addBroadTypes
,
exonsToTranscripts
,
filterGtfOverlap
,
removeSameExon
,
reorderExonNumbers
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] exons.altName <- exons exons.altName$transcript_id <- paste(exons.altName$transcript_id, "duplicated", sep="_") exons.duplicated <- c(exons, exons.altName) length(exons.duplicated) exons.deduplicated <- removeDuplicateTranscripts(exons.duplicated) length(exons.deduplicated)
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] exons.altName <- exons exons.altName$transcript_id <- paste(exons.altName$transcript_id, "duplicated", sep="_") exons.duplicated <- c(exons, exons.altName) length(exons.duplicated) exons.deduplicated <- removeDuplicateTranscripts(exons.duplicated) length(exons.deduplicated)
Removes structural duplicates of exons in a GRanges object
removeSameExon(exons)
removeSameExon(exons)
exons |
GRanges object with exons |
GRanges object with unique exons
Beth Signal
Other gtf manipulation: UTR2UTR53
,
addBroadTypes
,
exonsToTranscripts
,
filterGtfOverlap
,
removeDuplicateTranscripts
,
reorderExonNumbers
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] exons.duplicated <- c(exons[1:4], exons[1:4]) length(exons.duplicated) exons.deduplicated <- removeSameExon(exons.duplicated) length(exons.deduplicated)
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] exons.duplicated <- c(exons[1:4], exons[1:4]) length(exons.duplicated) exons.deduplicated <- removeSameExon(exons.duplicated) length(exons.deduplicated)
Remove version number from ensembl gene/transcript ids
removeVersion(ids)
removeVersion(ids)
ids |
vector of ensembl ids |
vector of ensembl ids without the version number
Beth Signal
removeVersion("ENSMUSG00000001017.15")
removeVersion("ENSMUSG00000001017.15")
Reorder the exon numbers in a gtf annotation
reorderExonNumbers(exons, by = "transcript_id")
reorderExonNumbers(exons, by = "transcript_id")
exons |
GRanges object made from a GTF with ONLY exon annotations (no gene, transcript, CDS etc.) |
by |
what column are the transcripts grouped by? |
The same input GRanges, but with exon numbers reordered.
Beth Signal
Other gtf manipulation: UTR2UTR53
,
addBroadTypes
,
exonsToTranscripts
,
filterGtfOverlap
,
removeDuplicateTranscripts
,
removeSameExon
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] exons <- reorderExonNumbers(exons)
gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] exons <- reorderExonNumbers(exons)
Find transcripts containing/overlapping junctions and replace them with alternative junctions
replaceJunction(whippetDataSet, junctionPairs, exons, type = NA)
replaceJunction(whippetDataSet, junctionPairs, exons, type = NA)
whippetDataSet |
whippetDataSet generated from |
junctionPairs |
GRanges object with alternative Whippet junctions. Generated by findJunctionPairs() |
exons |
GRanges object made from a GTF containing exon coordinates |
type |
type of Whippet event (AA/AD/AF/AL). Note only one event type should be processed at a time. |
GRanges object with transcripts containing alternative junctions.
Beth Signal
Other whippet splicing isoform creation: addIntronInTranscript
,
findExonContainingTranscripts
,
findIntronContainingTranscripts
,
findJunctionPairs
,
skipExonInTranscript
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.altAce <- filterWhippetEvents(wds, eventTypes="AA") jncPairs.altAce <- findJunctionPairs(wds.altAce, type="AA") transcripts.altAce <- replaceJunction(wds.altAce, jncPairs.altAce, exons, type="AA") wds.altDon <- filterWhippetEvents(wds, eventTypes="AD") jncPairs.altDon <- findJunctionPairs(wds.altDon, type="AD") transcripts.altDon <- replaceJunction(wds.altDon, jncPairs.altDon, exons, type="AD") wds.altFirst <- filterWhippetEvents(wds, eventTypes="AF", psiDelta=0.2) jncPairs.altFirst <- findJunctionPairs(wds.altFirst, type="AF") transcripts.altFirst <- replaceJunction(wds.altFirst, jncPairs.altFirst, exons, type="AF") wds.altLast <- filterWhippetEvents(wds, eventTypes="AL", psiDelta=0.2) jncPairs.altLast <- findJunctionPairs(wds.altLast, type="AL") transcripts.altLast <- replaceJunction(wds.altLast, jncPairs.altLast, exons, type="AL")
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.altAce <- filterWhippetEvents(wds, eventTypes="AA") jncPairs.altAce <- findJunctionPairs(wds.altAce, type="AA") transcripts.altAce <- replaceJunction(wds.altAce, jncPairs.altAce, exons, type="AA") wds.altDon <- filterWhippetEvents(wds, eventTypes="AD") jncPairs.altDon <- findJunctionPairs(wds.altDon, type="AD") transcripts.altDon <- replaceJunction(wds.altDon, jncPairs.altDon, exons, type="AD") wds.altFirst <- filterWhippetEvents(wds, eventTypes="AF", psiDelta=0.2) jncPairs.altFirst <- findJunctionPairs(wds.altFirst, type="AF") transcripts.altFirst <- replaceJunction(wds.altFirst, jncPairs.altFirst, exons, type="AF") wds.altLast <- filterWhippetEvents(wds, eventTypes="AL", psiDelta=0.2) jncPairs.altLast <- findJunctionPairs(wds.altLast, type="AL") transcripts.altLast <- replaceJunction(wds.altLast, jncPairs.altLast, exons, type="AL")
Remove and include a skipped exon from the transcripts it overlaps
skipExonInTranscript(skippedExons, exons, glueExons = TRUE, whippetDataSet = NULL, match = "exact")
skipExonInTranscript(skippedExons, exons, glueExons = TRUE, whippetDataSet = NULL, match = "exact")
skippedExons |
data.frame generataed by findExonContainingTranscripts() |
exons |
GRanges object made from a GTF with ONLY exon annotations (no gene, transcript, CDS etc.) |
glueExons |
Join together exons that are not seperated by exons? |
whippetDataSet |
whippetDataSet generated from |
match |
what type of match replacement should be done? exact: exact matches to the skipped event only, also removes any intron overlaps skip: keep non-exact exon match coordinates in included sets, and skip them in skipped sets replace: replace non-exact exon match coordinates with event coordinates in included sets, and skip them in skipped sets |
GRanges with transcripts skipping exons
Beth Signal
Other whippet splicing isoform creation: addIntronInTranscript
,
findExonContainingTranscripts
,
findIntronContainingTranscripts
,
findJunctionPairs
,
replaceJunction
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) exonFromGRanges <- exons[exons$exon_id == "ENSMUSE00001271768.1"] exons.exonSkip <- findExonContainingTranscripts(exonFromGRanges, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, match="skip")
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] transcripts <- gtf[gtf$type=="transcript"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) exonFromGRanges <- exons[exons$exon_id == "ENSMUSE00001271768.1"] exons.exonSkip <- findExonContainingTranscripts(exonFromGRanges, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, match="skip")
Summarise exon biotypes to broader categories
summariseExonTypes(types)
summariseExonTypes(types)
types |
vector of exon biotypes |
vector of broader exon biotypes
Beth Signal
Other DEXSeq processing methods: DEXSeqIdsToGeneIds
,
findDEXexonType
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") DEXSeqGtfFile <- system.file("extdata","gencode.vM14.dexseq.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) DEXSeqGtf <- rtracklayer::import(DEXSeqGtfFile) findDEXexonType("ENSMUSG00000032366.15:E028", DEXSeqGtf, gtf) DEXSeqResultsFile <- system.file("extdata","dexseq_results_significant.txt", package = "GeneStructureTools") DEXSeqResults <- read.table(DEXSeqResultsFile, sep="\t") types <- findDEXexonType(rownames(DEXSeqResults), DEXSeqGtf, gtf) summarisedTypes <- summariseExonTypes(types) table(types, summarisedTypes)
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") DEXSeqGtfFile <- system.file("extdata","gencode.vM14.dexseq.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) DEXSeqGtf <- rtracklayer::import(DEXSeqGtfFile) findDEXexonType("ENSMUSG00000032366.15:E028", DEXSeqGtf, gtf) DEXSeqResultsFile <- system.file("extdata","dexseq_results_significant.txt", package = "GeneStructureTools") DEXSeqResults <- read.table(DEXSeqResultsFile, sep="\t") types <- findDEXexonType(rownames(DEXSeqResults), DEXSeqGtf, gtf) summarisedTypes <- summariseExonTypes(types) table(types, summarisedTypes)
Compare open reading frames for two sets of paired transcripts
transcriptChangeSummary(transcriptsX, transcriptsY, BSgenome, exons, NMD = FALSE, NMDModel = NULL, compareBy = "gene", orfPrediction = "allFrames", compareToGene = FALSE, whippetDataSet = NULL, exportGTF = NULL)
transcriptChangeSummary(transcriptsX, transcriptsY, BSgenome, exons, NMD = FALSE, NMDModel = NULL, compareBy = "gene", orfPrediction = "allFrames", compareToGene = FALSE, whippetDataSet = NULL, exportGTF = NULL)
transcriptsX |
GRanges object with exon annotations for all transcripts to be compared for the 'normal' condition |
transcriptsY |
GRanges object with exon annotations for all transcripts to be compared for the 'alternative' condition |
BSgenome |
BSGenome object containing the genome for the species analysed |
exons |
GRanges object made from a GTF containing exon coordinates |
NMD |
Use NMD predictions? (Note: notNMD must be installed to use this feature) |
NMDModel |
Use the "base" or "lncRNA" NMD model? |
compareBy |
compare isoforms by 'transcript' id, or aggregate all changes occuring by 'gene' |
orfPrediction |
What type of orf predictions to return. default= |
compareToGene |
compare alternative isoforms to all normal gene isoforms (in exons) |
whippetDataSet |
whippetDataSet generated from |
exportGTF |
file name to export alternative isoform GTFs (default= |
Summarised ORF changes data.frame
Beth Signal
Other transcript isoform comparisons: attrChangeAltSpliced
,
orfDiff
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) transcriptChangeSummary(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="included_exon"], ExonSkippingTranscripts[ExonSkippingTranscripts$set=="skipped_exon"], BSgenome=g,exons)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) exons <- gtf[gtf$type=="exon"] g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 wds.exonSkip <- filterWhippetEvents(wds, eventTypes="CE",psiDelta = 0.2) exons.exonSkip <- findExonContainingTranscripts(wds.exonSkip, exons, variableWidth=0, findIntrons=FALSE, transcripts) ExonSkippingTranscripts <- skipExonInTranscript(exons.exonSkip, exons, whippetDataSet=wds.exonSkip) transcriptChangeSummary(ExonSkippingTranscripts[ExonSkippingTranscripts$set=="included_exon"], ExonSkippingTranscripts[ExonSkippingTranscripts$set=="skipped_exon"], BSgenome=g,exons)
Annotate UTRs from Gencode GTF as 5' or 3'
UTR2UTR53(gtf)
UTR2UTR53(gtf)
gtf |
GRanges object of the GTF |
gtf annotation GRanges object
Beth Signal
Other gtf manipulation: addBroadTypes
,
exonsToTranscripts
,
filterGtfOverlap
,
removeDuplicateTranscripts
,
removeSameExon
,
reorderExonNumbers
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) table(gtf$type)
gtfFile <- system.file("extdata","example_gtf.gtf", package = "GeneStructureTools") gtf <- rtracklayer::import(gtfFile) gtf <- UTR2UTR53(gtf) table(gtf$type)
Class whippetDataSet
contains information read from whippet output files
Compare open reading frames for whippet differentially spliced events
whippetTranscriptChangeSummary(whippetDataSet, gtf.all = NULL, BSgenome, eventTypes = "all", exons = NULL, transcripts = NULL, NMD = FALSE, exportGTF = NULL)
whippetTranscriptChangeSummary(whippetDataSet, gtf.all = NULL, BSgenome, eventTypes = "all", exons = NULL, transcripts = NULL, NMD = FALSE, exportGTF = NULL)
whippetDataSet |
whippetDataSet generated from |
gtf.all |
GRanges gtf annotation (can be used instead of specifying exons and transcripts) |
BSgenome |
BSGenome object containing the genome for the species analysed |
eventTypes |
which event type to filter for? default = "all" |
exons |
GRanges gtf annotation of exons |
transcripts |
GRanges gtf annotation of transcripts |
NMD |
Use NMD predictions? (Note: notNMD must be installed to use this feature) |
exportGTF |
file name to export alternative isoform GTFs (default=NULL) |
data.frame containing signficant whippet diff data and ORF change summaries
Beth Signal
Other whippet data processing: coordinates
,
diffSplicingResults
,
filterWhippetEvents
,
formatWhippetEvents
,
junctions
, readCounts
,
readWhippetDIFFfiles
,
readWhippetDataSet
,
readWhippetJNCfiles
,
readWhippetPSIfiles
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 whippetTranscriptChangeSummary(wds, gtf.all=gtf,BSgenome = g)
whippetFiles <- system.file("extdata","whippet/", package = "GeneStructureTools") wds <- readWhippetDataSet(whippetFiles) wds <- filterWhippetEvents(wds) gtf <- rtracklayer::import(system.file("extdata","example_gtf.gtf", package = "GeneStructureTools")) g <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10 whippetTranscriptChangeSummary(wds, gtf.all=gtf,BSgenome = g)