Recent breakthrough in bioinformatics now allows us to accurately reconstruct and quantify full-length gene isoforms from RNA-sequencing data (via tools such as Cufflinks, Kallisto and Salmon). These tools made it possible to start analyzing alternative isoform usage, but unfortunately RNA-sequencing data is still underutilized since such analyses are both hard to make and therfore only rarely done.
To solve this problem we developed IsoformSwitchAnalyzeR. IsoformSwitchAnalyzeR is an easy to use R package that facilitates statistical identification of isoform switching from RNA-seq derived quantification of novel and/or annotated full-length isoforms. IsoformSwitchAnalyzeR furthermore facilitate integration of many sources of annotation including features such as Open Reading Frame (ORF), protein domains (via Pfam), signal peptides (via SignalP), coding potential (via CPAT) as well as sensitivity to Non-sense Mediated Decay (NMD). The combination of identified isoform switches and their annotation also enables IsoformSwitchAnalyzeR to predict potential functional consequence of the identified isoform switches - such as loss of protein domains or coding potential - thereby identifying isoform switches of particular interest. Lastly, IsoformSwitchAnalyzeR provide article ready visualization of isoform switches as well as multiple layers of summary statistics describing the genome wide occurence of isoform switches and their consequences.
In summary IsoformSwitchAnalyzeR enables analysis of RNA-seq data with isoform resolution with a focus on isoform switching (with predicted consequences) thereby expanding the usability of RNA-seq data.
What To Cite (please remember)
[Quick Start]
Frequently Asked Questions and Problems
The combination of alternative Transcription Start sites (aTSS), Alternative Splicing (AS) and alternative Transcription Termination Sites (aTTS) is often referred to as alternative transcription and is considered the major factors in modifying the pre-RNA and contributing to the complexity of higher organisms. Alternative transcription is widely used as recently demonstrated by The ENCODE Consortium, which found that an average of 6.3 different transcripts were generated per gene, although the individual number of transcripts from a single gene have been reported anywhere from one to thousands.
The importance of analyzing isoforms instead of genes has been highlighted by many examples showing functionally important changes that cannot be detected at gene level. One of these examples is the pyruvate kinase. In normal adult homeostasis, cells use the adult isoform (M1), which supports oxidative phosphorylation. But almost all cancers use the embryonic isoform (M2), which promotes aerobic glycolysis, one of the hallmarks of cancer. Such a shift in isoform usage has been termed isoform switching and frequently will not be detected at the gene level.
In 2010 a breakthrough in bioinformatics with the emergence of tools such as Cufflinks, which allows researchers to reconstruct and quantify full length transcripts from RNA-seq data. Tools for fast transcript quantification such as Salmon and Kallisto were the next breakthrough making it very fast to perform isoform quantification. Such data has the potential to facilitate both genome wide analysis of alternative isoform usage and identification of isoform switching - but unfortunately these types of analysis are still only rarely done.
We hypothesis that there are multiple reasons why RNA-seq data is not used to its full potential:
To solve these problems we developed IsoformSwitchAnalyzeR.
IsoformSwitchAnalyzeR is an easy to use R package that enables the user to import the (novel) full-length derived isoforms from an RNA-seq experiment into R. If annotated transcripts are analyzed, IsoformSwitchAnalyzeR offers integration with the multi-layer information stored in a GTF file including the annotated coding sequences (CDS). If transcript structure were predicted, de-novo or guided, IsoformSwitchAnalyzeR offers a highly accurate tool for identifying the dominant ORF of the isoforms. The knowledge of isoform positions for the CDS/ORF furthermore allows for prediction of sensitivity to Nonsense Mediated Decay (NMD) - the mRNA quality control machinery that degrades isoforms with pre-mature termination codons (PTC).
Next, IsoformSwitchAnalyzeR enables identification of isoform switches via newly developed statistical methods that test each individual isoform for differential usage and thereby identifies the exact isoforms involved in isoform switch.
Since we know the exon structure of the full-length isoform, IsoformSwitchAnalyzeR can extract the underlying nucleotide sequence from a reference genome. This enables integration with the Coding Potential Assessment Tool (CPAT) which predicts the coding potential of an isoform and can furthermore be used to increase accuracy of ORF predictions. By combining the CDS/ORF isoform positions with the nucleotide sequence, we can also extract the most likely amino acid (AA) sequence of the CDS/ORF. The AA sequence enables integration of analysis of protein domains (via Pfam) and signal peptides (via SignalP) - both of which are supported by IsoformSwitchAnalyzeR. Lastly, since the structures of all expressed isoforms from a given gene are known, one can also annotate intron retentions (via spliceR).
Combined, IsoformSwitchAnalyzeR enables annotation of isoforms with intron retentions, ORF, NMD sensitivity, coding potential, protein domains as well as signal peptides, resulting in the identification of important functional consequences of the isoform switches.
IsoformSwitchAnalyzeR contains tools that allow the user to create article ready visualization of both individual isoform switches as well as general common consequences of isoform switches. These visualizations are easy to understand and integrate all the information gathered throughout the workflow. An example of visualization can be found here Examples of visualization.
Lastly IsoformSwitchAnalyzeR is based on standard Bioconductor classes such as GRanges and BSgenome, whereby it supports all species and annotation versions facilitated in the Bioconductor annotation packages.
Back to Table of Content.
IsoformSwitchAnalyzeR is part of the Bioconductor repository and community which means it is distributed with, and dependent on, Bioconductor. Installation of IsoformSwitchAnalyzeR is very easy and can be done from within the R terminal. If it is the first time you use Bioconductor, simply copy-paste the following into your R session to install the basic bioconductor packages:
source("http://bioconductor.org/biocLite.R")
biocLite()
If you already have installed Bioconductor, running these two commands will check whether updates for installed packages are available.
After you have installed the basic bioconductor packages you can install IsoformSwitchAnalyzeR by copy pasting the following into your R session:
source("http://bioconductor.org/biocLite.R")
biocLite("IsoformSwitchAnalyzeR")
This will install the IsoformSwitchAnalyzeR package as well as other R packages that are needed for IsoformSwitchAnalyzeR to work.
This R package comes with a lot of documentation. Much information can be found in the R help files (which can easily be accessed by running the following command in R “?functionName”, for example “?isoformSwitchTest”). Furthermore this vignette contains a lot of information, make sure to read both sources carefully as it will contain the answer to the most Frequently Asked Questions and Problems.
If you have unanswered questions or comments regarding IsoformSwitchAnalyzeR please post them on the associted google group: https://groups.google.com/forum/#!forum/isoformswitchanalyzer (after making sure the question have not already been answered there).
If you want to report a bug (found in the newest version of the R package) please make an issue with a reproducible example at github https://github.com/kvittingseerup/IsoformSwitchAnalyzeR - remember to add the appropriate label.
If you have suggestions for improvements also put them on github (https://github.com/kvittingseerup/IsoformSwitchAnalyzeR) this will allow other people to upvote you idea by reactions thereby showing us there is wide support of implementing your idea.
Back to Table of Content.
The IsoformSwitchAnalyzeR tool is only made possible by a string of other tools and scientific discoveries - please read this section thoroughly and cite the appropriate articles. Note that due to the references being divided into sections some references appear more than once.
If you are using the
Refrences:
The idea behind IsoformSwitchAnalyzeR is to make it easy to do advanced post analysis of full length RNA-seq derived transcripts with a focus on finding, annotating and visualizing isoform switches with functional consequences. IsoformSwitchAnalyzeR therefore performs 3 specific tasks:
A normal workflow for identification and analysis of isoform switches with functional consequences can be divide into two parts (also illustrated below in Figure 1).
1) Extract Isoform Switches and Their Sequences. This part includes importing the data into R, identifying isoform swithces, annotating those switches with open reading frames (ORF) and extract both the nucleotide and peptide sequence. The later step enables the usage of external sequence analysis tools such as
All of the above steps is performed by the high level function:
isoformSwitchAnalysisPart1()
See below for example of usage, and Detailed Workflow for details on the individual steps.
2) Plot All Isoform Switches and Their annotation. This part involves importing and incorporating the results of the external sequence analyssi, identifying intron retentions, predicting functional consequences and lastly plotting all genes with isoform switches as well as summarizing general consequences of switching.
All of this can be done using the function:
isoformSwitchAnalysisPart2()
See below for usage example, and Detailed Workflow for details on the individual steps.
Alternatively if one does not plan to incorporate external sequence analysis, it is possible to run the full workflow using:
isoformSwitchAnalysisCombined()
Which correspond to running isoformSwitchAnalysisPart1() and isoformSwitchAnalysisPart2() without adding the external results.
Before you start on this section, we recommend you to read though the [Quick Start] section, as it gives the large overview, introduces some basic concepts and give a couple of tips not repeated here.
Compared to the workflow presented in [Quick Start], a full workflow for analyzing isoform switches naturally have a lot of sub-steps which each can be customized/optimized. In this section, we will go more into depth with each of the steps as well as provide tips and shortcuts for working with IsoformSwitchAnalyzeR. Specifically each of the main functions behind these steps will be described and illustrated. For a more comprehensive and detailed description of each individual function please refer to the individual function documentation (via ?functionName).
The first section contains IsoformSwitchAnalyzeR Background Information, followed by a detailed isoform switch analysis workflow and lastly an overview of useful Other Tools in IsoformSwitchAnalyzeR is provided.
Back to Table of Content
The detailed workflow consists of the following steps (illustrated in Figure 2) which, just like before, can be divided into two parts:
Part 1) Extract Isoform Switches and Their Sequences.
This corresponds to running the following functions in sequential order (which incidentally is just what the isoformSwitchAnalysisPart1() function does):
### Import data and create SwitchAnalyzeRlist
# either:
mySwitchList <- importCufflinksData()
# or:
mySwitchList <- importIsoformExpression()
mySwitchList <- importRdata()
### Run analysis
mySwitchList <- preFilter(mySwitchList)
mySwitchList <- isoformSwitchTestDRIMSeq(mySwitchList) # OR
mySwitchList <- isoformSwitchTest(mySwitchList)
mySwitchList <- analyzeORF(mySwitchList)
mySwitchList <- extractSequence(mySwitchList)
Part 2) Plot All Isoform Switches and Their annotation.
This corresponds to running the following functions in sequential order (just like the isoformSwitchAnalysisPart2() function does):
### Add annotation
mySwitchList <- analyzeCPAT(mySwitchList)
mySwitchList <- analyzePFAM(mySwitchList)
mySwitchList <- analyzeSignalP(mySwitchList)
mySwitchList <- analyzeIntronRetention(mySwitchList)
### Analyse consequences
mySwitchList <- analyzeSwitchConsequences(mySwitchList)
### visualize results
switchPlotTopSwitches(mySwitchList)
extractConsequenceSummary(mySwitchList)
The combined workflow is visualized here:
The switchAnalyzeRlist object is created to specifically contain and summarize all relevant information about the isoforms involved in isoform switching. The switchAnalyzeRlist is a named list, meaning each entry in the list can be accessed by its name via the ‘$’ symbol or by using “[[‘entryName’]]”. A newly created switchAnalyzeRlist object contains 6 entries, and as the isoforms are gradually annotated and analyzed more entries are added.
data("exampleSwitchList") # A newly created switchAnalyzeRlist + switch analysis
names(exampleSwitchList)
#> [1] "isoformFeatures" "exons" "conditions"
#> [4] "designMatrix" "isoformCountMatrix" "sourceId"
#> [7] "isoformSwitchAnalysis"
data("exampleSwitchListAnalyzed") # A fully analyzed switchAnalyzeRlist
names(exampleSwitchListAnalyzed)
#> [1] "isoformFeatures" "exons"
#> [3] "conditions" "designMatrix"
#> [5] "isoformCountMatrix" "sourceId"
#> [7] "isoformSwitchAnalysis" "orfAnalysis"
#> [9] "ntSequence" "aaSequence"
#> [11] "signalPeptideAnalysis" "domainAnalysis"
#> [13] "intronRetentionAnalysis" "switchConsequence"
The first entry ‘isoformFeatures’ is a data.frame where all relevant data about each comparison of an isoform (between conditions), as well as the analysis performed and annotation added via IsoformSwitchAnalyzeR is stored. Amongst the default information is isoform and gene id, gene and isoform expression as well as the isoform_switch_q_value and isoform_switch_q_value where the result of the differential isoform analysis is stored. The comparisons made can be identified from ‘condition_1’ to ‘condition_2’, meaning ‘condition_1’ is considered the ground state and ‘condition_2’ the changed state. This also means that a positive dIF value indicates the isoform usage is increased in ‘condition_2’ compared to ‘condition_1’. Since the ‘isoformFeatures’ entry is the most relevant part of the switchAnalyzeRlist object, the most used standard methods have also been implemented to works directly on isoformFeatures.
### Preview
head(exampleSwitchList, 2)
#> iso_ref gene_ref isoform_id gene_id
#> 19 isoComp_00000007 geneComp_00000005 TCONS_00000007 XLOC_000005
#> 22 isoComp_00000008 geneComp_00000005 TCONS_00000008 XLOC_000005
#> condition_1 condition_2 gene_name nearest_ref_id class_code
#> 19 hESC Fibroblasts <NA> uc009vjk.2 =
#> 22 hESC Fibroblasts <NA> uc001aau.2 =
#> TSS_group_id CDS_id length locus gene_status gene_value_1
#> 19 TSS2 NA 2750 chr1:322036-328580 OK 696.704
#> 22 TSS3 NA 4369 chr1:322036-328580 OK 696.704
#> gene_value_2 gene_stderr_1 gene_stderr_2 gene_log2_fold_change
#> 19 48.0566 3.592857 2.307488 -3.85774
#> 22 48.0566 3.592857 2.307488 -3.85774
#> gene_p_value gene_q_value gene_significant iso_status iso_value_1
#> 19 2.66665e-09 3.20379e-08 yes OK 358.383
#> 22 2.66665e-09 3.20379e-08 yes OK 338.308
#> iso_value_2 iso_stderr_1 iso_stderr_2 iso_log2_fold_change iso_p_value
#> 19 29.28480 2.091049 17.489950 -3.61328 4.83698e-05
#> 22 5.01291 1.322809 6.999295 -6.07655 2.64331e-03
#> iso_q_value iso_significant IF1 IF2 dIF
#> 19 0.000548629 yes 0.5143978 0.6093814 0.09498364
#> 22 0.015898000 yes 0.4855835 0.1043126 -0.38127092
#> isoform_switch_q_value gene_switch_q_value
#> 19 NA 1
#> 22 NA 1
identical(
head(exampleSwitchList), head(exampleSwitchList$isoformFeatures)
)
#> [1] TRUE
identical(
tail(exampleSwitchList), tail(exampleSwitchList$isoformFeatures)
)
#> [1] TRUE
### Dimentions
dim(exampleSwitchList$isoformFeatures)
#> [1] 259 36
nrow(exampleSwitchList)
#> [1] 259
ncol(exampleSwitchList)
#> [1] 36
dim(exampleSwitchList)
#> [1] 259 36
A very useful functionality implemented in IsoformSwitchAnalyzeR is subsetSwitchAnalyzeRlist() function, which allows for removal of isoforms and all their associated information across all entries in a switchAnalyzeRlist. The function subsets the switchAnalyzeRlist based on a vector of logicals matching the isoformFeatures entry of the list.
exampleSwitchListAnalyzed
#> This switchAnalyzeRlist list contains:
#> 490 isoforms from 120 genes
#> 3 comparison from 3 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 hESC vs Fibroblasts 148 79
#> 2 iPS vs Fibroblasts 135 69
#> 3 iPS vs hESC 135 76
#> 4 combined 264 120
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification, ORFs, ntSequence, aaSequence, Signal Peptides, Protein Domains, Intron Retentions, Switch Consequences, Coding Potential"
### Subset
subsetSwitchAnalyzeRlist(
exampleSwitchListAnalyzed,
exampleSwitchListAnalyzed$isoformFeatures$gene_name == 'ARHGEF19'
)
#> This switchAnalyzeRlist list contains:
#> 3 isoforms from 1 genes
#> 3 comparison from 3 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 hESC vs Fibroblasts 2 1
#> 2 iPS vs Fibroblasts 3 1
#> 3 iPS vs hESC 2 1
#> 4 combined 3 1
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification, ORFs, ntSequence, aaSequence, Signal Peptides, Protein Domains, Intron Retentions, Switch Consequences, Coding Potential"
Transcript structure information is stored in the exon entry of the switchAnalyzeRlist and contains the genomic coordinates for each exon in each isoform, as well as a column indicating which isoform it originates from. This information is stored as GenomicRanges (GRanges), which is very useful for overlapping genomic features and interacting with other Bioconductor packages.
head(exampleSwitchList$exons,2)
#> GRanges object with 2 ranges and 2 metadata columns:
#> seqnames ranges strand | isoform_id gene_id
#> <Rle> <IRanges> <Rle> | <character> <character>
#> [1] chr1 [322037, 322228] + | TCONS_00000007 XLOC_000005
#> [2] chr1 [324288, 324345] + | TCONS_00000007 XLOC_000005
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths
A full description of the initial switchAnalyzeRlist can be found via ?switchAnalyzeRlist and each additional entry added is described in the “value” part of the documentation for the specific function adding the entry.
Back to Table of Content
In this section we will give a brief overview of the functions implemented in IsoformSwitchAnalyzeR. With a few exceptions, these functions can be divided into 5 groups:
The exceptions to these categories are the following:
createSwitchAnalyzeRlist()
preFilter()
isoformSwitchTest()
isoformSwitchTestDRIMSeq()
Where:
As well as some addtional tools not essetial for the workflow (more about them later).
Note that we have tried to be systematic with function argument names for cutoff, meaning that if the argument name contains “cutoff” the value supplied is not included whereas if it contains “min” or “max” the given value is included.
Back to Table of Content
IsoformSwitchAnalyzeR can analyze the output from any tool that quantify (de-novo/guided deconvoluted) full-length isoforms. All it requires is the following 3 sets of data:
And furthermore it is highly recommended to also obtain the estimated replicate abundances and also use them in the construction of the switchAnalyzeRlist.
From these data, a minimum switchAnalyzeRlist (see The switchAnalyzeRlist for description) can be constructed. To facilitate the usage of IsoformSwitchAnalyzeR, several dedicated wrappers for constructing the switchAnalyzeRlist from different data sources are included in switchAnalyzeRlist. See the following sections for details about the specific import methods:
The data from Cufflinks/Cuffdiff is of special interest because Cuffdiff is amongst the most advanced tools for analyzing RNA-seq data with isoform resolution. Furthermore, Cuffdiff also supports identification of isoform switching by: a) Cuffdiff have a test of isoform switches amongst isoforms with shared promoter (isoforms from same pre-mRNA), b) Cuffdiff gives better estimates of the expression uncertainties (variance) than estimated from the raw data, which can be utilized with one of the isoform switch test (at isoform level) implemented (see Identifying Isoform Switches).
There are two ways of obtaining a switchAnalyzeRlist from cufflinks.
importCufflinksCummeRbund()
importCufflinksFiles()
Creating a switchAnalyzeRlist via cummeRbund is a two-step process. First use the cummeRbund function readCufflinks() to create the SQL backend (which in the following example is done by prepareCuffExample() on the cummeRbund example data).
cuffDB <- prepareCuffExample()
cuffDB
#> CuffSet instance with:
#> 3 samples
#> 400 genes
#> 1203 isoforms
#> 662 TSS
#> 906 CDS
#> 1062 promoters
#> 1986 splicing
#> 990 relCDS
This SQL backend is also extremely useful since it, via other functions implemented in cummeRbund, allows for a lot of (necessary!) quality control (QC) and exploratory data analysis (EDA). See section 5 in the cummeRbund manual for more information).
Now we have the sql connection, simply use importCufflinksCummeRbund() to create the switchAnalyzeRlist.
aSwitchList <- importCufflinksCummeRbund(cuffDB)
#> Reading cuffDB, isoforms...
#> Reading cuffDB, exons...
#> Analyzing cufflinks annotation problem...
#> Fixing cufflinks annotation problem...
#> Cufflinks annotation problem was fixed for 65 Cuff_genes
#> Extracting analysis of alternative splicing...
#> Making IsoformSwitchAanalyzeRlist...
#> Done
aSwitchList
#> This switchAnalyzeRlist list contains:
#> 1203 isoforms from 468 genes
#> 3 comparison from 3 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 hESC vs Fibroblasts 0 0
#> 2 iPS vs Fibroblasts 0 0
#> 3 iPS vs hESC 0 0
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification"
An alternative way of importing data from Cufflinks/Cuffdiff is to directly use the files outputted by Cuffdiff. This can be done with the importCufflinksFiles() function. In the example below we link to the example files in the cummeRbund package via the system.file() function. To analyze your own data simply supply the (full) path to the required files as a string (aka do not use the system.file() function).
testSwitchList <- importCufflinksFiles(
pathToGTF = system.file('extdata/chr1_snippet.gtf', package = "cummeRbund"),
pathToGeneDEanalysis = system.file('extdata/gene_exp.diff', package = "cummeRbund"),
pathToIsoformDEanalysis = system.file('extdata/isoform_exp.diff', package = "cummeRbund"),
pathToGeneFPKMtracking = system.file('extdata/genes.fpkm_tracking', package = "cummeRbund"),
pathToIsoformFPKMtracking = system.file('extdata/isoforms.fpkm_tracking', package = "cummeRbund"),
pathToIsoformReadGroupTracking = system.file('extdata/isoforms.read_group_tracking', package = "cummeRbund"),
pathToSplicingAnalysis = system.file('extdata/splicing.diff', package = "cummeRbund"),
pathToReadGroups = system.file('extdata/read_groups.info', package = "cummeRbund"),
pathToRunInfo = system.file('extdata/run.info', package = "cummeRbund"),
fixCufflinksAnnotationProblem=TRUE,
quiet=TRUE
)
testSwitchList
#> This switchAnalyzeRlist list contains:
#> 1203 isoforms from 468 genes
#> 3 comparison from 3 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 hESC vs Fibroblasts 0 0
#> 2 iPS vs Fibroblasts 0 0
#> 3 iPS vs hESC 0 0
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification"
Note that both cufflinks import functions per default:
The resulting switchAnalyzeRlist can then be used with the rest of the isoform switch analysis pipeline. The next step is typically Filtering.
While RSEM is a well-established and robust isoform quantification tool, a new generation of tools often refered to as pseudo/quasi-aligners, such as Salmon and Kallisto, are on the rise due to their speed and increased accuracy. We therefore naturally support analysis of quantifications made with all three tools by allowing for easy import of the isoform replicate count matrix via the importIsoformExpression() function.
This function is a wrapper for the tximport R package with extra functionalities. By using txtimport it allows the isoform counts to be estimated from the abundance estimates - which is recomended as it will incooperate the bias correction performed by the tools into the switch identification. The importIsoformExpression can furthermore perform a inter-library normalization of the abundance estimats which is nessesary for a fair comparison of libraries (for in deapth discussion of normalization and why it is needed see section 2.7 of the edgeR user guide).
Using importIsoformExpression results in the creation of a list which contains both a count and an abundance matrix can then easily be used to generate a switchAnalyzeRlist with the importRdata() function (as described in the next section).
As described above, all you need to create a switchAnalyzeRlist is:
Note that:
Once you have imported the isoform expression estimates (via importIsoformExpression() or otherwise) as well as the isoform structre and the isoform annotation into R an switchAnalyzeRlist can be constructed using the general-purpose wrapper:
importRdata()
All you need to supply to the importRdata() is the estimated replicate isoform count expression, optionally the replicate isoform abundance estimates, the isoform annotation and a design matrix indicating which samples are from the same condition. The importRdata() function will then:
Note that:
### Construct data needed from example data in cummeRbund package
### (The recomended way of analyzing Cufflinks/Cuffdiff datat is via importCufflinksCummeRbund()
### This is jus an easy way to get some example data ).
cuffDB <- prepareCuffExample()
isoRepCount <- repCountMatrix(isoforms(cuffDB))
isoRepCount$isoform_id <- rownames(isoRepCount)
### Design matrix
designMatrix <- cummeRbund::replicates(cuffDB)[,c('rep_name','sample_name')]
colnames(designMatrix) <- c('sampleID','condition')
localAnnotaion <- import(system.file("extdata/chr1_snippet.gtf", package="cummeRbund"))[,c('transcript_id','gene_id')]
colnames(localAnnotaion@elementMetadata)[1] <- 'isoform_id'
### Take a look at the data
head(isoRepCount, 2)
#> iPS_0 iPS_1 hESC_0 hESC_1 Fibroblasts_0
#> TCONS_00000001 0.0000 14.99310 0 0.000000 11.6729
#> TCONS_00000002 16.0613 8.06806 0 0.654207 0.0000
#> Fibroblasts_1 isoform_id
#> TCONS_00000001 12.1679 TCONS_00000001
#> TCONS_00000002 0.0000 TCONS_00000002
designMatrix
#> sampleID condition
#> 1 iPS_0 iPS
#> 2 iPS_1 iPS
#> 3 hESC_0 hESC
#> 4 hESC_1 hESC
#> 5 Fibroblasts_0 Fibroblasts
#> 6 Fibroblasts_1 Fibroblasts
head(localAnnotaion, 3)
#> GRanges object with 3 ranges and 2 metadata columns:
#> seqnames ranges strand | isoform_id gene_id
#> <Rle> <IRanges> <Rle> | <character> <character>
#> [1] chr1 [11874, 12227] + | TCONS_00000003 XLOC_000001
#> [2] chr1 [12646, 12697] + | TCONS_00000003 XLOC_000001
#> [3] chr1 [13221, 14409] + | TCONS_00000003 XLOC_000001
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths
### Create switchAnalyzeRlist
aSwitchList <- importRdata(
isoformCountMatrix=isoRepCount,
designMatrix=designMatrix,
isoformExonAnnoation=localAnnotaion,
showProgress=FALSE
)
#> Step 1 of 5: Obtaining annotation...
#> Step 2 of 5: Calculating gene expression...
#> Step 3 of 5: Merging gene and isoform expression...
#> Step 4 of 5: Making comparisons...
#> Step 5 of 5: Making switchAnalyzeRlist object...
#> Done
aSwitchList
#> This switchAnalyzeRlist list contains:
#> 1203 isoforms from 400 genes
#> 3 comparison from 3 conditions
This switchAnalyzeRlist can then be used with the rest of the pipeline. The next step in a typical workflow Filtering.
Alternatives to using importRdata() are either importGTF() or createSwitchAnalyzeRlist().
If the importGTF() function is used it will generate a switchAnalyzeRlist with dummy variables:
aSwitchList <- importGTF(pathToGTF = system.file("extdata/chr1_snippet.gtf", package = "cummeRbund"))
aSwitchList
#> This switchAnalyzeRlist list contains:
#> 1203 isoforms from 400 genes
#> 1 comparison from 2 conditions
head(aSwitchList,2)
#> iso_ref gene_ref isoform_id gene_id
#> 3 isoComp_00000001 geneComp_00000001 TCONS_00000001 XLOC_000001
#> 2 isoComp_00000002 geneComp_00000001 TCONS_00000002 XLOC_000001
#> condition_1 condition_2 gene_name class_code gene_value_1 gene_value_2
#> 3 plaseholder1 plaseholder2 <NA> NA 0 0
#> 2 plaseholder1 plaseholder2 <NA> NA 0 0
#> gene_stderr_1 gene_stderr_2 gene_log2_fold_change gene_p_value
#> 3 NA NA 0 1
#> 2 NA NA 0 1
#> gene_q_value iso_value_1 iso_value_2 iso_stderr_1 iso_stderr_2
#> 3 1 0 0 NA NA
#> 2 1 0 0 NA NA
#> iso_log2_fold_change iso_p_value iso_q_value IF1 IF2 dIF
#> 3 0 1 1 NA NA NA
#> 2 0 1 1 NA NA NA
#> isoform_switch_q_value gene_switch_q_value
#> 3 NA NA
#> 2 NA NA
head(aSwitchList$conditions,2)
#> condition nrReplicates
#> 1 plaseholder1 1
#> 2 plaseholder2 1
From above, it is observed that dummy variables have been inserted in both isoformFeatures and conditions entries of the switchAnalyzeRlist. This approach is well suited if you just want to annotate a transcriptome and are not interested in expression. If you are interested in expression estimates, it is probably easier to use importRdata.
Back to Table of Content
Once you have a switchAnalyzeRlist, there is a good chance that it contains a lot of genes/isoforms that are irrelevant for an analysis of isoform swhiches. Examples of such could be single isoform gene or non-expressed isoforms. These extra genes/isoforms will make the downstream analysis take (much) longer than necessary. Therfore we have implemented a pre-filtering step to remove thse features before continuing with the analysus. Furthermore, filtering can enhance the reliability of the downstream analysis as described below.
By using preFilter() it is possible to remove genes and isoforms from all aspects of the switchAnalyzeRlist via filtering on:
Removal of single isoform genes is the default setting in preFilter() since these genes per definition cannot have changes in isoform usage.
Filtering on isoforms expression allows removal of non-used isoforms that only appear in the switchAnalyzeRlist because they were in the isoform/gene annotation used. Furthermore, the expression filtering allows removal oflowly expressed isoforms where the expression levels might be untrustworthy. The filtering on gene expression allows for removal of lowly expressed genes which causes the calculations of the Isoform Fractions (IF) to become untrustworthy.
The filter on Isoform Fraction allows removal of isoforms that only contribute minimally to the gene expression thereby speeding up and simplifying the rest of the downstream analysis.
Filtering on unwanted isoform classes is only an option available if using data form Cufflinks/Cuffdiff since the Tuxedo workflow includes classification of transcript types. This allows for removal of for example transcripts classified as “Possible polymerase run-on fragment” or “Repeat”. See full description of Cufflinks/Cuffdiff transcript classification here
The preFilter() function is used as follows:
data("exampleSwitchList")
exampleSwitchList
#> This switchAnalyzeRlist list contains:
#> 259 isoforms from 84 genes
#> 1 comparison from 2 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 hESC vs Fibroblasts 0 0
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification"
exampleSwitchListFiltered <- preFilter(exampleSwitchList, geneExpressionCutoff = 1, isoformExpressionCutoff = 0, removeSingleIsoformGenes = TRUE)
#> The fitering removed 53 ( 20.46% of ) transcripts. There is now 206 isoforms left
exampleSwitchListFilteredStrict <- preFilter(exampleSwitchList, geneExpressionCutoff = 10, isoformExpressionCutoff = 3, removeSingleIsoformGenes = TRUE)
#> The fitering removed 64 ( 24.71% of ) transcripts. There is now 195 isoforms left
Back to Table of Content
As identification of isoform switches are essential for IsoformSwitchAnalyzeR, multipe ways of identifying isoform switches are supported. Currently IsoformSwitchAnalyzeR directly supports three different approachs:
However results of other tools for which can test for differential isoform usage at gene and/or isoform level can also be used. See Testing Isoform Switches with other Tools.
No matter which method is used, all the downstream functionality in IsoformSwitchAnalyzeR (e.g. not only the test) uses tow paramters to defines a significant isoform switch:
This combination is used since a q-value is only a measure of the statistical certainty of the difference between two groups and thereby not reflects the effect size which however is measured by the dIF values.
IsoformSwitchAnalyzeR supports tests for differential isoform usage performed both at isoform and gene resolution. Note however that for gene level analysis you reliy on cutoffs on dIF values for identifying the isoforms involved in a switch - something that, when compared to the isoform level analysis, could give false positives. We therefore recommend the use of isoform level analysis whenever possible.
A recent breakthrough in the analysis of isoform usage was when the DRIMSeq package (By Nowicka et al, see What To Cite - please remember to cite the tool) was updated to not only analyze genes for differential isoform usage but provide a statistical test for each isoform analyzed thereby allowing for identification of the excat isoforms involved in a switch.
IsoformSwitchAnalyzeR supports identification of isoform switches with DRIMSeq via the wrapper isoformSwitchTestDRIMSeq() which performs the full analysis of all isoforms and comparisons in switchAnalyzeRlist. Afterwards, the results are integrated back into the switchAnalyzeRlist which is returned to the user as usual. Note that we support both the use of the gene-level and isoform-level analysis perfromed by DRIMSeq as controlled by the ‘testIntegration’ argument. Using this argument, it is also possible to be even more stringent by requiring that both the gene- and isoform-level analysis must be significant.
The isoformSwitchTestDRIMSeq() function is used as follows:
# Load example data and prefilter
data("exampleSwitchList")
exampleSwitchList <- preFilter(exampleSwitchList) # preFilter for fast runtime
#> The fitering removed 53 ( 20.46% of ) transcripts. There is now 206 isoforms left
# Perfom test (takes ~10 sec)
exampleSwitchListAnalyzed <- isoformSwitchTestDRIMSeq(
switchAnalyzeRlist = exampleSwitchList,
testIntegration='isoform_only',
reduceToSwitchingGenes=TRUE
)
#> Step 1 of 6: Creating DM data object...
#> Step 2 of 6: Filtering DM data object...
#> Step 3 of 6: Estimating precision paramter (this may take a while)...
#> Step 4 of 6: Fitting linear models (this may take a while)...
#> Step 5 of 6: Testing pairwise comparison(s)...
#> Step 6 of 6: Preparing output...
#> Result added switchAnalyzeRlist
#> An isoform switch analysis was performed for 54 gene comparisons (90%).
#> Done
# Summarize swiching geatures
extractSwitchSummary(exampleSwitchListAnalyzed)
#> Comparison nrIsoforms nrGenes
#> 1 hESC vs Fibroblasts 85 41
An important argument in isoformSwitchTestDRIMSeq() is the ‘reduceToSwitchingGenes’. This is simply an option to reduce the switchAnalyzeRlist to only genes that contains at least one isoform passing the supplied alpha and dIFcutoff cutoffs. This option ensures the rest of the workflow goes significantly faster since isoforms from genes without isoform switching is not analyzed.
Please note:
The DRIMSeq approach utilizes Bayesian information sharing across genes/isoforms. We therfore recommend most users to use the DRIMSeq approach as it is much more sensitive when having a small number of samples. The exception to this recommendation is if a large number of samples (per condition) are analyzed. For such cases we recommend the statistical test implemented in IsoformSwitchAnalyzeR as it fast and the gain of the information sharing in DRIMSeq decreases as the number of samples increases.
IsoformSwitchAnalyzeR contains an implementation of a test for differential isoform usage (in the function isoformSwitchTest()) described in Vitting-Seerup et al (2017) (see What To Cite). Here, we will introduce two sections describing the test for differential isoform usage implemented in isoformSwitchTest(). First, the idea behind the test will be introduced and afterwards we will show how to use it (jump to Usage of The Statistical Test).
One of the better ways of descripting isoform usage is as Isoform Fractions (IF). These are simply values describing how large a fraction of the expression of a given gene originate from a given isoform. So the Isoform Fraction of isoform x (\(IF_x\)) is given by: \(IF_x = i_x / g\) where \(i_x\) is the expression of isoform x and g is the parent gene expression. Note that the parent gene expression is the sum of the expression of all the associated isoforms. The difference between two conditions, condition 1 and 2, is then given by \(dIF= IF_2 - IF_1\).
The Idea Behind The Statistical Test implemented in IsoformSwitchAnalyzeR is:
More information can be found in the documentation page of isoformSwitchTest().
The test described above in The Idea Behind The Statistical Test section is implemented in IsoformSwitchAnalyzeR as the function:
isoformSwitchTest()
### Show arguments of function
args(isoformSwitchTest)
#> function (switchAnalyzeRlist, alpha = 0.05, dIFcutoff = 0.1,
#> reduceToSwitchingGenes = TRUE, calibratePvalues = TRUE, showProgress = FALSE,
#> quiet = FALSE)
#> NULL
If a single isoform from a given gene is termed significant, we define this gene as having an isoform switch, since if there is a change in the isoform usage of one isoform there must be changes in another (or multiple) isoforms compensating for the identified change. To extrapolate this to gene level, we define the q-value at the gene level (in the ‘gene_switch_q_value’ column) as equivalent to the smallest q-value of the associated isoforms (in the ’isoform_switch_q_value column).
Another important argument in isoformSwitchTestDRIMSeq() is the ‘reduceToSwitchingGenes’. This is simply an option to reduce the switchAnalyzeRlist to only genes that contains at least one isoform passing the given alpha and dIFcutoff cutoffs. This option ensures that rest of the workflow runs significantly faster, since isoforms from genes without isoform switching are not analyzed.
Note that if a switchAnalyzeRlist already contains the result of a switch test this will be completely overwritten by isoformSwitchTestDRIMSeq().
The DRIMSeq approach utilizes Bayesian information sharing across genes/isoforms. We therfore recommend most users to use the DRIMSeq approach as it is much more sensitive when having a small number of samples. The exception to this recommendation is if a large number of samples (per condition) are analyzed. For such cases we recommend the statistical test implemented in IsoformSwitchAnalyzeR as it fast and the gain of the information sharing in DRIMSeq decreases as the number of samples increases.
# Load example data and prefilter
data("exampleSwitchList")
exampleSwitchList <- preFilter(exampleSwitchList) # preFilter for fast runtime
#> The fitering removed 53 ( 20.46% of ) transcripts. There is now 206 isoforms left
# Perfom test
exampleSwitchListAnalyzed <- isoformSwitchTest(exampleSwitchList)
#> Step 1 of 3: Filtering for eligible data...
#> Found 76 ( 36.89 %) of isoform comparisons eligable for switch analysis
#> Step 2 of 3: Testing isoform usage...
#> Step 3 of 3: Preparing output...
# Summarize swiching geatures
extractSwitchSummary(exampleSwitchListAnalyzed)
#> Comparison nrIsoforms nrGenes
#> 1 hESC vs Fibroblasts 1 1
As suggested by the Ferguson et al (2014) the p-value calibration is not always performed (even though calibratePvalues=TRUE, see details in ?isoformSwitchTest). It is therefore useful to check in which comparison the p-value calibration was performed. To this end extractCalibrationStatus() have been implemented.
extractCalibrationStatus(exampleSwitchListAnalyzed)
#> Comparison pvalsAreCalibrated
#> 1 hESC vs Fibroblasts TRUE
From which we see that we should be careful with comparing the two conditions since they have been treated differently. To fix this problem we can simply turn the calibration off:
# Perfom test
exampleSwitchListAnalyzed <- isoformSwitchTest(exampleSwitchList, calibratePvalues = FALSE)
#> Step 1 of 3: Filtering for eligible data...
#> Found 76 ( 36.89 %) of isoform comparisons eligable for switch analysis
#> Step 2 of 3: Testing isoform usage...
#> Step 3 of 3: Preparing output...
# Summarize swiching geatures
extractSwitchSummary(exampleSwitchListAnalyzed)
#> Comparison nrIsoforms nrGenes
#> 1 hESC vs Fibroblasts 1 1
# check callibration status
extractCalibrationStatus(exampleSwitchListAnalyzed)
#> Comparison pvalsAreCalibrated
#> 1 hESC vs Fibroblasts FALSE
Back to Table of Content
IsoformSwichtAnalyzeR also supports the analysis of isoform switches found via other tools. All you need to do is to generate a switchAnalyzeRlist with the corresponding data via one of the import options (see Importing Data Into R) and then fill in either the isoform_switch_q_value or gene_switch_q_value columns in the isoformFeatures entry of the switchAnalyzeRlist with the multiple testing corrected p-values (q-values) from the testing tool.
If the external tool used has isoform resolution (one test per isoform), the q-values should be added to the isoform_switch_q_value column and the smallest q-value of a given gene (in a specific comparison) should be added to all the isoforms for that gene in the gene_switch_q_value column. If the external tool has lower resolution (pre-mRNA or gene level resolution) the q-values should only be added to the gene_switch_q_value column (and the isoform_switch_q_value should be set to NA).
In this way, IsoformSwitchAnalyzeR can also support non-isoform resolution testing - note however that:
Once the isoform switches have been found, the next step is to annotate the isoforms involved in the isoform switches. The first step in such annotation is to obtain Open Reading Frames (ORFs).
If known annotated isoforms were quantified (aka you did not perform a (guided) de novo isoform reconstruction), you should use the annotated CDS (Coding Sequence) form the annotation - which can be imported and integrated though one of the implemented methods, see Importing Data Into R.
If you performed (guided) de-novo isoform reconstruction (isoform deconvolution), prediction of ORFs can be done with high accuracy from the transcript sequence alone (see supplementary data in Vitting-Seerup et al 2017). Such prediction can be don with:
analyzeORF()
This function utilizes that we know the genomic coordinates of each transcript to extract the transcript nucleotide sequence from a reference genome (supplied via the genomeObject argument). In analyzeORF four different methods for predicting the ORF are implemented, suitable for different purposes and circumstances. The four methods are:
One important argument is the minORFlength, which will ensure that only ORFs longer than minORFlength nucleotides are annotated. Besides predicting the ORF, information about the most likely stop codon also allows for prediction of Pre-mature Termination Codon (PTC) and thereby sensitivity to degradation via the Nonsense Mediated Decay (NMD) machinery. This analysis is also implemented in the analyzeORF() function and is controlled by the PTCDistance argument.
Note that the ORF prediction can be integrated with both the CPAT results (via the removeNoncodinORFs paramter, see ?analyzeCPAT) as well as Pfam results (see Augmenting ORF Predictions with Pfam Results)
The analyzeORF() function can be used as follows:
### This example relies on the example data from the 'Usage of The Statistical Test' section above
### analyzeORF needs the genomic sequence to predict ORFs.
# These are readily advailable from Biocindoctor as BSgenome orbjects:
# http://bioconductor.org/packages/release/BiocViews.html#___BSgenome
# Here we use Hg19 - which can be download by copy/pasting the following two lines into the R termminal:
# source("https://bioconductor.org/biocLite.R")
# biocLite("BSgenome.Hsapiens.UCSC.hg19")
library(BSgenome.Hsapiens.UCSC.hg19)
exampleSwitchListAnalyzed <- analyzeORF(exampleSwitchListAnalyzed, genomeObject = Hsapiens, showProgress=FALSE)
#> Step 1 of 3 : Extracting transcript sequences...
#> Step 2 of 3 : Locating ORFs of interest...
#> Step 3 of 3 : Scanning for PTCs...
#> 3 putative ORFs were identified and analyzed
#> Done
head(exampleSwitchListAnalyzed$orfAnalysis, 3)
#> isoform_id orfTransciptStart orfTransciptEnd orfTransciptLength
#> 1 TCONS_00003919 381 641 261
#> 2 TCONS_00003920 139 660 522
#> 3 TCONS_00003921 139 1335 1197
#> orfStarExon orfEndExon orfStartGenomic orfEndGenomic
#> 1 3 3 1354877 1354617
#> 2 2 3 1356385 1354617
#> 3 2 4 1356385 1354483
#> stopDistanceToLastJunction stopIndex PTC
#> 1 -313 0 FALSE
#> 2 -313 0 FALSE
#> 3 -447 0 FALSE
As seen above, including genomic and transcript coordinates of ORF start and stop as well as PTC status, is added to the ‘orfAnalysis’ entry of the switchAnalyzeRlist.
If the user wants to use the longestAnnotated or mostUpstreamAnnoated methods, the analyzeORF() function requires known CDS to be supplied as described above. The CDS must be stored as a CDSSeq (see ?CDSSet) and a wrapper for downloading the CDS of the most frequently used datasets from UCSC genome browser is available via:
getCDS()
Now that we know the ORF of a transcript, we can obtain the corresponding protein amino acid sequence of the ORF simply by translating the nucleotide sequence of the ORF into amino acids. This opens the possibility for also performing both internal and external sequence analysis to annotate the isoforms involved in isoform switches further. To facilitate this we have implemented:
extractSequence()
Which allows for the extraction of both nucleotide and amino acid sequences from the switchAnalyzeRlist. To facilitate external sequence analysis extractSequence can output fasta files (one per sequence type) and to facilitate internal sequence analysis the sequences can be added to the switchAnalyzeRlist. An example of the internal sequence analysis is a pairwise comparison of the ORFs in two switching isoforms (see Predicting Switch Consequences below).
### This example relies on the example data from the 'Analyzing Open Reading Frames' section above
exampleSwitchListAnalyzed <- extractSequence(
exampleSwitchListAnalyzed,
genomeObject = Hsapiens,
pathToOutput = '<insert_path>',
writeToFile=FALSE # to avoid output when running this example data
)
#> Step 1 of 3 : Extracting transcript nucleotide sequences...
#> Step 2 of 3 : Extracting ORF AA sequences...
#> Step 3 of 3 : Preparing output...
#> Done
head(exampleSwitchListAnalyzed$ntSequence,2)
#> A DNAStringSet instance of length 2
#> width seq names
#> [1] 1456 CTTTTCAGCAGCAGACACTCC...AATAAAGCCCTTCCTTCTAC TCONS_00003919
#> [2] 1475 CTTTTCAGCAGCAGACACTCC...AATAAAGCCCTTCCTTCTAC TCONS_00003920
head(exampleSwitchListAnalyzed$aaSequence,2)
#> A AAStringSet instance of length 2
#> width seq names
#> [1] 87 MAALRCTGLPPEDTCLPSSCW...ERPPYTWLQSEGMGLPWGFC TCONS_00003919
#> [2] 174 MDSQRPEPREEEEEEQELRWM...ERPPYTWLQSEGMGLPWGFC TCONS_00003920
Back to Table of Content
The two fasta files outputted by extractSequence() (if writeToFile=TRUE) can be used as input to amongst others:
These three tools are the once currently supported but if you have additional ideas please do not hesitate to contact us as described in How To Get Help.
We suggest the external sequence analysis tools are run as follows:
Back to Table of Content
After the external sequence analysis with CPAT, Pfam, and SignalP have been performed (please remember to cite as describe in What To Cite), the results can be extracted and added to the switchAnalyzeRlist via respectively
analyzeCPAT()
analyzePFAM()
analyzeSignalP()
The functions are used as:
### Load test data (maching the external sequence analysis results)
data("exampleSwitchListIntermediary")
exampleSwitchListIntermediary
#> This switchAnalyzeRlist list contains:
#> 162 isoforms from 40 genes
#> 1 comparison from 2 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 hESC vs Fibroblasts 83 40
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification, ORFs, ntSequence, aaSequence"
### Add PFAM analysis
exampleSwitchListAnalyzed <- analyzePFAM(
switchAnalyzeRlist = exampleSwitchListIntermediary,
pathToPFAMresultFile = system.file("extdata/pfam_results.txt", package = "IsoformSwitchAnalyzeR"),
showProgress=FALSE
)
#> Converting AA coordinats to transcript and genomic coordinats...
#> Added domain information to 140 (86.42%) transcripts
### Add SignalP analysis
exampleSwitchListAnalyzed <- analyzeSignalP(
switchAnalyzeRlist = exampleSwitchListAnalyzed,
pathToSignalPresultFile = system.file("extdata/signalP_results.txt", package = "IsoformSwitchAnalyzeR")
)
#> Added signal peptide information to 17 (10.49%) transcripts
### Add CPAT analysis
exampleSwitchListAnalyzed <- analyzeCPAT(
switchAnalyzeRlist = exampleSwitchListAnalyzed,
pathToCPATresultFile = system.file("extdata/cpat_results.txt", package = "IsoformSwitchAnalyzeR"),
codingCutoff = 0.725, # the coding potential cutoff we suggested for human
removeNoncodinORFs = TRUE # because ORF was predicted de novo
)
#> Added coding potential to 162 (100%) transcripts
exampleSwitchListAnalyzed
#> This switchAnalyzeRlist list contains:
#> 162 isoforms from 40 genes
#> 1 comparison from 2 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 hESC vs Fibroblasts 83 40
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification, ORFs, ntSequence, aaSequence, Protein Domains, Signal Peptides, Coding Potential"
Where the pathTo<tool_name>resultFile points to the result files constructed as suggested in Advise for Running External Sequence Analysis Tools and in the detailed documentation of each of the individual functions.
Of particular interest is the removeNoncodinORFs argument in analyzeCPAT() since it allows to integrate the CPAT and ORF analysis by removing ORFs from isoforms not predicted to be coding by CPAT. This can be particular useful if isoforms and ORFs have been predicted de novo (both guided or non-guided). Note that if enabled (by setting to TRUE) it will affect all downstream analysis and plots as both analysis of domains and signal peptides requires that ORFs are annotated (e.g. analyzeSwitchConsequences will for example not consider the domains (potentially) found by Pfam if the ORF have been removed).
Back to Table of Content
Another annotation we easily can obtain, since we know the exon structure of all isoforms in a given gene (with isoform switching), is intron retentions. This can be done via the spliceR R package via the wrapper implemented in analyzeIntronRetention().
### This example relies on the example data from the 'Importing External Sequences Analysis' section above
exampleSwitchListAnalyzed <- analyzeIntronRetention(exampleSwitchListAnalyzed, quiet=TRUE)
### overview of number of intron retentions (IR)
table(exampleSwitchListAnalyzed$isoformFeatures$IR)
#>
#> 0 1 2
#> 135 25 2
Meaning 25 isoforms each with one intron retention (IR) and 2 isoforms with 2 IRs was identified.
Back to Table of Content
If an isoform has a significant change in its contribution to the gene expression, there must per definition be a reciprocal changes in one (or more) isoforms in the opposite direction compensating for the initial change. We utilize this by extracting the isoforms that are significantly differentially used and compare them to the isoforms that are compensating. By comparing all the information gathered during the workflow described above for the isoform used more (positive dIF) to the isoforms used less (negative dIF) we can identify differences in the functional annotation of the isoforms and thereby identify potential function consequences of the isoform switch.
Specifically, IsoformSwitchAnalyzeR contains a function analyzeSwitchConsequences() which extract the isoforms with significant changes in their isoform usage (defined by the alpha and dIFcutoff parameters, see Identifying Isoform Switches for details) as well as the isoforms, with a lage opposite change in isoform usage (also controlled via the dIFcutoff parameters), that compensate for the changes. Note that if an isoform level test was not used the gene is require to be significant (defined by the alpha paramter) but isoforms are then selected purely based on their changes in dIF values.
These isoforms are then divided into the isoforms that increases their contribution to gene expression (positive dIF values larger than dIFcutoff) and the isoforms that decrease their contribution (negative dIF values smaller than -dIFcutoff). The isoforms with increased contribution are then (in a pairwise manner) compared to the isoform with decreasing contribution. In each of these comparisons the isoforms compared are analyzed for differences in their annotation (controlled by the consequencesToAnalyze parameter). Currently 22 different features of the isoforms can be compared, which include features such as intron retention, coding potential, NMD status, protein domains and the sequence similarity of the amino acid sequence of the annotated ORFs. For a full list see under details at the manual page for analyzeSwitchConsequences() (?analyzeSwitchConsequences).
A more strict analysis can be performed by enabling the onlySigIsoforms argument, which causes analyzeSwitchConsequences() to only consider significant isoforms (defined by the alpha and dIFcutoff parameters) meaning the compensatory changes in isoform usage are ignored unless they themselves are significant.
The analysis of consequences can be performed as follows:
### This example relies on the example data from the 'Predicting Intron Retentions' section above
# the consequences highlighted in the text above
consequencesOfInterest <- c('intron_retention','coding_potential','NMD_status','domains_identified','ORF_seq_similarity')
exampleSwitchListAnalyzed <- analyzeSwitchConsequences(
exampleSwitchListAnalyzed,
consequencesToAnalyze = consequencesOfInterest,
dIFcutoff = 0.4, # very high cutoff for fast runtimes
showProgress=FALSE
)
#> Step 1 of 4: Extracting genes with isoform switches...
#> Step 2 of 4: Analyzing 5 pairwise isoforms comparisons...
#> Step 3 of 4: Massaging isoforms comparisons results...
#> Step 4 of 4: Preparing output...
#> Identified genes with containing isoforms switching with functional consequences...
extractSwitchSummary(exampleSwitchListAnalyzed, dIFcutoff = 0.4, filterForConsequences = FALSE)
#> Comparison nrIsoforms nrGenes
#> 1 hESC vs Fibroblasts 12 7
extractSwitchSummary(exampleSwitchListAnalyzed, dIFcutoff = 0.4, filterForConsequences = TRUE)
#> Comparison nrIsoforms nrGenes
#> 1 hESC vs Fibroblasts 10 5
Meaning that ~5/7 genes with switches (with dIF > 0.4) have isoform switches with functional consequences. For a larger dataset, it looks more like:
### Load already analyzed data
data("exampleSwitchListAnalyzed")
extractSwitchSummary(exampleSwitchListAnalyzed, filterForConsequences = FALSE)
#> Comparison nrIsoforms nrGenes
#> 1 hESC vs Fibroblasts 148 79
#> 2 iPS vs Fibroblasts 135 69
#> 3 iPS vs hESC 135 76
#> 4 combined 264 120
extractSwitchSummary(exampleSwitchListAnalyzed, filterForConsequences = TRUE)
#> Comparison nrIsoforms nrGenes
#> 1 hESC vs Fibroblasts 117 58
#> 2 iPS vs Fibroblasts 107 53
#> 3 iPS vs hESC 109 58
#> 4 combined 215 97
Where the combined row indicates the number of unique features across all comparisons.
Please note that the analyzeSwitchConsequences() function contains a lot of parameters and cutoffs for deciding when a difference in the annoatation is sufficiently different (e.g. supplying ‘ntCutoff=50’ ensures that differences in nucleotide lengths of two ORF are larger than 50 nucleotides (aka not just 3 nucleotide)).
Back to Table of Content
Now that we have performed a genome-wide analysis of isoform switches with potential consequences, there are two types of post-analysis we can perform:
IsoformSwitchAnalyzeR contains tools for both, as illustrated in each of the following sections - starting with the analysis of individual genes.
When analyzing the individual genes with isoform switches, the genes/isoforms with the largest changes in isoform usage (aka “most” switching genes/isoforms) are of particular interest. IsoformSwitchAnalyzeR can help you obtain these, either by sorting for the smallest q-values (getting the most significant genes) or the largest absolute dIF values (getting the largest effect sizes (switches) that are still significant). Both methods are implemented for both genes and isoforms in the extractTopSwitches() function and are controlled via the sortByQvals argument.
### Load already analyzed data
data("exampleSwitchListAnalyzed")
### Lets reduce the switchAnalyzeRlist to only one condition
exampleSwitchListAnalyzedSubset <- subsetSwitchAnalyzeRlist(
exampleSwitchListAnalyzed,
exampleSwitchListAnalyzed$isoformFeatures$condition_1 == 'iPS' &
exampleSwitchListAnalyzed$isoformFeatures$condition_2 == 'hESC'
)
exampleSwitchListAnalyzedSubset
#> This switchAnalyzeRlist list contains:
#> 321 isoforms from 76 genes
#> 1 comparison from 2 conditions
#>
#> Switching features:
#> Comparison switchingIsoforms switchingGenes
#> 1 iPS vs hESC 135 76
#>
#> Feature analyzed:
#> [1] "Isoform Swich Identification, ORFs, ntSequence, aaSequence, Signal Peptides, Protein Domains, Intron Retentions, Switch Consequences, Coding Potential"
### Extract top 2 switching genes (by dIF values)
extractTopSwitches(
exampleSwitchListAnalyzedSubset,
filterForConsequences = TRUE,
n = 2,
sortByQvals = FALSE
)
#> gene_ref gene_id gene_name condition_1 condition_2
#> 1 geneComp_00001060 XLOC_000108 TNFRSF1B iPS hESC
#> 2 geneComp_00000963 XLOC_000027 SCNN1D iPS hESC
#> gene_switch_q_value combinedDIF switchConsequencesGene
#> 1 3.942358e-19 1.999947 TRUE
#> 2 1.189401e-16 1.587289 TRUE
### Extract top 2 switching genes (by q-value)
extractTopSwitches(
exampleSwitchListAnalyzedSubset,
filterForConsequences = TRUE,
n = 2,
sortByQvals = TRUE
)
#> gene_ref gene_id gene_name condition_1 condition_2
#> 1 geneComp_00001378 XLOC_001391 USP48 iPS hESC
#> 2 geneComp_00001109 XLOC_000150:CROCC CROCC iPS hESC
#> gene_switch_q_value switchConsequencesGene
#> 1 2.225470e-27 TRUE
#> 2 1.058416e-26 TRUE
Let us take a look at the switching isoforms in the TBC1D22B gene which plays a role in protecting from apoptosis:
### Extract data.frame with all switching isoforms
switchingIso <- extractTopSwitches(
exampleSwitchListAnalyzedSubset,
filterForConsequences = TRUE,
n = NA, # n=NA: all features are returned
extractGenes = FALSE, # when FALSE isoforms are returned
sortByQvals = TRUE
)
subset(switchingIso, gene_name == 'TNFRSF1B')
#> iso_ref gene_ref isoform_id gene_id gene_name
#> 6 isoComp_00002713 geneComp_00001060 TCONS_00000307 XLOC_000108 TNFRSF1B
#> 20 isoComp_00002711 geneComp_00001060 TCONS_00000305 XLOC_000108 TNFRSF1B
#> 67 isoComp_00002712 geneComp_00001060 TCONS_00000306 XLOC_000108 TNFRSF1B
#> condition_1 condition_2 iso_significant IF1 IF2 dIF
#> 6 iPS hESC no 1 0.000 -1.000
#> 20 iPS hESC yes 0 0.832 0.832
#> 67 iPS hESC no 0 0.168 0.168
#> isoform_switch_q_value switchConsequencesGene comparison
#> 6 3.942358e-19 TRUE iPS_vs_hESC
#> 20 1.003206e-09 TRUE iPS_vs_hESC
#> 67 2.125305e-03 TRUE iPS_vs_hESC
The isoform switch can also be visually analyzed by the Isoform Switch Analysis Plot. Note that since there is only one comparison in the switchAnalyzeRlist (after the subset), it is not necessary to specify the conditions.
switchPlot(exampleSwitchListAnalyzedSubset, gene = 'TNFRSF1B')
From this plot, it is clear to see that in the iPS cells it is mainly the isoform ‘TCONS_00000307’ that is used. The isoform is predicted to be NMD sensitive (due to the inclusion of a premature termination codon (PTC) in exon 7). The identified isoform switch does however result in a switch away from the NMD sensitive isoform whereby only the protein coding isoforms are used in normal human Embryonic Stemm Cells (hESC) and not in iPS cells. Note that no differential expression was identified by Cufflinks/Cuffdiff meaning that a gene-level analysis might have missed this change.
One advantage of the Isoform Switch Analysis Plot is that it contains all information needed to judge the potential impact of an isoform switch. This also means a systematic screening of the top isoform switches can be made by generating this plot for each of those genes.
To facilitate such screening we have implemented switchPlotTopSwitches() which will extract the top n genes (controlled by the n argument) with significant switches (as defined by alpha and dIFcutoff) and output a pdf or png version of the corresponding Isoform Switch Analysis Plot to the directory given by the outputDestination argument. The function furthermore automatically ranks (by p-value or switch size) the switches and supports to either filter for isoform switches with predicted functional consequences or to output those with and those without consequences to separate folders.
switchPlotTopSwitches(
switchAnalyzeRlist = exampleSwitchListAnalyzed,
n= 10,
filterForConsequences = FALSE,
splitFunctionalConsequences = TRUE
)
If you want to output a plot for all genes simply set “n=NA”.
Back to Table of Content
A genome wide analysis is both useful for getting an overview of the extent of isoform switching as well as discovering larger patterns.
The extractConsequenceSummary() function can globally summarize isoform switches with consequences as shown here:
### Load the large example dataset
data("exampleSwitchListAnalyzed")
### Extract summary
consequenceSummary <- extractConsequenceSummary(
exampleSwitchListAnalyzed,
returnResult = TRUE,
plotGenes = TRUE
)
subset(consequenceSummary, featureCompared=='Intron retention')
#> Comparison featureCompared switchConsequence
#> 12 hESC vs Fibroblasts Intron retention Intron retention gain
#> 13 iPS vs Fibroblasts Intron retention Intron retention gain
#> 14 iPS vs hESC Intron retention Intron retention gain
#> 15 hESC vs Fibroblasts Intron retention Intron retention loss
#> 16 iPS vs Fibroblasts Intron retention Intron retention loss
#> 17 iPS vs hESC Intron retention Intron retention loss
#> 18 hESC vs Fibroblasts Intron retention Intron retention switch
#> nrGenesWithConsequences nrIsoWithConsequences
#> 12 4 8
#> 13 2 5
#> 14 5 8
#> 15 13 25
#> 16 7 12
#> 17 11 21
#> 18 1 1
Illustrating that both a tabular and a visual summary can be obtained thereby providing a general overview of the isoform switches.
From this summary plot, we can see for large fraction of isoform switches where the upregulated isoforms more frequently have intron retentions.
Note furthermore that the focus of the plot generated is on the number of genes with isoform switches with functional consequences (controlled by plotGenes = TRUE and asFractionTotal=FALSE), whereas in Short Example Workflow the y-axis was the fraction of isoforms involved in isoform switching with functional consequences.
If the difference between conditions are substantial such effects could indicate a more systematic change on genome wide scale. This can also be annalyzed with IsoformSwitchAnalyzeR by simultaniously analyzing all isoforms with a specific feature (fx intron retention) for changes in isoform usage:
symmaryStatistics <- extractGenomeWideAnalysis(
switchAnalyzeRlist = exampleSwitchListAnalyzed,
featureToExtract = 'isoformUsage', # default - alternatives are 'isoformExp', 'geneExp' and 'all'
plot=TRUE,
returnResult = TRUE
)
(Note that the 3 dots in each violin plot correspond to the 25th, 50th (median) and 75th percentiles).
Which shows that, although many isoform switches are identified, no global changes are identified (in this toy data).
Back to Table of Content
Apart from the analysis presented above, IsoformSwitchAnalyzeR also contains a couple of other tools which will be presented in this section
The central visualization is the Isoform Switch Analysis Plot, created with switchPlot() as shown above. This plot is made from 4 subplots, which can be created individually using respectively:
switchPlotTranscript() # Visualizes the transcripts and their annotation
switchPlotGeneExp() # Visualizes the gene exression
switchPlotIsoExp() # Visualizes the isoform exression
switchPlotIsoUsage() # Visualizes the isoform usage
As illustrated here:
### Load already analyzed data
data("exampleSwitchListAnalyzed")
switchPlotTranscript(exampleSwitchListAnalyzedSubset, gene = 'TNFRSF1B')
switchPlotGeneExp (exampleSwitchListAnalyzedSubset, gene = 'TNFRSF1B', condition1='iPS', condition2='hESC')
switchPlotIsoExp(exampleSwitchListAnalyzedSubset, gene = 'TNFRSF1B', condition1='iPS', condition2='hESC')
switchPlotIsoUsage(exampleSwitchListAnalyzedSubset, gene = 'TNFRSF1B', condition1='iPS', condition2='hESC')
Another usefull tool implemented is the isoformToGeneExp() function which takes a data.frame of isoform expression and summarizes it up to gene level expression values. This function is quite usefull in conjunction with the importIsoformExpression() function since the importIsoformExpression() function imports isoform-level data - which can then be summarized to gene level data with isoformToGeneExp(). This is quite analogous to what tximport can do - except it allows for incooperation of the additional features, such as inter-library normalization, implemented in importIsoformExpression().
The last tool currently build into IsoformSwitchAnalyzeR is the extractExpressionMatrix() function. The expression information stored in the switchAnalyzeRlist’s isoformFeatures is ideal for comparing multiple annotation in a specific comparison of two conditions, but is not well suited for the comparison of multiple conditions. The extractExpressionMatrix() function solves this by converting the average expression information (gene expression, isoform expression or Isoform Fraction values as controlled via the feature argument) into a matrix format as illustrated here:
data("exampleSwitchListIntermediary")
ifMatrix <- extractExpressionMatrix(exampleSwitchListAnalyzed)
head(ifMatrix)
#> hESC iPS Fibroblasts
#> TCONS_00000007 5.143978e-01 0.5318571 0.60938144
#> TCONS_00000008 4.855835e-01 0.2155097 0.10431262
#> TCONS_00000009 1.832227e-05 0.2526343 0.28630615
#> TCONS_00000018 4.470920e-01 0.0000000 0.20569781
#> TCONS_00000019 0.000000e+00 0.3351936 0.05855237
#> TCONS_00000020 5.529069e-01 0.6648072 0.73574777
dim(ifMatrix)
#> [1] 269 3
Such a matrix can be used for global comparisons of multiple condtions and potential analysis are sample correlations:
# correlation plot
ggplot(melt(cor(ifMatrix)), aes(x=Var1, y=Var2, fill=value)) +
geom_tile() +
scale_fill_continuous('Correlation') +
labs(x='Condition', y='Condition') +
theme_minimal()
Or expression (via a heatmap):
library(pheatmap)
pheatmap(t(ifMatrix), show_colnames=FALSE)
Back to Table of Content
The workflow described here is an extension of the workflow enabled in the analyzeCPAT() function via the “removeNoncodinORFs” argument which if TRUE remove ORF information from the isoforms where the CPAT analysis classifies them as non-coding. Here we will “rescue”" the ORF of isoforms predicted to be noncoding by CPAT but which have a predicted protein domain. This workflow is an alternative apporach to setting “removeNoncodinORFs=TRUE” and should be used instead of that option.
Note that this is only recommended if ORFs were predicted (aka not imported from a GTF file). After this procedure, isoforms with ORFs will be isoforms with an ORF longer than minORFlength (if specified in analyzeORF) which are predicted to be coding by CPAT OR have a predicted protein domain (by Pfam).
Since the ORF information is stored in the ‘orfAnalysis’ analysis entry of the switchList we can remove it (by replacing it with NAs) as follows:
### Test data
data("exampleSwitchListAnalyzed")
exampleSwitchListAnalyzed
### Extract coding isoforms
nonCodingIsoforms <- unique(exampleSwitchListAnalyzed$isoformFeatures$isoform_id[
which( ! exampleSwitchListAnalyzed$isoformFeatures$codingPotential )
])
### Rescue those with protein domains
nonCodingIsoformsRescued <- setdiff(nonCodingIsoforms, exampleSwitchListAnalyzed$domainAnalysis$isoform_id)
# nr rescued
length(nonCodingIsoforms) - length(nonCodingIsoformsRescued)
### Remove noncoding isoforms ORF annotation
sum(is.na(exampleSwitchListAnalyzed$orfAnalysis$orfTransciptStart))
exampleSwitchListAnalyzed$orfAnalysis[
which( exampleSwitchListAnalyzed$orfAnalysis$isoform_id %in% nonCodingIsoformsRescued), 2:ncol(exampleSwitchListAnalyzed$orfAnalysis)
] <- NA
exampleSwitchListAnalyzed$isoformFeatures$PTC[which(exampleSwitchListAnalyzed$isoformFeatures$isoform_id %in% nonCodingIsoforms)] <- NA
sum(is.na(exampleSwitchListAnalyzed$orfAnalysis$orfTransciptStart))
Back to Table of Content
Recent research suggests that small upstream ORFs are far more frequent that previously assumed. It is therefore of particular interesting to start analyzing these and here we have indirectly presented a tool which can do just so: analyzeORF().
Here we show how one could start such an analysis of small upstream ORFs.
# run ORF analysis on longest ORF
exampleSwitchListAnalyzed <- analyzeORF(exampleSwitchListAnalyzed, genomeObject = Hsapiens, method='longest')
mean(exampleSwitchListAnalyzed$orfAnalysis$orfTransciptLength)
# run ORF analysis on most upstream ORF
exampleSwitchListAnalyzed2 <- analyzeORF(exampleSwitchListAnalyzed, genomeObject = Hsapiens, orfMethod = 'mostUpstream', minORFlength = 50)
mean(exampleSwitchListAnalyzed2$orfAnalysis$orfTransciptLength)
# calculate pairwise difference
summary(
exampleSwitchListAnalyzed2$orfAnalysis$orfTransciptLength -
exampleSwitchListAnalyzed$orfAnalysis$orfTransciptLength[
match(exampleSwitchListAnalyzed2$orfAnalysis$isoform_id ,exampleSwitchListAnalyzed$orfAnalysis$isoform_id)
]
)
Back to Table of Content
The sequences stored in the SwitchAnalyzeRlist are not needed after consequences have been predicted and can thereby be removed thereby reducing the object size as well as loading/saving times. This has for example been done for the example dataset ‘exampleSwitchListAnalyzed’. This is simply done as follows
summary(exampleSwitchListIntermediary)
exampleSwitchListIntermediary$ntSequence <- NULL
exampleSwitchListIntermediary$aaSequence <- NULL
summary(exampleSwitchListIntermediary)
Back to Table of Content
There has been quite a bit of debate about whether the default parameters for the codingPotential calculated for CPAT are to lenient (aka to low). This will always be a problem by having a single cutoff. One possible solution is to introduce an “unknown” class with medium size coding potential which we can then disregard.
Lets start by looking at the distribution of coding potential values
data("exampleSwitchListAnalyzed")
hist(exampleSwitchListAnalyzed$isoformFeatures$codingPotentialValue)
These coding potential values are summarized by the cutoff supplied in the ‘codingPotential’ column - which is what IsoformSwitchAnalyzeR uses in the downstream analysis
### These are summarized by the cutoff supplied in the 'codingPotential' column - which is what is used by IsoformSwitchAnalyzeR in the downstream analysis
table(exampleSwitchListAnalyzed$isoformFeatures$codingPotential, exclude = NULL)
By simply setting the mid-range values to NA it will cause IsoformSwitchAnalyzeR to ignore them thereby removing the isoforms with “unknown” coding potential. This can be done as follows:
exampleSwitchListAnalyzed$isoformFeatures$codingPotential <- NA
exampleSwitchListAnalyzed$isoformFeatures$codingPotential[which(exampleSwitchListAnalyzed$isoformFeatures$codingPotentialValue > 0.75)] <- TRUE
exampleSwitchListAnalyzed$isoformFeatures$codingPotential[which(exampleSwitchListAnalyzed$isoformFeatures$codingPotentialValue < 0.25)] <- FALSE
table(exampleSwitchListAnalyzed$isoformFeatures$codingPotential, exclude = NULL)
Back to Table of Content
As we have shown there are quite a lot of problems with known CDS annotation (see supplementary data in Vitting-Seerup et al 2017) we have implemented two ways to ensure a high quality of the CDS imported from a GTF annotation file:
The difference between the isoforms involved in an isoform switch can arise by changes in three distinct biological mechanisms:
Since we how the structure of the isoforms involved in a isoform switch we can also analyze which (combination) of these biological mechanisms gives rise to the difference between the two isoforms involved in a isoform switches.
This is simply done by, in addition to running analyzeSwitchConsequences with the consequences you find interesting, you make a separate consequence analysis of consequences (also with analyzeSwitchConsequences) where the consequences you analyze (supplied to the consequencesToAnalyze argument) are:
Then we can simply compare the result of this analysis to the isoform switches with consequences we already have identified to be of interest and thereby identify which biological mechanisms gives rise to the isoform switches with consequence you are interested in. One suggestion for such an analysis are illustrated here:
### Load example data
data("exampleSwitchListAnalyzed")
### Reduce datasize for fast runtime
randomGenes <- sample(unique(exampleSwitchListAnalyzed$isoformFeatures$gene_id), size = 40)
exampleSwitchListAnalyzedSubset <- subsetSwitchAnalyzeRlist(
exampleSwitchListAnalyzed,
exampleSwitchListAnalyzed$isoformFeatures$gene_id %in% randomGenes
)
### analyze the biological mechanismes
bioMechanismeAnalysis <- analyzeSwitchConsequences(
exampleSwitchListAnalyzedSubset,
consequencesToAnalyze = c('tss','tts','intron_structure'),
showProgress = FALSE
)$switchConsequence # only the consequences are interesting here
### subset to those with differences
bioMechanismeAnalysis <- bioMechanismeAnalysis[which(bioMechanismeAnalysis$isoformsDifferent),]
### extract the consequences of interest alerady stored in the switchAnalyzeRlist
myConsequences <- exampleSwitchListAnalyzedSubset$switchConsequence
myConsequences <- myConsequences[which(myConsequences$isoformsDifferent),]
myConsequences$isoPair <- paste(myConsequences$isoformUpregulated, myConsequences$isoformDownregulated) # id for specific iso comparison
### Obtain the mechanisms of the isoform switches with consequences
bioMechanismeAnalysis$isoPair <- paste(bioMechanismeAnalysis$isoformUpregulated, bioMechanismeAnalysis$isoformDownregulated)
bioMechanismeAnalysis <- bioMechanismeAnalysis[which(bioMechanismeAnalysis$isoPair %in% myConsequences$isoPair),] # id for specific iso comparison
This result is best summarized in a Venn diagram:
### Create list with the isoPair ids for each consequencee
AS <- bioMechanismeAnalysis$isoPair[ which( bioMechanismeAnalysis$featureCompared == 'intron_structure')]
aTSS <- bioMechanismeAnalysis$isoPair[ which( bioMechanismeAnalysis$featureCompared == 'tss' )]
aTTS <- bioMechanismeAnalysis$isoPair[ which( bioMechanismeAnalysis$featureCompared == 'tts' )]
mechList <- list(
AS=AS,
aTSS=aTSS,
aTTS=aTTS
)
### Create venn diagram
library(VennDiagram)
myVenn <- venn.diagram(mechList, col='transparent', alpha=0.4, fill=brewer.pal(n=3,name='Dark2'), filename=NULL)
### Plot the venn diagram
grid.newpage() ; grid.draw(myVenn)
From which the relative importance of each of the three mechanisms, as well as the combination of these, can be seen.
Back to Table of Content
None yet but we will continously add the most common problems and questions raised on github or in the google group (see How To Get Help).
Back to Table of Content
With this vignette, we hope to provide a thorough introduction to IsoformSwitchAnalyzeR as well as give some examples of what IsoformSwitchAnalyzeR can be used for.
We aim to continuously keep IsoformSwitchAnalyzeR up to date and update it. The update aspect includes the integration of new tools as they are developed (isoform quantification, isoform switch identification, statistical ORF detection, external sequence analysis etc.) so please feel free to suggest new tools to us (see the How To Get Help section for info of how to get in contact). Tools should preferably either be light weight (runnable on an old laptop) and/or available via web-services as that will allow a larger audience to use it.
Back to Table of Content
sessionInfo()
#> R version 3.4.2 (2017-09-28)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 16.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.6-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.6-bioc/R/lib/libRlapack.so
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> attached base packages:
#> [1] grid stats4 parallel stats graphics grDevices utils
#> [8] datasets methods base
#>
#> other attached packages:
#> [1] pheatmap_1.0.8 BSgenome.Hsapiens.UCSC.hg19_1.4.0
#> [3] BSgenome_1.46.0 Biostrings_2.46.0
#> [5] XVector_0.18.0 IsoformSwitchAnalyzeR_1.0.0
#> [7] spliceR_1.20.0 plyr_1.8.4
#> [9] RColorBrewer_1.1-2 VennDiagram_1.6.17
#> [11] futile.logger_1.4.3 cummeRbund_2.20.0
#> [13] Gviz_1.22.0 rtracklayer_1.38.0
#> [15] GenomicRanges_1.30.0 GenomeInfoDb_1.14.0
#> [17] IRanges_2.12.0 S4Vectors_0.16.0
#> [19] fastcluster_1.1.24 reshape2_1.4.2
#> [21] ggplot2_2.2.1 RSQLite_2.0
#> [23] BiocGenerics_0.24.0
#>
#> loaded via a namespace (and not attached):
#> [1] ProtGenerics_1.10.0 bitops_1.0-6
#> [3] matrixStats_0.52.2 bit64_0.9-7
#> [5] progress_1.1.2 httr_1.3.1
#> [7] rprojroot_1.2 tools_3.4.2
#> [9] backports_1.1.1 R6_2.2.2
#> [11] rpart_4.1-11 Hmisc_4.0-3
#> [13] DBI_0.7 lazyeval_0.2.1
#> [15] colorspace_1.3-2 nnet_7.3-12
#> [17] gridExtra_2.3 prettyunits_1.0.2
#> [19] RMySQL_0.10.13 bit_1.1-12
#> [21] curl_3.0 compiler_3.4.2
#> [23] Biobase_2.38.0 htmlTable_1.9
#> [25] DelayedArray_0.4.0 labeling_0.3
#> [27] scales_0.5.0 checkmate_1.8.5
#> [29] stringr_1.2.0 digest_0.6.12
#> [31] Rsamtools_1.30.0 foreign_0.8-69
#> [33] rmarkdown_1.6 pkgconfig_2.0.1
#> [35] base64enc_0.1-3 dichromat_2.0-0
#> [37] htmltools_0.3.6 limma_3.34.0
#> [39] ensembldb_2.2.0 htmlwidgets_0.9
#> [41] rlang_0.1.2 BiocInstaller_1.28.0
#> [43] shiny_1.0.5 BiocParallel_1.12.0
#> [45] acepack_1.4.1 VariantAnnotation_1.24.0
#> [47] RCurl_1.95-4.8 magrittr_1.5
#> [49] GenomeInfoDbData_0.99.1 Formula_1.2-2
#> [51] Matrix_1.2-11 Rcpp_0.12.13
#> [53] munsell_0.4.3 edgeR_3.20.0
#> [55] stringi_1.1.5 yaml_2.1.14
#> [57] MASS_7.3-47 SummarizedExperiment_1.8.0
#> [59] zlibbioc_1.24.0 AnnotationHub_2.10.0
#> [61] blob_1.1.0 lattice_0.20-35
#> [63] DRIMSeq_1.6.0 splines_3.4.2
#> [65] GenomicFeatures_1.30.0 locfit_1.5-9.1
#> [67] knitr_1.17 biomaRt_2.34.0
#> [69] futile.options_1.0.0 XML_3.98-1.9
#> [71] evaluate_0.10.1 biovizBase_1.26.0
#> [73] latticeExtra_0.6-28 lambda.r_1.2
#> [75] data.table_1.10.4-3 httpuv_1.3.5
#> [77] gtable_0.2.0 assertthat_0.2.0
#> [79] mime_0.5 xtable_1.8-2
#> [81] AnnotationFilter_1.2.0 survival_2.41-3
#> [83] tibble_1.3.4 GenomicAlignments_1.14.0
#> [85] AnnotationDbi_1.40.0 memoise_1.1.0
#> [87] tximport_1.6.0 cluster_2.0.6
#> [89] interactiveDisplayBase_1.16.0
Back to Table of Content