regionReport 1.18.2
R is an open-source statistical environment which can be easily modified to enhance its functionality via packages. regionReport is a R package available via the Bioconductor repository for packages. R can be installed on any operating system from CRAN after which you can install regionReport by using the following commands in your R session:
install.packages("BiocManager")
BiocManager::install("regionReport")
## Check that you have a valid Bioconductor installation
BiocManager::valid()regionReport is based on many other packages and in particular in those that have implemented the infrastructure needed for dealing with RNA-seq data. That is, packages like GenomicFeatures that allow you to import the data and DESeq2 for generating differential expression results. A regionReport user is not expected to deal with those packages directly.
If you are asking yourself the question “Where do I start using Bioconductor?” you might be interested in this blog post.
As package developers, we try to explain clearly how to use our packages and in which order to use the functions. But R and Bioconductor have a steep learning curve so it is critical to learn where to ask for help. The blog post quoted above mentions some but we would like to highlight the Bioconductor support site as the main resource for getting help: remember to use the regionReport tag and check the older posts. Other alternatives are available such as creating GitHub issues and tweeting. However, please note that if you want to receive help you should adhere to the posting guidelines. It is particularly critical that you provide a small reproducible example and your session information so package developers can track down the source of the error.
We hope that regionReport will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
## Citation info
citation('regionReport')## 
## Collado-Torres L, Jaffe AE, Leek JT (2016). "regionReport:
## Interactive reports for region-level and feature-level genomic
## analyses [version2; referees: 2 approved, 1 approved with
## reservations]." _F1000Research_, *4*, 105. doi:
## 10.12688/f1000research.6379.2 (URL:
## https://doi.org/10.12688/f1000research.6379.2), <URL:
## http://f1000research.com/articles/4-105/v2>.
## 
## Collado-Torres L, Nellore A, Frazee AC, Wilks C, Love MI, Langmead
## B, Irizarry RA, Leek JT, Jaffe AE (2017). "Flexible expressed region
## analysis for RNA-seq with derfinder." _Nucl. Acids Res._. doi:
## 10.1093/nar/gkw852 (URL: https://doi.org/10.1093/nar/gkw852), <URL:
## http://nar.oxfordjournals.org/content/early/2016/09/29/nar.gkw852>.
## 
## Collado-Torres L, Jaffe AE, Leek JT (2017). _regionReport: Generate
## HTML or PDF reports for a set of genomic regions or DESeq2/edgeR
## results_. doi: 10.18129/B9.bioc.regionReport (URL:
## https://doi.org/10.18129/B9.bioc.regionReport),
## https://github.com/leekgroup/regionReport - R package version
## 1.18.2, <URL: http://www.bioconductor.org/packages/regionReport>.
## 
## To see these entries in BibTeX format, use 'print(<citation>,
## bibtex=TRUE)', 'toBibtex(.)', or set
## 'options(citation.bibtex.max=999)'.regionReport (Collado-Torres, Jaffe, and Leek, 2016) creates HTML or PDF reports for a set of genomic regions such as derfinder (Collado-Torres, Nellore, Frazee, Wilks, et al., 2017) results or for feature-level analyses performed with DESeq2 (Love, Huber, and Anders, 2014) or edgeR (Robinson, McCarthy, and Smyth, 2010; McCarthy, Chen, and Smyth, 2012; Chen, Lun, and Smyth, 2014). The HTML reports are styled with rmarkdown by default but can optionally be styled with knitrBootstrap (Hester, 2018).
This package includes a basic exploration for a general set of genomic regions which can be easily customized to include the appropriate conclusions and/or further exploration of the results. Such a report can be generated using renderReport(). regionReport has a separate template for running a basic exploration analysis of derfinder
results by using derfinderReport(). derfinderReport() is specific to single base-level approach derfinder results. A third template is included for exploring DESeq2 or edgeR differential expression results.
All reports are written in R Markdown
format and include all the code for making the plots and explorations in the report itself. For all templates, regionReport relies on
knitr (Xie, 2014), rmarkdown
, DT (Xie, Cheng, and Tan, 2019) and optionally knitrBootstrap
(Hester, 2018) for generating the report. The reports can be either in HTML or PDF format and can be easily customized.
The plots in regionReport for exploring DESeq2 are powered by ggplot2 (Wickham, 2016) and pheatmap (Kolde, 2019).
The regionReport supplementary website regionReportSupp has examples of using regionReport with DESeq2 results. In particular, please look at DESeq2.html which has the code for generating some DESeq2 results based on the DESeq2 vignette. Then it uses those results to create HTML and PDF versions of the same report. The resulting reports are available in the following locations:
Note that in both examples we changed the ggplot2 theme to theme_bw(). Also, in the PDF version we used the option device = 'pdf' instead of the default device = 'png' in DESeq2Report() since PDF figures are more appropriate for PDF reports: they look better than PNG figures.
If you want to create a similar HTML report as the one linked in this section, simply run example('DESeq2Report', 'regionReport', ask=FALSE). The only difference will be the ggplot2 theme for the plots.
regionReport has the edgeReport() function that takes as input a DGEList object and the results from the differential expression analysis using edgeR. edgeReport() internally uses DEFormats to convert the results to DESeq2’s format and then uses DESeqReport() to create the final report. The report looks nearly the same whether you performed the differential expression analysis with DESeq2 or edgeR in order to make more homogenous the exploratory data analysis step.
The regionReport supplementary website regionReportSupp has examples of using regionReport with edgeR results. In particular, please look at edgeR.html which has the code for generating some random data with DEFormats and performing the differential expression analysis with edgeR. Then it uses those results to create HTML and PDF versions of the same report. The resulting reports are available in the following locations:
Note that in both examples we changed the ggplot2 theme to theme_linedraw(). Also, in the PDF version we used the option device = 'pdf' instead of the default device = 'png' in edgeReport() since PDF figures are more appropriate for PDF reports: they look better than PNG figures.
If you want to create a similar HTML report as the one linked in this section, simply run example('edgeReport', 'regionReport', ask=FALSE). The only difference will be the ggplot2 theme for the plots and the amount of data simulated with DEFormats.
The plots in regionReport for region reports are powered by derfinderPlot (Collado-Torres, Jaffe, and Leek, 2017), ggbio (Yin, Cook, and Lawrence, 2012), and ggplot2 (Wickham, 2016).
The regionReport supplementary website regionReportSupp has examples of using regionReport with results from DiffBind and derfinder. Included as a vignette, this package also has an example using a small data set derived from bumphunter. These represent different uses of regionReport for results from ChIP-seq, methylation, and RNA-seq data. In particular, the DiffBind example illustrates how to expand a basic report created with renderReport().
For a general use case, you first have to identify a set of genomic regions of interest and store it as a GRanges object. In a typical workflow you will have some variables measured for each of the regions, such as p-values and scores. renderReport() uses the set of regions and three main arguments:
pvalueVars: this is a character vector (named optionally) with the names of the variables that are bound between 0 and 1, such as p-values. For each of these variables, renderReport() explores the distribution by chromosome, the overall distribution, and makes a table with commonly used cutoffs.densityVars: is another character vector (named optionally) with another set of variables you wish to explore by making density graphs. This is commonly used for scores and other similar numerical variables.significantVar: is a logical vector separating the regions into by whether they are statistically significant. For example, this information is used to explore the width of all the regions and compare it the significant ones.Other parameters control the name of the report, where it’ll be located, the transcripts database used to annotate the nearest genes, graphical parameters, etc.
Here is a short example of how to use renderReport(). Note that we are using regions produced by derfinder just for convenience sake. You can also run this example by using example('renderReport', 'regionReport', ask=FALSE).
## Load derfinder
library('derfinder')
regions <- genomeRegions$regions
## Assign chr length
library('GenomicRanges')
seqlengths(regions) <- c('chr21' = 48129895)
## The output will be saved in the 'derfinderReport-example' directory
dir.create('renderReport-example', showWarnings = FALSE, recursive = TRUE)
## Generate the HTML report
report <- renderReport(regions, 'Example run', pvalueVars = c(
    'Q-values' = 'qvalues', 'P-values' = 'pvalues'), densityVars = c(
    'Area' = 'area', 'Mean coverage' = 'meanCoverage'), 
    significantVar = regions$qvalues <= 0.05, nBestRegions = 20,
    outdir = 'renderReport-example')For derfinder results created via the expressed regions-level approach you can use renderReport() to explore the results. If you use DESeq2 to perform the differential expression analysis of the expressed regions you can then use DESeq2Report().
Prior to using regionReport::derfinderReport() you must use derfinder to analyze a specific data set. While there are many ways to do so, we recommend using analyzeChr() with the same prefix argument. Then merging the results with mergeResults(). This is the recommended pipeline for the single base-level approach.
Below, we run derfinder for the example data included in the package. The steps are:
## Load derfinder
library('derfinder')
## The output will be saved in the 'report' directory
dir.create('report', showWarnings = FALSE, recursive = TRUE)The following code runs derfinder.
## Save the current path
initialPath <- getwd()
setwd(file.path(initialPath, 'report'))
## Generate output from derfinder
## Collapse the coverage information
collapsedFull <- collapseFullCoverage(list(genomeData$coverage), 
verbose=TRUE)
## Calculate library size adjustments
sampleDepths <- sampleDepth(collapsedFull, probs=c(0.5), nonzero=TRUE, 
verbose=TRUE)
## Build the models
group <- genomeInfo$pop
adjustvars <- data.frame(genomeInfo$gender)
models <- makeModels(sampleDepths, testvars=group, adjustvars=adjustvars)
## Analyze chromosome 21
analysis <- analyzeChr(chr='21', coverageInfo=genomeData, models=models, 
cutoffFstat=1, cutoffType='manual', seeds=20140330, groupInfo=group, 
mc.cores=1, writeOutput=TRUE, returnOutput=TRUE)
## Save the stats options for later
optionsStats <- analysis$optionsStats
## Change the directory back to the original one
setwd(initialPath)For convenience, we have included the derfinder results as part of regionReport. Note that the above functions are routinely checked as part of derfinder.
## Copy previous results
file.copy(system.file(file.path('extdata', 'chr21'), package='derfinder', 
mustWork=TRUE), 'report', recursive=TRUE)## [1] TRUENext, proceed to merging the results.
## Merge the results from the different chromosomes. In this case, there's 
## only one: chr21
mergeResults(chrs = 'chr21', prefix = 'report',
    genomicState = genomicState$fullGenome)## 2019-06-18 21:06:38 mergeResults: Saving options used## 2019-06-18 21:06:38 Loading chromosome chr21## Neither 'cutoffFstatUsed' nor 'optionsStats' were supplied, so the FWER calculation step will be skipped.## 2019-06-18 21:06:38 mergeResults: Saving fullNullSummary## 2019-06-18 21:06:38 mergeResults: Re-calculating the p-values## 2019-06-18 21:06:38 mergeResults: Saving fullRegions## 2019-06-18 21:06:38 mergeResults: assigning genomic states## 2019-06-18 21:06:38 annotateRegions: counting## 2019-06-18 21:06:38 annotateRegions: annotating## 2019-06-18 21:06:38 mergeResults: Saving fullAnnotatedRegions## 2019-06-18 21:06:38 mergeResults: Saving fullFstats## 2019-06-18 21:06:38 mergeResults: Saving fullTime## Load optionsStats
load(file.path('report', 'chr21', 'optionsStats.Rdata'))Once the derfinder output has been generated and merged, use derfinderReport() to create the HTML report.
## Load derfindeReport
library('regionReport')## Generate the HTML report
report <- derfinderReport(prefix='report', browse=FALSE,
    nBestRegions=15, makeBestClusters=TRUE, outdir='html',
    fullCov=list('21'=genomeDataRaw$coverage), optionsStats=optionsStats)## Writing 10 Bibtex entries ... OK
## Results written to file 'report/html/basicExploration.bib'
## extendedMapSeqlevels: sequence names mapped from NCBI to UCSC for species homo_sapiens
## processing file: basicExploration.Rmd
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
## Parsing transcripts...
## No transcripts found at this region.
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## Parsing transcripts...
## Parsing exons...
## Parsing cds...
## Parsing utrs...
## ------exons...
## ------cdss...
## ------introns...
## ------utr...
## aggregating...
## Done
## Constructing graphics...
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## Coordinate system already present. Adding new coordinate system, which will replace the existing one.
## output file: basicExploration.knit.md## /usr/bin/pandoc +RTS -K512m -RTS basicExploration.utf8.md --to html4 --from markdown+autolink_bare_uris+ascii_identifiers+tex_math_single_backslash+smart --output basicExploration.html --email-obfuscation none --self-contained --wrap preserve --standalone --section-divs --table-of-contents --toc-depth 3 --variable toc_float=1 --variable toc_selectors=h1,h2,h3 --variable toc_collapsed=1 --variable toc_smooth_scroll=1 --variable toc_print=1 --template /tmp/Rtmp2Eb4on/BiocStyle/template.html --no-highlight --variable highlightjs=1 --number-sections --css /home/biocbuild/bbs-3.9-bioc/R/library/BiocStyle/resources/html/bioconductor.css --variable 'theme:bootstrap' --include-in-header /tmp/Rtmp2Eb4on/rmarkdown-str2e111bab5510.html --mathjax --variable 'mathjax-url:https://mathjax.rstudio.com/latest/MathJax.js?config=TeX-AMS-MML_HTMLorMML' --metadata pagetitle=basicExploration.utf8.md --variable code_folding=hide --variable code_menu=1## 
## Output created: basicExploration.htmlOnce the output is generated, you can browse the report from R using
browseURL() as shown below.
## Browse the report
browseURL(report)Note that the reports require an active Internet connection to render correctly.
The report is self-explanatory and will change some of the text depending on the input options.
If the report is taking too long to compile (say more than 3 hours), you might
want to consider setting nBestCluters to a small number or even set
makeBestClusters to FALSE.
This package was made possible thanks to:
Code for creating the vignette
## Create the vignette
library('rmarkdown')
system.time(render('regionReport.Rmd', 'BiocStyle::html_document'))
## Extract the R code
library('knitr')
knit('regionReport.Rmd', tangle = TRUE)## Clean up
file.remove('regionReportRef.bib')## [1] TRUEunlink('report', recursive = TRUE)Date the vignette was generated.
## [1] "2019-06-18 21:07:52 EDT"Wallclock time spent generating the vignette.
## Time difference of 1.283 minsR session information.
## ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
##  setting  value                       
##  version  R version 3.6.0 (2019-04-26)
##  os       Ubuntu 18.04.2 LTS          
##  system   x86_64, linux-gnu           
##  ui       X11                         
##  language (EN)                        
##  collate  C                           
##  ctype    en_US.UTF-8                 
##  tz       America/New_York            
##  date     2019-06-18                  
## 
## ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
##  package                           * version   date       lib source        
##  acepack                             1.4.1     2016-10-29 [2] CRAN (R 3.6.0)
##  annotate                            1.62.0    2019-06-18 [2] Bioconductor  
##  AnnotationDbi                     * 1.46.0    2019-06-18 [2] Bioconductor  
##  AnnotationFilter                    1.8.0     2019-06-18 [2] Bioconductor  
##  assertthat                          0.2.1     2019-03-21 [2] CRAN (R 3.6.0)
##  backports                           1.1.4     2019-04-10 [2] CRAN (R 3.6.0)
##  base64enc                           0.1-3     2015-07-28 [2] CRAN (R 3.6.0)
##  bibtex                              0.4.2     2017-06-30 [2] CRAN (R 3.6.0)
##  Biobase                           * 2.44.0    2019-06-18 [2] Bioconductor  
##  BiocGenerics                      * 0.30.0    2019-06-18 [2] Bioconductor  
##  BiocManager                         1.30.4    2018-11-13 [2] CRAN (R 3.6.0)
##  BiocParallel                        1.18.0    2019-06-18 [2] Bioconductor  
##  BiocStyle                         * 2.12.0    2019-06-18 [2] Bioconductor  
##  biomaRt                             2.40.0    2019-06-18 [2] Bioconductor  
##  Biostrings                          2.52.0    2019-06-18 [2] Bioconductor  
##  biovizBase                        * 1.32.0    2019-06-18 [2] Bioconductor  
##  bit                                 1.1-14    2018-05-29 [2] CRAN (R 3.6.0)
##  bit64                               0.9-7     2017-05-08 [2] CRAN (R 3.6.0)
##  bitops                              1.0-6     2013-08-17 [2] CRAN (R 3.6.0)
##  blob                                1.1.1     2018-03-25 [2] CRAN (R 3.6.0)
##  bookdown                            0.11      2019-05-28 [2] CRAN (R 3.6.0)
##  BSgenome                            1.52.0    2019-06-18 [2] Bioconductor  
##  bumphunter                        * 1.26.0    2019-06-18 [2] Bioconductor  
##  checkmate                           1.9.3     2019-05-03 [2] CRAN (R 3.6.0)
##  cli                                 1.1.0     2019-03-19 [2] CRAN (R 3.6.0)
##  cluster                             2.0.9     2019-05-01 [2] CRAN (R 3.6.0)
##  codetools                           0.2-16    2018-12-24 [2] CRAN (R 3.6.0)
##  colorspace                          1.4-1     2019-03-18 [2] CRAN (R 3.6.0)
##  crayon                              1.3.4     2017-09-16 [2] CRAN (R 3.6.0)
##  crosstalk                           1.0.0     2016-12-21 [2] CRAN (R 3.6.0)
##  curl                                3.3       2019-01-10 [2] CRAN (R 3.6.0)
##  data.table                          1.12.2    2019-04-07 [2] CRAN (R 3.6.0)
##  DBI                                 1.0.0     2018-05-02 [2] CRAN (R 3.6.0)
##  DEFormats                           1.12.0    2019-06-18 [2] Bioconductor  
##  DelayedArray                        0.10.0    2019-06-18 [2] Bioconductor  
##  derfinder                         * 1.18.3    2019-06-18 [2] Bioconductor  
##  derfinderHelper                     1.18.1    2019-06-18 [2] Bioconductor  
##  derfinderPlot                     * 1.18.1    2019-06-18 [2] Bioconductor  
##  DESeq2                              1.24.0    2019-06-18 [2] Bioconductor  
##  dichromat                           2.0-0     2013-01-24 [2] CRAN (R 3.6.0)
##  digest                              0.6.19    2019-05-20 [2] CRAN (R 3.6.0)
##  doRNG                               1.7.1     2018-06-22 [2] CRAN (R 3.6.0)
##  dplyr                               0.8.1     2019-05-14 [2] CRAN (R 3.6.0)
##  DT                                * 0.7       2019-06-11 [2] CRAN (R 3.6.0)
##  edgeR                               3.26.4    2019-06-18 [2] Bioconductor  
##  ensembldb                           2.8.0     2019-06-18 [2] Bioconductor  
##  evaluate                            0.14      2019-05-28 [2] CRAN (R 3.6.0)
##  foreach                           * 1.4.4     2017-12-12 [2] CRAN (R 3.6.0)
##  foreign                             0.8-71    2018-07-20 [2] CRAN (R 3.6.0)
##  Formula                             1.2-3     2018-05-03 [2] CRAN (R 3.6.0)
##  genefilter                          1.66.0    2019-06-18 [2] Bioconductor  
##  geneplotter                         1.62.0    2019-06-18 [2] Bioconductor  
##  GenomeInfoDb                      * 1.20.0    2019-06-18 [2] Bioconductor  
##  GenomeInfoDbData                    1.2.1     2019-04-26 [2] Bioconductor  
##  GenomicAlignments                   1.20.1    2019-06-18 [2] Bioconductor  
##  GenomicFeatures                   * 1.36.1    2019-06-18 [2] Bioconductor  
##  GenomicFiles                        1.20.0    2019-06-18 [2] Bioconductor  
##  GenomicRanges                     * 1.36.0    2019-06-18 [2] Bioconductor  
##  GGally                              1.4.0     2018-05-17 [2] CRAN (R 3.6.0)
##  ggbio                             * 1.32.0    2019-06-18 [2] Bioconductor  
##  ggplot2                           * 3.2.0     2019-06-16 [2] CRAN (R 3.6.0)
##  glue                                1.3.1     2019-03-12 [2] CRAN (R 3.6.0)
##  graph                               1.62.0    2019-06-18 [2] Bioconductor  
##  gridExtra                         * 2.3       2017-09-09 [2] CRAN (R 3.6.0)
##  gtable                              0.3.0     2019-03-25 [2] CRAN (R 3.6.0)
##  highr                               0.8       2019-03-20 [2] CRAN (R 3.6.0)
##  Hmisc                               4.2-0     2019-01-26 [2] CRAN (R 3.6.0)
##  hms                                 0.4.2     2018-03-10 [2] CRAN (R 3.6.0)
##  htmlTable                           1.13.1    2019-01-07 [2] CRAN (R 3.6.0)
##  htmltools                           0.3.6     2017-04-28 [2] CRAN (R 3.6.0)
##  htmlwidgets                         1.3       2018-09-30 [2] CRAN (R 3.6.0)
##  httpuv                              1.5.1     2019-04-05 [2] CRAN (R 3.6.0)
##  httr                                1.4.0     2018-12-11 [2] CRAN (R 3.6.0)
##  IRanges                           * 2.18.1    2019-06-18 [2] Bioconductor  
##  iterators                         * 1.0.10    2018-07-13 [2] CRAN (R 3.6.0)
##  jsonlite                            1.6       2018-12-07 [2] CRAN (R 3.6.0)
##  knitcitations                     * 1.0.8     2017-07-04 [2] CRAN (R 3.6.0)
##  knitr                             * 1.23      2019-05-18 [2] CRAN (R 3.6.0)
##  knitrBootstrap                      1.0.2     2018-05-24 [2] CRAN (R 3.6.0)
##  labeling                            0.3       2014-08-23 [2] CRAN (R 3.6.0)
##  later                               0.8.0     2019-02-11 [2] CRAN (R 3.6.0)
##  lattice                             0.20-38   2018-11-04 [2] CRAN (R 3.6.0)
##  latticeExtra                        0.6-28    2016-02-09 [2] CRAN (R 3.6.0)
##  lazyeval                            0.2.2     2019-03-15 [2] CRAN (R 3.6.0)
##  limma                               3.40.2    2019-06-18 [2] Bioconductor  
##  locfit                            * 1.5-9.1   2013-04-20 [2] CRAN (R 3.6.0)
##  lubridate                           1.7.4     2018-04-11 [2] CRAN (R 3.6.0)
##  magrittr                            1.5       2014-11-22 [2] CRAN (R 3.6.0)
##  markdown                            1.0       2019-06-07 [2] CRAN (R 3.6.0)
##  Matrix                              1.2-17    2019-03-22 [2] CRAN (R 3.6.0)
##  matrixStats                         0.54.0    2018-07-23 [2] CRAN (R 3.6.0)
##  memoise                             1.1.0     2017-04-21 [2] CRAN (R 3.6.0)
##  mgcv                              * 1.8-28    2019-03-21 [2] CRAN (R 3.6.0)
##  mime                                0.7       2019-06-11 [2] CRAN (R 3.6.0)
##  munsell                             0.5.0     2018-06-12 [2] CRAN (R 3.6.0)
##  nlme                              * 3.1-140   2019-05-12 [2] CRAN (R 3.6.0)
##  nnet                                7.3-12    2016-02-02 [2] CRAN (R 3.6.0)
##  org.Hs.eg.db                      * 3.8.2     2019-06-18 [2] Bioconductor  
##  OrganismDbi                         1.26.0    2019-06-18 [2] Bioconductor  
##  pillar                              1.4.1     2019-05-28 [2] CRAN (R 3.6.0)
##  pkgconfig                           2.0.2     2018-08-16 [2] CRAN (R 3.6.0)
##  pkgmaker                            0.27      2018-05-25 [2] CRAN (R 3.6.0)
##  plyr                                1.8.4     2016-06-08 [2] CRAN (R 3.6.0)
##  prettyunits                         1.0.2     2015-07-13 [2] CRAN (R 3.6.0)
##  progress                            1.2.2     2019-05-16 [2] CRAN (R 3.6.0)
##  promises                            1.0.1     2018-04-13 [2] CRAN (R 3.6.0)
##  ProtGenerics                        1.16.0    2019-06-18 [2] Bioconductor  
##  purrr                               0.3.2     2019-03-15 [2] CRAN (R 3.6.0)
##  qvalue                              2.16.0    2019-06-18 [2] Bioconductor  
##  R6                                  2.4.0     2019-02-14 [2] CRAN (R 3.6.0)
##  RBGL                                1.60.0    2019-06-18 [2] Bioconductor  
##  RColorBrewer                      * 1.1-2     2014-12-07 [2] CRAN (R 3.6.0)
##  Rcpp                                1.0.1     2019-03-17 [2] CRAN (R 3.6.0)
##  RCurl                               1.95-4.12 2019-03-04 [2] CRAN (R 3.6.0)
##  RefManageR                          1.2.12    2019-04-03 [2] CRAN (R 3.6.0)
##  regionReport                      * 1.18.2    2019-06-18 [1] Bioconductor  
##  registry                            0.5-1     2019-03-05 [2] CRAN (R 3.6.0)
##  reshape                             0.8.8     2018-10-23 [2] CRAN (R 3.6.0)
##  reshape2                            1.4.3     2017-12-11 [2] CRAN (R 3.6.0)
##  rlang                               0.3.4     2019-04-07 [2] CRAN (R 3.6.0)
##  rmarkdown                           1.13      2019-05-22 [2] CRAN (R 3.6.0)
##  rngtools                            1.3.1.1   2019-04-26 [2] CRAN (R 3.6.0)
##  rpart                               4.1-15    2019-04-12 [2] CRAN (R 3.6.0)
##  Rsamtools                           2.0.0     2019-06-18 [2] Bioconductor  
##  RSQLite                             2.1.1     2018-05-06 [2] CRAN (R 3.6.0)
##  rstudioapi                          0.10      2019-03-19 [2] CRAN (R 3.6.0)
##  rtracklayer                         1.44.0    2019-06-18 [2] Bioconductor  
##  S4Vectors                         * 0.22.0    2019-06-18 [2] Bioconductor  
##  scales                              1.0.0     2018-08-09 [2] CRAN (R 3.6.0)
##  sessioninfo                       * 1.1.1     2018-11-05 [2] CRAN (R 3.6.0)
##  shiny                               1.3.2     2019-04-22 [2] CRAN (R 3.6.0)
##  stringi                             1.4.3     2019-03-12 [2] CRAN (R 3.6.0)
##  stringr                             1.4.0     2019-02-10 [2] CRAN (R 3.6.0)
##  SummarizedExperiment                1.14.0    2019-06-18 [2] Bioconductor  
##  survival                            2.44-1.1  2019-04-01 [2] CRAN (R 3.6.0)
##  tibble                              2.1.3     2019-06-06 [2] CRAN (R 3.6.0)
##  tidyselect                          0.2.5     2018-10-11 [2] CRAN (R 3.6.0)
##  TxDb.Hsapiens.UCSC.hg19.knownGene * 3.2.2     2019-04-26 [2] Bioconductor  
##  VariantAnnotation                   1.30.1    2019-06-18 [2] Bioconductor  
##  whisker                           * 0.3-2     2013-04-28 [2] CRAN (R 3.6.0)
##  withr                               2.1.2     2018-03-15 [2] CRAN (R 3.6.0)
##  xfun                                0.7       2019-05-14 [2] CRAN (R 3.6.0)
##  XML                                 3.98-1.20 2019-06-06 [2] CRAN (R 3.6.0)
##  xml2                                1.2.0     2018-01-24 [2] CRAN (R 3.6.0)
##  xtable                              1.8-4     2019-04-21 [2] CRAN (R 3.6.0)
##  XVector                             0.24.0    2019-06-18 [2] Bioconductor  
##  yaml                                2.2.0     2018-07-25 [2] CRAN (R 3.6.0)
##  zlibbioc                            1.30.0    2019-06-18 [2] Bioconductor  
## 
## [1] /tmp/RtmpkOavtq/Rinst17bd6aad58ea
## [2] /home/biocbuild/bbs-3.9-bioc/R/libraryThis vignette was generated using BiocStyle (Oleś, Morgan, and Huber, 2019) with knitr (Xie, 2014) and rmarkdown running behind the scenes.
Citations made with knitcitations (Boettiger, 2017).
[1] S. Arora, M. Morgan, M. Carlson, and H. Pagès. GenomeInfoDb: Utilities for manipulating chromosome and other ‘seqname’ identifiers. 2017. DOI: 10.18129/B9.bioc.GenomeInfoDb.
[1] B. Auguie. gridExtra: Miscellaneous Functions for “Grid” Graphics. R package version 2.3. 2017. URL: https://CRAN.R-project.org/package=gridExtra.
[1] C. Boettiger. knitcitations: Citations for ‘Knitr’ Markdown Files. R package version 1.0.8. 2017. URL: https://CRAN.R-project.org/package=knitcitations.
[1] M. Carlson and B. P. Maintainer. TxDb.Hsapiens.UCSC.hg19.knownGene: Annotation package for TxDb object(s). R package version 3.2.2. 2015.
[1] Y. Chen, A. T. L. Lun, and G. K. Smyth. “Differential expression analysis of complex RNA-seq experiments using edgeR”. In: Statistical Analysis of Next Generation Sequencing Data. Ed. by S. Datta and D. Nettleton. New York: Springer, 2014, pp. 51-74.
[1] L. Collado-Torres, A. E. Jaffe, and J. T. Leek. derfinderPlot: Plotting functions for derfinder. https://github.com/leekgroup/derfinderPlot - R package version 1.18.1. 2017. DOI: 10.18129/B9.bioc.derfinderPlot. URL: http://www.bioconductor.org/packages/derfinderPlot.
[1] G. Csárdi, R. core, H. Wickham, W. Chang, et al. sessioninfo: R Session Information. R package version 1.1.1. 2018. URL: https://CRAN.R-project.org/package=sessioninfo.
[1] J. Hester. knitrBootstrap: ‘knitr’ Bootstrap Framework. R package version 1.0.2. 2018. URL: https://CRAN.R-project.org/package=knitrBootstrap.
[1] R. Kolde. pheatmap: Pretty Heatmaps. R package version 1.0.12. 2019. URL: https://CRAN.R-project.org/package=pheatmap.
[1] E. Neuwirth. RColorBrewer: ColorBrewer Palettes. R package version 1.1-2. 2014. URL: https://CRAN.R-project.org/package=RColorBrewer.
[1] A. Oleś. DEFormats: Differential gene expression data formats converter. R package version 1.12.0. 2019. URL: https://github.com/aoles/DEFormats.
[1] A. Oleś, M. Morgan, and W. Huber. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.12.0. 2019. URL: https://github.com/Bioconductor/BiocStyle.
[1] R Core Team. R: A Language and Environment for Statistical Computing. R Foundation for Statistical Computing. Vienna, Austria, 2019. URL: https://www.R-project.org/.
[1] R Core Team. R: A Language and Environment for Statistical Computing. R Foundation for Statistical Computing. Vienna, Austria, 2019. URL: https://www.R-project.org/.
[1] H. Wickham. ggplot2: Elegant Graphics for Data Analysis. Springer-Verlag New York, 2016. ISBN: 978-3-319-24277-4. URL: https://ggplot2.tidyverse.org.
[1] Y. Xie. “knitr: A Comprehensive Tool for Reproducible Research in R”. In: Implementing Reproducible Computational Research. Ed. by V. Stodden, F. Leisch and R. D. Peng. ISBN 978-1466561595. Chapman and Hall/CRC, 2014. URL: http://www.crcpress.com/product/isbn/9781466561595.
[1] Y. Xie, J. Cheng, and X. Tan. DT: A Wrapper of the JavaScript Library ‘DataTables’. R package version 0.7. 2019. URL: https://CRAN.R-project.org/package=DT.
[1] T. Yin, D. Cook, and M. Lawrence. “ggbio: an R package for extending the grammar of graphics for genomic data”. In: Genome Biology 13.8 (2012), p. R77.
[1] T. Yin, M. Lawrence, and D. Cook. biovizBase: Basic graphic utilities for visualization of genomic data. R package version 1.32.0. 2019.