Package version: InPAS 1.6.0

Contents

1 Introduction

Alternative polyadenylation (APA) is one of the most important post-transcriptional regulation mechanisms which occurs in most human genes. APA in a gene can result in altered expression of the gene, which can lead pathological effect to the cells, such as uncontrolled cell cycle and growth. However, there are only limited ways to identify and quantify APA in genes, and most of them suffers from complicated process for library construction and requires large amount of RNAs.

RNA-seq has become one of the most commonly used methods for quantifying genome-wide gene expression. There are massive RNA-seq datasets available publicly 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

To use InPAS, BSgenome and TxDb object have to be loaded before run.

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

Users can prepare annotaiton by utr3Annotation with a TxDb and org annotation. The 3UTR annotation prepared by utr3Annotation includes the gaps to next annotation 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 load mm10 and hg19 annotation from pre-prepared data. Here we loaded the prepared mm10 3UTR annotation file.

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

The coverage is calculated from BEDGraph file. The RNA-seq BAM files could be converted to BEDGraph files by bedtools genomecov tool (parameter: -bg -split). PWM and a classifier of polyA signal can be used for adjusting CP sites prediction.

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.01273315                    0.6582754         0.02317380
## UM15         0.01280700                    0.6621169         0.02337027
##      UTR3.expressed.gene.subset.coverage.rate
## Baf3                                0.8102631
## UM15                                0.8171327
##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.6737
##      uc009dhz.2_19|Atg7|utr3       74244        Atg7     FALSE 27383.4125
##      uc009dmb.2_4|Lrtm2|utr3      211187       Lrtm2     FALSE   168.5509
##   uc009eet.1_3|BC035044|utr3      232406    BC035044     FALSE   128.8275
##                              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              128849128            128846244
##                                        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 [128846244, 128850081]      - |             unknown
##               transcript        gene      symbol truncated  fit_value
##              <character> <character> <character> <logical>  <numeric>
##   uc009daz.2  uc009daz.2       17342        Mitf     FALSE 12594.6737
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE 27383.4125
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE   168.5509
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE   128.8275
##              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              128849128            128846244   novel distal
##              utr3start   utr3end Baf3_short.form.usage
##              <numeric> <numeric>             <numeric>
##   uc009daz.2  98018276  98021358             33.531338
##   uc009dhz.2 114859443 114860614            520.914790
##   uc009dmb.2 119316955 119315133              8.857865
##   uc009eet.1 128849981 128848044             22.671784
##              UM15_short.form.usage Baf3_long.form.usage
##                          <numeric>            <numeric>
##   uc009daz.2             1.0521497           282.824024
##   uc009dhz.2           172.1523492           208.024604
##   uc009dmb.2            49.5356160             8.535131
##   uc009eet.1             0.7334097             7.876603
##              UM15_long.form.usage Baf3_PDUI UM15_PDUI short.mean.gp1
##                         <numeric> <numeric> <numeric>      <numeric>
##   uc009daz.2            189.14112 0.8940074 0.9944680      33.531338
##   uc009dhz.2            456.54462 0.2853798 0.7261760     520.914790
##   uc009dmb.2             70.81263 0.4907223 0.5883977       8.857865
##   uc009eet.1             21.96014 0.2578402 0.9676820      22.671784
##              long.mean.gp1 short.mean.gp2 long.mean.gp2  PDUI.gp1
##                  <numeric>      <numeric>     <numeric> <numeric>
##   uc009daz.2    282.824024      1.0521497     189.14112 0.8940074
##   uc009dhz.2    208.024604    172.1523492     456.54462 0.2853798
##   uc009dmb.2      8.535131     49.5356160      70.81263 0.4907223
##   uc009eet.1      7.876603      0.7334097      21.96014 0.2578402
##               PDUI.gp2       dPDUI      logFC      P.Value    adj.P.Val
##              <numeric>   <numeric>  <numeric>    <numeric>    <numeric>
##   uc009daz.2 0.9944680 -0.10046063 -0.1536382 1.586980e-06 2.115973e-06
##   uc009dhz.2 0.7261760 -0.44079612 -1.3474358 1.709476e-60 6.837904e-60
##   uc009dmb.2 0.5883977 -0.09767539 -0.2618847 5.930309e-01 5.930309e-01
##   uc009eet.1 0.9676820 -0.70984179 -1.9080557 2.632879e-08 5.265757e-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 [128846244, 128850081]      - |             unknown
##               transcript        gene      symbol truncated  fit_value
##              <character> <character> <character> <logical>  <numeric>
##   uc009daz.2  uc009daz.2       17342        Mitf     FALSE 12594.6737
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE 27383.4125
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE   168.5509
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE   128.8275
##              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              128849128            128846244   novel distal
##              utr3start   utr3end Baf3_short.form.usage
##              <numeric> <numeric>             <numeric>
##   uc009daz.2  98018276  98021358             33.531338
##   uc009dhz.2 114859443 114860614            520.914790
##   uc009dmb.2 119316955 119315133              8.857865
##   uc009eet.1 128849981 128848044             22.671784
##              UM15_short.form.usage Baf3_long.form.usage
##                          <numeric>            <numeric>
##   uc009daz.2             1.0521497           282.824024
##   uc009dhz.2           172.1523492           208.024604
##   uc009dmb.2            49.5356160             8.535131
##   uc009eet.1             0.7334097             7.876603
##              UM15_long.form.usage Baf3_PDUI UM15_PDUI short.mean.gp1
##                         <numeric> <numeric> <numeric>      <numeric>
##   uc009daz.2            189.14112 0.8940074 0.9944680      33.531338
##   uc009dhz.2            456.54462 0.2853798 0.7261760     520.914790
##   uc009dmb.2             70.81263 0.4907223 0.5883977       8.857865
##   uc009eet.1             21.96014 0.2578402 0.9676820      22.671784
##              long.mean.gp1 short.mean.gp2 long.mean.gp2  PDUI.gp1
##                  <numeric>      <numeric>     <numeric> <numeric>
##   uc009daz.2    282.824024      1.0521497     189.14112 0.8940074
##   uc009dhz.2    208.024604    172.1523492     456.54462 0.2853798
##   uc009dmb.2      8.535131     49.5356160      70.81263 0.4907223
##   uc009eet.1      7.876603      0.7334097      21.96014 0.2578402
##               PDUI.gp2       dPDUI      logFC      P.Value    adj.P.Val
##              <numeric>   <numeric>  <numeric>    <numeric>    <numeric>
##   uc009daz.2 0.9944680 -0.10046063 -0.1536382 1.586980e-06 2.115973e-06
##   uc009dhz.2 0.7261760 -0.44079612 -1.3474358 1.709476e-60 6.837904e-60
##   uc009dmb.2 0.5883977 -0.09767539 -0.2618847 5.930309e-01 5.930309e-01
##   uc009eet.1 0.9676820 -0.70984179 -1.9080557 2.632879e-08 5.265757e-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 process described above can be done in one step.

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 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 14 with logLik: -724.6964 
## converged at iteration 1 with logLik: -997.3022 
## converged at iteration 6 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.70792
##   uc009dhz.2  uc009dhz.2       74244        Atg7     FALSE   7630.49809
##   uc009die.2  uc009die.2      232334       Vgll4     FALSE  10704.69891
##   uc009dmb.2  uc009dmb.2      211187       Lrtm2     FALSE     18.63507
##   uc009eet.1  uc009eet.1      232406    BC035044     FALSE    227.54183
##   uc009eol.1  uc009eol.1       11569       Aebp2      TRUE 204539.69842
##              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.660257
##   uc009dhz.2 114859343 114860614    184,184,184,...              3.725449
##   uc009die.2 114862164 114862092       76,76,76,...            121.888138
##   uc009dmb.2 119317055 119315133       11,12,12,...              9.226098
##   uc009eet.1 128850081 128848044       22,22,22,...             22.044952
##   uc009eol.1 140651362 140651622 1408,1412,1384,...           1065.980134
##              Baf3_long.form.usage Baf3_PDUI       short.mean
##                         <numeric> <numeric>           <list>
##   uc009daz.2           283.389388 0.9052538 29.6602572856635
##   uc009dhz.2           209.978806 0.9825673  3.7254488495449
##   uc009die.2           174.122892 0.5882311   121.8881378455
##   uc009dmb.2             9.117988 0.4970533 9.22609762692124
##   uc009eet.1            10.327422 0.3190196 22.0449523788087
##   uc009eol.1           221.456140 0.1720133 1065.98013415893
##                     long.mean              PDUI P.Value adj.P.Val
##                        <list>            <list>  <list>    <list>
##   uc009daz.2 283.389388104407 0.905253822444981       1         1
##   uc009dhz.2 209.978806469604 0.982567268751943       1         1
##   uc009die.2 174.122891566265 0.588231093659866       1         1
##   uc009dmb.2 9.11798839458414 0.497053294663731       1         1
##   uc009eet.1 10.3274224192527 0.319019611124456       1         1
##   uc009eol.1 221.456140350877 0.172013283092553       1         1
##                  dPDUI      PASS
##              <numeric> <logical>
##   uc009daz.2      <NA>     FALSE
##   uc009dhz.2      <NA>     FALSE
##   uc009die.2      <NA>     FALSE
##   uc009dmb.2      <NA>     FALSE
##   uc009eet.1      <NA>     FALSE
##   uc009eol.1      <NA>     FALSE
##   -------
##   seqinfo: 1 sequence from an unspecified genome; no seqlengths

3 Session Info

sessionInfo()

R version 3.3.1 (2016-06-21) Platform: x86_64-pc-linux-gnu (64-bit) Running under: Ubuntu 16.04.1 LTS

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.12.0
[2] e1071_1.6-7
[3] seqinr_3.3-3
[4] BSgenome.Drerio.UCSC.danRer7_1.4.0
[5] org.Hs.eg.db_3.4.0
[6] TxDb.Mmusculus.UCSC.mm10.knownGene_3.4.0 [7] BSgenome.Mmusculus.UCSC.mm10_1.4.0
[8] BSgenome_1.42.0
[9] rtracklayer_1.34.0
[10] Biostrings_2.42.0
[11] XVector_0.14.0
[12] InPAS_1.6.0
[13] GenomicFeatures_1.26.0
[14] AnnotationDbi_1.36.0
[15] GenomicRanges_1.26.0
[16] GenomeInfoDb_1.10.0
[17] IRanges_2.8.0
[18] S4Vectors_0.12.0
[19] Biobase_2.34.0
[20] BiocGenerics_0.20.0
[21] BiocStyle_2.2.0

loaded via a namespace (and not attached): [1] httr_1.2.1 AnnotationHub_2.6.0
[3] splines_3.3.1 Formula_1.2-1
[5] shiny_0.14.1 assertthat_0.1
[7] interactiveDisplayBase_1.12.0 latticeExtra_0.6-28
[9] Rsamtools_1.26.0 yaml_2.1.13
[11] RSQLite_1.0.0 lattice_0.20-34
[13] biovizBase_1.22.0 limma_3.30.0
[15] chron_2.3-47 digest_0.6.10
[17] RColorBrewer_1.1-2 colorspace_1.2-7
[19] preprocessCore_1.36.0 htmltools_0.3.5
[21] httpuv_1.3.3 Matrix_1.2-7.1
[23] plyr_1.8.4 XML_3.98-1.4
[25] biomaRt_2.30.0 zlibbioc_1.20.0
[27] xtable_1.8-2 scales_0.4.0
[29] BiocParallel_1.8.0 tibble_1.2
[31] ggplot2_2.1.0 SummarizedExperiment_1.4.0
[33] nnet_7.3-12 Gviz_1.18.0
[35] survival_2.39-5 magrittr_1.5
[37] mime_0.5 evaluate_0.10
[39] MASS_7.3-45 truncnorm_1.0-7
[41] foreign_0.8-67 class_7.3-14
[43] BiocInstaller_1.24.0 tools_3.3.1
[45] data.table_1.9.6 formatR_1.4
[47] matrixStats_0.51.0 stringr_1.1.0
[49] depmixS4_1.3-3 munsell_0.4.3
[51] cluster_2.0.5 ensembldb_1.6.0
[53] ade4_1.7-4 grid_3.3.1
[55] RCurl_1.95-4.8 dichromat_2.0-0
[57] VariantAnnotation_1.20.0 Rsolnp_1.16
[59] bitops_1.0-6 rmarkdown_1.1
[61] gtable_0.2.0 DBI_0.5-1
[63] R6_2.2.0 GenomicAlignments_1.10.0
[65] gridExtra_2.2.1 knitr_1.14
[67] Hmisc_3.17-4 stringi_1.1.2
[69] Rcpp_0.12.7 rpart_4.1-10
[71] acepack_1.3-3.3