Contents

1 About APAlyzer

APAlyzer is a toolkit for bioinformatic analysis of alternative polyadenylation (APA) events using RNA sequencing data. Our main approach is comparison of sequencing reads in regions demarcated by high quality polyadenylation sites (PASs) annotated in the PolyA_DB database (https://exon.apps.wistar.org/PolyA_DB/v3/) (Wang et al. 2017, 2018). The current version uses RNA-seq data to examine APA events in 3’ untranslated regions (3’UTRs) and in introns. The coding regions are used for gene expression calculation.

2 Program installation

APAlyzer should be installed using BiocManager:

if (!"BiocManager" %in% rownames(installed.packages()))
     install.packages("BiocManager")
BiocManager::install("APAlyzer")

Alternatively, it can also be installed as follows:

R CMD INSTALL APAlyzer.tar.gz

In additions, user can also install development version of APAlyzer directly from GitHub:

BiocManager::install('RJWANGbioinfo/APAlyzer')

After installation, APAlyzer can be used by:

library(APAlyzer)

3 Sample data and PAS references

3.1 RNA-seq BAM files

The package reads BAM file(s) to obtain read coverage information in different genomic regions. The following example shows that we first specify paths to example BAM files in the TBX20BamSubset (Bindreither 2018) data package. In this example, BAM files correspond to mouse RNA-seq data, (mapped to mm9).

suppressMessages(library("TBX20BamSubset"))
suppressMessages(library("Rsamtools"))
flsall = getBamFileList()
flsall
##                                                                           SRR316184 
## "/home/biocbuild/bbs-3.17-bioc/R/site-library/TBX20BamSubset/extdata/SRR316184.bam" 
##                                                                           SRR316185 
## "/home/biocbuild/bbs-3.17-bioc/R/site-library/TBX20BamSubset/extdata/SRR316185.bam" 
##                                                                           SRR316186 
## "/home/biocbuild/bbs-3.17-bioc/R/site-library/TBX20BamSubset/extdata/SRR316186.bam" 
##                                                                           SRR316187 
## "/home/biocbuild/bbs-3.17-bioc/R/site-library/TBX20BamSubset/extdata/SRR316187.bam" 
##                                                                           SRR316188 
## "/home/biocbuild/bbs-3.17-bioc/R/site-library/TBX20BamSubset/extdata/SRR316188.bam" 
##                                                                           SRR316189 
## "/home/biocbuild/bbs-3.17-bioc/R/site-library/TBX20BamSubset/extdata/SRR316189.bam"

3.2 Genomic reference

PAS references in the genome (both 3’UTRs and introns) are required by our package. We have pre-built a reference file for the mouse genome (mm9), which can be loaded from extdata:

library(repmis)
URL="https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/"
file="mm9_REF.RData"
source_data(paste0(URL,file,"?raw=True"))
## Downloading data from: https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/mm9_REF.RData?raw=True
## SHA-1 hash of the downloaded data file is:
## fa117a77976a6931e4755a3182544df1d1de7eb6
## [1] "refUTRraw" "dfIPA"     "dfLE"

This extdata covers 3’UTR APA regions (refUTRraw), IPA regions (dfIPA), and 3’-most exon regions (dfLE). The refUTRraw is a data frame containing 6 columns for genomic information of 3’UTR PASs:

head(refUTRraw,2)
##   gene_symbol Chrom Strand   First    Last  cdsend
## 1        Xkr4  chr1      - 3204561 3189814 3206103
## 2      Lypla1  chr1      + 4835237 4836816 4835097

dfIPA is a data frame containing 8 columns for Intronic PASs; ‘upstreamSS’ means the closest 5’ or 3’ splice site to IPA, ‘downstreamSS’ means closest 3’ splice site:

head(dfIPA,2)
##              PASid gene_symbol Chrom Strand      Pos upstreamSS downstreamSS
## 5  chr4:32818455:+        Mdn1  chr4      + 32818455   32818305     32818889
## 10 chr4:32837613:+        Mdn1  chr4      + 32837613   32837452     32837887

dfLE is a data frame containing 5 colmuns for 3’ least exon; ‘LEstart’ means the start genomic position of last 3’ exon.

head(dfLE,2)
##     gene_symbol Chrom Strand  LEstart      TES
## 1 0610009O20Rik chr18      + 38421534 38425118
## 2 0610010F05Rik chr11      - 23465136 23462079

In additions to mouse mm9, our package has also a pre-build version for mouse mm10, human hg38 and human hg19 genome:

URL="https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/"
file="hg19_REF.RData"
source_data(paste0(URL,file,"?raw=True"))
## Downloading data from: https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/hg19_REF.RData?raw=True
## SHA-1 hash of the downloaded data file is:
## 7ba8e91beba626318e48736692d6442480d1ef06
## [1] "refUTRraw_hg19" "dfIPA_hg19"     "dfLE_hg19"

More pre-build refercence can be found at the reference and testing data repo: https://github.com/RJWANGbioinfo/PAS_reference_RData_and_testing_data

3.3 Building 3’UTR and intronic PAS reference region at once

To quantify the relative expression of PAS, we will need to build the reference regions for them, although this can be build separately in previous version. We also provide a new fouction REF4PAS starting from APAlyzer 1.2.1 to build these regions at once:

refUTRraw=refUTRraw[which(refUTRraw$Chrom=='chr19'),]
dfIPAraw=dfIPA[which(dfIPA$Chrom=='chr19'),]
dfLEraw=dfLE[which(dfLE$Chrom=='chr19'),]   
PASREF=REF4PAS(refUTRraw,dfIPAraw,dfLEraw)
UTRdbraw=PASREF$UTRdbraw
dfIPA=PASREF$dfIPA
dfLE=PASREF$dfLE    

In this case, UTRdbraw is the reference region used for 3’UTR APA analysis, while dfIPA and dfLE are needed in the intronic APA analysis.

3.4 Building 3’UTR PAS and IPA reference using GTF files

Although we highly suggest user use references regions genrated from PolyA_DB. Start from APAlyzer 1.2.1, we also provide a new fouction that can help users to build their reference directly from gene annotation GTF files, we hope this can help the species which are not covered by the PolyA_DB yet:

## build Reference ranges for 3'UTR PASs in mouse
download.file(url='ftp://ftp.ensembl.org/pub/release-99/gtf/mus_musculus/Mus_musculus.GRCm38.99.gtf.gz',
          destfile='Mus_musculus.GRCm38.99.gtf.gz')           
GTFfile="Mus_musculus.GRCm38.99.gtf.gz" 
PASREFraw=PAS2GEF(GTFfile)  
refUTRraw=PASREFraw$refUTRraw
dfIPAraw=PASREFraw$dfIPA
dfLEraw=PASREFraw$dfLE
PASREF=REF4PAS(refUTRraw,dfIPAraw,dfLEraw)

Start from APAlyzer 1.9.1, two annotation methods can be set through AnnoMethod= to annotate the PAS, ‘legacy’ or ‘V2’ (default), ‘legacy’ is the old method used within previous versions, while ‘V2’ is the updated and improved version.

4 Analysis of APA in 3’UTRs

4.1 Building aUTR and cUTR references

To calculate 3’UTR APA relative expression (RE), we first need to define the refence regions of aUTR and cUTR using refUTRraw. Since the sample data only contains mapping information on chr19, we can zoom into reference regions on chr19 only:

refUTRraw=refUTRraw[which(refUTRraw$Chrom=='chr19'),]
UTRdbraw=REF3UTR(refUTRraw)

The REF3UTR function returns a genomic range containing aUTR(pPAS to dPAS) and cUTR(cdsend to pPAS) regions for each gene:

head(UTRdbraw,2)
## GRangesList object of length 2:
## $`1700020D05Rik__aUTR`
## GRanges object with 1 range and 2 metadata columns:
##         seqnames          ranges strand |   gene_symbol           gene_name
##            <Rle>       <IRanges>  <Rle> |      <factor>         <character>
##   10082    chr19 5492231-5495277      - | 1700020D05Rik 1700020D05Rik__aUTR
##   -------
##   seqinfo: 21 sequences from an unspecified genome; no seqlengths
## 
## $`1700020D05Rik__cUTR`
## GRanges object with 1 range and 2 metadata columns:
##          seqnames          ranges strand |   gene_symbol           gene_name
##             <Rle>       <IRanges>  <Rle> |      <factor>         <character>
##   100821    chr19 5495277-5502867      - | 1700020D05Rik 1700020D05Rik__cUTR
##   -------
##   seqinfo: 21 sequences from an unspecified genome; no seqlengths

4.2 Calculation of relative expression

Once cUTR and aUTR regions are defined, the RE of 3’UTR APA of each gene can be calculated by PASEXP_3UTR:

DFUTRraw=PASEXP_3UTR(UTRdbraw, flsall, Strandtype="forward")
## [1] "SRR316184, Strand: forward, finished"
## [1] "SRR316185, Strand: forward, finished"
## [1] "SRR316186, Strand: forward, finished"
## [1] "SRR316187, Strand: forward, finished"
## [1] "SRR316188, Strand: forward, finished"
## [1] "SRR316189, Strand: forward, finished"

The PASEXP_3UTR 3UTR requires two inputs: 1) aUTR and cUTR reference regions, and 2) path of BAM file(s). In additions to input, one can also define the strand of the sequencing using Strandtype. The detailed usage can also be obtained by the command ?PASEXP_3UTR.The output data frame covers reads count (in aUTR or cUTR), RPKM (in aUTR or cUTR) and RE (log2(aUTR/cUTR)) for each gene:

head(DFUTRraw,2)
##     gene_symbol SRR316184_areads SRR316184_aRPKM SRR316184_creads
## 1 1700020D05Rik                0               0               31
## 2 1810009A15Rik                0               0                0
##   SRR316184_cRPKM SRR316184_3UTR_RE SRR316185_areads SRR316185_aRPKM
## 1        118.4358              -Inf                0               0
## 2          0.0000               NaN                0               0
##   SRR316185_creads SRR316185_cRPKM SRR316185_3UTR_RE SRR316186_areads
## 1               20        74.34672              -Inf                0
## 2                0         0.00000               NaN                0
##   SRR316186_aRPKM SRR316186_creads SRR316186_cRPKM SRR316186_3UTR_RE
## 1               0               23        87.90994              -Inf
## 2               0                0         0.00000               NaN
##   SRR316187_areads SRR316187_aRPKM SRR316187_creads SRR316187_cRPKM
## 1                0               0               33        88.85001
## 2                0               0                0         0.00000
##   SRR316187_3UTR_RE SRR316188_areads SRR316188_aRPKM SRR316188_creads
## 1              -Inf                0               0               33
## 2               NaN                0               0                0
##   SRR316188_cRPKM SRR316188_3UTR_RE SRR316189_areads SRR316189_aRPKM
## 1        77.68085              -Inf                0               0
## 2         0.00000               NaN                0               0
##   SRR316189_creads SRR316189_cRPKM SRR316189_3UTR_RE
## 1               11        31.76981              -Inf
## 2                0         0.00000               NaN

5 Analysis of APA in introns

5.1 Building intronic polyA references

Analysis of IPA requires two genomic regions: IPA regions and 3’-most exons. As mentioned above, these regions in mouse and human genomes have been pre-built in the package:

#mouse(mm9):
URL="https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/"
file="mm9_REF.RData"
source_data(paste0(URL,file,"?raw=True"))
## Downloading data from: https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/mm9_REF.RData?raw=True
## SHA-1 hash of the downloaded data file is:
## fa117a77976a6931e4755a3182544df1d1de7eb6
## [1] "refUTRraw" "dfIPA"     "dfLE"
#human(hg19):
URL="https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/"
file="hg19_REF.RData"
source_data(paste0(URL,file,"?raw=True"))
## Downloading data from: https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/hg19_REF.RData?raw=True
## SHA-1 hash of the downloaded data file is:
## 7ba8e91beba626318e48736692d6442480d1ef06
## [1] "refUTRraw_hg19" "dfIPA_hg19"     "dfLE_hg19"

5.2 Calculation of relative expression

Similar to 3’UTR APA, RE of IPAs can be calculated using PASEXP_IPA: PASEXP_IPA:

dfIPA=dfIPA[which(dfIPA$Chrom=='chr19'),]
dfLE=dfLE[which(dfLE$Chrom=='chr19'),]
IPA_OUTraw=PASEXP_IPA(dfIPA, dfLE, flsall, Strandtype="forward", nts=1)

Note that, as a specific feature for IPA, one can set more threads using ‘nts=’(the parameter passed to Rsubread::featureCounts, check ?Rsubread::featureCounts for details) to increase calculation speed. The detailed usage can be obtained by the command ?PASEXP_IPA.

The output data frame contains read count and read density IPA upstream (a), IPA downstream (b) and 3’-most exon region (c). The RE of IPA is calculated as log2((a - b)/c).

head(IPA_OUTraw,2)
##     gene_symbol                         IPAID            PASid
## 1 9130011E15Rik chr19:45894071:-9130011E15Rik chr19:45894071:-
## 2 9130011E15Rik chr19:45952090:-9130011E15Rik chr19:45952090:-
##   SRR316184_IPA_UPreads SRR316184_IPA_UPRPK SRR316184_IPA_DNreads
## 1                     1            1.445087                     0
## 2                   104           12.453598                    26
##   SRR316184_IPA_DNRPK SRR316184_LEreads SRR316184_LERPK SRR316184_IPA_RE
## 1           0.0000000                96        59.47955        -5.363166
## 2           0.8276828                96        59.47955        -2.355049
##   SRR316185_IPA_UPreads SRR316185_IPA_UPRPK SRR316185_IPA_DNreads
## 1                     0             0.00000                     1
## 2                   121            14.48928                    20
##   SRR316185_IPA_DNRPK SRR316185_LEreads SRR316185_LERPK SRR316185_IPA_RE
## 1           1.6611296                96        59.47955         0.000000
## 2           0.6366791                96        59.47955        -2.102237
##   SRR316186_IPA_UPreads SRR316186_IPA_UPRPK SRR316186_IPA_DNreads
## 1                     0             0.00000                     0
## 2                   108            12.93258                    13
##   SRR316186_IPA_DNRPK SRR316186_LEreads SRR316186_LERPK SRR316186_IPA_RE
## 1           0.0000000                69        42.75093         0.000000
## 2           0.4138414                69        42.75093        -1.771866
##   SRR316187_IPA_UPreads SRR316187_IPA_UPRPK SRR316187_IPA_DNreads
## 1                     0             0.00000                     0
## 2                   199            23.82948                    22
##   SRR316187_IPA_DNRPK SRR316187_LEreads SRR316187_LERPK SRR316187_IPA_RE
## 1            0.000000                40        24.78315       0.00000000
## 2            0.700347                40        24.78315      -0.09964814
##   SRR316188_IPA_UPreads SRR316188_IPA_UPRPK SRR316188_IPA_DNreads
## 1                     1            1.445087                     0
## 2                   207           24.787451                    15
##   SRR316188_IPA_DNRPK SRR316188_LEreads SRR316188_LERPK SRR316188_IPA_RE
## 1           0.0000000                55        34.07683       -4.5595631
## 2           0.4775093                55        34.07683       -0.4872446
##   SRR316189_IPA_UPreads SRR316189_IPA_UPRPK SRR316189_IPA_DNreads
## 1                     1            1.445087                     0
## 2                   226           27.062627                    15
##   SRR316189_IPA_DNRPK SRR316189_LEreads SRR316189_LERPK SRR316189_IPA_RE
## 1           0.0000000                46        28.50062       -4.3017653
## 2           0.4775093                46        28.50062       -0.1003744

6 Significance analysis of APA events

Once the read coverage information is obtained for each sample, one can compare APA regulation difference between two different groups. In this analysis, there are two types of experimental design: 1) without replicates; 2) with replicates. A sample table will be generated according to the design:

# Build the sample table with replicates
sampleTable1 = data.frame(samplename = c(names(flsall)),
                    condition = c(rep("NT",3),rep("KD",3)))
sampleTable1
##   samplename condition
## 1  SRR316184        NT
## 2  SRR316185        NT
## 3  SRR316186        NT
## 4  SRR316187        KD
## 5  SRR316188        KD
## 6  SRR316189        KD
# Build the sample table without replicates
sampleTable2 = data.frame(samplename = c("SRR316184","SRR316187"),
                    condition = c("NT","KD")) 
sampleTable2
##   samplename condition
## 1  SRR316184        NT
## 2  SRR316187        KD

6.1 Significantly regulated APA in 3’UTRs

To analyze 3’UTR APA between samples (KD and NT groups in the example) without replicates, sampleTable2 is used. The function used here is called APAdiff (detailed information can be obtained by the command ?APAdiff). It will fist to go through the sample table to determine whether it is a replicate design or non-replicate design. Then the APA compassion will be performed.

# Analysis 3'UTR APA between KD and NT group using non-repilicate design
test_3UTRsing=APAdiff(sampleTable2,DFUTRraw, 
                        conKET='NT',
                        trtKEY='KD',
                        PAS='3UTR',
                        CUTreads=0,
                        p_adjust_methods="fdr")

The APAdiff function requires two inputs: 1) A sample table defining groups/conditions of the samples, and 2) read coverage information of aUTRs and cUTRs, which can be obtained by PASEXP_3UTR from the previous step. The group name, i.e., treatment or control, can be defined b trtKEY= and conKET=; the PAS type analyzed should be defined by PAS=; the read cutoff used for aUTR and cUTR is defined by CUTreads= with the default value being 0; and the method used for p value correction is defined through p_adjust_methods= with the default value being ‘fdr’. In the non-replicate design, the APA pattern will be compared between two samples and output will be shown in a data frame:

head(test_3UTRsing,2)
##     gene_symbol        RED    pvalue     p_adj APAreg
## 3 1810055G02Rik -0.1195816 1.0000000 1.0000000     NC
## 4 2700081O15Rik  0.3862662 0.3292784 0.8007581     NC
table(test_3UTRsing$APAreg)
## 
##  DN  NC  UP 
##  31 246   8

The output contains 4 columns: ‘gene symbol’ describes gene information; ‘RED’ is relative expression difference between two groups; ‘pvalue’ is statistical significance based on the Fisher’s exact test; ‘p_adj’ is FDR adjusted pvalue and ‘APAreg’ is 3’UTR APA regulation pattern in the gene. We define 3 types in ‘APAreg’, ‘UP’ means aUTR abundance in the treatment group (‘KD’ in this case) is at least 5% higher than that in control (‘NT’ in this case), and ‘pvalue’<0.05; ‘DN’ means aUTR abundance is 5% lower in treatment than that in control and p-value<0.05; ‘NC’ are the remaining genes. With respect to 3’UTR size changes, ‘UP’ means 3’UTR lengthening, and ‘DN’ 3’UTR shortening.

For the replicate design, we use t-test for significance analysis. However, other tools based on negative binomial data distribution, such as
DEXSeq (Anders, Reyes, and W. 2012) might also be used.

# Analysis 3'UTR APA between KD and NT group using multi-repilicate design
test_3UTRmuti=APAdiff(sampleTable1,
                        DFUTRraw, 
                        conKET='NT',
                        trtKEY='KD',
                        PAS='3UTR',
                        CUTreads=0,
                        p_adjust_methods="fdr",
                        MultiTest='unpaired t-test')
head(test_3UTRmuti,2)
##     gene_symbol        RED     pvalue     p_adj APAreg
## 3 1810055G02Rik -0.7631689 0.08590526 0.4961374     NC
## 4 2700081O15Rik  0.1872116 0.07691742 0.4961374     NC
table(test_3UTRmuti$APAreg)
## 
##  DN  NC  UP 
##  24 226   8

In the replicate design, additional parameter MultiTest can be used to set the statistical method. In the current version, there are 3 methods can be set: “unpaired t-test” (default), “paired t-test”, “ANOVA”. After the running of APAdiff, a few columns will be generated: ‘RED’ is difference of averaged relative expression between two groups; ‘pvalue’ is the p-value calculated through the method defined by MultiTest. ‘APAreg’ is the predicted APA regulation direction; ‘UP’ is defined as ‘RED’ >0 and ‘pvalue’ <0.05; while ‘DN’ is the opposite; and ‘NC’ is the remaining genes.

6.2 Significantly regulated APA in introns

IPA comparison is similar to 3’UTR APA using APAdiff, except that it (1) uses IPA expression as input, and (2) ‘PAS=’ needs to be defined as ‘IPA’, and (3) the analysis is performed on each IPA. Note that, the direction of IPA regulation is similar to that of 3’UTR APA. This means ‘UP’ is defined as up-regulation of IPA (RED > 0); ‘DN’ is the opposite; and ‘NC’ is the remaining genes.

Analysis of IPA between KD and NT groups without replicates is shown here:

test_IPAsing=APAdiff(sampleTable2,
                        IPA_OUTraw, 
                        conKET='NT',
                        trtKEY='KD',
                        PAS='IPA',
                        CUTreads=0,
                        p_adjust_methods="fdr") 
head(test_IPAsing,2)
##     gene_symbol            PASid      RED       pvalue        p_adj APAreg
## 2 9130011E15Rik chr19:45952090:- 2.255401 1.790010e-12 3.615821e-10     UP
## 5        Ablim1 chr19:57151655:- 0.626637 5.392007e-01 1.000000e+00     NC

Analysis of IPA between KD and NT groups using replicate data is shown here:

test_IPAmuti=APAdiff(sampleTable1,
                        IPA_OUTraw, 
                        conKET='NT',
                        trtKEY='KD',
                        PAS='IPA',
                        CUTreads=0,
                        p_adjust_methods="fdr",
                        MultiTest='unpaired t-test')
head(test_IPAmuti,2)
##     gene_symbol            PASid         RED      pvalue     p_adj APAreg
## 2 9130011E15Rik chr19:45952090:-  1.84729507 0.001294432 0.2808917     UP
## 5        Ablim1 chr19:57151655:- -0.09466715 0.900027947 0.9480877     NC

7 Visualization of analysis results

Start from APAlyzer 1.2.1, we provides two new fouction called APAVolcano and APABox for users to plot their RED results using volcano plot and box plot. In the volcano plot, users can also label the top genes using top= or a set of specific gene using markergenes=, for example:

APAVolcano(test_3UTRsing, PAS='3UTR', Pcol = "pvalue", top=5, main='3UTR APA')

In the box plot, RED is ploted on ‘UP’, ‘DN’, and ‘NC’ genes:

APABox(test_3UTRsing, xlab = "APAreg", ylab = "RED", plot_title = NULL)

In addtion to volcano and box plots, APA comparison result can be also plotted using either boxplots or violin plots or CDF curves. For the previous 3’UTR APA and IPA comparison outputs, one needs to first build the plotting data frame:

test_3UTRmuti$APA="3'UTR"
test_IPAmuti$APA="IPA"
dfplot=rbind(test_3UTRmuti[,c('RED','APA')],test_IPAmuti[,c('RED','APA')])

To make violin plots and CDF curves using ggplot2:

library(ggplot2)
###violin
ggplot(dfplot, aes(x = APA, y = RED)) + 
    geom_violin(trim = FALSE) + 
    geom_boxplot(width = 0.2)+ theme_bw() + 
    geom_hline(yintercept=0, linetype="dashed", color = "red")

###CDF
ggplot(dfplot, aes( x = RED, color = APA)) + 
    stat_ecdf(geom = "step") +
    ylab("cumulative fraction")+ 
    geom_vline(xintercept=0, linetype="dashed", color = "gray")+ theme_bw() + 
    geom_hline(yintercept=0.5, linetype="dashed", color = "gray")

8 Gene expression analysis using coding regions

APA is frequently involved in gene expression regulation. To compare gene expression vs. APA in different samples, our package provides a simple function to assess the expression changes using RNA-seq reads mapped to coding sequences.

8.1 Building coding region references

suppressMessages(library("GenomicFeatures"))
suppressMessages(library("org.Mm.eg.db"))
extpath = system.file("extdata", "mm9.chr19.refGene.R.DB", package="APAlyzer")
txdb=loadDb(extpath, packageName='GenomicFeatures')
IDDB = org.Mm.eg.db
CDSdbraw=REFCDS(txdb,IDDB)
## 'select()' returned many:1 mapping between keys and columns

8.2 Calculation of expression

DFGENEraw=GENEXP_CDS(CDSdbraw, flsall, Strandtype="forward")
## [1] "SRR316184, Strand: forward, finished"
## [1] "SRR316185, Strand: forward, finished"
## [1] "SRR316186, Strand: forward, finished"
## [1] "SRR316187, Strand: forward, finished"
## [1] "SRR316188, Strand: forward, finished"
## [1] "SRR316189, Strand: forward, finished"

9 Extract 3’most alignment from the pair-end bam file

Start from APAlyzer 1.5.5, we provides a new fouction called ThreeMostPairBam for users to extract three prime most alignment from the paired-end bam file and save it into a new bam file. For example:

## Extract 3 prime most alignment of a paired-end 
## bam file and saved into a new bam file
Bamfile='/path/to/inputdir/input.bam'
Outdir='/path/to/outdir/'  
StrandType="forward-reverse"    ## "forward-reverse",  or "reverse-forward" or "NONE"   
ThreeMostPairBam (BamfilePath=Bamfile, 
                    OutDirPath=Outdir, 
                    StrandType='forward-reverse')

The output bamfile (A bam file named as *.3most.bam located in Outdir) can be used as input file in PASEXP_3UTR. To use this bam file for PASEXP_IPA, the user need to set the SeqType to “ThreeMostPairEnd”, for example:

ThreemostBamfile="test.3most.bam"
IPA_OUTraw=PASEXP_IPA(dfIPA, dfLE, ThreemostBamfile, Strandtype="forward", SeqType="ThreeMostPairEnd")

10 Complete Analysis Example: APA analysis in mouse testis versus heart

10.1 About this dataset

To provide a complete tutorial of APA analysis using our package, we have now prepared a testing dataset through down sampling of mouse RNA-Seq data in heart (GSM900193) and testis (GSM900199):

Sample ID SRRID Sample Name Down sampling reads
GSM900199 SRR453175 Heart_Rep1 5 Million
GSM900199 SRR453174 Heart_Rep2 5 Million
GSM900199 SRR453173 Heart_Rep3 5 Million
GSM900199 SRR453172 Heart_Rep4 5 Million
GSM900193 SRR453143 Testis_rep1 5 Million
GSM900193 SRR453142 Testis_rep2 5 Million
GSM900193 SRR453141 Testis_rep3 5 Million
GSM900193 SRR453140 Testis_rep4 5 Million

10.2 Download the bam files

download_testbam()
flsall <- dir(getwd(),".bam")
flsall<-paste0(getwd(),'/',flsall)
names(flsall)<-gsub('.bam','',dir(getwd(),".bam"))

10.3 Build the PAS reference regions

library(repmis)
URL="https://github.com/RJWANGbioinfo/PAS_reference_RData/blob/master/"
file="mm9_REF.RData"
source_data(paste0(URL,file,"?raw=True"))
PASREF=REF4PAS(refUTRraw,dfIPA,dfLE)
UTRdbraw=PASREF$UTRdbraw
dfIPA=PASREF$dfIPA
dfLE=PASREF$dfLE

10.4 Calculation of relative expression of 3’UTR APA and IPA

UTR_APA_OUT=PASEXP_3UTR(UTRdbraw, flsall, Strandtype="invert")
IPA_OUT=PASEXP_IPA(dfIPA, dfLE, flsall, Strandtype="invert", nts=1)

10.5 Significantly regulated APA in 3’UTRs

############# 3utr APA #################
sampleTable = data.frame(samplename = c('Heart_rep1',
                                        'Heart_rep2',
                                        'Heart_rep3',
                                        'Heart_rep4',
                                        'Testis_rep1',
                                        'Testis_rep2',
                                        'Testis_rep3',
                                        'Testis_rep4'),
                        condition = c(rep("Heart",4),
                                    rep("Testis",4)))
                        
test_3UTRAPA=APAdiff(sampleTable,UTR_APA_OUT, 
                        conKET='Heart',
                        trtKEY='Testis',
                        PAS='3UTR',
                        CUTreads=5)
table(test_3UTRAPA$APAreg)                      
## 
##  DN  NC  UP 
## 490 719  80
APAVolcano(test_3UTRAPA, PAS='3UTR', Pcol = "pvalue", plot_title='3UTR APA')

APABox(test_3UTRAPA, xlab = "APAreg", ylab = "RED", plot_title = NULL)

10.6 Significantly regulated APA in Intron

############# IPA #################
test_IPA=APAdiff(sampleTable,
                        IPA_OUT, 
                        conKET='Heart',
                        trtKEY='Testis',
                        PAS='IPA',
                        CUTreads=5)
table(test_IPA$APAreg)  
## 
##  DN  NC  UP 
## 125 924 312
APAVolcano(test_IPA, PAS='IPA', Pcol = "pvalue", plot_title='IPA')

APABox(test_IPA, xlab = "APAreg", ylab = "RED", plot_title = NULL)

11 FAQs

How to generate a BAM file list for analysis?

A BAM file list containing both BAM file names and paths of the files. Let’s say all the BAM files are stored in bamdir, then BAM file lists can be obtained through:

flsall = dir(bamdir,".bam")
flsall=paste0(bamdir,flsall)
names(flsall)=dir(bamdir,".bam")

Why am I getting error messages when I try to get txdb using makeTxDbFromUCSC?

You can try either upgrade your Bioconductor, or load the genome annotation using GTF, or load the prebuild genome annotation using ‘.R.DB’ file, e.g., mm9.refGene.R.DB.

12 Session Information

The session information records all the package versions used in generation of the present document.

sessionInfo()
## R version 4.3.0 RC (2023-04-13 r84269)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.17-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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] org.Mm.eg.db_3.17.0    GenomicFeatures_1.52.0 AnnotationDbi_1.62.0  
##  [4] Biobase_2.60.0         ggplot2_3.4.2          knitr_1.42            
##  [7] repmis_0.5             TBX20BamSubset_1.35.0  Rsamtools_2.16.0      
## [10] Biostrings_2.68.0      XVector_0.40.0         GenomicRanges_1.52.0  
## [13] GenomeInfoDb_1.36.0    IRanges_2.34.0         S4Vectors_0.38.0      
## [16] BiocGenerics_0.46.0    APAlyzer_1.14.0        BiocStyle_2.28.0      
## 
## loaded via a namespace (and not attached):
##  [1] DBI_1.1.3                   bitops_1.0-7               
##  [3] biomaRt_2.56.0              rlang_1.1.0                
##  [5] magrittr_2.0.3              matrixStats_0.63.0         
##  [7] compiler_4.3.0              RSQLite_2.3.1              
##  [9] png_0.1-8                   vctrs_0.6.2                
## [11] stringr_1.5.0               pkgconfig_2.0.3            
## [13] crayon_1.5.2                Rsubread_2.14.0            
## [15] fastmap_1.1.1               magick_2.7.4               
## [17] dbplyr_2.3.2                labeling_0.4.2             
## [19] utf8_1.2.3                  rmarkdown_2.21             
## [21] purrr_1.0.1                 bit_4.0.5                  
## [23] xfun_0.39                   zlibbioc_1.46.0            
## [25] cachem_1.0.7                jsonlite_1.8.4             
## [27] progress_1.2.2              blob_1.2.4                 
## [29] highr_0.10                  DelayedArray_0.26.0        
## [31] BiocParallel_1.34.0         parallel_4.3.0             
## [33] HybridMTest_1.44.0          prettyunits_1.1.1          
## [35] VariantAnnotation_1.46.0    R6_2.5.1                   
## [37] bslib_0.4.2                 stringi_1.7.12             
## [39] rtracklayer_1.60.0          jquerylib_0.1.4            
## [41] Rcpp_1.0.10                 bookdown_0.33              
## [43] SummarizedExperiment_1.30.0 R.utils_2.12.2             
## [45] Matrix_1.5-4                R.cache_0.16.0             
## [47] tidyselect_1.2.0            yaml_2.3.7                 
## [49] codetools_0.2-19            curl_5.0.0                 
## [51] lattice_0.21-8              tibble_3.2.1               
## [53] plyr_1.8.8                  withr_2.5.0                
## [55] KEGGREST_1.40.0             evaluate_0.20              
## [57] BiocFileCache_2.8.0         xml2_1.3.3                 
## [59] pillar_1.9.0                BiocManager_1.30.20        
## [61] filelock_1.0.2              MatrixGenerics_1.12.0      
## [63] generics_0.1.3              RCurl_1.98-1.12            
## [65] hms_1.1.3                   munsell_0.5.0              
## [67] scales_1.2.1                glue_1.6.2                 
## [69] tools_4.3.0                 BiocIO_1.10.0              
## [71] data.table_1.14.8           BSgenome_1.68.0            
## [73] locfit_1.5-9.7              GenomicAlignments_1.36.0   
## [75] XML_3.99-0.14               grid_4.3.0                 
## [77] tidyr_1.3.0                 colorspace_2.1-0           
## [79] GenomeInfoDbData_1.2.10     restfulr_0.0.15            
## [81] cli_3.6.1                   rappdirs_0.3.3             
## [83] fansi_1.0.4                 dplyr_1.1.2                
## [85] gtable_0.3.3                R.methodsS3_1.8.2          
## [87] DESeq2_1.40.0               sass_0.4.5                 
## [89] digest_0.6.31               ggrepel_0.9.3              
## [91] farver_2.1.1                rjson_0.2.21               
## [93] memoise_2.0.1               htmltools_0.5.5            
## [95] R.oo_1.25.0                 lifecycle_1.0.3            
## [97] httr_1.4.5                  bit64_4.0.5

13 Acknowledgements

We thank members of the Bin Tian lab for helpful discussions and package testing.

References

Anders, S., A. Reyes, and Huber W. 2012. “Detecting Differential Usage of Exons from Rna-Seq Data.” Genome Research 22: 4025. https://doi.org/10.1101/gr.133744.111.

Bindreither, D. 2018. TBX20BamSubset: Subset of Bam Files from the "Tbx20" Experiment.

Wang, R., R. Nambiar, D. Zheng, and Tian B. 2017. “PolyA_DB 3 Catalogs Cleavage and Polyadenylation Sites Identified by Deep Sequencing in Multiple Genomes.” Nucleic Acids Research 46 (D1): D315–D319. https://doi.org/10.1093/nar/gkx1000.

Wang, R., D. Zheng, G. Yehia, and Tian B. 2018. “A Compendium of Conserved Cleavage and Polyadenylation Events in Mammalian Genes.” Genome Research 28 (10): 1427–41. https://doi.org/10.1101/gr.237826.118.