Contents

1 Introduction

Alternative polyadenylation (APA) is one of the most important post-transcriptional regulation mechanisms which is prevalent in human genes. Like alternative splicing, APA can create proteome diversity. In addition, it defines 3’ UTR and results in altered expression of the gene. It is a tightly controlled process and mis-regulation of APA can lead to pathological effects to the cells, such as uncontrolled cell cycle and growth. Although several high throughput sequencing methods have been developed, there are still limited data available.

RNA-seq has become one of the most commonly used methods for quantifying genome-wide gene expression. There are massive RNA-seq datasets available in public databases such as GEO and TCGA. For this reason, we develop the InPAS algorithm for identifying APA from conventional RNA-seq data.

The workflow for InPAS is:

2 How to

First, load the required libraries InPAS, species specific BSgenome, and TxDb as follows.

library(InPAS)
library(BSgenome.Mmusculus.UCSC.mm10)
library(TxDb.Mmusculus.UCSC.mm10.knownGene)
path <- file.path(find.package("InPAS"), "extdata")

Next, prepare annotation using the function utr3Annotation with a TxDb and org annotation. Please note that the 3’ UTR annotation prepared by utr3Annotation includes the gaps to the next annotated region.

library(org.Hs.eg.db)
samplefile <- system.file("extdata", "hg19_knownGene_sample.sqlite",
            package="GenomicFeatures")
txdb <- loadDb(samplefile)
utr3.sample.anno <- utr3Annotation(txdb=txdb, 
                   orgDbSYMBOL="org.Hs.egSYMBOL")
utr3.sample.anno
## GRanges object with 155 ranges and 7 metadata columns:
##                             seqnames              ranges strand |     feature
##                                <Rle>           <IRanges>  <Rle> | <character>
##      uc001bum.2_5|IQCC|utr3     chr1   32673684-32674288      + |        utr3
##    uc001fbq.3_3|S100A9|utr3     chr1 153333315-153333503      + |        utr3
##    uc001gde.2_2|LRRC52|utr3     chr1 165533062-165533185      + |        utr3
##   uc001hfg.3_15|PFKFB2|utr3     chr1 207245717-207251162      + |        utr3
##   uc001hfh.3_15|PFKFB2|utr3     chr1 207252365-207254368      + |        utr3
##                         ...      ...                 ...    ... .         ...
##      uc004dsv.3_19|PHF8|CDS     chrX   53964414-53964492      - |         CDS
##      uc004dsx.3_15|PHF8|CDS     chrX   53969797-53969835      - |         CDS
##     uc004ehz.1_5|ARMCX3|CDS     chrX 100879970-100881109      + |         CDS
##    uc004elw.3_6|FAM199X|CDS     chrX 103434289-103434459      + |         CDS
##      uc004fmj.1_10|GAB3|CDS     chrX 153906455-153906571      - |         CDS
##                             annotatedProximalCP          exon  transcript
##                                     <character>   <character> <character>
##      uc001bum.2_5|IQCC|utr3             unknown  uc001bum.2_5  uc001bum.2
##    uc001fbq.3_3|S100A9|utr3             unknown  uc001fbq.3_3  uc001fbq.3
##    uc001gde.2_2|LRRC52|utr3             unknown  uc001gde.2_2  uc001gde.2
##   uc001hfg.3_15|PFKFB2|utr3             unknown uc001hfg.3_15  uc001hfg.3
##   uc001hfh.3_15|PFKFB2|utr3             unknown uc001hfh.3_15  uc001hfh.3
##                         ...                 ...           ...         ...
##      uc004dsv.3_19|PHF8|CDS             unknown uc004dsv.3_19  uc004dsv.3
##      uc004dsx.3_15|PHF8|CDS             unknown uc004dsx.3_15  uc004dsx.3
##     uc004ehz.1_5|ARMCX3|CDS             unknown  uc004ehz.1_5  uc004ehz.1
##    uc004elw.3_6|FAM199X|CDS             unknown  uc004elw.3_6  uc004elw.3
##      uc004fmj.1_10|GAB3|CDS             unknown uc004fmj.1_10  uc004fmj.1
##                                    gene      symbol truncated
##                             <character> <character> <logical>
##      uc001bum.2_5|IQCC|utr3       55721        IQCC     FALSE
##    uc001fbq.3_3|S100A9|utr3        6280      S100A9     FALSE
##    uc001gde.2_2|LRRC52|utr3      440699      LRRC52     FALSE
##   uc001hfg.3_15|PFKFB2|utr3        5208      PFKFB2     FALSE
##   uc001hfh.3_15|PFKFB2|utr3        5208      PFKFB2     FALSE
##                         ...         ...         ...       ...
##      uc004dsv.3_19|PHF8|CDS       23133        PHF8     FALSE
##      uc004dsx.3_15|PHF8|CDS       23133        PHF8     FALSE
##     uc004ehz.1_5|ARMCX3|CDS       51566      ARMCX3     FALSE
##    uc004elw.3_6|FAM199X|CDS      139231     FAM199X     FALSE
##      uc004fmj.1_10|GAB3|CDS      139716        GAB3     FALSE
##   -------
##   seqinfo: 27 sequences from hg19 genome; no seqlengths

Users can also directly load the 3’ UTR annotation included in this package for mm10 and hg19. Here we show how to load the pre-built mm10 3’ UTR annotation file.

##step1 annotation
# load from dataset
data(utr3.mm10)

For coverage calculation, alignment files in BEDGraph format are required, which can be generated from BAM files using the genomecov tool in bedtools with parameter: -bg -split. Potential pA sites identified from the coverage data can be filtered/adjusted using the classifier provided by cleanUpdTseq. The following scripts illustrate the function calls needed to perform the complete analysis using InPAS.

load(file.path(path, "polyA.rds"))
library(cleanUpdTSeq)
data(classifier)
bedgraphs <- c(file.path(path, "Baf3.extract.bedgraph"), 
               file.path(path, "UM15.extract.bedgraph"))
hugeData <- FALSE
##step1 Calculate coverage
coverage <- coverageFromBedGraph(bedgraphs, 
                                 tags=c("Baf3", "UM15"),
                                 genome=BSgenome.Mmusculus.UCSC.mm10,
                                 hugeData=hugeData)

## we hope the coverage rate of should be greater than 0.75 in the expressed gene.
## which is used because the demo data is a subset of genome.
coverageRate(coverage=coverage,
             txdb=TxDb.Mmusculus.UCSC.mm10.knownGene,
             genome=BSgenome.Mmusculus.UCSC.mm10,
             which = GRanges("chr6", ranges=IRanges(98013000, 140678000)))
## strand information will be ignored.
##      gene.coverage.rate expressed.gene.coverage.rate UTR3.coverage.rate
## Baf3         0.01031381                    0.6300834         0.01830495
## UM15         0.01020966                    0.6236834         0.01846496
##      UTR3.expressed.gene.subset.coverage.rate
## Baf3                                0.8082201
## UM15                                0.8152853
##step2 Predict cleavage sites
CPs <- CPsites(coverage=coverage, 
               genome=BSgenome.Mmusculus.UCSC.mm10,
               utr3=utr3.mm10, 
               search_point_START=200, 
               cutEnd=.2, 
               long_coverage_threshold=3,
               background="10K",
               txdb=TxDb.Mmusculus.UCSC.mm10.knownGene,
               PolyA_PWM=pwm, 
               classifier=classifier,
               shift_range=50,
               step=10)
head(CPs)
## GRanges object with 4 ranges and 12 metadata columns:
##                              seqnames              ranges strand |
##                                 <Rle>           <IRanges>  <Rle> |
##      uc009daz.2_10|Mitf|utr3     chr6   98018176-98021358      + |
##      uc009dhz.2_19|Atg7|utr3     chr6 114859343-114860614      + |
##      uc009dmb.2_4|Lrtm2|utr3     chr6 119315133-119317055      - |
##   uc009eet.1_3|BC035044|utr3     chr6 128848044-128850081      - |
##                              annotatedProximalCP          exon  transcript
##                                      <character>   <character> <character>
##      uc009daz.2_10|Mitf|utr3             unknown uc009daz.2_10  uc009daz.2
##      uc009dhz.2_19|Atg7|utr3             unknown uc009dhz.2_19  uc009dhz.2
##      uc009dmb.2_4|Lrtm2|utr3             unknown  uc009dmb.2_4  uc009dmb.2
##   uc009eet.1_3|BC035044|utr3             unknown  uc009eet.1_3  uc009eet.1
##                                     gene      symbol truncated fit_value
##                              <character> <character> <logical> <numeric>
##      uc009daz.2_10|Mitf|utr3       17342        Mitf     FALSE 12594.674
##      uc009dhz.2_19|Atg7|utr3       74244        Atg7     FALSE 27383.413
##      uc009dmb.2_4|Lrtm2|utr3      211187       Lrtm2     FALSE   168.551
##   uc009eet.1_3|BC035044|utr3      232406    BC035044     FALSE   123.974
##                              Predicted_Proximal_APA Predicted_Distal_APA
##                                           <numeric>            <numeric>
##      uc009daz.2_10|Mitf|utr3               98018978             98021358
##      uc009dhz.2_19|Atg7|utr3              114859674            114862071
##      uc009dmb.2_4|Lrtm2|utr3              119316541            119315133
##   uc009eet.1_3|BC035044|utr3              128849177            128846744
##                                        type utr3start   utr3end
##                                 <character> <numeric> <numeric>
##      uc009daz.2_10|Mitf|utr3 novel proximal  98018276  98021358
##      uc009dhz.2_19|Atg7|utr3   novel distal 114859443 114860614
##      uc009dmb.2_4|Lrtm2|utr3 novel proximal 119316955 119315133
##   uc009eet.1_3|BC035044|utr3   novel distal 128849981 128848044
##   -------
##   seqinfo: 42 sequences from mm10 genome; no seqlengths
##step3 Estimate 3UTR usage
res <- testUsage(CPsites=CPs, 
                 coverage=coverage, 
                 genome=BSgenome.Mmusculus.UCSC.mm10,
                 utr3=utr3.mm10, 
                 method="fisher.exact",
                 gp1="Baf3", gp2="UM15")

##step4 view the results
as(res, "GRanges") 
## GRanges object with 4 ranges and 27 metadata columns:
##              seqnames              ranges strand | annotatedProximalCP
##                 <Rle>           <IRanges>  <Rle> |         <character>
##   uc009daz.2     chr6   98018176-98021358      + |             unknown
##   uc009dhz.2     chr6 114859343-114862071      + |             unknown
##   uc009dmb.2     chr6 119315133-119317055      - |             unknown
##   uc009eet.1     chr6 128846744-128850081      - |             unknown
##               transcript        gene      symbol truncated fit_value
##              <character> <character> <character> <logical> <numeric>
##   uc009daz.2  uc009daz.2       17342        Mitf     FALSE 12594.674
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE 27383.413
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE   168.551
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE   123.974
##              Predicted_Proximal_APA Predicted_Distal_APA           type
##                           <numeric>            <numeric>    <character>
##   uc009daz.2               98018978             98021358 novel proximal
##   uc009dhz.2              114859674            114862071   novel distal
##   uc009dmb.2              119316541            119315133 novel proximal
##   uc009eet.1              128849177            128846744   novel distal
##              utr3start   utr3end Baf3_short.form.usage UM15_short.form.usage
##              <numeric> <numeric>             <numeric>             <numeric>
##   uc009daz.2  98018276  98021358              33.53134               1.05215
##   uc009dhz.2 114859443 114860614             520.91479             172.15235
##   uc009dmb.2 119316955 119315133               8.85786              49.53562
##   uc009eet.1 128849981 128848044              21.64617               0.00000
##              Baf3_long.form.usage UM15_long.form.usage Baf3_PDUI UM15_PDUI
##                         <numeric>            <numeric> <numeric> <numeric>
##   uc009daz.2            282.82402             189.1411  0.894007  0.994468
##   uc009dhz.2            208.02460             456.5446  0.285380  0.726176
##   uc009dmb.2              8.53513              70.8126  0.490722  0.588398
##   uc009eet.1              8.90222              23.2913  0.291414  1.000000
##              short.mean.gp1 long.mean.gp1 short.mean.gp2 long.mean.gp2
##                   <numeric>     <numeric>      <numeric>     <numeric>
##   uc009daz.2       33.53134     282.82402        1.05215      189.1411
##   uc009dhz.2      520.91479     208.02460      172.15235      456.5446
##   uc009dmb.2        8.85786       8.53513       49.53562       70.8126
##   uc009eet.1       21.64617       8.90222        0.00000       23.2913
##               PDUI.gp1  PDUI.gp2      dPDUI     logFC     P.Value   adj.P.Val
##              <numeric> <numeric>  <numeric> <numeric>   <numeric>   <numeric>
##   uc009daz.2  0.894007  0.994468 -0.1004606 -0.153638 1.58698e-06 2.11597e-06
##   uc009dhz.2  0.285380  0.726176 -0.4407961 -1.347436 1.70948e-60 6.83790e-60
##   uc009dmb.2  0.490722  0.588398 -0.0976754 -0.261885 5.93031e-01 5.93031e-01
##   uc009eet.1  0.291414  1.000000 -0.7085863 -1.778859 4.13501e-08 8.27003e-08
##   -------
##   seqinfo: 1 sequence from an unspecified genome; no seqlengths
filterRes(res, gp1="Baf3", gp2="UM15",
          background_coverage_threshold=3,
          adj.P.Val_cutoff=0.05, 
          dPDUI_cutoff=0.3, 
          PDUI_logFC_cutoff=0.59)
## GRanges object with 4 ranges and 28 metadata columns:
##              seqnames              ranges strand | annotatedProximalCP
##                 <Rle>           <IRanges>  <Rle> |         <character>
##   uc009daz.2     chr6   98018176-98021358      + |             unknown
##   uc009dhz.2     chr6 114859343-114862071      + |             unknown
##   uc009dmb.2     chr6 119315133-119317055      - |             unknown
##   uc009eet.1     chr6 128846744-128850081      - |             unknown
##               transcript        gene      symbol truncated fit_value
##              <character> <character> <character> <logical> <numeric>
##   uc009daz.2  uc009daz.2       17342        Mitf     FALSE 12594.674
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE 27383.413
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE   168.551
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE   123.974
##              Predicted_Proximal_APA Predicted_Distal_APA           type
##                           <numeric>            <numeric>    <character>
##   uc009daz.2               98018978             98021358 novel proximal
##   uc009dhz.2              114859674            114862071   novel distal
##   uc009dmb.2              119316541            119315133 novel proximal
##   uc009eet.1              128849177            128846744   novel distal
##              utr3start   utr3end Baf3_short.form.usage UM15_short.form.usage
##              <numeric> <numeric>             <numeric>             <numeric>
##   uc009daz.2  98018276  98021358              33.53134               1.05215
##   uc009dhz.2 114859443 114860614             520.91479             172.15235
##   uc009dmb.2 119316955 119315133               8.85786              49.53562
##   uc009eet.1 128849981 128848044              21.64617               0.00000
##              Baf3_long.form.usage UM15_long.form.usage Baf3_PDUI UM15_PDUI
##                         <numeric>            <numeric> <numeric> <numeric>
##   uc009daz.2            282.82402             189.1411  0.894007  0.994468
##   uc009dhz.2            208.02460             456.5446  0.285380  0.726176
##   uc009dmb.2              8.53513              70.8126  0.490722  0.588398
##   uc009eet.1              8.90222              23.2913  0.291414  1.000000
##              short.mean.gp1 long.mean.gp1 short.mean.gp2 long.mean.gp2
##                   <numeric>     <numeric>      <numeric>     <numeric>
##   uc009daz.2       33.53134     282.82402        1.05215      189.1411
##   uc009dhz.2      520.91479     208.02460      172.15235      456.5446
##   uc009dmb.2        8.85786       8.53513       49.53562       70.8126
##   uc009eet.1       21.64617       8.90222        0.00000       23.2913
##               PDUI.gp1  PDUI.gp2      dPDUI     logFC     P.Value   adj.P.Val
##              <numeric> <numeric>  <numeric> <numeric>   <numeric>   <numeric>
##   uc009daz.2  0.894007  0.994468 -0.1004606 -0.153638 1.58698e-06 2.11597e-06
##   uc009dhz.2  0.285380  0.726176 -0.4407961 -1.347436 1.70948e-60 6.83790e-60
##   uc009dmb.2  0.490722  0.588398 -0.0976754 -0.261885 5.93031e-01 5.93031e-01
##   uc009eet.1  0.291414  1.000000 -0.7085863 -1.778859 4.13501e-08 8.27003e-08
##                   PASS
##              <logical>
##   uc009daz.2     FALSE
##   uc009dhz.2      TRUE
##   uc009dmb.2     FALSE
##   uc009eet.1      TRUE
##   -------
##   seqinfo: 1 sequence from an unspecified genome; no seqlengths

The steps described above can be done in one function call.

if(interactive()){
    res <- inPAS(bedgraphs=bedgraphs, tags=c("Baf3", "UM15"), 
              genome=BSgenome.Mmusculus.UCSC.mm10, 
              utr3=utr3.mm10, gp1="Baf3", gp2="UM15",
              txdb=TxDb.Mmusculus.UCSC.mm10.knownGene,
              search_point_START=200,
              short_coverage_threshold=15,
              long_coverage_threshold=3, 
              cutStart=0, cutEnd=.2, 
              hugeData=FALSE,
              shift_range=50,
              PolyA_PWM=pwm, classifier=classifier,
              method="fisher.exact",
              adj.P.Val_cutoff=0.05,
              dPDUI_cutoff=0.3, 
              PDUI_logFC_cutoff=0.59)
}

InPAS can also handle single group data.

inPAS(bedgraphs=bedgraphs[1], tags=c("Baf3"), 
      genome=BSgenome.Mmusculus.UCSC.mm10, 
      utr3=utr3.mm10, gp1="Baf3", gp2=NULL,
      txdb=TxDb.Mmusculus.UCSC.mm10.knownGene,
      search_point_START=200,
      short_coverage_threshold=15,
      long_coverage_threshold=3, 
      cutStart=0, cutEnd=.2, 
      hugeData=FALSE, 
      PolyA_PWM=pwm, classifier=classifier,
      method="singleSample",
      adj.P.Val_cutoff=0.05,
      step=10)
## converged at iteration 1 with logLik: -1835.501 
## converged at iteration 5 with logLik: -838.8306 
## converged at iteration 1 with logLik: -1496.738 
## converged at iteration 13 with logLik: -724.6964 
## converged at iteration 1 with logLik: -997.3022 
## converged at iteration 7 with logLik: -555.5656 
## converged at iteration 1 with logLik: -188.938 
## converged at iteration 23 with logLik: -152.6663 
## converged at iteration 1 with logLik: -462.3804 
## converged at iteration 10 with logLik: -214.1651
## dPDUI is calculated by gp2 - gp1.
## GRanges object with 6 ranges and 22 metadata columns:
##              seqnames              ranges strand | annotatedProximalCP
##                 <Rle>           <IRanges>  <Rle> |         <character>
##   uc009daz.2     chr6   98018176-98021358      + |             unknown
##   uc009dhz.2     chr6 114859343-114862075      + |             unknown
##   uc009die.2     chr6 114860617-114862164      - |             unknown
##   uc009dmb.2     chr6 119315133-119317055      - |             unknown
##   uc009eet.1     chr6 128847265-128850081      - |             unknown
##   uc009eol.1     chr6 140651362-140651622      + |             unknown
##               transcript        gene      symbol truncated   fit_value
##              <character> <character> <character> <logical>   <numeric>
##   uc009daz.2  uc009daz.2       17342        Mitf     FALSE  17843.7079
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE   7630.4981
##   uc009die.2  uc009die.2      232334       Vgll4     FALSE  10704.6989
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE     18.6351
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE    227.5418
##   uc009eol.1  uc009eol.1       11569       Aebp2      TRUE 204539.6984
##              Predicted_Proximal_APA Predicted_Distal_APA           type
##                           <numeric>            <numeric>    <character>
##   uc009daz.2               98019022             98021358 novel proximal
##   uc009dhz.2              114860283            114862075   novel distal
##   uc009die.2              114861446            114860617   novel distal
##   uc009dmb.2              119316683            119315133 novel proximal
##   uc009eet.1              128848843            128847265   novel distal
##   uc009eol.1              140651566            140651622 novel proximal
##              utr3start   utr3end              data2 Baf3_short.form.usage
##              <numeric> <numeric>             <list>             <numeric>
##   uc009daz.2  98018176  98021358    473,460,457,...              29.66026
##   uc009dhz.2 114859343 114860614    184,184,184,...               3.72545
##   uc009die.2 114862164 114862092       76,76,76,...             121.88814
##   uc009dmb.2 119317055 119315133       11,12,12,...               9.22610
##   uc009eet.1 128850081 128848044       22,22,22,...              22.04495
##   uc009eol.1 140651362 140651622 1408,1412,1384,...            1065.98013
##              Baf3_long.form.usage Baf3_PDUI short.mean long.mean     PDUI
##                         <numeric> <numeric>     <list>    <list>   <list>
##   uc009daz.2            283.38939  0.905254    29.6603   283.389 0.905254
##   uc009dhz.2            209.97881  0.982567    3.72545   209.979 0.982567
##   uc009die.2            174.12289  0.588231    121.888   174.123 0.588231
##   uc009dmb.2              9.11799  0.497053     9.2261   9.11799 0.497053
##   uc009eet.1             10.32742  0.319020     22.045   10.3274  0.31902
##   uc009eol.1            221.45614  0.172013    1065.98   221.456 0.172013
##              P.Value adj.P.Val     dPDUI      PASS
##               <list>    <list> <numeric> <logical>
##   uc009daz.2       1         1       NaN     FALSE
##   uc009dhz.2       1         1       NaN     FALSE
##   uc009die.2       1         1       NaN     FALSE
##   uc009dmb.2       1         1       NaN     FALSE
##   uc009eet.1       1         1       NaN     FALSE
##   uc009eol.1       1         1       NaN     FALSE
##   -------
##   seqinfo: 1 sequence from an unspecified genome; no seqlengths

3 Session Info

sessionInfo()

R version 4.0.3 (2020-10-10) Platform: x86_64-pc-linux-gnu (64-bit) Running under: Ubuntu 18.04.5 LTS

Matrix products: default BLAS: /home/biocbuild/bbs-3.12-bioc/R/lib/libRblas.so LAPACK: /home/biocbuild/bbs-3.12-bioc/R/lib/libRlapack.so

locale: [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
[5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
[7] LC_PAPER=en_US.UTF-8 LC_NAME=C
[9] LC_ADDRESS=C LC_TELEPHONE=C
[11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C

attached base packages: [1] stats4 parallel stats graphics grDevices utils datasets [8] methods base

other attached packages: [1] cleanUpdTSeq_1.28.0
[2] org.Hs.eg.db_3.12.0
[3] TxDb.Mmusculus.UCSC.mm10.knownGene_3.10.0 [4] BSgenome.Mmusculus.UCSC.mm10_1.4.0
[5] BSgenome_1.58.0
[6] rtracklayer_1.50.0
[7] Biostrings_2.58.0
[8] XVector_0.30.0
[9] InPAS_1.22.0
[10] GenomicFeatures_1.42.0
[11] AnnotationDbi_1.52.0
[12] GenomicRanges_1.42.0
[13] GenomeInfoDb_1.26.0
[14] IRanges_2.24.0
[15] S4Vectors_0.28.0
[16] Biobase_2.50.0
[17] BiocGenerics_0.36.0
[18] BiocStyle_2.18.0

loaded via a namespace (and not attached): [1] colorspace_1.4-1 seqinr_4.2-4
[3] ellipsis_0.3.1 class_7.3-17
[5] depmixS4_1.4-2 biovizBase_1.38.0
[7] htmlTable_2.1.0 base64enc_0.1-3
[9] dichromat_2.0-0 rstudioapi_0.11
[11] bit64_4.0.5 xml2_1.3.2
[13] splines_4.0.3 knitr_1.30
[15] ade4_1.7-15 Formula_1.2-4
[17] Rsamtools_2.6.0 cluster_2.1.0
[19] dbplyr_1.4.4 png_0.1-7
[21] BiocManager_1.30.10 compiler_4.0.3
[23] httr_1.4.2 backports_1.1.10
[25] assertthat_0.2.1 Matrix_1.2-18
[27] lazyeval_0.2.2 limma_3.46.0
[29] htmltools_0.5.0 prettyunits_1.1.1
[31] tools_4.0.3 gtable_0.3.0
[33] glue_1.4.2 GenomeInfoDbData_1.2.4
[35] dplyr_1.0.2 rappdirs_0.3.1
[37] Rcpp_1.0.5 vctrs_0.3.4
[39] preprocessCore_1.52.0 nlme_3.1-150
[41] xfun_0.18 stringr_1.4.0
[43] lifecycle_0.2.0 ensembldb_2.14.0
[45] XML_3.99-0.5 zlibbioc_1.36.0
[47] MASS_7.3-53 scales_1.1.1
[49] VariantAnnotation_1.36.0 hms_0.5.3
[51] MatrixGenerics_1.2.0 ProtGenerics_1.22.0
[53] SummarizedExperiment_1.20.0 AnnotationFilter_1.14.0
[55] RColorBrewer_1.1-2 yaml_2.2.1
[57] curl_4.3 memoise_1.1.0
[59] gridExtra_2.3 ggplot2_3.3.2
[61] biomaRt_2.46.0 rpart_4.1-15
[63] latticeExtra_0.6-29 stringi_1.5.3
[65] RSQLite_2.2.1 e1071_1.7-4
[67] checkmate_2.0.0 BiocParallel_1.24.0
[69] truncnorm_1.0-8 rlang_0.4.8
[71] pkgconfig_2.0.3 matrixStats_0.57.0
[73] bitops_1.0-6 Rsolnp_1.16
[75] evaluate_0.14 lattice_0.20-41
[77] purrr_0.3.4 GenomicAlignments_1.26.0
[79] htmlwidgets_1.5.2 bit_4.0.4
[81] tidyselect_1.1.0 magrittr_1.5
[83] bookdown_0.21 R6_2.4.1
[85] generics_0.0.2 Hmisc_4.4-1
[87] DelayedArray_0.16.0 DBI_1.1.0
[89] pillar_1.4.6 foreign_0.8-80
[91] BSgenome.Drerio.UCSC.danRer7_1.4.0 survival_3.2-7
[93] RCurl_1.98-1.2 nnet_7.3-14
[95] tibble_3.0.4 crayon_1.3.4
[97] BiocFileCache_1.14.0 rmarkdown_2.5
[99] jpeg_0.1-8.1 progress_1.2.2
[101] grid_4.0.3 data.table_1.13.2
[103] blob_1.2.1 digest_0.6.27
[105] openssl_1.4.3 munsell_0.5.0
[107] Gviz_1.34.0 askpass_1.1

1. Sheppard, S., Lawson, N. D. & Zhu, L. J. Accurate identification of polyadenylation sites from 3′ end deep sequencing using a naive bayes classifier. Bioinformatics 29, 2564–2571 (2013).