EventPointer is an R package to identify alternative splicing events that involve either simple (case-control experiment) or complex experimental designs such as time course experiments and studies including paired-samples. The algorithm can be used to analyze data from either junction arrays (Affymetrix Arrays) or sequencing data (RNA-Seq).
The algorithm can generate a series of files to visualize the detected alternative splicing events in IGV. This eases the interpretation of results and the design of primers for standard PCR validation. The software returns a data.frame with the detected alternative splicing events: gene name, type of event (cassette, alternative 3’,…,etc), genomic position, statistical significance and increment of the percent spliced in (\(\Delta \Psi\)) for all the events.
EventPointer can be installed from Bioconductor using the BiocManager package:
library(BiocManager)
if (!requireNamespace("BiocManager", quietly=TRUE))
install.packages("BiocManager")
BiocManager::install("EventPointer")
EventPointer R package provides users a simplified way to identify, classify and visualize alternative splicing events. EventPointer can work with arrays and bulk-RNA-Seq data. For arrays there is a unique pipeline. For RNA-Seq there are two main pipelines: one able to detect novel events and the other one relying on annotated events. The steps required by the algorithm are almost the same for all approaches:
Splicing Graph Creation The Splicing Graph (SG) is a directed graph used to represent the structure of a gene. EventPointer relies on SGSeq package to obtain the corresponding SGs for every gene. For arrays, the SG is built according to the reference transcriptome selected by the user and for RNA-Seq, it is created by predicting the features in the .bam files provided by the user or by the reference transcriptome.
Event Detection Once the graphs are created, EventPointer analyzes each SG in order to find the alternative splicing events. The definition of splicing events by EventPointer consists in a triplet of subgraphs (P1,P2 and PR) i.e. for a cassette exon, PR correspond to the flanking exons, P1 the junctions and exon, and P2 the junction that skips the exon. This is depicted in Figure 1.
Figure 1. Definition of alternative splicing events for EventPointer
Figure 2 shows each of the main steps in a graphical layout.
EventPointer was originally developed for array-based analysis and has since been adapted for bulk RNA-Seq. This vignette is structured into two main sections: RNA-Seq analysis and array analysis. The RNA-Seq section is divided into two parts:
EventPointerBAM: Alignment-based analysis, which detects novel events specific to the samples (novel version).
EventPonterST: Analysis of events annotated in a reference transcriptome.
To cite EventPointer:
For EventPonterST:
For EventPointer for arryas:
Romero, Juan P., et al. “EventPointer: an effective identification of alternative splicing events using junction arrays.” BMC genomics 17.1 (2016): 467. doi:10.1186/s12864-016-2816-x
Romero, Juan P., et al. “Comparison of RNA-seq and microarray platforms for splice event detection using a cross-platform algorithm.” BMC genomics 19.1 (2018): 703. doi:10.1186/s12864-018-5082-2
Figure 2. Overview of EventPointer
EventPointer has two pipelines for RNA-Seq analysis:
EventPointerBAM (EP_BAM): The detection, classification and quantification of splicing events is carried out using the BAM files generated from RNA-seq sequence alignment. It should be noted that the alignment must be performed using an aligner that takes splicing into account (recommended tool: STAR). The main advantage of this methodology is its ability to detect de novo events based on reads alignment.
EventPointerST (EP_ST): The detection and classification of events are determined by the corresponding transcriptome annotation in GTF format. In this case, the quantification files of transcripts come from pseudoaligners (examples: salmon or kallisto).
As mentioned above, EventPointer consists of several core steps: (i) building the splicing graph, (ii) detecting and classifying events, (iii) quantification, and (iv) performing statistical analysis. Finally, users can visualize the results with IGV by files generated by EventPointer.
The Figure 3 shows the pipeline of EventPointerBAM and EventPointerST:
Figure 3 (A) Worflow and advantages and disadvantages of EventPointerBAM. (B) Worflow and advantages and disadvantages of EventPointerST.
The aim of EventPointerBAM is to perform detection, classification,
quantification and differential analysis between conditions, using as
input data sequence alignments of RNA-seq data performed using aligners
sensitive to junction reads.
There are two requirements for the BAM files in order
to get EventPointerBAM working correctly:
The BAM files should include the XS-flag, the flag can be included in the files when running the alignment. Most of the available software has parameters to incude the flag. For example, in the case of STAR the flag –outSAMattributes XS must be included when mapping.
All files to be analyzed should have the corresponding index files (.bai) in the same directory as the BAM files. Create the index before running EventPointer.
To demonstrate the use of EventPointer, we utilise synthetically generated FASTQ files containing genes with predefined alternative splicing events, including alternative 3’ sites, alternative 5’ sites, and cassette exons. These FASTQ files are first processed with the STAR aligner to generate the BAM files required as input for EventPointer. In order to reduce computation time in this example, the data have been limited to genes from chromosome 3 with positive strand (strand +).
The function EventsDetection_BAM() is responsible for the detection of events from the input alignment files, their classification and the calculation of the associated PSI. The inputs for the function are described below:
Bams and GTF file parameters
PathSamplesAbundance: Path to bam and bai files or path to folder with bam and bai files.
PathTranscriptomeGTF: Path to file containing the regions to be analysed from the bam files in GTF format.
Detection parameters
Rsamtools library:path_bam_files = system.file("extdata/bams", package = "EventPointer")
path_file_1 = dir(path_bam_files,pattern = "*.bam$",full.names = TRUE)[1]
file_tmp <- BamFile(path_file_1)
si <- seqinfo(file_tmp)
sl <- rep(seqlevels(si),2)
st <- c(rep(c("+"), length(si)),rep(c("-"), length(si)))
ranges_info <- GRanges(sl, IRanges(1, seqlengths(si)[sl]), st)
| index | seqnames | start | end | width | strand |
|---|---|---|---|---|---|
| 1 | 1 | 1 | 249250621 | 249250621 | p |
| 2 | 10 | 1 | 135534747 | 135534747 | p |
| 3 | 11 | 1 | 135006516 | 135006516 | p |
| 4 | 12 | 1 | 133851895 | 133851895 | p |
| 5 | 13 | 1 | 115169878 | 115169878 | p |
| 6 | 14 | 1 | 107349540 | 107349540 | p |
| 7 | 15 | 1 | 102531392 | 102531392 | p |
| 8 | 16 | 1 | 90354753 | 90354753 | p |
| 9 | 17 | 1 | 81195210 | 81195210 | p |
| 10 | 18 | 1 | 78077248 | 78077248 | p |
| 11 | 19 | 1 | 59128983 | 59128983 | p |
| 12 | 2 | 1 | 243199373 | 243199373 | p |
| 13 | 20 | 1 | 63025520 | 63025520 | p |
| 14 | 21 | 1 | 48129895 | 48129895 | p |
| 15 | 22 | 1 | 51304566 | 51304566 | p |
| 16 | 3 | 1 | 198022430 | 198022430 | p |
| 17 | 4 | 1 | 191154276 | 191154276 | p |
| 18 | 5 | 1 | 180915260 | 180915260 | p |
| 19 | 6 | 1 | 171115067 | 171115067 | p |
| 20 | 7 | 1 | 159138663 | 159138663 | p |
| 21 | 8 | 1 | 146364022 | 146364022 | p |
| 22 | 9 | 1 | 141213431 | 141213431 | p |
| 23 | MT | 1 | 16569 | 16569 | p |
| 24 | X | 1 | 155270560 | 155270560 | p |
| 25 | Y | 1 | 59373566 | 59373566 | p |
| 26 | GL000192.1 | 1 | 547496 | 547496 | p |
| 27 | GL000225.1 | 1 | 211173 | 211173 | p |
| 28 | GL000194.1 | 1 | 191469 | 191469 | p |
| 29 | GL000193.1 | 1 | 189789 | 189789 | p |
| 30 | GL000200.1 | 1 | 187035 | 187035 | p |
| 31 | GL000222.1 | 1 | 186861 | 186861 | p |
| 32 | GL000212.1 | 1 | 186858 | 186858 | p |
| 33 | GL000195.1 | 1 | 182896 | 182896 | p |
| 34 | GL000223.1 | 1 | 180455 | 180455 | p |
| 35 | GL000224.1 | 1 | 179693 | 179693 | p |
| 36 | GL000219.1 | 1 | 179198 | 179198 | p |
| 37 | GL000205.1 | 1 | 174588 | 174588 | p |
| 38 | GL000215.1 | 1 | 172545 | 172545 | p |
| 39 | GL000216.1 | 1 | 172294 | 172294 | p |
| 40 | GL000217.1 | 1 | 172149 | 172149 | p |
| 41 | GL000199.1 | 1 | 169874 | 169874 | p |
| 42 | GL000211.1 | 1 | 166566 | 166566 | p |
| 43 | GL000213.1 | 1 | 164239 | 164239 | p |
| 44 | GL000220.1 | 1 | 161802 | 161802 | p |
| 45 | GL000218.1 | 1 | 161147 | 161147 | p |
| 46 | GL000209.1 | 1 | 159169 | 159169 | p |
| 47 | GL000221.1 | 1 | 155397 | 155397 | p |
| 48 | GL000214.1 | 1 | 137718 | 137718 | p |
| 49 | GL000228.1 | 1 | 129120 | 129120 | p |
| 50 | GL000227.1 | 1 | 128374 | 128374 | p |
| 51 | GL000191.1 | 1 | 106433 | 106433 | p |
| 52 | GL000208.1 | 1 | 92689 | 92689 | p |
| 53 | GL000198.1 | 1 | 90085 | 90085 | p |
| 54 | GL000204.1 | 1 | 81310 | 81310 | p |
| 55 | GL000233.1 | 1 | 45941 | 45941 | p |
| 56 | GL000237.1 | 1 | 45867 | 45867 | p |
| 57 | GL000230.1 | 1 | 43691 | 43691 | p |
| 58 | GL000242.1 | 1 | 43523 | 43523 | p |
| 59 | GL000243.1 | 1 | 43341 | 43341 | p |
| 60 | GL000241.1 | 1 | 42152 | 42152 | p |
| 61 | GL000236.1 | 1 | 41934 | 41934 | p |
| 62 | GL000240.1 | 1 | 41933 | 41933 | p |
| 63 | GL000206.1 | 1 | 41001 | 41001 | p |
| 64 | GL000232.1 | 1 | 40652 | 40652 | p |
| 65 | GL000234.1 | 1 | 40531 | 40531 | p |
| 66 | GL000202.1 | 1 | 40103 | 40103 | p |
| 67 | GL000238.1 | 1 | 39939 | 39939 | p |
| 68 | GL000244.1 | 1 | 39929 | 39929 | p |
| 69 | GL000248.1 | 1 | 39786 | 39786 | p |
| 70 | GL000196.1 | 1 | 38914 | 38914 | p |
| 71 | GL000249.1 | 1 | 38502 | 38502 | p |
| 72 | GL000246.1 | 1 | 38154 | 38154 | p |
| 73 | GL000203.1 | 1 | 37498 | 37498 | p |
| 74 | GL000197.1 | 1 | 37175 | 37175 | p |
| 75 | GL000245.1 | 1 | 36651 | 36651 | p |
| 76 | GL000247.1 | 1 | 36422 | 36422 | p |
| 77 | GL000201.1 | 1 | 36148 | 36148 | p |
| 78 | GL000235.1 | 1 | 34474 | 34474 | p |
| 79 | GL000239.1 | 1 | 33824 | 33824 | p |
| 80 | GL000210.1 | 1 | 27682 | 27682 | p |
| 81 | GL000231.1 | 1 | 27386 | 27386 | p |
| 82 | GL000229.1 | 1 | 19913 | 19913 | p |
| 83 | GL000226.1 | 1 | 15008 | 15008 | p |
| 84 | GL000207.1 | 1 | 4262 | 4262 | p |
| 85 | 1 | 1 | 249250621 | 249250621 | n |
| 86 | 10 | 1 | 135534747 | 135534747 | n |
| 87 | 11 | 1 | 135006516 | 135006516 | n |
| 88 | 12 | 1 | 133851895 | 133851895 | n |
| 89 | 13 | 1 | 115169878 | 115169878 | n |
| 90 | 14 | 1 | 107349540 | 107349540 | n |
| 91 | 15 | 1 | 102531392 | 102531392 | n |
| 92 | 16 | 1 | 90354753 | 90354753 | n |
| 93 | 17 | 1 | 81195210 | 81195210 | n |
| 94 | 18 | 1 | 78077248 | 78077248 | n |
| 95 | 19 | 1 | 59128983 | 59128983 | n |
| 96 | 2 | 1 | 243199373 | 243199373 | n |
| 97 | 20 | 1 | 63025520 | 63025520 | n |
| 98 | 21 | 1 | 48129895 | 48129895 | n |
| 99 | 22 | 1 | 51304566 | 51304566 | n |
| 100 | 3 | 1 | 198022430 | 198022430 | n |
| 101 | 4 | 1 | 191154276 | 191154276 | n |
| 102 | 5 | 1 | 180915260 | 180915260 | n |
| 103 | 6 | 1 | 171115067 | 171115067 | n |
| 104 | 7 | 1 | 159138663 | 159138663 | n |
| 105 | 8 | 1 | 146364022 | 146364022 | n |
| 106 | 9 | 1 | 141213431 | 141213431 | n |
| 107 | MT | 1 | 16569 | 16569 | n |
| 108 | X | 1 | 155270560 | 155270560 | n |
| 109 | Y | 1 | 59373566 | 59373566 | n |
| 110 | GL000192.1 | 1 | 547496 | 547496 | n |
| 111 | GL000225.1 | 1 | 211173 | 211173 | n |
| 112 | GL000194.1 | 1 | 191469 | 191469 | n |
| 113 | GL000193.1 | 1 | 189789 | 189789 | n |
| 114 | GL000200.1 | 1 | 187035 | 187035 | n |
| 115 | GL000222.1 | 1 | 186861 | 186861 | n |
| 116 | GL000212.1 | 1 | 186858 | 186858 | n |
| 117 | GL000195.1 | 1 | 182896 | 182896 | n |
| 118 | GL000223.1 | 1 | 180455 | 180455 | n |
| 119 | GL000224.1 | 1 | 179693 | 179693 | n |
| 120 | GL000219.1 | 1 | 179198 | 179198 | n |
| 121 | GL000205.1 | 1 | 174588 | 174588 | n |
| 122 | GL000215.1 | 1 | 172545 | 172545 | n |
| 123 | GL000216.1 | 1 | 172294 | 172294 | n |
| 124 | GL000217.1 | 1 | 172149 | 172149 | n |
| 125 | GL000199.1 | 1 | 169874 | 169874 | n |
| 126 | GL000211.1 | 1 | 166566 | 166566 | n |
| 127 | GL000213.1 | 1 | 164239 | 164239 | n |
| 128 | GL000220.1 | 1 | 161802 | 161802 | n |
| 129 | GL000218.1 | 1 | 161147 | 161147 | n |
| 130 | GL000209.1 | 1 | 159169 | 159169 | n |
| 131 | GL000221.1 | 1 | 155397 | 155397 | n |
| 132 | GL000214.1 | 1 | 137718 | 137718 | n |
| 133 | GL000228.1 | 1 | 129120 | 129120 | n |
| 134 | GL000227.1 | 1 | 128374 | 128374 | n |
| 135 | GL000191.1 | 1 | 106433 | 106433 | n |
| 136 | GL000208.1 | 1 | 92689 | 92689 | n |
| 137 | GL000198.1 | 1 | 90085 | 90085 | n |
| 138 | GL000204.1 | 1 | 81310 | 81310 | n |
| 139 | GL000233.1 | 1 | 45941 | 45941 | n |
| 140 | GL000237.1 | 1 | 45867 | 45867 | n |
| 141 | GL000230.1 | 1 | 43691 | 43691 | n |
| 142 | GL000242.1 | 1 | 43523 | 43523 | n |
| 143 | GL000243.1 | 1 | 43341 | 43341 | n |
| 144 | GL000241.1 | 1 | 42152 | 42152 | n |
| 145 | GL000236.1 | 1 | 41934 | 41934 | n |
| 146 | GL000240.1 | 1 | 41933 | 41933 | n |
| 147 | GL000206.1 | 1 | 41001 | 41001 | n |
| 148 | GL000232.1 | 1 | 40652 | 40652 | n |
| 149 | GL000234.1 | 1 | 40531 | 40531 | n |
| 150 | GL000202.1 | 1 | 40103 | 40103 | n |
| 151 | GL000238.1 | 1 | 39939 | 39939 | n |
| 152 | GL000244.1 | 1 | 39929 | 39929 | n |
| 153 | GL000248.1 | 1 | 39786 | 39786 | n |
| 154 | GL000196.1 | 1 | 38914 | 38914 | n |
| 155 | GL000249.1 | 1 | 38502 | 38502 | n |
| 156 | GL000246.1 | 1 | 38154 | 38154 | n |
| 157 | GL000203.1 | 1 | 37498 | 37498 | n |
| 158 | GL000197.1 | 1 | 37175 | 37175 | n |
| 159 | GL000245.1 | 1 | 36651 | 36651 | n |
| 160 | GL000247.1 | 1 | 36422 | 36422 | n |
| 161 | GL000201.1 | 1 | 36148 | 36148 | n |
| 162 | GL000235.1 | 1 | 34474 | 34474 | n |
| 163 | GL000239.1 | 1 | 33824 | 33824 | n |
| 164 | GL000210.1 | 1 | 27682 | 27682 | n |
| 165 | GL000231.1 | 1 | 27386 | 27386 | n |
| 166 | GL000229.1 | 1 | 19913 | 19913 | n |
| 167 | GL000226.1 | 1 | 15008 | 15008 | n |
| 168 | GL000207.1 | 1 | 4262 | 4262 | n |
In this case, as we are interested in the chr3 with + starnd, we want the
region 16. Therefore, we will set the parameter region to 16.
min_junction_count: Minimum number of junctions detected in the alignment to be considered in the splicing graph. Default value is 2.
max_complexity: Maximum allowed complexity. If a locus exceeds this threshold, it is skipped, resulting in a warning. Complexity is defined as the maximum number of unique predicted splice junctions overlapping a given position. High complexity regions are often due to spurious read alignments and can slow down processing. Default value is 30.
min_n_sample: Minimum number of samples that a junction must have to be considered. Default value is 3.
min_anchor: Minimum number of aligned bases at one end of an exon to consider a junction. Default value is 6.
PSI calculation parameters
nboot: Number of resamples of the quantification of the samples to perform bootstrap. Default value is 20.
lambda:
Performance parameters
cores: Nº of cores used for the function. Default value is 1.
verbose: Logical indicating whether to show warnings and messages. If FALSE, warnings from internal functions (e.g., makeTxDbFromGRanges) will be suppressed. Default is FALSE.
Results parameters
PathSGResult: Path where to save the following 4 files:
The following code chunk describe how to use this function:
PathSamplesAbundance <- system.file("extdata/bams", package = "EventPointer")
PathTranscriptomeGTF <- list.files(PathSamplesAbundance,"*.gtf",full.names = T)
PathSGResult <- tempdir()
EventsDetection_BAM(PathSamplesAbundance,
PathTranscriptomeGTF,
region = 16,
min_junction_count = 2,
max_complexity = 30,
min_n_sample = NULL,
min_anchor = 6,
nboot = 20,
lambda = NULL,
cores = 1,
PathSGResult = PathSGResult)
The 4 files created by this function are:
A further explanation of these files are depicted in the following subsecction “step-by-step”
As mentioned above, the function EventsDetection_BAM unifies three main steps.
As each task has a long duration and it may be of interest
to a user to perform any of these
tasks separately, the following section explains how to execute this
function step by step.
get Information for the Splicing graph
The first step of EventsDetection_BAM is to get information to build the splicing graph.
For this, the function PrepareBam_EP is called. The parameters of this
functions are depicted above.
PathSamplesAbundance <- system.file("extdata/bams", package = "EventPointer")
PathTranscriptomeGTF <- list.files(PathSamplesAbundance,"*.gtf",full.names = T)
PathSGResult <- tempdir()
SgFC <- PrepareBam_EP(PathSamplesAbundance = PathSamplesAbundance,
PathTranscriptomeGTF = PathTranscriptomeGTF,
region = 16,
min_junction_count = 2,
max_complexity = 30,
min_n_sample = NULL,
min_anchor = 6,
cores = 1,
PathSGResult = PathSGResult,
verbose = FALSE)
This functions is the most demading regarding time computing. The SgFC object can be laoded from the EventPointer package:
load(system.file("extdata/SgFC.RData", package = "EventPointer"))
SgFC
## class: SGFeatureCounts
## dim: 1186 9
## metadata(0):
## assays(2): counts FPKM
## rownames: NULL
## rowData names(0):
## colnames(9):
## SRR1513329_switched_A5_A3Aligned.sortedByCoord.out_chr3.bam
## SRR1513329_switched_SEAligned.sortedByCoord.out_chr3.bam ...
## SRR1513331_switched_SEAligned.sortedByCoord.out_chr3.bam
## SRR1513331Aligned.sortedByCoord.out_chr3.bam
## colData names(6): sample_name file_bam ... frag_length lib_size
The SgFC object is a SGFeatureCounts object. for more information check
SGSeq
Vignette.
Briefly the SGFeaturesCounts contains a GRanges object with all the
required elements to create the different splicing graphs found in the
given samples. It also contains the number of counts associated with
each element of the splicing graph.
Event Detection
The second step of EventsDetection_BAM takes the SgFC object to create
the splicing graph and then detect and classify splicing events.
This is done by calling the function EventDetection.
This function only requieres two parameters: the SgFC object, and the
number of cores:
EventsDetection_EPBAM<-EventDetection(SgFC, cores=1)
The function retunes a list called EventsDetection_EPBAM that
corresponds to the EventsDetection_EPBAM.RData generated by
EventsDetection_BAM function.
The output is a list vector with
one element for each gene. Each gene (element of the list) is itself a
list with a length equal to the number of events in that gene. Finally,
each event is also a list containing the event information. In order to
see the values of this object, the follownig code chunk show an example
of how to see the information of the first event of the first gene.
index_gene <- 1
index_event <- 1
length(EventsDetection_EPBAM) #number of genes
length(EventsDetection_EPBAM[[index_gene]]) #number of events in the corresponding gene
names(EventsDetection_EPBAM[[index_gene]][[index_event]]) #the different information stored in the object.
In order to replicate the .csv file stored by the
EventsDetection_BAM function, you can use the following code chunk:
event_list <- list()
index <- 1
for (geneEvent in EventsDetection_EPBAM) {
for (event in geneEvent) {
event_list[[index]] <- event$Info
index <- index + 1
}
}
totalEventTable <- do.call(rbind, event_list)
write.csv(totalEventTable, file = paste0(PathSGResult, "/TotalEventsFound.csv"), row.names = FALSE)
Get PSI
Finally, the EventsDetection_BAM function compute the PSI values by
calling the function getPSI_RNASeq_boot which needs the following
parameters: lambda, and nboot (depicted above).
cores <- 1
lambda <- NULL
nboot <- 20
PSI_boots <- getPSI_RNASeq_boot(EventsDetection_EPBAM, lambda = lambda,
cores = cores, nboot = nboot)
PSI_boots is an array of number_of_events x (number_of_boots+1) x number_of_samples.
In the second dimension, we have the PSI values corresponding to the maximum likelihood estimation of the expression and the PSI corresponding to the nboot bootstraps.
The main function to perform the statistical analysis is EventPointerStats_BAM.
data(PSI_boots)
Design <- cbind(rep(1, 9), rep(c(1, 0, 0), 3), rep(c(0, 1, 0), 3))
Contrast <- cbind(c(0, 1, 0), c(0, 0, 1))
PathSGResult <- tempdir()
result = EventPointerStats_BAM(PSI_boots,
Design,
Contrast,
Threshold = 0,
nbootstraps = 1000,
cores = 1,
ram = 0.1,
pathResult = PathSGResult)
##
## The Contrast and Design Matrices have been correctly designed.
##
##
## The Contrast and Design Matrices have been correctly designed.
##
## Number of chunks: 1
##
## Performing resampling of chunk number 1 , this may take a while...
##
## ****Analysis 100 % Completed****
##
## Preparing output...
## finish chunk 1
##
## The program has succesfully ended.
EventPointerStats_BAM returns a list with the reulsts of the statistical analysis.
The list contains the summary of the Bootstrap analysis: deltaPSI, Pvalues, iqr_events,
median_events, and LocalFDR. The first for are data.frame with the corresponding information
for each event (row) and contrast. LocalFDR is a list of length equal to the number
of contrast. Each elemnt has information of the local flase dicovery rate.
Brielfy:
EventPointerStats_BAM also return files with the results of the statistical analysis
for each contrast in a .csv file.
Alternative splicing events detected by EventPointerBAM can be visualized using the EventPointerBAM_IGV() function as shown in Figure 4. This function displays the paths that make up the splicing graph of each event in a gtf format that can be viewed using the IGV program. Each event is composed of 3 paths: path1, path2, and Reference, which will be indicated along with the event id for identification as “path1_”, “path2_”, “Ref_1”, and “Ref_2”.
To create the GTF files, the algorithm uses the EventPointerBAM_IGV() function that
has the following parameters:
Events and reference transcriptome parameters
EventsDetection_BAM with the events to be included in the GTF file.Results parameters
EventsTxt<-paste(system.file("extdata",package="EventPointer"),"/TotalEventsFound.csv",sep="")
load(paste(system.file("extdata",package="EventPointer"),"/SgFC.RData",sep=""))
SG_RNASeq <- SgFC
PathGTF <- tempdir()
EventPointerBAM_IGV(SG_RNASeq = SgFC,
EventsCSV = EventsTxt,
PathGTF = "./")
The GTF file can be loaded in IGV and will display as following:
Figure 4. EventPointer visualization using IGV.
This approach of EventPointer can be run with tow main functions: EventsDetection_BAM() and EventPointerStats_BAM(). The former functin can be run step by step by using the function:
PrepareBam_EP, EventDetection, and getPSI_RNASeq_boot.
When running both EventsDetection_BAM() and EventPointerStats_BAM(), the resulting files at the defined output folder will be those shown in Figure 5.
Figure 5. Output EventPointerBAM functions.
Finally, results can be visualized using the function EventPointerBAM_IGV.
In this pipeline alternative splicng events are detected from a reference transcriptome without finding novel events.
The events quantification relies on isoform expression estimate from pseudo-alignment process such as Kallisto or Salmon. Besides, we provide a function to leverage the bootstrap data from kallisto or salmon. Preveous statistical analysis have also been adapted to this data. Further, Primers design for PCR validation and protein domain enrichment analysis can be performed. Figure 6 shows an overview of this branch of EventPointer.
Figure 6. Analysis of RNA-Seq data with kallisto output
We use EventDetection_transcriptome to identify and classify
alternative splicing events of a given reference transcriptome. The
required parameters are:
Following code shows how to apply EventDetection_transcriptome
function using the Gencode 24 transcriptome (GRCH 38) as reference
annotation. In this example we do not use the full reference but only
two genes (ENSG00000185252.17 and ENSG00000254709.7).
# Set input variables
PathFiles<-system.file("extdata",package="EventPointer")
inputFile <- paste(PathFiles,"/gencode.v24.ann_2genes.gtf",sep="")
Pathtxt <- tempdir()
# Run the function
EventXtrans <- EventDetection_transcriptome(inputFile = inputFile,
Pathtxt=Pathtxt,
cores=1)
This function takes a couple of hours to complete (depending on the computer), and creates a file called EventsFound.txt. This file contains the information of the events. This .txt file will be used as an input in further functions.
The object returned by the function is a list containing five elements: three sparce matrices that relate which isoforms build up the paths (path1,path2 and pathRef) of each event. The fourth element contains the name of the reference annotation: only appear the name of the transcript. The final element is SG_List: a list with the information of the graph of each gene, this variable is necesary for Primers design step.
names(EventXtrans)
NOTE: as the events are detected based on the reference transcriptome, the function
EventDetection_transcriptome can be applied only once for each reference
transcriptome.
Given (I) the output of EventDetection_transcriptome, which contains
events information and (II) the result of pseudo-alignment, EventPointer
is able to get the expression of the path of each event and also to
compute the PSI value.
To first step is to load the data from pseudo-alignment. EventPointer works with both Kallisto and Salmon methods.
In this example we have applied kallisto to the samples ERR315326, ERR315400, ERR315407 and ERR315494 using the transcriptome reference used before. These samples have been download from EMBL-EBI dataset.
EventPointer offers the option to take advantage of bootstraps that
return both Kallisto and Salmon. In case we want to leverage the
Bootstrap data we have to use the function getbootstrapdata:
PathSamples<-system.file("extdata",package="EventPointer")
PathSamples <- paste0(PathSamples,"/output")
PathSamples <- dir(PathSamples,full.names = TRUE)
data_exp <- getbootstrapdata(PathSamples = PathSamples,type = "kallisto")
## 1 2 3 4
getbootstrapdata returns a list containing the quantification data
with the bootstrap information. The list is of length equal to the
number of samples. Each element of the list storage a matrix where the
number of rows is equal to the number of transcripts and the number of
columns is equal to the number of bootstrats + 1. (the first column
corresponds to the maximum likelihood expression and the rest to the
bootstrap data).
In case you don’t want to use the bootstrap returned by the pseudo-alignmnet step (or is no available), you should load a matrix with the expression of the transcripts. The following code chunk shown an example of how to do this:
# this code chunk is an example of how to load the data from kallisto output.
# the expression of the isoforms are counts
PathFiles <- system.file("extdata",package="EventPointer")
filesnames <- dir(paste0(PathFiles,"/output"))
PathFiles <- dir(paste0(PathFiles,"/output"),full.names = TRUE)
dirtoload <- paste0(PathFiles,"/","abundance.tsv")
RNASeq <- read.delim(dirtoload[1],sep = "\t", colClasses = c(NA,"NULL","NULL",NA,"NULL"))
for (n in 2:length(dirtoload)){
RNASeq[,n+1] <- read.delim(dirtoload[n],sep = '\t', colClasses = c('NULL','NULL','NULL',NA,'NULL'))
}
rownames(RNASeq)<-RNASeq[,1]
RNASeq<-RNASeq[,-1]
colnames(RNASeq) <- filesnames
Once we have the expression of the isoforms loaded, the next step is to
compute the PSI value. For both cases (with and withuot bootstrap data
from psudoaligment) EventPointer provides the function
GetPSI_FromTranRef for this purpose.
The function GetPSI_FromTranRef returns the values of \(\Psi\) and the
expression of the paths of each events. This functions requires the
following inputs:
This function also requires that the same reference transcriptome have been used in both the pseudo-alignment and event detections steps. Besides, in the variable Samples the type of annotation used must be equal to the one used in the variable EventXtrans.
Get psi if NO bootstrap from pseudo-alignmnet is used
First step is to verify if same annotation is used in RNASeq and EventXtrans:
rownames(RNASeq)[1:5]
## [1] "ENST00000400451.6|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319648.2|ZNF74-005|ZNF74|3159|protein_coding|"
## [2] "ENST00000611540.4|ENSG00000185252.17|OTTHUMG00000150687.4|-|ZNF74-202|ZNF74|3714|protein_coding|"
## [3] "ENST00000403682.7|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319645.2|ZNF74-002|ZNF74|2789|protein_coding|"
## [4] "ENST00000357502.5|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319622.3|ZNF74-001|ZNF74|2788|protein_coding|"
## [5] "ENST00000493734.5|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319646.2|ZNF74-003|ZNF74|3758|retained_intron|"
EventXtrans$transcritnames[1:5]
## [1] "ENST00000611540.4" "ENST00000400451.6" "ENST00000403682.7"
## [4] "ENST00000357502.5" "ENST00000493734.5"
We need to change the rownames of RNASeq variable before applying
GetPSI_FromTranRef function:
#change rownames of RNASeq variable
rownames(RNASeq) <- sapply(strsplit(rownames(RNASeq),"\\|"),function(X) return(X[1]))
RNASeq<-as.matrix(RNASeq) #must be a matrix variable
PSIss_nb <- GetPSI_FromTranRef(PathsxTranscript = EventXtrans,
Samples = RNASeq,
Bootstrap = FALSE,
Filter = FALSE)
PSI <- PSIss_nb$PSI
Expression <- PSIss_nb$ExpEvs
The output of the function is a list containing two elements: a matrix with the \(\Psi\) values, and a list containing as many matrices as number of events. In each matrix is stored the expression of the different paths of the correspondig event along the samples.
If bootstrap data from pseudo-alignment IS required
As done above, we need to check if same annotation is used in data_exp and EventXtrans. The difference from above is that data_exp is a list and not a matrix as RNASeq. Therefore, we are going to used the rownames of the first list of data_exp to check if same annotations is used:
rownames(data_exp[[1]])[1:5]
## [1] "ENST00000400451.6|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319648.2|ZNF74-005|ZNF74|3159|protein_coding|"
## [2] "ENST00000611540.4|ENSG00000185252.17|OTTHUMG00000150687.4|-|ZNF74-202|ZNF74|3714|protein_coding|"
## [3] "ENST00000403682.7|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319645.2|ZNF74-002|ZNF74|2789|protein_coding|"
## [4] "ENST00000357502.5|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319622.3|ZNF74-001|ZNF74|2788|protein_coding|"
## [5] "ENST00000493734.5|ENSG00000185252.17|OTTHUMG00000150687.4|OTTHUMT00000319646.2|ZNF74-003|ZNF74|3758|retained_intron|"
EventXtrans$transcritnames[1:5]
## [1] "ENST00000611540.4" "ENST00000400451.6" "ENST00000403682.7"
## [4] "ENST00000357502.5" "ENST00000493734.5"
We need to change the rownames of the first element of the list
data_exp before applying GetPSI_FromTranRef function:
#change rownames of the first element of teh list data_exp
rownames(data_exp[[1]]) <- sapply(strsplit(rownames(data_exp[[1]]),"\\|"),function(X) return(X[1]))
PSIss <- GetPSI_FromTranRef(PathsxTranscript = EventXtrans,
Samples = data_exp,
Bootstrap = TRUE,
Filter = FALSE)
PSI <- PSIss$PSI
Expression <- PSIss$ExpEvs
The output of the function is a list containing two elements: an array with the \(\Psi\) values with dimension = c(number of events, number of bootstrap+1,number of samples), and a list containing as many matrices as number of events. In each matrix is stored the expression of the different paths of the corresponding event along the samples.
We have implemented a bootstrap test to evaluate the difference in \(\Psi\) among the conditions under study. This test can be done both for PSI values that contain information from the pseudo-alignment bootstrap and for cases in which we do not have this information.
The function EventPointer_Bootstraps proceed to calculate the increase
in PSI and its p.values by performing bootstrap statistic test. This
function can be parallelized. Users may set the number of Cores (by
default 1) and the ram available (by default 1 Gb).
Dmatrix <- cbind(1,rep(c(0,1),each=2))
Cmatrix <- matrix(c(0,1),nrow=2)
Fit_bots <- EventPointer_Bootstraps(PSI = PSI,
Design = Dmatrix,
Contrast = Cmatrix,
cores = 1,
ram = 1,
nbootstraps = 10,
UsePseudoAligBootstrap = TRUE)
##
## The Contrast and Design Matrices have been correctly designed.
##
## Number of chunks: 1
##
## Performing resampling of chunk number 1 , this may take a while...
##
## ****Analysis 100 % Completed****
##
## Preparing output...
## finish chunk 1
##
## The program has succesfully ended.
You can extract a table of the top-ranked events from the results with the function ResulTable
ResulTable(EP_Result = Fit_bots,coef = 1,number = 5)
## deltaPSI pvalue lfdr qvalues
## ENSG00000185252.17_5 -0.0001254954 0.00828407 0.2716446 0.09940884
## ENSG00000254709.7_3 0.0453065788 0.06340562 0.9251835 0.35061198
## ENSG00000185252.17_22 0.0688425350 0.08765300 0.9671472 0.35061198
## ENSG00000185252.17_6 0.0032836810 0.16558769 1.0000000 0.49670313
## ENSG00000185252.17_23 0.0688425350 0.20695964 1.0000000 0.49670313
The statistical analysis of the alternative splicing events is done in exactly the same way as for junction arrays (section 4).
The function EventPointer_RNASeq_TranRef perform this statistical
analysis and requires the following parameters:
# Design and contrast matrix:
Design <- matrix(c(1,1,1,1,0,0,1,1),nrow=4)
Contrast <- matrix(c(0,1),nrow=2)
# Statistical analysis:
Fit <- EventPointer_RNASeq_TranRef(Count_Matrix = Expression,Statistic = "LogFC",Design = Design, Contrast = Contrast)
## Done
## Analysis Finished
## Done
##
## 17:30:25 Analysis Completed
The output of this function is a data.frame with the information of the names of the event, its p.values and the corresponding z.value. If there is more than one contrast, the function returns as many data.frames as number of contrast and all these data.frame are sotred in an unique list.
| Event_ID | Pvalue | Zvalue |
|---|---|---|
| ENSG00000254709.7_3 | 0.99939 | 0.00076 |
| ENSG00000254709.7_1 | 0.96698 | 0.04140 |
| ENSG00000185252.17_22 | 0.97798 | 0.02760 |
| ENSG00000185252.17_23 | 0.97798 | 0.02760 |
| ENSG00000185252.17_6 | 0.96521 | 0.04361 |
| ENSG00000185252.17_7 | 0.70442 | 0.37936 |
| ENSG00000185252.17_8 | 0.70442 | 0.37936 |
| ENSG00000185252.17_9 | 0.70442 | 0.37936 |
| ENSG00000185252.17_25 | 0.70442 | 0.37936 |
| ENSG00000185252.17_5 | 0.63875 | 0.46945 |
| ENSG00000185252.17_21 | 0.57791 | 0.55644 |
| ENSG00000185252.17_24 | 0.44775 | -0.75917 |
EventPointer creates one GTF file that can be loaded into IGV to visualize the alternative splicing events. Figure XX displays an example result showed in IGV (5th ranked event in Table 2). Also, in the figure a reference transcriptome is displayed (blue track).
To create the GTF files, the algorithm uses the
EventPointer_RNASeq_TranRef_IGV functions with the following
parameters:
EventDetection_transcriptome function.EventDetection_transcriptome that contains the information of each
event, or table with specific events that the user want to load into
IGV to visualize.As EventDetection_transcriptome function returns the splicing graph
information, this function does not need to create the splicing graph.
SG_List <- EventXtrans$SG_List
PathEventsTxt<-system.file('extdata',package='EventPointer')
PathEventsTxt <- paste0(PathEventsTxt,"/EventsFound_Gencode24_2genes.txt")
PathGTF <- tempdir()
EventPointer_RNASeq_TranRef_IGV(SG_List = SG_List,pathtoeventstable = PathEventsTxt,PathGTF = PathGTF)
##
## Creating .GTF ...
|
| | 0%
|
|====== | 8%
|
|============ | 17%
|
|================== | 25%
|
|======================= | 33%
|
|============================= | 42%
|
|=================================== | 50%
|
|========================================= | 58%
|
|=============================================== | 67%
|
|==================================================== | 75%
|
|========================================================== | 83%
|
|================================================================ | 92%
|
|======================================================================| 100%
##
## .GTF created
The main idea of protein domain analysis is to analyze whether the
presence of a protein domain increases or decreases in the condition
under study. To do this, we have implemented the function
Protein_Domain_Enrichment. This function requires the following
parameters:
As mentioned, it is first necessary to know the structure of protein
domains within the transcriptome or in other words: to relate the
isoforms to protein domains. This can be done using bioMart R
package. The input needed by EP to perform this analysis is a matrix
that relates transcripts to protein domains (TxD matrix). Therefore, you
can work either with Pfam, Interpro or superfamily annotation.
The following conde chunck shows an example of how to get this TxD matrix:
library(biomaRt)
mart <- useMart(biomart = "ENSEMBL_MART_ENSEMBL", host = "mar2016.archive.ensembl.org")
mart<-useDataset("hsapiens_gene_ensembl",mart)
mistranscritos <- EventXtrans$transcritnames
head(mistranscritos)
## we need to remove the ".x":
mistranscritos <- gsub("\\..*","",mistranscritos)
Dominios <- getBM(attributes = c("ensembl_transcript_id","interpro","interpro_description"),
filters = "ensembl_transcript_id",
values = mistranscritos,
mart=mart)
#we build the isoform x protein domain matrix
library(Matrix)
ii <- match(Dominios$ensembl_transcript_id,mistranscritos)
misDominios <- unique(Dominios[,2])
jj <- match(Dominios[,2],misDominios)
TxD <- sparseMatrix(i=ii,j = jj,dims = c(length(mistranscritos),length(misDominios)),x = 1)
rownames(TxD) <- mistranscritos
colnames(TxD) <- misDominios
The above code chunk build the required Transcript X Domain (TxD) sparse matrix where \(TxD_{ij}\) = 1 means that Domain \(j\) match to transcript \(i\).
The next step is to perform the enrichment analysis with the function Protein_Domain_Enrichment Internally, given the TxD matrix Protein_Domain_Enrichment relates protein domains to events and performs the enrichment analysis:
NOTE: Check that EventXtrans$transcritnames annotation type match to the rownames of the TxD matrix.
data("TxD")
#check same annotation for transcripts:
EventXtrans$transcritnames[1]
# [1] "ENST00000611540.4"
rownames(TxD)[1]
# [1] "ENST00000611540"
## as is not the same, we change EventXtrans$transcritnames annotation
transcriptnames <- EventXtrans$transcritnames
transcriptnames <- gsub("\\..*","",transcriptnames)
EventXtrans$transcritnames <- transcriptnames
Result_PDEA <- Protein_Domain_Enrichment(PathsxTranscript = EventXtrans,
TxD = TxD,
Diff_PSI = Fit_bots$deltaPSI)
The function returns a list containing the results of the protein domain enrichment analysis. This list contains 3 matrices in which the rows indicate the protein domains and the columns the number of contrasts. The 3 matrices are the following:
EventPointer can also be used to design primers for PCR validation. The
aim of FindPrimers function is the design of PCR primers and TaqMan
probes for detection and quantification of alternative splicing.
Depending on the assay we want to carry out the the algorithm will
design the primers for a conventional PCR or the primers and TaqMan
probes if we are performing a TaqMan assay. In the case of a
conventional PCR we will be able to detect the alternative splicing
event. Besides, the algorithm gives as an output the length of the PCR
bands that are going to appear. In the case of a TaqMan assay, we will
not only detect but also quantify alternative splicing.
The
Primers Design step has been developed with Primer3 software and works
with versions >= 2.3.6. In order to use this option you need to install
this sofware. You can download from
sourceforge
(for Windows, Mac OSX or Unix/Linux) or from
github
(for Mac OSX or Unix/Linux).
To work with Primer Design option
you also need to add Primer3 to your environment PATH. For Windows you
can add environment variables in:
Control Panel -> System and
Security -> System -> Advanced system settings -> Environment
Variables…
For Unix/Linux or Max you need to add the path to
Primer3 to the variable PATH using terminal.
Then, as is
explained below, one of the input variables for the main function is
Primer3Path which is a string variable with the complete path where
primer3_core.exe is placed. As we have added this path to the
environment variable we can set down this variable with the following
command as is shown later in the example:
Primer3Path <- Sys.which("primer3_core")
The main function of this step is called Find_Primers and requires the
following main parameters:
Other required parameters that have a default value but can be modified by user:
The following example shows how to design primers and TaqMan probes for a specific alternative splicing event:
data("EventXtrans")
#From the output of EventsGTFfromTranscriptomeGTF we take the splicing graph information
SG_list <- EventXtrans$SG_List
#SG_list contains the information of the splicing graphs for each gene
#Let's supone we want to design primers for the event 1 of the gene ENSG00000254709.7
#We take the splicing graph information of the required gene
SG <- SG_list$ENSG00000254709.7
#We point the event number
EventNum <- 1
#Define rest of variables:
Primer3Path <- Sys.which("primer3_core")
Dir <- "C:\\PROGRA~2\\primer3\\"
MyPrimers_taqman <- FindPrimers(SG = SG,
EventNum = EventNum,
Primer3Path = Primer3Path,
Dir = Dir,
mygenomesequence = BSgenome.Hsapiens.UCSC.hg38::Hsapiens,
taqman = 1,
nProbes=1,
nPrimerstwo=4,
ncommonForward=4,
ncommonReverse=4,
nExons=10,
nPrimers =5,
maxLength = 1200)
FindPrimers return a data.frame. Either for conventional PCR and taqman option the data.frame has the following columns:
| For1Seq | For2Seq | Rev1Seq | Rev2Seq | For1Exon | For2Exon | Rev1Exon | Rev2Exon | FINALvalue | DistPath1 | DistPath2 | DistNoPath |
|---|---|---|---|---|---|---|---|---|---|---|---|
| CTGAAGGCCAATGAGACCCA | NA | CCCAGTTCCGAAGACATAACAC | NA | 2.a | NA | 4.a | NA | 2600.8 | 319 | 316 | |
| TAGGGACAGGGACCAGAGC | NA | TGAGGCTCAGACCAAAACCC | GTTTGGAGGGTTTGGTGGTC | 2.a | NA | 4.a | 5.a | 4316.8 | 435 - 662 | 432 - 659 | 540 |
| GACCAGAGCCAGTCCAGG | NA | CCCAGTTCCGAAGACATAACAC | GTTTGGAGGGTTTGGTGGTC | 2.a | NA | 4.a | 5.a | 4334.8 | 453 - 652 | 450 - 649 | 530 |
| AGAGCCAGTCCAGGGAGAG | NA | CTTGGTCCCAGTTCCGAAGA | ATTCTGTAGGGGCCACTGTC | 2.a | NA | 4.a | 5.a | 4336.8 | 455 - 780 | 452 - 777 | 658 |
| AGAGCCAGTCCAGGGAGAG | NA | TGACCTTGGTCCCAGTTCC | GTTTGGAGGGTTTGGTGGTC | 2.a | NA | 4.a | 5.a | 4340.8 | 459 - 648 | 456 - 645 | 526 |
For the tacman option, the data.frame contains three additional columns:
| For1Seq | For2Seq | Rev1Seq | Rev2Seq | For1Exon | For2Exon | Rev1Exon | Rev2Exon | FINALvalue | DistPath1 | DistPath2 | DistNoPath | SeqProbeRef | SeqProbeP1 | SeqProbeP2 |
|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|
| CTGAAGGCCAATGAGACCCA | NA | CCCAGTTCCGAAGACATAACAC | NA | 2.a | NA | 4.a | NA | 1038 | 319 | 316 | TTGGTCTGAGCCTCAGTCAC | NA | NA | |
| TAGGGACAGGGACCAGAGC | NA | TGAGGCTCAGACCAAAACCC | GTTTGGAGGGTTTGGTGGTC | 2.a | NA | 4.a | 5.a | 2870 | 435 - 662 | 432 - 659 | 540 | GAGAGCAGACCCCAGGTG | NA | NA |
| GACCAGAGCCAGTCCAGG | NA | CCCAGTTCCGAAGACATAACAC | GTTTGGAGGGTTTGGTGGTC | 2.a | NA | 4.a | 5.a | 2906 | 453 - 652 | 450 - 649 | 530 | TTGGTCTGAGCCTCAGTCAC | NA | NA |
| AGAGCCAGTCCAGGGAGAG | NA | CTTGGTCCCAGTTCCGAAGA | ATTCTGTAGGGGCCACTGTC | 2.a | NA | 4.a | 5.a | 2910 | 455 - 780 | 452 - 777 | 658 | TCTGAGCCTCAGTCACTGTG | NA | NA |
| AGAGCCAGTCCAGGGAGAG | NA | TGACCTTGGTCCCAGTTCC | GTTTGGAGGGTTTGGTGGTC | 2.a | NA | 4.a | 5.a | 2918 | 459 - 648 | 456 - 645 | 526 | TCTGAGCCTCAGTCACTGTG | NA | NA |
We can visualize the result in IGV using the Find Motif tool given the output of FindPrimers function.
Forward primers are represented in blue whilst reverse primers are represented in red as long as the gene is transcribed from the forward strand. Otherwise the colors we will be swapped. TaqMan probes will always have the color of the forward primers.
Figure 7. Transcriptome, paths and primers visualization using IGV.
EventPointer is prepared to work with different Affymetrix arrays, such as: HTA 2.0, Clariom-D, RTA and MTA.
To build the CDF file (to work under the aroma.affymetrix framework), EventPointer requires a GTF file with the reference transcriptome information. In case not provided, the algorithm automatically downloads the required information from different databases such as ENSEMBL or UCSC. As the probes for the HTA 2.0 array are mapped to the HG19 genome, the latests versions of the ensembl and ucsc genome, mapped to the hg19 version, will be downloaded. For the other arrays, the following genomes are used: ClariomD = GRCh38, MTA = mm10 and RTA = rn6.
The required files are:
Files 1 and 2 contain probe information for the array. These files and the corresponding CDF files can be downloaded from the following URL: EventPointer Dropbox
The format of these files is briefly explained in the following paragraphs:
For the Exon Probes, the file corresponds to a tab separated .txt file composed of 11 columns that include: Probe Id, X coordinate in the array, Y coordinate in the array, Transcript Cluster Id, Probeset Id, Probeset name, Probe sequence, chromosome, start, end and strand.
The Junction probes file is also a tab separated .txt composed of 10 columns: Probe Id, X coordinate in the array, Y coordinate in the array, Transcript Cluster Id, Probeset Id, Probeset name, Probe sequence, chromosome and probe alignments.
For the GTF the standard format is used. (For more information see GTF)
This vignette uses the probes and annotation file for the DONSON gene in order to run the examples in a short amount of time. To generate the corresponding CDF file for the whole genome, the files from the Dropbox link must be used.
Note: It is advisable to work with reference GTF files, as the probes are annotated to them. However, if other database is used, the algorithm will only include probes that are mapped to such database.
This step can be skipped if the corresponding CDF file is doownloaded from the Dropbox link. The available CDF files were generated using the GTF files for each of the arrays, if users want to generate a CDF for other databases (Ensembl or UCSC), this step should be used.
The function CDFfromGTF generates the CDF file used afterwards in the aroma.affymetrix pre-processing pipeline. Internally, it calls flat2cdf function written by Elizabeth Purdom. More information about it can be seen in the following webpage: Create CDF from scratch
The required input for the function is described below:
This function takes a couple of hours to complete (depending on the computer), and creates the following files:
The following code chunks show examples on how to run the function using part of the Affymetrix GTF and the example data included in the package. This data corresponds to the Affymetrix HTA 2.0 GTF representing only the DONSON gene and the probes that are mapped to it.
Using Affymetrix GTF as reference transcriptome
# Set input variables
PathFiles<-system.file("extdata",package="EventPointer")
DONSON_GTF<-paste(PathFiles,"/DONSON.gtf",sep="")
PSRProbes<-paste(PathFiles,"/PSR_Probes.txt",sep="")
JunctionProbes<-paste(PathFiles,"/Junction_Probes.txt",sep="")
Directory<-tempdir()
array<-"HTA-2_0"
# Run the function
CDFfromGTF(input="AffyGTF",inputFile=DONSON_GTF,
PSR=PSRProbes,Junc=JunctionProbes,
PathCDF=Directory,microarray=array)
Note: Both the .flat and .CDF file take large amounts of space in the hard drive, it is recommended that the directory has at least 1.5 GB of free space.
Figure 3 shows a screen shot with the output of the CDFfromGTF function
Figure 9. Output of CDFfromGTF
Once the file is created, the name of the cdf file can be changed. Internally, the algorithm gives the name as HTA-2_0, but the actual name of the file can be different. In the rest of the vignette, we have renamed our CDF file as EP_HTA-2_0.
Once the CDF file has been created, it can be used for the analysis of different experiments.
For microarray data, a pre-processing pipeline must be applied to the .cel files of the experiment. Usually this pre-processing is done using the aroma.affymetrix R package. This procedure normalizes and summarizes the expression of the different probesets into single values.
The aroma.affymetrix R package provides users different functions to work with affymetrix arrays. The functions are used to extract all the information contained in the .cel files and to do all the required pre-processing steps such as background correction, normalization and summarization. The package requires a certain folder structure in order to work correctly. For more information about aroma.affymetrix visit the webpage:aroma project
The following code chunk displays the pipeline used to get the results required by EventPointer after the pre-processing using aroma.affymetrix. The code is not intended to be a runable example, but just to show users the settings and functions that can be used. In order for users to have a runable example, the corrrespoding folder structure and files are required.
# Simple example of Aroma.Affymetrix Preprocessing Pipeline
verbose <- Arguments$getVerbose(-8);
timestampOn(verbose);
projectName <- "Experiment"
cdfGFile <- "EP_HTA-2_0,r"
cdfG <- AffymetrixCdfFile$byChipType(cdfGFile)
cs <- AffymetrixCelSet$byName(projectName, cdf=cdfG)
bc <- NormExpBackgroundCorrection(cs, method="mle", tag=c("*","r11"));
csBC <- process(bc,verbose=verbose,ram=20);
qn <- QuantileNormalization(csBC, typesToUpdate="pm");
csN <- process(qn,verbose=verbose,ram=20);
plmEx <- ExonRmaPlm(csN, mergeGroups=FALSE)
fit(plmEx, verbose=verbose)
cesEx <- getChipEffectSet(plmEx)
ExFit <- extractDataFrame(cesEx, addNames = TRUE)
EventPointer is the main function of the algorithm
The function requires the following parameters:
If the Filter variable is TRUE, the following is performed:
The summarized expression value of all the reference paths is obtained
and the threshold is set depending on the Qn value used.
An event is considered if at least one sample in all paths is expressed above the threshold.
The rest of the events are not shown unless the Filter variable is set
to FALSE
data(ArraysData)
Dmatrix<-matrix(c(1,1,1,1,0,0,1,1),nrow=4,ncol=2,byrow=FALSE)
Cmatrix<-t(t(c(0,1)))
EventsFound<-paste(system.file("extdata",package="EventPointer"),"/EventsFound.txt",sep="")
Events<-EventPointer(Design=Dmatrix,
Contrast=Cmatrix,
ExFit=ArraysData,
Eventstxt=EventsFound,
Filter=FALSE,
Qn=0.25,
Statistic="LogFC",
PSI=TRUE)
## 17:30:26 Running EventPointer:
## ** Statistical Analysis: Logarithm of the fold change of both isoforms
## ** Delta PSI will be calculated
## ** No expression filter
## ----------------------------------------------------------------
## ** Calculating PSI...Done
## ** Running Statistical Analysis...Done
##
## 17:30:26 Analysis Completed
Table 1 displays the output of EventPointer function
| Gene name | Event Type | Genomic Position | Splicing Z Value | Splicing Pvalue | Delta PSI | |
|---|---|---|---|---|---|---|
| TC21001058.hg_8 | TC21001058.hg | Alternative 3’ Splice Site | 21:34957032-34958284 | 6.86631 | 0.0000 | 0.00000 |
| TC21001058.hg_6 | TC21001058.hg | Complex Event | 21:34950750-34953608 | 6.33338 | 0.0000 | -0.09861 |
| TC21001058.hg_9 | TC21001058.hg | Alternative Last Exon | 21:34951868-34956896 | 6.08929 | 0.0000 | -0.44545 |
| TC21001058.hg_10 | TC21001058.hg | Complex Event | 21:34955972-34958284 | -5.03597 | 0.0000 | 0.04857 |
| TC21001058.hg_4 | TC21001058.hg | Complex Event | 21:34955972-34958284 | 1.43180 | 0.1522 | 0.00000 |
The columns of the data.frame are:
EventPointer creates two different GTF files to visualize the alternative splicing events. Figure 4 displays the cassette exon for the COPS7A gene (4th ranked event in Table 1). In the IGV visualization, the reference path is shown in gray color, the path 1 in red and path 2 in green. Below the transcripts, the different probes that are measuring each of the paths are represented in the same colors.
Figure 10. EventPointer visualization using IGV
To create the GTF files, the algorithm uses the EventPointer_IGV functions with the following parameters:
The inital process of the function is slow, as the splicing graphs must be created every time the function is executed. A progress bar is displayed to the user to inform about the progress of the function.
Once the process is completed two GTF files are generated in the specified directory:
Figure 11. GTF files generated with EventPointer_IGV
# Set Input Variables
DONSON_GTF<-paste(PathFiles,"/DONSON.gtf",sep="")
PSRProbes<-paste(PathFiles,"/PSR_Probes.txt",sep="")
JunctionProbes<-paste(PathFiles,"/Junction_Probes.txt",sep="")
Directory<-tempdir()
EventsFound<-paste(system.file("extdata",package="EventPointer"),"/EventsFound.txt",sep="")
array<-"HTA-2_0"
# Generate Visualization files
EventPointer_IGV(Events[1,,drop=FALSE],"AffyGTF",DONSON_GTF,PSRProbes,JunctionProbes,Directory,EventsFound,array)
EventPointer can also identify Multi-Path events. The multi-path events are composed of more elements than a triplet where the concentration of the reference path should be equal to the sum of the concentration of the rest of paths. This is depicted in Figure 8.
Figure 12. Multi-Path alternative splicing events for EventPointer
EventPointer identify Multi-Path events and estimate the percent spliced in value \(\Psi\).
The function CDFfromGTF_Multipath generates the CDF file used afterwards in the aroma.affymetrix preprocesing pipeline. This function is equivalent to the function CDFfromGTF.
The function requires the following parameters:
The function CDFfromGTF_Multipath detects events with number of paths from 2 to variable paths. This function classifies the events with two paths as the CDFfromGTF function does and the events with more than two paths as “Multi-Path”.
This function takes a couple of hours to complete (depending on the computer), and creates the same files as the function CDFfromGTF.
The following code chunks show examples on how to run the function using part of the Affymetrix GTF and the example data included in the package. his data corresponds to the Affymetrix HTA 2.0 GTF representing only the DONSON gene and the probes that are mapped to it.
Using Affymetrix GTF as reference transcriptome
# Set input variables
PathFiles<-system.file("extdata",package="EventPointer")
DONSON_GTF<-paste(PathFiles,"/DONSON.gtf",sep="")
PSRProbes<-paste(PathFiles,"/PSR_Probes.txt",sep="")
JunctionProbes<-paste(PathFiles,"/Junction_Probes.txt",sep="")
Directory<-tempdir()
array<-"HTA-2_0"
# Run the function
CDFfromGTF_Multipath(input="AffyGTF",inputFile=DONSON_GTF,
PSR=PSRProbes,Junc=JunctionProbes,
PathCDF=Directory,microarray=array,paths=3)
##RNA-Seq (Event detection for Multi-Path)
As for two-paths events, the only required files are the corresponding .BAM files from the experiment. After the BAM Preparation step explained in section 4.2 we have all the elements needed to analyze each of the splicing graphs. To detect Multi-Path events, the function EventDetectionMultipath is used.
This function requires the following parameters:
The function EventDetectionMultipath detects events with number of paths from 2 to variable paths. This function classifies the events with two paths as the EventDetection function does and the events with more than two paths as “Multi-Path”.
# Run EventDetection function
data(SG_RNASeq)
TxtPath<-tempdir()
AllEvents_RNASeq_MP<-EventDetectionMultipath(SG_RNASeq,cores=1,Path=TxtPath,paths=3)
This function retireves a list with all the events found for all the genes present in the experiment. It also generates a file called EventsFound_RNASeq.txt with the information for every detected event.
The list is organized in the following way:
Events[[i]][[j]]
The list will have as many \(i\) values as genes and \(j\) values as many events detected for the \(i_{th}\) gene. In other words, the command above will display the \(j_{th}\) event detected for the \(i_{th}\) gene.
As described above, EventPointer is able to detect splicing events that
cannot be classified as canonical splicing events. EventPointer
classifies these events as “Complex Event”. Within this category one can
observe simple events such as ‘Multiple Skipping Exons’ and events whose
structure is very complex. Thus, EP is classifying in the same category
events that have nothing to do with each other.
EventPointer
includes a function called Events_ReClassification to reclassify
complex events. The input to this function is:
Thus, EventPointer reclassifies the complex events according to how
similar the event is to the struture canonical events. The same complex
event can have several types (see). Further, EP adds a new type of
event: “multiple skipping exon”. These events are characterized by
presenting several exons in a row as alternative exons. I.e., If there
is only one alternative exon we would be talking about a “Casstte
Exon”.
The following code shows how this function works for a subset of 5 splicing events:
#load splicing graph
data("SG_reclassify")
#load table with info of the events
PathFiles<-system.file("extdata",package="EventPointer")
inputFile <- paste(PathFiles,"/Events_found_class.txt",sep="")
EventTable <- read.delim(file=inputFile)
#this table has the information of 5 complex events.
EventTable_new <- Events_ReClassification(EventTable = EventTable,
SplicingGraph = SG_reclassify)
##
## Starting reclassification
##
|
| | 0%
|
|============== | 20%
|
|============================ | 40%
|
|========================================== | 60%
|
|======================================================== | 80%
|
|======================================================================| 100%
##
## Reclassification Finished
The following figures show the structure of these events and their corresponding gene:
Example 1
EventPointer provides three different statistical tests that can be used to determine the statistical significance to the alternative splicing events.
In Table 5, the most relevatn coefficients from the statistical and events information point of view are shown
Table 6. Coefficients of the linear model used to detect alternative splicing
There are a number of alternatives to detect AS using these coefficients. Either of \(\beta_4\) , \(\beta_5\) , \(\beta_4\) + \(\beta_5\) is theoretically able to detect AS events. Some of them are more sensitive than others depending on the relative concentrations of the isoforms. For example, if isoform 2 is much more highly expressed than isoform 1 in both conditions, \(\beta_4\) will be more sensitive than \(\beta_4\) + \(\beta_5\) since in the latter case, the numerator and denominator of the logarithms of both terms are similar, and hence their logs are close to zero.
A contrast on \(\beta_5\) would seem to be more sensitive than a contrast on \(\beta_4\) or \(\beta_4\) + \(\beta_5\) ; however, in practice, this contrast proved to be “too sensitive” and led to many false positives especially in weakly expressed isoforms. If one of the paths is not expressed in any condition, its signal will be similar in either condition (the background level) and a change in the expression of the other isoform will drive to a false positive detection. This contrast can be used only if the signals are filtered to guarantee that they are above the background.
In the PCR validation, the contrast that provided the best performance was the combination of the fold changes of both isoforms (i.e. \(\beta_3\) + \(\beta_4\) and \(\beta_3\)+\(\beta_4\) + \(\beta_5\) in Table 4) plus the requirement that the fold-changes have opposite directions, i.e. if isoform 1 significantly increases its expression, isoform 2 must significantly decrease its expression and vice versa. Therefore, if this test requires that the detected AS events show a significant change of the expression both paths and this change must be in opposite direction.
In order to compute this contrast, we summed up the p-values (one-tailed) for both contrasts. If the null hypothesis holds, the expected null distribution is triangular from 0 to 2 with the peak at 1, and the summation of the p-values must be close to 0 or close to 2 for genes with differential AS. Using this triangular distribution, it is possible to assign an overall p-value to their sum. We preferred this combination rather than the classical Fisher method since in the latter a single good p-value yields a good summary p-value for the event. Using this approach, both p-values must be close to zero or one in order to generate a significant overall p-value.
The different options availabe in EventPointer are:
1. LogFC : Compute the contrast using \(\beta_3\) + \(\beta_4\) and \(\beta_3\)+\(\beta_4\) + \(\beta_5\)
2. Dif_LogFC : Compute the contrast using \(\beta_4\) and \(\beta_4\) + \(\beta_5\)
3. DRS: Compute the constast using \(\beta_5\)
Alpha is a parameter used by SGSeq R package to predict the features that are along the different bam files that are being analyzed. As stated in the help menu for the predictTxtFeatures function:
Alpha Minimum FPKM required for a splice junction to be included.
The user can change the value to be more or less restrictive when deciding if a feature is included or not. As the alpha value increases, the algorithm will slow down as the splicing graphs would became more complex.
EventPointer estimates the abundance of the isoforms mapped to each of the paths, in an splicing event, to obtain the PSI values. With this values, a simple linear model, using the provided design and contrast matrices, is solved and this increment is returned to the user in the data.frame (if PSI argument is TRUE).
It is possible to obtain the estimated PSI values using the internal functions getPSI,getPSImultipath for junction arrays, getPSI_RNASeq or getPSI_RNASeq_MultiPath for RNA-Seq data (data from the pipeline described in section 4.2: statistical analysis from de BAM files).
These functions not only calculate the value of \(\Psi\) but also the residuals of the simple linear model used to estimate the values of \(\Psi\).
# Microarrays (two paths)
data(ArraysData)
PSI_Arrays_list<-EventPointer:::getPSI(ArraysData)
PSI_Arrays <- PSI_Arrays_list$PSI
Residuals_Arrays <- PSI_Arrays_list$Residuals
# Microarrays (Multi-Path)
data(ArrayDatamultipath)
PSI_MP_Arrays_list <- EventPointer:::getPSImultipath(ArrayDatamultipath)
PSI_multipath_Arrays <- PSI_MP_Arrays_list$PSI
Residuals_multipath_Arrays <- PSI_MP_Arrays_list$Residuals
# RNASeq (two paths)
data(AllEvents_RNASeq)
PSI_RNASeq_list<-EventPointer:::getPSI_RNASeq(AllEvents_RNASeq)
PSI_RNASeq <- PSI_RNASeq_list$PSI
Residuals_RNASeq <- PSI_RNASeq_list$Residuals
# RNASeq (Multi-Path)
data(AllEvents_RNASeq_MP)
PSI_MP_RNASeq_list <- EventPointer:::getPSI_RNASeq_MultiPath(AllEvents_RNASeq_MP)
PSI_multipath_RNASeq <- PSI_MP_RNASeq_list$PSI
Residuals_multipath_RNASeq <- PSI_MP_RNASeq_list$Residuals
We can apply the function PSI_Statistic to the values of \(\Psi\) (only for two-paths events). This function takes as input the values of PSI and performs a statistical analysis based on permutation test.
Dmatrix<-matrix(c(1,1,1,1,0,0,1,1),nrow=4,ncol=2,byrow=FALSE)
Cmatrix<-t(c(0,1))
table <- PSI_Statistic(PSI = PSI_Arrays,Design = Dmatrix,Contrast = Cmatrix,nboot = 20)
The residual obtained in the linear model used to estimate the values of \(\Psi\) must be independent from the Design matrix of the experiment. The data of the residuals that returns the internal functions getPSI or getPSI_RNASeq are useful to validate if this is true. The next code shows an example of how to perform a statistical analysis of the residuals.
Dmatrix<-matrix(c(1,1,1,1,0,0,1,1),nrow=4,ncol=2,byrow=FALSE)
Cmatrix<-t(t(c(0,1)))
Ress <- vector("list", length = ncol(Cmatrix))
fitresiduals <- limma::lmFit(Residuals_Arrays,design = Dmatrix)
fitresiduals2 <- limma::contrasts.fit(fitresiduals, Cmatrix)
fitresiduals2 <- limma::eBayes(fitresiduals2)
# repeated if there is more than one contrast
for (jj in 1:ncol(Cmatrix)) {
TopPSI <- limma::topTable(fitresiduals2, coef = jj, number = Inf)[, 1, drop = FALSE]
colnames(TopPSI)<-"Residuals"
Ress[[jj]] <- TopPSI
}
sessionInfo()
## R Under development (unstable) (2025-12-22 r89219)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.23-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB 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
##
## time zone: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] kableExtra_1.4.0 dplyr_1.1.4
## [3] EventPointer_3.19.1 Matrix_1.7-4
## [5] SGSeq_1.45.0 SummarizedExperiment_1.41.0
## [7] Biobase_2.71.0 MatrixGenerics_1.23.0
## [9] matrixStats_1.5.0 Rsamtools_2.27.0
## [11] Biostrings_2.79.4 XVector_0.51.0
## [13] GenomicRanges_1.63.1 Seqinfo_1.1.0
## [15] IRanges_2.45.0 S4Vectors_0.49.0
## [17] BiocGenerics_0.57.0 generics_0.1.4
## [19] knitr_1.51 BiocStyle_2.39.0
##
## loaded via a namespace (and not attached):
## [1] RColorBrewer_1.1-3 rstudioapi_0.17.1 jsonlite_2.0.0
## [4] shape_1.4.6.1 tximport_1.39.1 magrittr_2.0.4
## [7] GenomicFeatures_1.63.1 farver_2.1.2 rmarkdown_2.30
## [10] BiocIO_1.21.0 vctrs_0.6.5 memoise_2.0.1
## [13] RCurl_1.98-1.17 progress_1.2.3 htmltools_0.5.9
## [16] S4Arrays_1.11.1 curl_7.0.0 Rhdf5lib_1.33.0
## [19] SparseArray_1.11.10 rhdf5_2.55.12 sass_0.4.10
## [22] parallelly_1.46.0 bslib_0.9.0 httr2_1.2.2
## [25] plyr_1.8.9 poibin_1.6 cachem_1.1.0
## [28] GenomicAlignments_1.47.0 igraph_2.2.1 lifecycle_1.0.4
## [31] iterators_1.0.14 pkgconfig_2.0.3 R6_2.6.1
## [34] fastmap_1.2.0 future_1.68.0 digest_0.6.39
## [37] AnnotationDbi_1.73.0 textshaping_1.0.4 RSQLite_2.4.5
## [40] filelock_1.0.3 nnls_1.6 httr_1.4.7
## [43] abind_1.4-8 compiler_4.6.0 bit64_4.6.0-1
## [46] doParallel_1.0.17 speedglm_0.3-5 S7_0.2.1
## [49] BiocParallel_1.45.0 DBI_1.2.3 biglm_0.9-3
## [52] R.utils_2.13.0 biomaRt_2.67.0 MASS_7.3-65
## [55] lava_1.8.2 quantreg_6.1 rappdirs_0.3.3
## [58] DelayedArray_0.37.0 rjson_0.2.23 affxparser_1.83.0
## [61] tools_4.6.0 otel_0.2.0 future.apply_1.20.1
## [64] R.oo_1.27.1 glue_1.8.0 restfulr_0.0.16
## [67] rhdf5filters_1.23.3 grid_4.6.0 reshape2_1.4.5
## [70] fgsea_1.37.4 lpSolve_5.6.23 gtable_0.3.6
## [73] BSgenome_1.79.1 R.methodsS3_1.8.2 hms_1.1.4
## [76] data.table_1.18.0 xml2_1.5.1 stringr_1.6.0
## [79] foreach_1.5.2 pillar_1.11.1 limma_3.67.0
## [82] splines_4.6.0 BiocFileCache_3.1.0 lattice_0.22-7
## [85] aroma.light_3.41.0 survival_3.8-3 rtracklayer_1.71.3
## [88] bit_4.6.0 SparseM_1.84-2 tidyselect_1.2.1
## [91] RBGL_1.87.0 bookdown_0.46 svglite_2.2.2
## [94] xfun_0.55 statmod_1.5.1 stringi_1.8.7
## [97] UCSC.utils_1.7.1 yaml_2.3.12 evaluate_1.0.5
## [100] codetools_0.2-20 cobs_1.3-9-1 cigarillo_1.1.0
## [103] tibble_3.3.0 qvalue_2.43.0 BiocManager_1.30.27
## [106] graph_1.89.1 cli_3.6.5 systemfonts_1.3.1
## [109] jquerylib_0.1.4 dichromat_2.0-0.1 Rcpp_1.1.0.8.2
## [112] GenomeInfoDb_1.47.2 globals_0.18.0 dbplyr_2.5.1
## [115] png_0.1-8 XML_3.99-0.20 RUnit_0.4.33.1
## [118] parallel_4.6.0 MatrixModels_0.5-4 ggplot2_4.0.1
## [121] blob_1.2.4 prettyunits_1.2.0 bitops_1.0-9
## [124] txdbmaker_1.7.3 listenv_0.10.0 glmnet_4.1-10
## [127] viridisLite_0.4.2 scales_1.4.0 prodlim_2025.04.28
## [130] crayon_1.5.3 rlang_1.1.6 cowplot_1.2.0
## [133] fastmatch_1.1-6 KEGGREST_1.51.1