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 |
##                                <Rle>           <IRanges>  <Rle> |
##      uc001bum.2_5|IQCC|utr3     chr1   32673684-32674288      + |
##    uc001fbq.3_3|S100A9|utr3     chr1 153333315-153333503      + |
##    uc001gde.2_2|LRRC52|utr3     chr1 165533062-165533185      + |
##   uc001hfg.3_15|PFKFB2|utr3     chr1 207245717-207251162      + |
##   uc001hfh.3_15|PFKFB2|utr3     chr1 207252365-207254368      + |
##                         ...      ...                 ...    ... .
##      uc004dsv.3_19|PHF8|CDS     chrX   53964414-53964492      - |
##      uc004dsx.3_15|PHF8|CDS     chrX   53969797-53969835      - |
##     uc004ehz.1_5|ARMCX3|CDS     chrX 100879970-100881109      + |
##    uc004elw.3_6|FAM199X|CDS     chrX 103434289-103434459      + |
##      uc004fmj.1_10|GAB3|CDS     chrX 153906455-153906571      - |
##                                 feature annotatedProximalCP          exon
##                             <character>         <character>   <character>
##      uc001bum.2_5|IQCC|utr3        utr3             unknown  uc001bum.2_5
##    uc001fbq.3_3|S100A9|utr3        utr3             unknown  uc001fbq.3_3
##    uc001gde.2_2|LRRC52|utr3        utr3             unknown  uc001gde.2_2
##   uc001hfg.3_15|PFKFB2|utr3        utr3             unknown uc001hfg.3_15
##   uc001hfh.3_15|PFKFB2|utr3        utr3             unknown uc001hfh.3_15
##                         ...         ...                 ...           ...
##      uc004dsv.3_19|PHF8|CDS         CDS             unknown uc004dsv.3_19
##      uc004dsx.3_15|PHF8|CDS         CDS             unknown uc004dsx.3_15
##     uc004ehz.1_5|ARMCX3|CDS         CDS             unknown  uc004ehz.1_5
##    uc004elw.3_6|FAM199X|CDS         CDS             unknown  uc004elw.3_6
##      uc004fmj.1_10|GAB3|CDS         CDS             unknown uc004fmj.1_10
##                              transcript        gene      symbol truncated
##                             <character> <character> <character> <logical>
##      uc001bum.2_5|IQCC|utr3  uc001bum.2       55721        IQCC     FALSE
##    uc001fbq.3_3|S100A9|utr3  uc001fbq.3        6280      S100A9     FALSE
##    uc001gde.2_2|LRRC52|utr3  uc001gde.2      440699      LRRC52     FALSE
##   uc001hfg.3_15|PFKFB2|utr3  uc001hfg.3        5208      PFKFB2     FALSE
##   uc001hfh.3_15|PFKFB2|utr3  uc001hfh.3        5208      PFKFB2     FALSE
##                         ...         ...         ...         ...       ...
##      uc004dsv.3_19|PHF8|CDS  uc004dsv.3       23133        PHF8     FALSE
##      uc004dsx.3_15|PHF8|CDS  uc004dsx.3       23133        PHF8     FALSE
##     uc004ehz.1_5|ARMCX3|CDS  uc004ehz.1       51566      ARMCX3     FALSE
##    uc004elw.3_6|FAM199X|CDS  uc004elw.3      139231     FAM199X     FALSE
##      uc004fmj.1_10|GAB3|CDS  uc004fmj.1      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
##                              <character> <character> <logical>
##      uc009daz.2_10|Mitf|utr3       17342        Mitf     FALSE
##      uc009dhz.2_19|Atg7|utr3       74244        Atg7     FALSE
##      uc009dmb.2_4|Lrtm2|utr3      211187       Lrtm2     FALSE
##   uc009eet.1_3|BC035044|utr3      232406    BC035044     FALSE
##                                     fit_value Predicted_Proximal_APA
##                                     <numeric>              <numeric>
##      uc009daz.2_10|Mitf|utr3 12594.6737481599               98018978
##      uc009dhz.2_19|Atg7|utr3 27383.4125399619              114859674
##      uc009dmb.2_4|Lrtm2|utr3 168.550899294257              119316541
##   uc009eet.1_3|BC035044|utr3 123.973723777853              128849177
##                              Predicted_Distal_APA           type utr3start
##                                         <numeric>    <character> <numeric>
##      uc009daz.2_10|Mitf|utr3             98021358 novel proximal  98018276
##      uc009dhz.2_19|Atg7|utr3            114862071   novel distal 114859443
##      uc009dmb.2_4|Lrtm2|utr3            119315133 novel proximal 119316955
##   uc009eet.1_3|BC035044|utr3            128846744   novel distal 128849981
##                                utr3end
##                              <numeric>
##      uc009daz.2_10|Mitf|utr3  98021358
##      uc009dhz.2_19|Atg7|utr3 114860614
##      uc009dmb.2_4|Lrtm2|utr3 119315133
##   uc009eet.1_3|BC035044|utr3 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.6737481599
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE 27383.4125399619
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE 168.550899294257
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE 123.973723777853
##              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.5313380764804      1.05214965526127
##   uc009dhz.2 114859443 114860614      520.914790102864      172.152349179872
##   uc009dmb.2 119316955 119315133      8.85786481015595      49.5356159541359
##   uc009eet.1 128849981 128848044      21.6461685265195                     0
##              Baf3_long.form.usage UM15_long.form.usage         Baf3_PDUI
##                         <numeric>            <numeric>         <numeric>
##   uc009daz.2      282.82402351953     189.141117177656 0.894007365933945
##   uc009dhz.2      208.02460383653     456.544620517098 0.285379834820432
##   uc009dmb.2     8.53513129879347     70.8126330731015 0.490722314046962
##   uc009eet.1     8.90221857025472     23.2912900575185 0.291413701877399
##                      UM15_PDUI   short.mean.gp1    long.mean.gp1
##                      <numeric>        <numeric>        <numeric>
##   uc009daz.2 0.994467997354577 33.5313380764804  282.82402351953
##   uc009dhz.2 0.726175952044354 520.914790102864  208.02460383653
##   uc009dmb.2 0.588397701216867 8.85786481015595 8.53513129879347
##   uc009eet.1                 1 21.6461685265195 8.90221857025472
##                short.mean.gp2    long.mean.gp2          PDUI.gp1
##                     <numeric>        <numeric>         <numeric>
##   uc009daz.2 1.05214965526127 189.141117177656 0.894007365933945
##   uc009dhz.2 172.152349179872 456.544620517098 0.285379834820432
##   uc009dmb.2 49.5356159541359 70.8126330731015 0.490722314046962
##   uc009eet.1                0 23.2912900575185 0.291413701877399
##                       PDUI.gp2               dPDUI              logFC
##                      <numeric>           <numeric>          <numeric>
##   uc009daz.2 0.994467997354577  -0.100460631420632 -0.153638226758031
##   uc009dhz.2 0.726175952044354  -0.440796117223922  -1.34743575844534
##   uc009dmb.2 0.588397701216867 -0.0976753871699045 -0.261884735309782
##   uc009eet.1                 1  -0.708586298122601  -1.77885938232419
##                           P.Value            adj.P.Val
##                         <numeric>            <numeric>
##   uc009daz.2 1.58697962749085e-06 2.11597283665446e-06
##   uc009dhz.2  1.7094761107206e-60  6.8379044428824e-60
##   uc009dmb.2    0.593030861268587    0.593030861268587
##   uc009eet.1 4.13501327860387e-08 8.27002655720775e-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.6737481599
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE 27383.4125399619
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE 168.550899294257
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE 123.973723777853
##              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.5313380764804      1.05214965526127
##   uc009dhz.2 114859443 114860614      520.914790102864      172.152349179872
##   uc009dmb.2 119316955 119315133      8.85786481015595      49.5356159541359
##   uc009eet.1 128849981 128848044      21.6461685265195                     0
##              Baf3_long.form.usage UM15_long.form.usage         Baf3_PDUI
##                         <numeric>            <numeric>         <numeric>
##   uc009daz.2      282.82402351953     189.141117177656 0.894007365933945
##   uc009dhz.2      208.02460383653     456.544620517098 0.285379834820432
##   uc009dmb.2     8.53513129879347     70.8126330731015 0.490722314046962
##   uc009eet.1     8.90221857025472     23.2912900575185 0.291413701877399
##                      UM15_PDUI   short.mean.gp1    long.mean.gp1
##                      <numeric>        <numeric>        <numeric>
##   uc009daz.2 0.994467997354577 33.5313380764804  282.82402351953
##   uc009dhz.2 0.726175952044354 520.914790102864  208.02460383653
##   uc009dmb.2 0.588397701216867 8.85786481015595 8.53513129879347
##   uc009eet.1                 1 21.6461685265195 8.90221857025472
##                short.mean.gp2    long.mean.gp2          PDUI.gp1
##                     <numeric>        <numeric>         <numeric>
##   uc009daz.2 1.05214965526127 189.141117177656 0.894007365933945
##   uc009dhz.2 172.152349179872 456.544620517098 0.285379834820432
##   uc009dmb.2 49.5356159541359 70.8126330731015 0.490722314046962
##   uc009eet.1                0 23.2912900575185 0.291413701877399
##                       PDUI.gp2               dPDUI              logFC
##                      <numeric>           <numeric>          <numeric>
##   uc009daz.2 0.994467997354577  -0.100460631420632 -0.153638226758031
##   uc009dhz.2 0.726175952044354  -0.440796117223922  -1.34743575844534
##   uc009dmb.2 0.588397701216867 -0.0976753871699045 -0.261884735309782
##   uc009eet.1                 1  -0.708586298122601  -1.77885938232419
##                           P.Value            adj.P.Val      PASS
##                         <numeric>            <numeric> <logical>
##   uc009daz.2 1.58697962749085e-06 2.11597283665446e-06     FALSE
##   uc009dhz.2  1.7094761107206e-60  6.8379044428824e-60      TRUE
##   uc009dmb.2    0.593030861268587    0.593030861268587     FALSE
##   uc009eet.1 4.13501327860387e-08 8.27002655720775e-08      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 15 with logLik: -553.6445 
## 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.7079181186
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE 7630.49809358993
##   uc009die.2  uc009die.2      232334       Vgll4     FALSE 10704.6989063837
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE  18.635065106255
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE  227.54183333935
##   uc009eol.1  uc009eol.1       11569       Aebp2      TRUE 204539.698420888
##              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
##              <numeric> <numeric>
##   uc009daz.2  98018176  98021358
##   uc009dhz.2 114859343 114860614
##   uc009die.2 114862164 114862092
##   uc009dmb.2 119317055 119315133
##   uc009eet.1 128850081 128848044
##   uc009eol.1 140651362 140651622
##                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                       data2
##                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                      <list>
##   uc009daz.2                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                        c(473, 460, 457, 450, 442, 436, 433, 435, 435, 431, 421, 420, 422, 423, 422, 422, 416, 417, 416, 417, 412, 416, 421, 424, 418, 408, 417, 406, 414, 415, 411, 402, 403, 393, 362, 348, 295, 266, 232, 230, 228, 224, 219, 216, 208, 206, 206, 205, 196, 196, 198, 197, 197, 199, 200, 197, 209, 213, 212, 211, 210, 213, 203, 204, 209, 207, 215, 216, 216, 215, 215, 215, 215, 206, 206, 206, 203, 205, 205, 206, 205, 203, 199, 199, 195, 193, 189, 196, 194, 194, 201, 203, 204, 202, 209, 211, 214, 220, 220, 220, \n221, 221, 221, 216, 215, 214, 214, 206, 205, 205, 208, 208, 203, 208, 208, 204, 204, 205, 205, 206, 204, 201, 198, 200, 198, 236, 224, 211, 203, 200, 203, 198, 195, 199, 199, 197, 196, 198, 198, 196, 196, 198, 200, 202, 201, 203, 214, 217, 217, 218, 218, 219, 220, 216, 215, 215, 202, 200, 200, 206, 206, 204, 201, 201, 196, 209, 201, 200, 204, 208, 208, 208, 208, 208, 208, 208, 208, 209, 209, 209, 224, 225, 225, 225, 226, 226, 226, 221, 221, 225, 218, 218, 218, 219, 213, 212, 213, 206, 209, 209, 213, \n235, 236, 234, 239, 246, 248, 262, 271, 282, 280, 282, 290, 284, 285, 285, 289, 299, 300, 301, 301, 301, 300, 296, 296, 254, 252, 252, 252, 259, 267, 266, 268, 268, 277, 282, 283, 289, 292, 291, 295, 294, 301, 303, 304, 306, 308, 307, 308, 324, 323, 342, 346, 347, 356, 355, 356, 358, 358, 352, 351, 351, 350, 350, 351, 338, 338, 345, 341, 337, 338, 340, 340, 340, 340, 339, 340, 336, 336, 334, 319, 317, 317, 315, 319, 319, 321, 319, 319, 315, 319, 316, 315, 314, 313, 319, 344, 360, 365, 363, 382, 366, \n366, 361, 362, 353, 356, 342, 334, 323, 322, 319, 313, 317, 325, 327, 323, 317, 324, 321, 321, 317, 321, 325, 325, 325, 329, 332, 328, 334, 323, 323, 322, 328, 320, 317, 318, 317, 317, 313, 323, 319, 313, 310, 318, 319, 306, 306, 306, 288, 288, 276, 273, 274, 268, 268, 269, 266, 266, 266, 266, 265, 266, 266, 269, 271, 270, 264, 264, 271, 270, 268, 274, 274, 287, 287, 286, 286, 286, 286, 286, 292, 293, 294, 291, 292, 290, 290, 291, 294, 292, 292, 294, 293, 300, 302, 275, 258, 250, 251, 229, 226, 225, \n225, 224, 225, 225, 229, 229, 231, 235, 235, 233, 230, 221, 220, 227, 221, 214, 214, 216, 215, 212, 206, 206, 206, 202, 198, 198, 197, 196, 204, 209, 202, 225, 226, 226, 219, 220, 220, 212, 223, 239, 242, 234, 234, 233, 233, 237, 241, 249, 252, 251, 257, 254, 255, 255, 266, 270, 282, 282, 283, 282, 282, 279, 279, 286, 298, 307, 310, 312, 328, 328, 334, 321, 323, 343, 349, 359, 364, 371, 365, 367, 381, 386, 387, 388, 388, 387, 384, 383, 387, 392, 390, 397, 391, 392, 392, 407, 407, 413, 411, 412, 413, \n409, 408, 408, 413, 417, 417, 414, 434, 438, 439, 440, 443, 435, 443, 457, 458, 462, 464, 462, 470, 470, 472, 472, 481, 482, 476, 490, 481, 481, 477, 452, 451, 451, 453, 450, 454, 453, 444, 428, 427, 427, 426, 430, 440, 442, 439, 434, 422, 425, 418, 425, 427, 439, 431, 429, 425, 428, 427, 431, 430, 436, 434, 431, 420, 411, 404, 406, 394, 390, 389, 389, 386, 366, 360, 352, 347, 342, 348, 345, 337, 332, 330, 337, 337, 338, 338, 338, 340, 357, 364, 357, 359, 392, 402, 388, 396, 388, 393, 401, 400, 400, \n400, 399, 393, 390, 394, 395, 381, 377, 376, 378, 380, 386, 382, 368, 367, 362, 361, 361, 355, 356, 373, 375, 368, 366, 375, 362, 361, 360, 366, 371, 385, 384, 383, 384, 379, 379, 377, 377, 376, 380, 403, 397, 391, 403, 404, 406, 417, 417, 425, 420, 417, 407, 413, 415, 428, 447, 451, 450, 453, 466, 467, 467, 467, 469, 466, 466, 460, 462, 457, 464, 468, 468, 469, 469, 479, 501, 502, 509, 502, 502, 504, 497, 498, 497, 499, 502, 502, 478, 472, 465, 459, 425, 418, 421, 413, 413, 407, 400, 405, 405, 405, \n408, 404, 403, 401, 403, 396, 397, 401, 398, 393, 388, 388, 386, 386, 385, 387, 386, 391, 394, 376, 374, 373, 376, 363, 363, 364, 361, 356, 353, 341, 356, 359, 359, 359, 362, 362, 369, 366, 363, 337, 339, 342, 326, 329, 328, 333, 331, 321, 323, 323, 319, 310, 307, 289, 271, 268, 272, 272, 252, 251, 247, 248, 246, 247, 243, 243, 241, 272, 265, 265, 269, 268, 267, 259, 234, 228, 221, 221, 219, 219, 219, 220, 222, 220, 216, 210, 210, 209, 211, 211, 213, 211, 211, 215, 215, 215, 215, 210, 210, 211, 207, \n207, 210, 208, 205, 206, 205, 205, 205, 205, 204, 203, 209, 209, 210, 208, 218, 232, 231, 236, 240, 262, 292, 298, 297, 321, 320, 322, 321, 316, 303, 301, 299, 304, 303, 304, 299, 316, 318, 318)
##   uc009dhz.2                                                                                                                                                                                                                          c(184, 184, 184, 184, 187, 187, 188, 190, 196, 201, 208, 212, 212, 216, 217, 218, 221, 221, 222, 222, 222, 228, 228, 228, 231, 235, 236, 236, 241, 241, 243, 244, 245, 245, 245, 245, 247, 242, 243, 245, 246, 243, 248, 251, 251, 250, 236, 235, 236, 232, 224, 223, 226, 228, 227, 228, 228, 228, 232, 235, 235, 235, 234, 234, 230, 230, 235, 236, 234, 234, 245, 245, 244, 239, 239, 239, 240, 242, 231, 227, 227, 223, 219, 215, 214, 214, 213, 189, 183, 178, 180, 177, 172, 169, 169, 174, 180, 179, 177, 175, \n177, 180, 180, 183, 180, 181, 181, 181, 176, 173, 168, 168, 170, 168, 168, 171, 172, 174, 174, 175, 181, 177, 177, 182, 181, 178, 178, 184, 187, 187, 185, 186, 185, 187, 201, 201, 198, 203, 202, 198, 196, 193, 188, 186, 188, 188, 190, 190, 189, 188, 186, 188, 188, 187, 194, 194, 197, 199, 199, 198, 198, 205, 212, 213, 211, 211, 205, 206, 206, 206, 195, 209, 207, 206, 207, 207, 208, 206, 206, 209, 209, 208, 212, 216, 234, 234, 232, 235, 238, 248, 246, 251, 256, 255, 261, 257, 248, 247, 246, 244, 242, \n239, 239, 236, 242, 241, 241, 239, 244, 247, 243, 240, 241, 247, 250, 250, 251, 257, 270, 271, 277, 282, 295, 289, 295, 295, 295, 289, 282, 284, 285, 283, 292, 292, 280, 281, 283, 281, 281, 282, 285, 289, 297, 300, 302, 303, 303, 305, 305, 306, 307, 302, 303, 302, 295, 302, 300, 298, 295, 293, 293, 286, 280, 279, 278, 279, 278, 276, 277, 277, 291, 278, 285, 290, 289, 291, 291, 305, 306, 303, 303, 306, 302, 299, 282, 283, 291, 288, 283, 274, 274, 269, 263, 267, 261, 266, 266, 264, 275, 278, 279, 279, \n282, 284, 278, 278, 278, 280, 281, 281, 293, 293, 293, 285, 281, 288, 283, 275, 262, 266, 256, 255, 248, 248, 242, 241, 244, 244, 249, 247, 252, 252, 244, 242, 249, 248, 247, 252, 252, 251, 250, 251, 244, 239, 235, 234, 233, 256, 259, 259, 260, 260, 255, 255, 259, 255, 255, 258, 261, 265, 269, 273, 271, 271, 277, 276, 279, 279, 280, 284, 286, 304, 298, 293, 292, 288, 286, 269, 268, 268, 268, 277, 277, 279, 285, 284, 278, 279, 283, 286, 286, 286, 292, 288, 292, 286, 286, 286, 279, 279, 280, 280, 282, \n280, 286, 289, 288, 292, 285, 286, 283, 287, 285, 285, 300, 289, 289, 290, 290, 289, 288, 290, 285, 287, 286, 287, 284, 287, 281, 283, 277, 278, 279, 281, 298, 300, 307, 302, 310, 315, 331, 330, 329, 329, 329, 335, 338, 313, 315, 318, 316, 316, 320, 320, 316, 314, 313, 312, 308, 304, 300, 303, 303, 303, 298, 298, 295, 295, 293, 289, 276, 257, 256, 256, 256, 256, 256, 257, 267, 267, 268, 256, 256, 255, 248, 247, 245, 244, 240, 236, 236, 237, 231, 231, 231, 231, 228, 228, 224, 218, 218, 220, 219, 219, \n217, 214, 214, 208, 208, 205, 196, 191, 190, 190, 175, 175, 177, 176, 176, 172, 170, 161, 160, 157, 158, 157, 157, 154, 154, 152, 152, 151, 149, 147, 121, 119, 111, 106, 98, 93, 75, 71, 71, 71, 69, 63, 61, 61, 56, 52, 52, 52, 48, 48, 48, 46, 46, 44, 44, 44, 44, 37, 37, 37, 36, 36, 36, 36, 36, 36, 33, 33, 33, 33, 33, 33, 33, 32, 22, 22, 21, 21, 21, 19, 19, 19, 19, 19, 19, 19, 19, 18, 18, 20, 20, 20, 21, 21, 21, 21, 19, 19, 17, 20, 16, 16, 25, 25, 25, 24, 25, 25, 27, 27, 27, 28, 27, 28, 28, 29, 33, \n35, 43, 45, 44, 44, 48, 56, 56, 58, 60, 77, 79, 79, 80, 82, 90, 88, 91, 21, 84, 92, 92, 92, 94, 95, 100, 101, 101, 101, 99, 101, 101, 102, 102, 104, 104, 104, 105, 109, 109, 112, 115, 115, 116, 117, 118, 122, 122, 122, 122, 125, 131, 132, 134, 141, 142, 142, 152, 157, 158, 164, 164, 173, 173, 177, 188, 188, 192, 192, 193, 196, 200, 202, 204, 203, 213, 214, 218, 226, 227, 227, 229, 226, 224, 224, 227, 220, 225, 236, 237, 246, 253, 256, 271, 298, 302, 313, 313, 317, 328, 324, 325, 316, 329, 329, 332, \n328, 320, 322, 321, 322, 307, 306, 311, 327, 317, 317, 322, 331, 333, 332, 332, 334, 341, 346, 347, 347, 348, 348, 354, 354, 360, 374, 378, 376, 377, 379, 379, 383, 391, 388, 388, 392, 392, 391, 395, 391, 395, 400, 402, 403, 411, 411, 409, 408, 416, 435, 446, 457, 459, 460, 461, 452, 453, 450, 439, 439, 435, 435, 434, 431, 427, 425, 421, 420, 410, 408, 404, 396, 395, 395, 393, 392, 391, 391, 388, 386, 381, 370, 369, 359, 352, 347, 332, 305, 300, 288, 287, 283, 271, 269, 267, 267, 254, 254, 251, 251, \n251, 248, 247, 245, 243, 241, 237, 223, 221, 221, 215, 206, 204, 204, 204, 201, 193, 187, 180, 180, 179, 179, 173, 173, 167, 152, 148, 148, 147, 145, 144, 136, 128, 128, 126, 122, 121, 121, 116, 116, 112, 107, 105, 101, 87, 86, 86, 80, 69, 53, 32, 16, 11, 4, 6, 16, 17, 20, 22, 22, 24, 24, 29, 36, 36, 38, 42, 45, 45, 45, 49, 49, 50, 57, 59, 59, 59, 60, 60, 60, 60, 60, 60, 61, 62, 63, 70, 70, 70, 70, 73, 73, 73, 73, 79, 83, 85, 86, 99, 99, 99, 100, 100, 101, 109, 109, 109, 113, 114, 119, 124, 128)
##   uc009die.2                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             c(76, 76, 76, 76, 76, 86, 86, 90, 90, 91, 94, 97, 102, 103, 106, 109, 125, 126, 128, 128, 128, 128, 129, 132, 132, 132, 133, 133, 134, 136, 143, 143, 143, 144, 156, 156, 158, 159, 169, 169, 172, 172, 172, 172, 173, 173, 175, 175, 175, 175, 178, 183, 183, 186, 186, 189, 190, 191, 191, 191, 191, 191, 191, 191, 191, 192, 193, 193, 195, 195, 195, 195, 195)
##   uc009dmb.2                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                  c(11, 12, 12, 11, 11, 11, 11, 11, 11, 11, 11, 11, 7, 7, 7, 7, 7, 7, 8, 8, 7, 7, 7, 7, 7, 7, 7, 7, 7, 9, 9, 9, 9, 9, 9, 9, 11, 11, 11, 11, 13, 14, 17, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 21, 21, 18, 18, 18, 18, 18, 18, 18, 18, 20, 22, 22, 22, 22, 22, 23, 23, 23, 23, 23, 23, 23, 23, 23, 22, 22, 22, 22, 22, 25, 25, 25, 25, 25, 25, 25, 25, 25, 25, 25, 25, 24, 24, 25, 25, 25, 25, 27, 27, 27, 27, 28, 28, 28, 28, 29, 29, 29, 28, 28, 27, 27, 27, 27, 27, 27, 28, 28, 28, 26, 26, \n26, 26, 26, 26, 26, 24, 24, 24, 24, 22, 21, 18, 17, 17, 17, 17, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 18, 15, 21, 21, 21, 21, 21, 21, 21, 21, 21, 19, 17, 17, 17, 17, 17, 16, 16, 16, 16, 16, 16, 16, 17, 17, 17, 17, 17, 17, 17, 14, 14, 14, 14, 16, 16, 16, 16, 16, 16, 16, 16, 16, 16, 15, 15, 15, 15, 13, 13, 13, 13, 12, 12, 12, 12, 11, 11, 11, 11, 11, 12, 12, 12, 12, 12, 12, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 13, 13, 13, 13, 13, 14, 13, 13, 13, 13, 13, 13, 19, 20, 20, 20, \n20, 20, 20, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 14, 14, 14, 14, 14, 14, 14, 17, 17, 17, 19, 17, 17, 17, 17, 17, 17, 18, 20, 20, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 22, 23, 23, 23, 23, 23, 24, 23, 23, 23, 24, 24, 25, 25, 25, 25, 25, 25, 25, 25, 25, 25, 25, 26, 26, 27, 27, 27, 25, 25, 25, 25, 26, 25, 25, 25, 25, 25, 27, 27, 21, 20, 22, 22, 22, 22, 22, 22, 22, 22, 23, 23, 23, 23, 23, 23, 23, 22, 22)
##   uc009eet.1 c(22, 22, 22, 22, 22, 21, 21, 21, 22, 22, 22, 22, 20, 20, 20, 19, 19, 19, 19, 19, 19, 21, 25, 25, 25, 24, 24, 24, 24, 24, 24, 22, 22, 22, 22, 21, 21, 21, 19, 19, 19, 19, 19, 19, 19, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 12, 12, 12, 12, 12, 12, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 11, 10, 10, 8, 8, 8, 7, 7, 7, 7, 7, 8, 8, 8, 8, 8, 8, 8, 8, 7, 7, 7, 7, 7, 7, 7, 7, 7, 7, 7, 7, 8, 6, 2, 2, 2, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, \n3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 4, 4, 4, 4, 4, 4, 4, 4, 4, 5, 5, 5, 5, 5, 5, 5, 5, 5, 4, 4, 4, 4, 4, 4, 4, 4, 4, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 4, 4, 4, 4, 4, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 2, 2, 2, 2, 6, 6, 6, 7, 9, 8, 8, 8, 8, 8, 8, 10, 10, 11, 16, 16, \n24, 24, 25, 26, 26, 27, 27, 26, 26, 26, 26, 30, 30, 30, 30, 30, 34, 34, 34, 40, 40, 40, 40, 41, 42, 42, 42, 42, 42, 42, 42, 42, 43, 44, 44, 44, 44, 44, 44, 44, 44, 44, 44, 44, 44, 44, 44, 44, 44, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 46, 47, 48, 49, 49, 49, 49, 49, 49, 49, 49, 52, 52, 53, 53, 53, 53, 53, 53, 53, 57, 57, 57, 57, 55, 55, 54, 51, 52, 44, 44, 43, 42, 45, 44, 44, 44, 44, 44, 44, 40, 40, 41, 42, 42, 39, 39, 39, 36, 37, 38, 38, 37, 43, 43, \n43, 43, 43, 44, 44, 44, 43, 42, 42, 42, 42, 47, 50, 50, 50, 50, 50, 50, 50, 50, 50, 50, 50, 48, 48, 53, 53, 57, 57, 57, 60, 60, 60, 62, 62, 62, 68, 68, 68, 68, 68, 70, 70, 70, 70, 70, 69, 68, 67, 67, 67, 68, 69, 69, 69, 69, 67, 67, 62, 62, 62, 61, 60, 60, 60, 56, 56, 56, 56, 56, 57, 57, 55, 54, 54, 55, 55, 55, 52, 52, 52, 52, 53, 53, 53, 56, 56, 55, 54, 54, 53, 55, 56, 53, 52, 51, 51, 51, 44, 45, 46, 46, 46, 48, 48, 48, 48, 48, 50, 51, 52, 47, 44, 44, 46, 46, 46, 46, 46, 46, 46, 46, 46, 49, 49, 44, \n46, 42, 42, 43, 42, 42, 42, 40, 40, 40, 35, 35, 35, 38, 38, 36, 36, 36, 36, 36, 36, 36, 36, 36, 36, 35, 34, 34, 36, 36, 44, 44, 44, 44, 44, 44, 46, 47, 47, 48, 48, 50, 50, 50, 49, 49, 49, 50, 50, 49, 49, 49, 49, 49, 49, 49, 48, 48, 48, 45, 47, 48, 48, 48, 48, 46, 45, 48, 48, 48, 48, 48, 48, 47, 46, 48, 48, 46, 46, 46, 46, 48, 46, 45, 44, 44, 44, 44, 42, 42, 42, 42, 42, 42, 42, 42, 42, 39, 39, 39, 37, 37, 37, 36, 34, 37, 37, 37, 37, 37, 36, 36, 36, 33, 33, 33, 33, 33, 33, 34, 34, 35, 35, 35, 35, 35, \n33, 35, 33, 33, 24, 24, 24, 24, 24, 24, 23, 21, 21, 20, 20, 18, 18, 18, 18, 18, 18, 17, 18, 18, 18, 18, 18, 18, 18, 18, 23, 23, 23, 23, 21, 20, 20, 20, 20, 20, 20, 17, 17, 17, 17, 17, 17, 17, 17, 15, 15, 14, 14, 14, 14, 12, 12, 12, 12, 12, 12, 12, 12, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 10, 10, 10, 10, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 12, 13, 12, 12, 12, 12, 13, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 14, 14, 14, 14, 14, 14, 14, 16, 16, 16, 16, 15, 15, 16, 16, \n16, 17, 17, 17, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 12, 12, 12, 13, 18, 20, 20, 20, 20, 20, 20, 20, 20, 20, 23, 23, 23, 23, 23, 23, 20, 20, 23, 24, 24, 24, 26, 26, 26, 26, 26, 25, 25, 25, 25, 26, 25, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 21, 21, 21, 21, 23, 24, 23, 23, 24, 23, 23, 23, 23, 23, 23, 23, 23, 23, 25, 25, 25, 25, 25, 25, 26, 26, 26, 26, 26, 25, 25, 25, 25, 28, \n28, 28, 28, 28, 28, 28, 31, 31, 32, 32, 32, 32, 32, 32, 31, 26, 24, 24, 24, 24, 34, 36, 39, 39, 39, 37, 38, 38, 38, 38, 39, 39, 39, 39, 44, 46, 46, 45, 45, 46, 46, 46, 46, 46, 46, 47, 48, 48, 48, 48, 48, 49, 49, 49, 50, 50, 50, 50, 50, 50, 50, 50, 50, 50, 50, 50, 52, 52, 52, 53, 51, 50, 54, 54, 53, 56, 56, 56, 56, 56, 56, 56, 56, 56, 56, 57, 57, 57, 57, 57, 56, 57, 57, 57, 57, 57, 57, 57, 57, 54, 54, 54, 54, 54, 54, 54, 51, 51, 53, 58, 59, 59, 59, 59, 59, 59, 59, 59, 59, 59, 49, 47, 46, 47, 47, 51, \n50, 50, 50, 53, 59, 60, 60, 57, 51, 52, 52, 51, 51, 50, 50, 50, 50, 50, 50, 53, 55, 58, 58, 58, 59, 58, 58, 58, 57, 57, 57, 57, 57, 57, 57, 57, 57, 61, 61, 62, 60, 60, 60, 59, 59, 59, 55, 55, 55, 52, 52, 52, 52, 52, 52, 53, 55, 56, 57, 56, 56, 57, 60, 60, 60, 59, 59, 59, 59, 62, 62, 62, 62, 66, 66, 66, 66, 66, 66, 66, 66, 66, 63, 58, 57, 57, 57, 57, 57, 57, 57, 62, 62, 63, 63, 63, 61, 62, 62, 57, 57, 58, 58, 56, 49, 48, 48, 49, 52, 49, 50, 50, 50, 51, 51, 51, 55, 55, 55, 51, 47, 45, 45, 45, 44, 44, \n45, 45, 45, 45, 45, 45, 45, 45, 45, 45, 45, 43, 43, 42, 43, 47, 50, 50, 50, 50, 50, 50, 51, 51, 51, 51, 51, 51, 51, 52, 50, 50, 47, 47, 47, 46, 43, 45, 45, 45, 45, 45, 45, 42, 42, 42, 42, 40, 41, 41, 41, 41, 41, 41)
##   uc009eol.1                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                c(1408, 1412, 1384, 1385, 1398, 1399, 1385, 1408, 1416, 1460, 1513, 1544, 1569, 1623, 1677, 1685, 1700, 1712, 1717, 1726, 1738, 1742, 1754, 1734, 1727, 1685, 1695, 1699, 1709, 1689, 1671, 1728, 1732, 1719, 1727, 1724, 1739, 1744, 1743, 1766, 1756, 1787, 1779, 1819, 1839, 1980, 1995, 1998, 1995, 1990, 1958, 1942, 1928, 1922, 1919, 1930, 1910, 1907, 1895, 1885, 1885, 1886, 1886, 1896, 1884, 1870, 1863, 1825, 1846, 1842, 1825, 1797, 1675, 1668, 1659, 1659, 1657, 1655, 1645, 1635, 1617, 1595, 1577, 1579, \n1578, 1562, 1557, 1536, 1506, 1495, 1486, 1482, 1480, 1488, 1454, 1455, 1435, 1426, 1407, 1387, 1352, 1357, 1352, 1389, 1369, 1363, 1377, 1365, 1376, 1377, 1322, 1315, 1288, 1239, 1174, 1163, 1151, 1137, 1130, 1119, 1128, 1128, 1128, 1141, 1142, 1143, 1155, 1187, 1174, 1182, 1180, 1129, 1125, 1150, 1139, 1129, 1106, 1100, 1096, 1075, 1052, 1017, 992, 957, 942, 765, 755, 761, 765, 767, 769, 771, 776, 787, 788, 773, 772, 767, 760, 760, 756, 755, 753, 743, 743, 742, 740, 742, 710, 707, 699, 694, 689, \n689, 687, 684, 684, 679, 679, 664, 657, 657, 648, 634, 634, 633, 632, 631, 629, 629, 629, 628, 628, 614, 614, 603, 601, 601, 594, 590, 576, 561, 545, 502)
##              Baf3_short.form.usage Baf3_long.form.usage         Baf3_PDUI
##                          <numeric>            <numeric>         <numeric>
##   uc009daz.2      29.6602572856635     283.389388104407 0.905253822444981
##   uc009dhz.2       3.7254488495449     209.978806469604 0.982567268751943
##   uc009die.2        121.8881378455     174.122891566265 0.588231093659866
##   uc009dmb.2      9.22609762692124     9.11798839458414 0.497053294663731
##   uc009eet.1      22.0449523788087     10.3274224192527 0.319019611124456
##   uc009eol.1      1065.98013415893     221.456140350877 0.172013283092553
##                    short.mean        long.mean              PDUI P.Value
##                        <list>           <list>            <list>  <list>
##   uc009daz.2 29.6602572856635 283.389388104407 0.905253822444981       1
##   uc009dhz.2  3.7254488495449 209.978806469604 0.982567268751943       1
##   uc009die.2   121.8881378455 174.122891566265 0.588231093659866       1
##   uc009dmb.2 9.22609762692124 9.11798839458414 0.497053294663731       1
##   uc009eet.1 22.0449523788087 10.3274224192527 0.319019611124456       1
##   uc009eol.1 1065.98013415893 221.456140350877 0.172013283092553       1
##              adj.P.Val     dPDUI      PASS
##                 <list> <numeric> <logical>
##   uc009daz.2         1       NaN     FALSE
##   uc009dhz.2         1       NaN     FALSE
##   uc009die.2         1       NaN     FALSE
##   uc009dmb.2         1       NaN     FALSE
##   uc009eet.1         1       NaN     FALSE
##   uc009eol.1         1       NaN     FALSE
##   -------
##   seqinfo: 1 sequence from an unspecified genome; no seqlengths

3 Session Info

sessionInfo()

R version 3.6.1 (2019-07-05) Platform: x86_64-pc-linux-gnu (64-bit) Running under: Ubuntu 18.04.3 LTS

Matrix products: default BLAS: /home/biocbuild/bbs-3.10-bioc/R/lib/libRblas.so LAPACK: /home/biocbuild/bbs-3.10-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.24.0
[2] org.Hs.eg.db_3.10.0
[3] TxDb.Mmusculus.UCSC.mm10.knownGene_3.10.0 [4] BSgenome.Mmusculus.UCSC.mm10_1.4.0
[5] BSgenome_1.54.0
[6] rtracklayer_1.46.0
[7] Biostrings_2.54.0
[8] XVector_0.26.0
[9] InPAS_1.18.0
[10] GenomicFeatures_1.38.0
[11] AnnotationDbi_1.48.0
[12] GenomicRanges_1.38.0
[13] GenomeInfoDb_1.22.0
[14] IRanges_2.20.0
[15] S4Vectors_0.24.0
[16] Biobase_2.46.0
[17] BiocGenerics_0.32.0
[18] BiocStyle_2.14.0

loaded via a namespace (and not attached): [1] colorspace_1.4-1 seqinr_3.6-1
[3] class_7.3-15 depmixS4_1.4-0
[5] biovizBase_1.34.0 htmlTable_1.13.2
[7] base64enc_0.1-3 dichromat_2.0-0
[9] rstudioapi_0.10 bit64_0.9-7
[11] splines_3.6.1 knitr_1.25
[13] zeallot_0.1.0 ade4_1.7-13
[15] Formula_1.2-3 Rsamtools_2.2.0
[17] cluster_2.1.0 dbplyr_1.4.2
[19] BiocManager_1.30.9 compiler_3.6.1
[21] httr_1.4.1 backports_1.1.5
[23] assertthat_0.2.1 Matrix_1.2-17
[25] lazyeval_0.2.2 limma_3.42.0
[27] acepack_1.4.1 htmltools_0.4.0
[29] prettyunits_1.0.2 tools_3.6.1
[31] gtable_0.3.0 glue_1.3.1
[33] GenomeInfoDbData_1.2.2 dplyr_0.8.3
[35] rappdirs_0.3.1 Rcpp_1.0.2
[37] vctrs_0.2.0 preprocessCore_1.48.0
[39] nlme_3.1-141 xfun_0.10
[41] stringr_1.4.0 ensembldb_2.10.0
[43] XML_3.98-1.20 zlibbioc_1.32.0
[45] MASS_7.3-51.4 scales_1.0.0
[47] VariantAnnotation_1.32.0 hms_0.5.1
[49] ProtGenerics_1.18.0 SummarizedExperiment_1.16.0
[51] AnnotationFilter_1.10.0 RColorBrewer_1.1-2
[53] yaml_2.2.0 curl_4.2
[55] memoise_1.1.0 gridExtra_2.3
[57] ggplot2_3.2.1 biomaRt_2.42.0
[59] rpart_4.1-15 latticeExtra_0.6-28
[61] stringi_1.4.3 RSQLite_2.1.2
[63] e1071_1.7-2 checkmate_1.9.4
[65] BiocParallel_1.20.0 truncnorm_1.0-8
[67] rlang_0.4.1 pkgconfig_2.0.3
[69] matrixStats_0.55.0 bitops_1.0-6
[71] Rsolnp_1.16 evaluate_0.14
[73] lattice_0.20-38 purrr_0.3.3
[75] GenomicAlignments_1.22.0 htmlwidgets_1.5.1
[77] bit_1.1-14 tidyselect_0.2.5
[79] magrittr_1.5 bookdown_0.14
[81] R6_2.4.0 Hmisc_4.2-0
[83] DelayedArray_0.12.0 DBI_1.0.0
[85] pillar_1.4.2 foreign_0.8-72
[87] BSgenome.Drerio.UCSC.danRer7_1.4.0 survival_2.44-1.1
[89] RCurl_1.95-4.12 nnet_7.3-12
[91] tibble_2.1.3 crayon_1.3.4
[93] BiocFileCache_1.10.0 rmarkdown_1.16
[95] progress_1.2.2 grid_3.6.1
[97] data.table_1.12.6 blob_1.2.0
[99] digest_0.6.22 openssl_1.4.1
[101] munsell_0.5.0 Gviz_1.30.0
[103] 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).