derfinder 1.18.9
R is an open-source statistical environment which can be easily modified to enhance its functionality via packages. derfinder 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 derfinder by using the following commands in your R session:
install.packages("BiocManager")
BiocManager::install("derfinder")
## Check that you have a valid Bioconductor installation
BiocManager::valid()derfinder 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 Rsamtools, GenomicAlignments and rtracklayer that allow you to import the data. A derfinder user is not expected to deal with those packages directly but will need to be familiar with GenomicRanges to understand the results derfinder generates. It might also prove to be highly beneficial to check the BiocParallel package for performing parallel computations.
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 derfinder 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 would like to highlight the derfinder user Jessica Hekman. She has used derfinder with non-human data, and in the process of doing so discovered some small bugs or sections of the documentation that were not clear.
We hope that derfinder will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
## Citation info
citation('derfinder')## 
## 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>.
## 
## Frazee AC, Sabunciyan S, Hansen KD, Irizarry RA, Leek JT (2014).
## "Differential expression analysis of RNA-seq data at single-base
## resolution." _Biostatistics_, *15 (3)*, 413-426. doi:
## 10.1093/biostatistics/kxt053 (URL:
## https://doi.org/10.1093/biostatistics/kxt053), <URL:
## http://biostatistics.oxfordjournals.org/content/15/3/413.long>.
## 
## Collado-Torres L, Jaffe AE, Leek JT (2017). _derfinder:
## Annotation-agnostic differential expression analysis of RNA-seq data
## at base-pair resolution via the DER Finder approach_. doi:
## 10.18129/B9.bioc.derfinder (URL:
## https://doi.org/10.18129/B9.bioc.derfinder),
## https://github.com/lcolladotor/derfinder - R package version 1.18.9,
## <URL: http://www.bioconductor.org/packages/derfinder>.
## 
## To see these entries in BibTeX format, use 'print(<citation>,
## bibtex=TRUE)', 'toBibtex(.)', or set
## 'options(citation.bibtex.max=999)'.Here is a very quick example of a DER Finder analysis. This analysis is explained in more detail later on in this document.
## Load libraries
library('derfinder')
library('derfinderData')
library('GenomicRanges')
## Determine the files to use and fix the names
files <- rawFiles(system.file('extdata', 'AMY', package = 'derfinderData'),
    samplepatt = 'bw', fileterm = NULL)
names(files) <- gsub('.bw', '', names(files))
## Load the data from disk -- On Windows you have to load data from Bam files
fullCov <- fullCoverage(files = files, chrs = 'chr21', verbose = FALSE)
## Get the region matrix of Expressed Regions (ERs)
regionMat <- regionMatrix(fullCov, cutoff = 30, L = 76, verbose = FALSE)
## Get pheno table
pheno <- subset(brainspanPheno, structure_acronym == 'AMY')
## Identify which ERs are differentially expressed, that is, find the DERs
library('DESeq2')
## Round matrix
counts <- round(regionMat$chr21$coverageMatrix)
## Round matrix and specify design
dse <- DESeqDataSetFromMatrix(counts, pheno, ~ group + gender)
## Perform DE analysis
dse <- DESeq(dse, test = 'LRT', reduced = ~ gender, fitType = 'local')
## Extract results
mcols(regionMat$chr21$regions) <- c(mcols(regionMat$chr21$regions), results(dse))
## Save info in an object with a shorter name
ers <- regionMat$chr21$regions
ersderfinder is an R package that implements the DER Finder approach (Frazee, Sabunciyan, Hansen, Irizarry, et al., 2014) for RNA-seq data. Briefly, this approach is an alternative to feature-counting and transcript assembly. The basic idea is to identify contiguous base-pairs in the genome that present differential expression signal. These base-pairs are grouped into _d_ifferentially _e_xpressed _r_regions (DERs). This approach is annotation-agnostic which is a feature you might be interested in. In particular, derfinder contains functions that allow you to identify DERs via two alternative methods. You can find more details on the full derfinder users guide.
This is a brief overview of what a DER Finder analysis looks like. In particular, here we will be identifying expressed regions (ERs) without relying on annotation. Next, we’ll identify candidate differentially expressed regions (DERs). Finally, we’ll compare the DERs with known annotation features.
We first load the required packages.
## Load libraries
library('derfinder')
library('derfinderData')
library('GenomicRanges')Next, we need to locate the chromosome 21 coverage files for a set of 12 samples. These samples are a small subset from the BrainSpan Atlas of the Human Brain (BrainSpan, 2011) publicly available data. The function rawFiles() helps us in locating these files.
## Determine the files to use and fix the names
files <- rawFiles(system.file('extdata', 'AMY', package = 'derfinderData'),
    samplepatt = 'bw', fileterm = NULL)
names(files) <- gsub('.bw', '', names(files))Next, we can load the full coverage data into memory using the fullCoverage() function. Note that the BrainSpan data is already normalized by the total number of mapped reads in each sample. However, that won’t be the case with most data sets in which case you might want to use the totalMapped and targetSize arguments. The function getTotalMapped() will be helpful to get this information.
## Load the data from disk
fullCov <- fullCoverage(files = files, chrs = 'chr21', verbose = FALSE,
    totalMapped = rep(1, length(files)), targetSize = 1)Now that we have the data, we can identify expressed regions (ERs) by using a cutoff of 30 on the base-level mean coverage from these 12 samples. Once the regions have been identified, we can calculate a coverage matrix with one row per ER and one column per sample (12 in this case). For doing this calculation we need to know the length of the sequence reads, which in this study were 76 bp long.
## Get the region matrix of Expressed Regions (ERs)
regionMat <- regionMatrix(fullCov, cutoff = 30, L = 76, verbose = FALSE)regionMatrix() returns a list of elements, each with three pieces of output. The actual ERs are arranged in a GRanges object named regions.
## regions output
regionMat$chr21$regions## GRanges object with 45 ranges and 6 metadata columns:
##      seqnames            ranges strand |            value             area
##         <Rle>         <IRanges>  <Rle> |        <numeric>        <numeric>
##    1    chr21   9827018-9827582      * | 313.671741938907 177224.534195483
##    2    chr21 15457301-15457438      * | 215.084607424943 29681.6758246422
##    3    chr21 20230140-20230192      * | 38.8325471608144 2058.12499952316
##    4    chr21 20230445-20230505      * | 41.3245082031834 2520.79500039419
##    5    chr21 27253318-27253543      * | 34.9131305372469  7890.3675014178
##   ..      ...               ...    ... .              ...              ...
##   41    chr21 33039644-33039688      * | 34.4705371008979 1551.17416954041
##   42    chr21 33040784-33040798      * |  32.134222178989 482.013332684835
##   43    chr21          33040890      * | 30.0925002098083 30.0925002098083
##   44    chr21 33040900-33040901      * | 30.1208333174388 60.2416666348775
##   45    chr21 48019401-48019414      * | 31.1489284609755 436.084998453657
##      indexStart  indexEnd cluster clusterL
##       <integer> <integer>   <Rle>    <Rle>
##    1          1       565       1      565
##    2        566       703       2      138
##    3        704       756       3      366
##    4        757       817       3      366
##    5        818      1043       4      765
##   ..        ...       ...     ...      ...
##   41       2180      2224      17       45
##   42       2225      2239      18      118
##   43       2240      2240      18      118
##   44       2241      2242      18      118
##   45       2243      2256      19       14
##   -------
##   seqinfo: 1 sequence from an unspecified genomebpCoverage is the base-level coverage list which can then be used for plotting.
## Base-level coverage matrices for each of the regions
## Useful for plotting
lapply(regionMat$chr21$bpCoverage[1:2], head, n = 2)## $`1`
##   HSB113 HSB123 HSB126 HSB130 HSB135 HSB136 HSB145 HSB153 HSB159 HSB178
## 1  93.20   3.32  28.22   5.62 185.17  98.34   5.88  16.71   3.52  15.71
## 2 124.76   7.25  63.68  11.32 374.85 199.28  10.39  30.53   5.83  29.35
##   HSB92 HSB97
## 1 47.40 36.54
## 2 65.04 51.42
## 
## $`2`
##     HSB113 HSB123 HSB126 HSB130 HSB135 HSB136 HSB145 HSB153 HSB159 HSB178
## 566  45.59   7.94  15.92  34.75 141.61 104.21  19.87  38.61   4.97   23.2
## 567  45.59   7.94  15.92  35.15 141.64 104.30  19.87  38.65   4.97   23.2
##     HSB92 HSB97
## 566 13.95 22.21
## 567 13.95 22.21## Check dimensions. First region is 565 long, second one is 138 bp long.
## The columns match the number of samples (12 in this case).
lapply(regionMat$chr21$bpCoverage[1:2], dim)## $`1`
## [1] 565  12
## 
## $`2`
## [1] 138  12The end result of the coverage matrix is shown below. Note that the coverage has been adjusted for read length. Because reads might not fully align inside a given region, the numbers are generally not integers but can be rounded if needed.
## Dimensions of the coverage matrix
dim(regionMat$chr21$coverageMatrix)## [1] 45 12## Coverage for each region. This matrix can then be used with limma or other pkgs
head(regionMat$chr21$coverageMatrix)##         HSB113     HSB123      HSB126      HSB130      HSB135      HSB136
## 1 3653.1093346 277.072106 1397.068687 1106.722895 8987.460401 5570.221054
## 2  333.3740816  99.987237  463.909476  267.354342 1198.713552 1162.313418
## 3   35.3828948  20.153553   30.725394   23.483947   16.786842   17.168947
## 4   42.3398681  29.931579   41.094474   24.724736   32.634080   19.309606
## 5   77.7402631 168.939342  115.059342  171.861974  180.638684   93.503158
## 6    0.7988158   1.770263    1.473421    2.231053    1.697368    1.007895
##        HSB145       HSB153     HSB159      HSB178       HSB92        HSB97
## 1 1330.158818 1461.2986829 297.939342 1407.288552 1168.519079 1325.9622371
## 2  257.114210  313.8513139  67.940131  193.695657  127.543553  200.7834228
## 3   22.895921   52.8756585  28.145395   33.127368   23.758816   20.4623685
## 4   33.802632   51.6146040  31.244343   33.576974   29.546183   28.2011836
## 5   90.950526   36.3046051  78.069605   97.151316  100.085790   35.5428946
## 6    1.171316    0.4221053   1.000132    1.139079    1.136447    0.3956579We can then use the coverage matrix and packages such as limma, DESeq2 or edgeR to identify which ERs are differentially expressed. Here we’ll use DESeq2 for which we need some phenotype data.
## Get pheno table
pheno <- subset(brainspanPheno, structure_acronym == 'AMY')Now we can identify the DERs using a rounded version of the coverage matrix.
## Identify which ERs are differentially expressed, that is, find the DERs
library('DESeq2')## Loading required package: SummarizedExperiment## Loading required package: Biobase## Welcome to Bioconductor
## 
##     Vignettes contain introductory material; view with
##     'browseVignettes()'. To cite Bioconductor, see
##     'citation("Biobase")', and for packages 'citation("pkgname")'.## Loading required package: DelayedArray## Loading required package: matrixStats## 
## Attaching package: 'matrixStats'## The following objects are masked from 'package:Biobase':
## 
##     anyMissing, rowMedians## Loading required package: BiocParallel## 
## Attaching package: 'DelayedArray'## The following objects are masked from 'package:matrixStats':
## 
##     colMaxs, colMins, colRanges, rowMaxs, rowMins, rowRanges## The following objects are masked from 'package:base':
## 
##     aperm, apply, rowsum## Round matrix
counts <- round(regionMat$chr21$coverageMatrix)
## Round matrix and specify design
dse <- DESeqDataSetFromMatrix(counts, pheno, ~ group + gender)## converting counts to integer mode## Perform DE analysis
dse <- DESeq(dse, test = 'LRT', reduced = ~ gender, fitType = 'local')## estimating size factors## estimating dispersions## gene-wise dispersion estimates## mean-dispersion relationship## final dispersion estimates## fitting model and testing## Extract results
mcols(regionMat$chr21$regions) <- c(mcols(regionMat$chr21$regions),
    results(dse))
## Save info in an object with a shorter name
ers <- regionMat$chr21$regions
ers## GRanges object with 45 ranges and 12 metadata columns:
##      seqnames            ranges strand |            value             area
##         <Rle>         <IRanges>  <Rle> |        <numeric>        <numeric>
##    1    chr21   9827018-9827582      * | 313.671741938907 177224.534195483
##    2    chr21 15457301-15457438      * | 215.084607424943 29681.6758246422
##    3    chr21 20230140-20230192      * | 38.8325471608144 2058.12499952316
##    4    chr21 20230445-20230505      * | 41.3245082031834 2520.79500039419
##    5    chr21 27253318-27253543      * | 34.9131305372469  7890.3675014178
##   ..      ...               ...    ... .              ...              ...
##   41    chr21 33039644-33039688      * | 34.4705371008979 1551.17416954041
##   42    chr21 33040784-33040798      * |  32.134222178989 482.013332684835
##   43    chr21          33040890      * | 30.0925002098083 30.0925002098083
##   44    chr21 33040900-33040901      * | 30.1208333174388 60.2416666348775
##   45    chr21 48019401-48019414      * | 31.1489284609755 436.084998453657
##      indexStart  indexEnd cluster clusterL          baseMean
##       <integer> <integer>   <Rle>    <Rle>         <numeric>
##    1          1       565       1      565  2846.28716155756
##    2        566       703       2      138  451.519643002767
##    3        704       756       3      366  29.5780895864663
##    4        757       817       3      366  36.0602688380417
##    5        818      1043       4      765  101.646763062388
##   ..        ...       ...     ...      ...               ...
##   41       2180      2224      17       45  20.7820346002786
##   42       2225      2239      18      118  6.41054227728292
##   43       2240      2240      18      118 0.129716561731455
##   44       2241      2242      18      118 0.702291488002413
##   45       2243      2256      19       14  5.29329263119445
##          log2FoldChange             lfcSE               stat
##               <numeric>         <numeric>          <numeric>
##    1  -1.69031816472548 0.831958703669481  0.215261662919431
##    2  -1.16404264400518 0.757489962635837   0.87112644422632
##    3 0.0461487746156371 0.458096633387309   3.13208182520646
##    4 -0.186620018986543 0.390920460570196   2.22570758632982
##    5 -0.138737745290957 0.320166382651628   3.95798720528119
##   ..                ...               ...                ...
##   41  -0.64205595270105   0.4276611764306  0.604781356173675
##   42 -0.634320864108627 0.512261897770782  0.545403917207132
##   43 -0.859548956942508  3.11653965011299 0.0206273206873746
##   44 -0.628284603274911  2.24737800886771  0.582510548110115
##   45  -1.69456340588639  1.25228952015329   5.78959103593262
##                  pvalue              padj
##               <numeric>         <numeric>
##    1  0.642674256598179 0.997155030203294
##    2  0.350643649949685 0.997155030203294
##    3 0.0767656530295047 0.863613596581927
##    4  0.135730461463207 0.997155030203294
##    5 0.0466494518014194   0.8620396310466
##   ..                ...               ...
##   41   0.43675951534617 0.997155030203294
##   42  0.460201756145957 0.997155030203294
##   43  0.885798852352319 0.997155030203294
##   44  0.445329945344273 0.997155030203294
##   45 0.0161213391338242 0.725460261022089
##   -------
##   seqinfo: 1 sequence from an unspecified genomeWe can then compare the DERs against known annotation to see which DERs overlap known exons, introns, or intergenic regions. A way to visualize this information is via a Venn diagram which we can create using vennRegions() from the derfinderPlot package as shown in Figure 1.
## Find overlaps between regions and summarized genomic annotation
annoRegs <- annotateRegions(ers, genomicState$fullGenome, verbose = FALSE)
library('derfinderPlot')## Registered S3 method overwritten by 'GGally':
##   method from   
##   +.gg   ggplot2venn <- vennRegions(annoRegs, counts.col = 'blue',
    main = 'Venn diagram using TxDb.Hsapiens.UCSC.hg19.knownGene annotation')
Figure 1: Venn diagram showing ERs by annotation class
We can also identify the nearest annotated feature. In this case, we’ll look for the nearest known gene from the UCSC hg19 annotation.
## Load database of interest
library('TxDb.Hsapiens.UCSC.hg19.knownGene')## Loading required package: GenomicFeatures## Loading required package: AnnotationDbitxdb <- keepSeqlevels(TxDb.Hsapiens.UCSC.hg19.knownGene, 'chr21')
## Find nearest feature
library('bumphunter')## Loading required package: foreach## Loading required package: iterators## Loading required package: locfit## locfit 1.5-9.1    2013-03-22genes <- annotateTranscripts(txdb)## No annotationPackage supplied. Trying org.Hs.eg.db.## Loading required package: org.Hs.eg.db## ## Getting TSS and TSE.## Getting CSS and CSE.## Getting exons.## Annotating genes.annotation <- matchGenes(ers, subject = genes)
## View annotation results
head(annotation)##       name
## 1     <NA>
## 2     LIPI
## 3 TMPRSS15
## 4 TMPRSS15
## 5      APP
## 6      APP
##                                                                                                                                                                                                                                                                    annotation
## 1                                                                                                                                                                                                                                                                        <NA>
## 2                                                                                                             NM_001302998 NM_001302999 NM_001303000 NM_001303001 NM_145317 NM_198996 NP_001289927 NP_001289928 NP_001289929 NP_001289930 NP_945347 XM_006723965 XP_006724028
## 3                                                                                             NM_002772 NP_002763 XM_011529654 XM_011529655 XM_011529656 XM_011529657 XM_011529658 XM_011529659 XP_011527956 XP_011527957 XP_011527958 XP_011527959 XP_011527960 XP_011527961
## 4                                                                                             NM_002772 NP_002763 XM_011529654 XM_011529655 XM_011529656 XM_011529657 XM_011529658 XM_011529659 XP_011527956 XP_011527957 XP_011527958 XP_011527959 XP_011527960 XP_011527961
## 5 NM_000484 NM_001136016 NM_001136129 NM_001136130 NM_001136131 NM_001204301 NM_001204302 NM_001204303 NM_201413 NM_201414 NP_000475 NP_001129488 NP_001129601 NP_001129602 NP_001129603 NP_001191230 NP_001191231 NP_001191232 NP_958816 NP_958817 XM_024452075 XP_024307843
## 6 NM_000484 NM_001136016 NM_001136129 NM_001136130 NM_001136131 NM_001204301 NM_001204302 NM_001204303 NM_201413 NM_201414 NP_000475 NP_001129488 NP_001129601 NP_001129602 NP_001129603 NP_001191230 NP_001191231 NP_001191232 NP_958816 NP_958817 XM_024452075 XP_024307843
##   description      region distance   subregion insideDistance exonnumber
## 1 close to 3' close to 3'      815        <NA>             NA         NA
## 2  downstream  downstream   125774        <NA>             NA         NA
## 3    upstream    upstream   454170        <NA>             NA         NA
## 4    upstream    upstream   454475        <NA>             NA         NA
## 5 inside exon      inside   289903 inside exon              0         16
## 6 inside exon      inside   289899 inside exon              0         16
##   nexons   UTR strand  geneL codingL    Geneid subjectHits
## 1      1  <NA>      +     60      NA 100500815          15
## 2      8  <NA>      - 102077  101886    149998         149
## 3     25  <NA>      - 134537  133653      5651         579
## 4     25  <NA>      - 134537  133653      5651         579
## 5     16 3'UTR      - 290585  230434       351         354
## 6     16 3'UTR      - 290585  230434       351         354## You can use derfinderPlot::plotOverview() to visualize this informationWe can check the base-level coverage information for some of our DERs. In this example we do so for the first 5 ERs (Figures 2, 3, 4, 5, 6).
## Extract the region coverage
regionCov <- regionMat$chr21$bpCoverage
plotRegionCoverage(regions = ers, regionCoverage = regionCov, 
    groupInfo = pheno$group, nearestAnnotation = annotation, 
    annotatedRegions = annoRegs, whichRegions = seq_len(5), txdb = txdb,
    scalefac = 1, ask = FALSE, verbose = FALSE)
Figure 2: Base-pair resolution plot of differentially expressed region 1
Figure 3: Base-pair resolution plot of differentially expressed region 2
Figure 4: Base-pair resolution plot of differentially expressed region 3
Figure 5: Base-pair resolution plot of differentially expressed region 4
Figure 6: Base-pair resolution plot of differentially expressed region 5
You can then use the regionReport package to generate interactive HTML reports exploring the results.
If you are interested in using derfinder we recommend checking the derfinder users guide and good luck with your analyses!
This package was made possible thanks to:
Code for creating the vignette
## Create the vignette
library('rmarkdown')
system.time(render('derfinder-quickstart.Rmd', 'BiocStyle::html_document'))
## Extract the R code
library('knitr')
knit('derfinder-quickstart.Rmd', tangle = TRUE)## Clean up
file.remove('quickstartRef.bib')## [1] TRUEDate the vignette was generated.
## [1] "2019-09-20 20:50:51 EDT"Wallclock time spent generating the vignette.
## Time difference of 1.075 minsR session information.
## ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
##  setting  value                       
##  version  R version 3.6.1 (2019-07-05)
##  os       Ubuntu 18.04.3 LTS          
##  system   x86_64, linux-gnu           
##  ui       X11                         
##  language (EN)                        
##  collate  C                           
##  ctype    en_US.UTF-8                 
##  tz       America/New_York            
##  date     2019-09-20                  
## 
## ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
##  package                           * version   date       lib source        
##  acepack                             1.4.1     2016-10-29 [2] CRAN (R 3.6.1)
##  annotate                            1.62.0    2019-09-20 [2] Bioconductor  
##  AnnotationDbi                     * 1.46.1    2019-09-20 [2] Bioconductor  
##  AnnotationFilter                    1.8.0     2019-09-20 [2] Bioconductor  
##  assertthat                          0.2.1     2019-03-21 [2] CRAN (R 3.6.1)
##  backports                           1.1.4     2019-04-10 [2] CRAN (R 3.6.1)
##  base64enc                           0.1-3     2015-07-28 [2] CRAN (R 3.6.1)
##  bibtex                              0.4.2     2017-06-30 [2] CRAN (R 3.6.1)
##  Biobase                           * 2.44.0    2019-09-20 [2] Bioconductor  
##  BiocGenerics                      * 0.30.0    2019-09-20 [2] Bioconductor  
##  BiocManager                         1.30.4    2018-11-13 [2] CRAN (R 3.6.1)
##  BiocParallel                      * 1.18.1    2019-09-20 [2] Bioconductor  
##  BiocStyle                         * 2.12.0    2019-09-20 [2] Bioconductor  
##  biomaRt                             2.40.4    2019-09-20 [2] Bioconductor  
##  Biostrings                          2.52.0    2019-09-20 [2] Bioconductor  
##  biovizBase                          1.32.0    2019-09-20 [2] Bioconductor  
##  bit                                 1.1-14    2018-05-29 [2] CRAN (R 3.6.1)
##  bit64                               0.9-7     2017-05-08 [2] CRAN (R 3.6.1)
##  bitops                              1.0-6     2013-08-17 [2] CRAN (R 3.6.1)
##  blob                                1.2.0     2019-07-09 [2] CRAN (R 3.6.1)
##  bookdown                            0.13      2019-08-21 [2] CRAN (R 3.6.1)
##  BSgenome                            1.52.0    2019-09-20 [2] Bioconductor  
##  bumphunter                        * 1.26.0    2019-09-20 [2] Bioconductor  
##  checkmate                           1.9.4     2019-07-04 [2] CRAN (R 3.6.1)
##  cli                                 1.1.0     2019-03-19 [2] CRAN (R 3.6.1)
##  cluster                             2.1.0     2019-06-19 [2] CRAN (R 3.6.1)
##  codetools                           0.2-16    2018-12-24 [2] CRAN (R 3.6.1)
##  colorspace                          1.4-1     2019-03-18 [2] CRAN (R 3.6.1)
##  crayon                              1.3.4     2017-09-16 [2] CRAN (R 3.6.1)
##  curl                                4.1       2019-09-16 [2] CRAN (R 3.6.1)
##  data.table                          1.12.2    2019-04-07 [2] CRAN (R 3.6.1)
##  DBI                                 1.0.0     2018-05-02 [2] CRAN (R 3.6.1)
##  DelayedArray                      * 0.10.0    2019-09-20 [2] Bioconductor  
##  derfinder                         * 1.18.9    2019-09-20 [1] Bioconductor  
##  derfinderData                     * 2.2.1     2019-09-19 [2] Bioconductor  
##  derfinderHelper                     1.18.1    2019-09-20 [2] Bioconductor  
##  derfinderPlot                     * 1.18.1    2019-09-20 [2] Bioconductor  
##  DESeq2                            * 1.24.0    2019-09-20 [2] Bioconductor  
##  dichromat                           2.0-0     2013-01-24 [2] CRAN (R 3.6.1)
##  digest                              0.6.20    2019-07-04 [2] CRAN (R 3.6.1)
##  doRNG                               1.7.1     2018-06-22 [2] CRAN (R 3.6.1)
##  dplyr                               0.8.3     2019-07-04 [2] CRAN (R 3.6.1)
##  ensembldb                           2.8.0     2019-09-20 [2] Bioconductor  
##  evaluate                            0.14      2019-05-28 [2] CRAN (R 3.6.1)
##  foreach                           * 1.4.7     2019-07-27 [2] CRAN (R 3.6.1)
##  foreign                             0.8-72    2019-08-02 [2] CRAN (R 3.6.1)
##  Formula                             1.2-3     2018-05-03 [2] CRAN (R 3.6.1)
##  genefilter                          1.66.0    2019-09-20 [2] Bioconductor  
##  geneplotter                         1.62.0    2019-09-20 [2] Bioconductor  
##  GenomeInfoDb                      * 1.20.0    2019-09-20 [2] Bioconductor  
##  GenomeInfoDbData                    1.2.1     2019-07-05 [2] Bioconductor  
##  GenomicAlignments                   1.20.1    2019-09-20 [2] Bioconductor  
##  GenomicFeatures                   * 1.36.4    2019-09-20 [2] Bioconductor  
##  GenomicFiles                        1.20.0    2019-09-20 [2] Bioconductor  
##  GenomicRanges                     * 1.36.1    2019-09-20 [2] Bioconductor  
##  GGally                              1.4.0     2018-05-17 [2] CRAN (R 3.6.1)
##  ggbio                               1.32.0    2019-09-20 [2] Bioconductor  
##  ggplot2                             3.2.1     2019-08-10 [2] CRAN (R 3.6.1)
##  glue                                1.3.1     2019-03-12 [2] CRAN (R 3.6.1)
##  graph                               1.62.0    2019-09-20 [2] Bioconductor  
##  gridExtra                           2.3       2017-09-09 [2] CRAN (R 3.6.1)
##  gtable                              0.3.0     2019-03-25 [2] CRAN (R 3.6.1)
##  highr                               0.8       2019-03-20 [2] CRAN (R 3.6.1)
##  Hmisc                               4.2-0     2019-01-26 [2] CRAN (R 3.6.1)
##  hms                                 0.5.1     2019-08-23 [2] CRAN (R 3.6.1)
##  htmlTable                           1.13.1    2019-01-07 [2] CRAN (R 3.6.1)
##  htmltools                           0.3.6     2017-04-28 [2] CRAN (R 3.6.1)
##  htmlwidgets                         1.3       2018-09-30 [2] CRAN (R 3.6.1)
##  httr                                1.4.1     2019-08-05 [2] CRAN (R 3.6.1)
##  IRanges                           * 2.18.2    2019-09-20 [2] Bioconductor  
##  iterators                         * 1.0.12    2019-07-26 [2] CRAN (R 3.6.1)
##  jsonlite                            1.6       2018-12-07 [2] CRAN (R 3.6.1)
##  knitcitations                     * 1.0.10    2019-09-15 [2] CRAN (R 3.6.1)
##  knitr                               1.25      2019-09-18 [2] CRAN (R 3.6.1)
##  lattice                             0.20-38   2018-11-04 [2] CRAN (R 3.6.1)
##  latticeExtra                        0.6-28    2016-02-09 [2] CRAN (R 3.6.1)
##  lazyeval                            0.2.2     2019-03-15 [2] CRAN (R 3.6.1)
##  limma                               3.40.6    2019-09-20 [2] Bioconductor  
##  locfit                            * 1.5-9.1   2013-04-20 [2] CRAN (R 3.6.1)
##  lubridate                           1.7.4     2018-04-11 [2] CRAN (R 3.6.1)
##  magrittr                            1.5       2014-11-22 [2] CRAN (R 3.6.1)
##  Matrix                              1.2-17    2019-03-22 [2] CRAN (R 3.6.1)
##  matrixStats                       * 0.55.0    2019-09-07 [2] CRAN (R 3.6.1)
##  memoise                             1.1.0     2017-04-21 [2] CRAN (R 3.6.1)
##  munsell                             0.5.0     2018-06-12 [2] CRAN (R 3.6.1)
##  nnet                                7.3-12    2016-02-02 [2] CRAN (R 3.6.1)
##  org.Hs.eg.db                      * 3.8.2     2019-07-24 [2] Bioconductor  
##  OrganismDbi                         1.26.0    2019-09-20 [2] Bioconductor  
##  pillar                              1.4.2     2019-06-29 [2] CRAN (R 3.6.1)
##  pkgconfig                           2.0.2     2018-08-16 [2] CRAN (R 3.6.1)
##  pkgmaker                            0.27      2018-05-25 [2] CRAN (R 3.6.1)
##  plyr                                1.8.4     2016-06-08 [2] CRAN (R 3.6.1)
##  prettyunits                         1.0.2     2015-07-13 [2] CRAN (R 3.6.1)
##  progress                            1.2.2     2019-05-16 [2] CRAN (R 3.6.1)
##  ProtGenerics                        1.16.0    2019-09-20 [2] Bioconductor  
##  purrr                               0.3.2     2019-03-15 [2] CRAN (R 3.6.1)
##  qvalue                              2.16.0    2019-09-20 [2] Bioconductor  
##  R6                                  2.4.0     2019-02-14 [2] CRAN (R 3.6.1)
##  RBGL                                1.60.0    2019-09-20 [2] Bioconductor  
##  RColorBrewer                        1.1-2     2014-12-07 [2] CRAN (R 3.6.1)
##  Rcpp                                1.0.2     2019-07-25 [2] CRAN (R 3.6.1)
##  RCurl                               1.95-4.12 2019-03-04 [2] CRAN (R 3.6.1)
##  RefManageR                          1.2.12    2019-04-03 [2] CRAN (R 3.6.1)
##  registry                            0.5-1     2019-03-05 [2] CRAN (R 3.6.1)
##  reshape                             0.8.8     2018-10-23 [2] CRAN (R 3.6.1)
##  reshape2                            1.4.3     2017-12-11 [2] CRAN (R 3.6.1)
##  rlang                               0.4.0     2019-06-25 [2] CRAN (R 3.6.1)
##  rmarkdown                           1.15      2019-08-21 [2] CRAN (R 3.6.1)
##  rngtools                            1.4       2019-07-01 [2] CRAN (R 3.6.1)
##  rpart                               4.1-15    2019-04-12 [2] CRAN (R 3.6.1)
##  Rsamtools                           2.0.1     2019-09-20 [2] Bioconductor  
##  RSQLite                             2.1.2     2019-07-24 [2] CRAN (R 3.6.1)
##  rstudioapi                          0.10      2019-03-19 [2] CRAN (R 3.6.1)
##  rtracklayer                         1.44.4    2019-09-20 [2] Bioconductor  
##  S4Vectors                         * 0.22.1    2019-09-20 [2] Bioconductor  
##  scales                              1.0.0     2018-08-09 [2] CRAN (R 3.6.1)
##  sessioninfo                       * 1.1.1     2018-11-05 [2] CRAN (R 3.6.1)
##  stringi                             1.4.3     2019-03-12 [2] CRAN (R 3.6.1)
##  stringr                             1.4.0     2019-02-10 [2] CRAN (R 3.6.1)
##  SummarizedExperiment              * 1.14.1    2019-09-20 [2] Bioconductor  
##  survival                            2.44-1.1  2019-04-01 [2] CRAN (R 3.6.1)
##  tibble                              2.1.3     2019-06-06 [2] CRAN (R 3.6.1)
##  tidyselect                          0.2.5     2018-10-11 [2] CRAN (R 3.6.1)
##  TxDb.Hsapiens.UCSC.hg19.knownGene * 3.2.2     2019-07-05 [2] Bioconductor  
##  VariantAnnotation                   1.30.1    2019-09-20 [2] Bioconductor  
##  vctrs                               0.2.0     2019-07-05 [2] CRAN (R 3.6.1)
##  withr                               2.1.2     2018-03-15 [2] CRAN (R 3.6.1)
##  xfun                                0.9       2019-08-21 [2] CRAN (R 3.6.1)
##  XML                                 3.98-1.20 2019-06-06 [2] CRAN (R 3.6.1)
##  xml2                                1.2.2     2019-08-09 [2] CRAN (R 3.6.1)
##  xtable                              1.8-4     2019-04-21 [2] CRAN (R 3.6.1)
##  XVector                             0.24.0    2019-09-20 [2] Bioconductor  
##  yaml                                2.2.0     2018-07-25 [2] CRAN (R 3.6.1)
##  zeallot                             0.1.0     2018-01-28 [2] CRAN (R 3.6.1)
##  zlibbioc                            1.30.0    2019-09-20 [2] Bioconductor  
## 
## [1] /tmp/RtmpGvZKFF/Rinst3a8b30ec392c
## [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, 2019).
[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.
[2] Bioconductor Package Maintainer, V. Obenchain, M. Love, L. Shepherd, et al. GenomicFiles: Distributed computing by file or by range. R package version 1.20.0. 2019.
[3] C. Boettiger. knitcitations: Citations for ‘Knitr’ Markdown Files. R package version 1.0.10. 2019. <URL: https://CRAN.R-project.org/package=knitcitations>.
[4] BrainSpan. “Atlas of the Developing Human Brain [Internet]. Funded by ARRA Awards 1RC2MH089921-01, 1RC2MH090047-01, and 1RC2MH089929-01.” 2011. <URL: http://www.brainspan.org/>.
[5] M. Carlson and B. P. Maintainer. TxDb.Hsapiens.UCSC.hg19.knownGene: Annotation package for TxDb object(s). R package version 3.2.2. 2015.
[6] L. Collado-Torres, A. Jaffe, and J. Leek. derfinderData: Processed BigWigs from BrainSpan for examples. R package version 2.2.1. 2019. <URL: https://github.com/leekgroup/derfinderData>.
[7] 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>.
[8] A. C. Frazee, S. Sabunciyan, K. D. Hansen, R. A. Irizarry, et al. “Differential expression analysis of RNA-seq data at single-base resolution”. In: Biostatistics 15 (3) (2014), pp. 413-426. DOI: 10.1093/biostatistics/kxt053. <URL: http://biostatistics.oxfordjournals.org/content/15/3/413.long>.
[9] F. E. Harrell Jr, w. c. f. C. Dupont, and m. others. Hmisc: Harrell Miscellaneous. R package version 4.2-0. 2019. <URL: https://CRAN.R-project.org/package=Hmisc>.
[10] A. E. Jaffe, P. Murakami, H. Lee, J. T. Leek, et al. “Bump hunting to identify differentially methylated regions in epigenetic epidemiology studies”. In: International journal of epidemiology 41.1 (2012), pp. 200-209. DOI: 10.1093/ije/dyr238.
[11] A. E. Jaffe, P. Murakami, H. Lee, J. T. Leek, et al. “Bump hunting to identify differentially methylated regions in epigenetic epidemiology studies”. In: International Journal of Epidemiology (2012).
[12] M. Lawrence, R. Gentleman, and V. Carey. “rtracklayer: an R package for interfacing with genome browsers”. In: Bioinformatics 25 (2009), pp. 1841-1842. DOI: 10.1093/bioinformatics/btp328. <URL: http://bioinformatics.oxfordjournals.org/content/25/14/1841.abstract>.
[13] M. Lawrence, W. Huber, H. Pagès, P. Aboyoun, et al. “Software for Computing and Annotating Genomic Ranges”. In: PLoS Computational Biology 9 (8 2013). DOI: 10.1371/journal.pcbi.1003118. <URL: http://www.ploscompbiol.org/article/info%3Adoi%2F10.1371%2Fjournal.pcbi.1003118}.>
[14] M. Lawrence, W. Huber, H. Pagès, P. Aboyoun, et al. “Software for Computing and Annotating Genomic Ranges”. In: PLoS Computational Biology 9 (8 2013). DOI: 10.1371/journal.pcbi.1003118. <URL: http://www.ploscompbiol.org/article/info%3Adoi%2F10.1371%2Fjournal.pcbi.1003118}.>
[15] M. Lawrence, W. Huber, H. Pagès, P. Aboyoun, et al. “Software for Computing and Annotating Genomic Ranges”. In: PLoS Computational Biology 9 (8 2013). DOI: 10.1371/journal.pcbi.1003118. <URL: http://www.ploscompbiol.org/article/info%3Adoi%2F10.1371%2Fjournal.pcbi.1003118}.>
[16] M. Lawrence, W. Huber, H. Pagès, P. Aboyoun, et al. “Software for Computing and Annotating Genomic Ranges”. In: PLoS Computational Biology 9 (8 2013). DOI: 10.1371/journal.pcbi.1003118. <URL: http://www.ploscompbiol.org/article/info%3Adoi%2F10.1371%2Fjournal.pcbi.1003118}.>
[17] M. Morgan, V. Obenchain, M. Lang, R. Thompson, et al. BiocParallel: Bioconductor facilities for parallel evaluation. R package version 1.18.1. 2019. <URL: https://github.com/Bioconductor/BiocParallel>.
[18] M. Morgan, H. Pagès, V. Obenchain, and N. Hayden. Rsamtools: Binary alignment (BAM), FASTA, variant call (BCF), and tabix file import. R package version 2.0.1. 2019. <URL: http://bioconductor.org/packages/Rsamtools>.
[19] 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>.
[20] H. Pagès, M. Carlson, S. Falcon, and N. Li. AnnotationDbi: Annotation Database Interface. 2017. DOI: 10.18129/B9.bioc.AnnotationDbi.
[21] H. Pagès, M. Lawrence, and P. Aboyoun. S4Vectors: S4 implementation of vector-like and list-like objects. 2017. DOI: 10.18129/B9.bioc.S4Vectors.
[22] 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/>.
[23] J. D. Storey, A. J. Bass, A. Dabney, and D. Robinson. qvalue: Q-value estimation for false discovery rate control. R package version 2.16.0. 2019. <URL: http://github.com/jdstorey/qvalue>.
[24] H. Wickham. ggplot2: Elegant Graphics for Data Analysis. Springer-Verlag New York, 2016. ISBN: 978-3-319-24277-4. <URL: https://ggplot2.tidyverse.org>.
[25] H. Wickham. “testthat: Get Started with Testing”. In: The R Journal 3 (2011), pp. 5-10. <URL: https://journal.r-project.org/archive/2011-1/RJournal_2011-1_Wickham.pdf>.
[26] 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>.
[27] T. Yin, M. Lawrence, and D. Cook. biovizBase: Basic graphic utilities for visualization of genomic data. R package version 1.32.0. 2019.