derfinderPlot 1.41.0
R
is an open-source statistical environment which can be easily modified to enhance its functionality via packages. derfinderPlot 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 derfinderPlot by using the following commands in your R
session:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("derfinderPlot")
## Check that you have a valid Bioconductor installation
BiocManager::valid()
derfinderPlot is based on many other packages and in particular in those that have implemented the infrastructure needed for dealing with RNA-seq data. A derfinderPlot user is not expected to deal with those packages directly but will need to be familiar with derfinder and for some plots with ggbio.
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
or derfinderPlot
tags 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 derfinderPlot will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
## Citation info
citation("derfinderPlot")
## To cite package 'derfinderPlot' in publications use:
##
## Collado-Torres L, Jaffe AE, Leek JT (2017). _derfinderPlot: Plotting
## functions for derfinder_. doi:10.18129/B9.bioc.derfinderPlot
## <https://doi.org/10.18129/B9.bioc.derfinderPlot>,
## https://github.com/leekgroup/derfinderPlot - R package version
## 1.41.0, <http://www.bioconductor.org/packages/derfinderPlot>.
##
## 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 <https://doi.org/10.1093/nar/gkw852>,
## <http://nar.oxfordjournals.org/content/early/2016/09/29/nar.gkw852>.
##
## To see these entries in BibTeX format, use 'print(<citation>,
## bibtex=TRUE)', 'toBibtex(.)', or set
## 'options(citation.bibtex.max=999)'.
derfinderPlot (Collado-Torres, Jaffe, and Leek, 2017) is an addon package for derfinder (Collado-Torres, Nellore, Frazee, Wilks, Love, Langmead, Irizarry, Leek, and Jaffe, 2017) with functions that allow you to visualize the results.
While the functions in derfinderPlot assume you generated the data with derfinder, they can be used with other GRanges
objects properly formatted.
The functions in derfinderPlot are:
plotCluster()
is a tailored ggbio (Yin, Cook, and Lawrence, 2012) plot that shows all the regions in a cluster (defined by distance). It shows the base-level coverage for each sample as well as the mean for each group. If these regions overlap any known gene, the gene and the transcript annotation is displayed.plotOverview()
is another tailored ggbio (Yin, Cook, and Lawrence, 2012) plot showing an overview of the whole genome. This plot can be useful to observe if the regions are clustered in a subset of a chromosome. It can also be used to check whether the regions match predominantly one part of the gene structure (for example, 3’ overlaps).plotRegionCoverage()
is a fast plotting function using R
base graphics that shows the base-level coverage for each sample inside a specific region of the genome. If the region overlaps any known gene or intron, the information is displayed. Optionally, it can display the known transcripts. This function is most likely the easiest to use with GRanges
objects from other packages.As an example, we will analyze a small subset of the samples from the BrainSpan Atlas of the Human Brain (BrainSpan, 2011) publicly available data.
We first load the required packages.
## Load libraries
suppressPackageStartupMessages(library("derfinder"))
library("derfinderData")
library("derfinderPlot")
## Registered S3 method overwritten by 'GGally':
## method from
## +.gg ggplot2
For this example, we created a small table with the relevant phenotype data for 12 samples: 6 from fetal samples and 6 from adult samples. We chose at random a brain region, in this case the primary auditory cortex (core) and for the example we will only look at data from chromosome 21. Other variables include the age in years and the gender. The data is shown below.
library("knitr")
## Get pheno table
pheno <- subset(brainspanPheno, structure_acronym == "A1C")
## Display the main information
p <- pheno[, -which(colnames(pheno) %in% c(
"structure_acronym",
"structure_name", "file"
))]
rownames(p) <- NULL
kable(p, format = "html", row.names = TRUE)
gender | lab | Age | group | |
---|---|---|---|---|
1 | M | HSB114.A1C | -0.5192308 | fetal |
2 | M | HSB103.A1C | -0.5192308 | fetal |
3 | M | HSB178.A1C | -0.4615385 | fetal |
4 | M | HSB154.A1C | -0.4615385 | fetal |
5 | F | HSB150.A1C | -0.5384615 | fetal |
6 | F | HSB149.A1C | -0.5192308 | fetal |
7 | F | HSB130.A1C | 21.0000000 | adult |
8 | M | HSB136.A1C | 23.0000000 | adult |
9 | F | HSB126.A1C | 30.0000000 | adult |
10 | M | HSB145.A1C | 36.0000000 | adult |
11 | M | HSB123.A1C | 37.0000000 | adult |
12 | F | HSB135.A1C | 40.0000000 | adult |
We can load the data from derfinderData (Collado-Torres, Jaffe, and Leek, 2024) by first identifying the paths to the BigWig files with derfinder::rawFiles()
and then loading the data with derfinder::fullCoverage()
.
## Determine the files to use and fix the names
files <- rawFiles(system.file("extdata", "A1C", package = "derfinderData"),
samplepatt = "bw", fileterm = NULL
)
names(files) <- gsub(".bw", "", names(files))
## Load the data from disk
system.time(fullCov <- fullCoverage(files = files, chrs = "chr21"))
## 2024-12-13 16:52:17.610884 fullCoverage: processing chromosome chr21
## 2024-12-13 16:52:17.63084 loadCoverage: finding chromosome lengths
## 2024-12-13 16:52:17.665901 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB103.bw
## 2024-12-13 16:52:17.956932 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB114.bw
## 2024-12-13 16:52:18.166625 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB123.bw
## 2024-12-13 16:52:18.317637 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB126.bw
## 2024-12-13 16:52:18.424821 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB130.bw
## 2024-12-13 16:52:18.556323 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB135.bw
## 2024-12-13 16:52:18.67727 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB136.bw
## 2024-12-13 16:52:18.797378 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB145.bw
## 2024-12-13 16:52:18.932188 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB149.bw
## 2024-12-13 16:52:19.068849 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB150.bw
## 2024-12-13 16:52:19.178278 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB154.bw
## 2024-12-13 16:52:19.336473 loadCoverage: loading BigWig file /home/biocbuild/bbs-3.21-bioc/R/site-library/derfinderData/extdata/A1C/HSB178.bw
## 2024-12-13 16:52:19.47859 loadCoverage: applying the cutoff to the merged data
## 2024-12-13 16:52:19.512034 filterData: originally there were 48129895 rows, now there are 48129895 rows. Meaning that 0 percent was filtered.
## user system elapsed
## 1.865 0.083 1.970
Alternatively, since the BigWig files are publicly available from BrainSpan (see here), we can extract the relevant coverage data using derfinder::fullCoverage()
. Note that as of rtracklayer 1.25.16 BigWig files are not supported on Windows: you can find the fullCov
object inside derfinderData to follow the examples.
## Determine the files to use and fix the names
files <- pheno$file
names(files) <- gsub(".A1C", "", pheno$lab)
## Load the data from the web
system.time(fullCov <- fullCoverage(files = files, chrs = "chr21"))
Once we have the base-level coverage data for all 12 samples, we can construct the models. In this case, we want to find differences between fetal and adult samples while adjusting for gender and a proxy of the library size.
## Get some idea of the library sizes
sampleDepths <- sampleDepth(collapseFullCoverage(fullCov), 1)
## 2024-12-13 16:52:19.618304 sampleDepth: Calculating sample quantiles
## 2024-12-13 16:52:19.768297 sampleDepth: Calculating sample adjustments
## Define models
models <- makeModels(sampleDepths,
testvars = pheno$group,
adjustvars = pheno[, c("gender")]
)
Next, we can find candidate differentially expressed regions (DERs) using as input the segments of the genome where at least one sample has coverage greater than 3. In this particular example, we chose a low theoretical F-statistic cutoff and used 20 permutations.
## Filter coverage
filteredCov <- lapply(fullCov, filterData, cutoff = 3)
## 2024-12-13 16:52:20.085433 filterData: originally there were 48129895 rows, now there are 90023 rows. Meaning that 99.81 percent was filtered.
## Perform differential expression analysis
suppressPackageStartupMessages(library("bumphunter"))
system.time(results <- analyzeChr(
chr = "chr21", filteredCov$chr21,
models, groupInfo = pheno$group, writeOutput = FALSE,
cutoffFstat = 5e-02, nPermute = 20, seeds = 20140923 + seq_len(20)
))
## 2024-12-13 16:52:21.246544 analyzeChr: Pre-processing the coverage data
## 2024-12-13 16:52:23.021404 analyzeChr: Calculating statistics
## 2024-12-13 16:52:23.025183 calculateStats: calculating the F-statistics
## 2024-12-13 16:52:24.196326 analyzeChr: Calculating pvalues
## 2024-12-13 16:52:24.197037 analyzeChr: Using the following theoretical cutoff for the F-statistics 5.31765507157871
## 2024-12-13 16:52:24.19846 calculatePvalues: identifying data segments
## 2024-12-13 16:52:24.206455 findRegions: segmenting information
## 2024-12-13 16:52:24.234482 findRegions: identifying candidate regions
## 2024-12-13 16:52:24.288839 findRegions: identifying region clusters
## 2024-12-13 16:52:24.442096 calculatePvalues: calculating F-statistics for permutation 1 and seed 20140924
## 2024-12-13 16:52:24.582569 findRegions: segmenting information
## 2024-12-13 16:52:24.606701 findRegions: identifying candidate regions
## 2024-12-13 16:52:24.669806 calculatePvalues: calculating F-statistics for permutation 2 and seed 20140925
## 2024-12-13 16:52:24.827055 findRegions: segmenting information
## 2024-12-13 16:52:24.851046 findRegions: identifying candidate regions
## 2024-12-13 16:52:24.896288 calculatePvalues: calculating F-statistics for permutation 3 and seed 20140926
## 2024-12-13 16:52:25.036293 findRegions: segmenting information
## 2024-12-13 16:52:25.060702 findRegions: identifying candidate regions
## 2024-12-13 16:52:25.106425 calculatePvalues: calculating F-statistics for permutation 4 and seed 20140927
## 2024-12-13 16:52:25.252255 findRegions: segmenting information
## 2024-12-13 16:52:25.285566 findRegions: identifying candidate regions
## 2024-12-13 16:52:25.330724 calculatePvalues: calculating F-statistics for permutation 5 and seed 20140928
## 2024-12-13 16:52:25.468399 findRegions: segmenting information
## 2024-12-13 16:52:25.4923 findRegions: identifying candidate regions
## 2024-12-13 16:52:25.536163 calculatePvalues: calculating F-statistics for permutation 6 and seed 20140929
## 2024-12-13 16:52:25.681477 findRegions: segmenting information
## 2024-12-13 16:52:25.705438 findRegions: identifying candidate regions
## 2024-12-13 16:52:25.75999 calculatePvalues: calculating F-statistics for permutation 7 and seed 20140930
## 2024-12-13 16:52:25.898296 findRegions: segmenting information
## 2024-12-13 16:52:25.922197 findRegions: identifying candidate regions
## 2024-12-13 16:52:25.966244 calculatePvalues: calculating F-statistics for permutation 8 and seed 20140931
## 2024-12-13 16:52:26.111007 findRegions: segmenting information
## 2024-12-13 16:52:26.134979 findRegions: identifying candidate regions
## 2024-12-13 16:52:26.179971 calculatePvalues: calculating F-statistics for permutation 9 and seed 20140932
## 2024-12-13 16:52:26.330541 findRegions: segmenting information
## 2024-12-13 16:52:26.354654 findRegions: identifying candidate regions
## 2024-12-13 16:52:26.400576 calculatePvalues: calculating F-statistics for permutation 10 and seed 20140933
## 2024-12-13 16:52:26.548942 findRegions: segmenting information
## 2024-12-13 16:52:26.573665 findRegions: identifying candidate regions
## 2024-12-13 16:52:26.620249 calculatePvalues: calculating F-statistics for permutation 11 and seed 20140934
## 2024-12-13 16:52:26.772172 findRegions: segmenting information
## 2024-12-13 16:52:26.79616 findRegions: identifying candidate regions
## 2024-12-13 16:52:26.841942 calculatePvalues: calculating F-statistics for permutation 12 and seed 20140935
## 2024-12-13 16:52:26.991589 findRegions: segmenting information
## 2024-12-13 16:52:27.024792 findRegions: identifying candidate regions
## 2024-12-13 16:52:27.069554 calculatePvalues: calculating F-statistics for permutation 13 and seed 20140936
## 2024-12-13 16:52:27.206728 findRegions: segmenting information
## 2024-12-13 16:52:27.230637 findRegions: identifying candidate regions
## 2024-12-13 16:52:27.274562 calculatePvalues: calculating F-statistics for permutation 14 and seed 20140937
## 2024-12-13 16:52:27.429873 findRegions: segmenting information
## 2024-12-13 16:52:27.453501 findRegions: identifying candidate regions
## 2024-12-13 16:52:27.497667 calculatePvalues: calculating F-statistics for permutation 15 and seed 20140938
## 2024-12-13 16:52:27.64746 findRegions: segmenting information
## 2024-12-13 16:52:27.67124 findRegions: identifying candidate regions
## 2024-12-13 16:52:27.715327 calculatePvalues: calculating F-statistics for permutation 16 and seed 20140939
## 2024-12-13 16:52:27.868175 findRegions: segmenting information
## 2024-12-13 16:52:27.892144 findRegions: identifying candidate regions
## 2024-12-13 16:52:27.936425 calculatePvalues: calculating F-statistics for permutation 17 and seed 20140940
## 2024-12-13 16:52:28.075556 findRegions: segmenting information
## 2024-12-13 16:52:28.099492 findRegions: identifying candidate regions
## 2024-12-13 16:52:28.144386 calculatePvalues: calculating F-statistics for permutation 18 and seed 20140941
## 2024-12-13 16:52:28.30785 findRegions: segmenting information
## 2024-12-13 16:52:28.33219 findRegions: identifying candidate regions
## 2024-12-13 16:52:28.378836 calculatePvalues: calculating F-statistics for permutation 19 and seed 20140942
## 2024-12-13 16:52:28.51874 findRegions: segmenting information
## 2024-12-13 16:52:28.542932 findRegions: identifying candidate regions
## 2024-12-13 16:52:28.588279 calculatePvalues: calculating F-statistics for permutation 20 and seed 20140943
## 2024-12-13 16:52:28.741893 findRegions: segmenting information
## 2024-12-13 16:52:28.766794 findRegions: identifying candidate regions
## 2024-12-13 16:52:28.83848 calculatePvalues: calculating the p-values
## 2024-12-13 16:52:28.90887 analyzeChr: Annotating regions
## No annotationPackage supplied. Trying org.Hs.eg.db.
## Loading required package: org.Hs.eg.db
## Loading required package: AnnotationDbi
## 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")'.
##
## Getting TSS and TSE.
## Getting CSS and CSE.
## Getting exons.
## Annotating genes.
## ...
## user system elapsed
## 57.828 0.977 58.834
## Quick access to the results
regions <- results$regions$regions
## Annotation database to use
suppressPackageStartupMessages(library("TxDb.Hsapiens.UCSC.hg19.knownGene"))
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
plotOverview()
Now that we have obtained the main results using derfinder, we can proceed to visualizing the results using derfinderPlot. The easiest to use of all the functions is plotOverview()
which takes a set of regions and annotation information produced by bumphunter::matchGenes()
.
Figure 1 shows the candidate DERs colored by whether their q-value was less than 0.10 or not.
## Q-values overview
plotOverview(regions = regions, annotation = results$annotation, type = "qval")
## 2024-12-13 16:53:20.257959 plotOverview: assigning chromosome lengths from hg19!
## 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.
Figure 2 shows the candidate DERs colored by the type of gene feature they are nearest too.
## Annotation overview
plotOverview(
regions = regions, annotation = results$annotation,
type = "annotation"
)
## 2024-12-13 16:53:22.134674 plotOverview: assigning chromosome lengths from hg19!
## Scale for x is already present.
## Adding another scale for x, which will replace the existing scale.
In this particular example, because we are only using data from one chromosome the above plot is not as informative as in a real case scenario. However, with this plot we can quickly observe that nearly all of the candidate DERs are inside an exon.
plotRegionCoverage()
The complete opposite of visualizing the candidate DERs at the genome-level is to visualize them one region at a time. plotRegionCoverage()
allows us to do this quickly for a large number of regions.
Before using this function, we need to process more detailed information using two derfinder functions: annotateRegions()
and getRegionCoverage()
as shown below.
## Get required information for the plots
annoRegs <- annotateRegions(regions, genomicState$fullGenome)
## 2024-12-13 16:53:23.165555 annotateRegions: counting
## 2024-12-13 16:53:23.25229 annotateRegions: annotating
regionCov <- getRegionCoverage(fullCov, regions)
## 2024-12-13 16:53:23.382745 getRegionCoverage: processing chr21
## 2024-12-13 16:53:23.446866 getRegionCoverage: done processing chr21
Once we have the relevant information we can proceed to plotting the first 10 regions. In this case, we will supply plotRegionCoverage()
with the information it needs to plot transcripts overlapping these 10 regions (Figures ??, ??, ??, ??, ??, ??, ??, ??, ??, ??).
## Plot top 10 regions
plotRegionCoverage(
regions = regions, regionCoverage = regionCov,
groupInfo = pheno$group, nearestAnnotation = results$annotation,
annotatedRegions = annoRegs, whichRegions = 1:10, txdb = txdb, scalefac = 1,
ask = FALSE, verbose = FALSE
)
The base-level coverage is shown in a log2 scale with any overlapping exons shown in dark blue and known introns in light blue.
plotCluster()
In this example, we noticed with the plotRegionCoverage()
plots that most of the candidate DERs are contained in known exons. Sometimes, the signal might be low or we might have used very stringent cutoffs in the derfinder analysis. One way we can observe this is by plotting clusters of regions where a cluster is defined as regions within 300 bp (default option) of each other.
To visualize the clusters, we can use plotCluster()
which takes similar input to plotOverview()
with the notable addition of the coverage information as well as the idx
argument. This argument specifies which region to focus on: it will be plotted with a red bar and will determine the cluster to display.
In Figure 13 we observe one large candidate DER with other nearby ones that do not have a q-value less than 0.10. In a real analysis, we would probably discard this region as the coverage is fairly low.
## First cluster
plotCluster(
idx = 1, regions = regions, annotation = results$annotation,
coverageInfo = fullCov$chr21, txdb = txdb, groupInfo = pheno$group,
titleUse = "pval"
)
## Parsing transcripts...
## Parsing exons...
## Parsing cds...
## Parsing utrs...
## ------exons...
## ------cdss...
## ------introns...
## ------utr...
## aggregating...
## Done
## Constructing graphics...
The second cluster (Figure 14) shows a larger number of potential DERs (again without q-values less than 0.10) in a segment of the genome where the coverage data is highly variable. This is a common occurrence with RNA-seq data.
## Second cluster
plotCluster(
idx = 2, regions = regions, annotation = results$annotation,
coverageInfo = fullCov$chr21, txdb = txdb, groupInfo = pheno$group,
titleUse = "pval"
)
## Parsing transcripts...
## Parsing exons...
## Parsing cds...
## Parsing utrs...
## ------exons...
## ------cdss...
## ------introns...
## ------utr...
## aggregating...
## Done
## Constructing graphics...
## Warning in !vapply(ggl, fixed, logical(1L)) & !vapply(PlotList, is, "Ideogram",
## : longer object length is not a multiple of shorter object length
## Warning in scale_y_continuous(trans = log2_trans()): log-2 transformation
## introduced infinite values.
These plots show an ideogram which helps quickly identify which region of the genome we are focusing on. Then, the base-level coverage information for each sample is displayed in log2. Next, the coverage group means are shown in the log2 scale. The plot is completed with the potential and candidate DERs as well as any known transcripts.
vennRegions
derfinder has functions for annotating regions given their genomic state. A typical visualization is to then view how many regions overlap known exons, introns, intergenic regions, none of them or several of these groups in a venn diagram. The function vennRegions()
makes this plot using the output from derfinder::annotateRegions()
as shown in Figure 15.
## Make venn diagram
venn <- vennRegions(annoRegs)
## It returns the actual venn counts information
venn
## exon intergenic intron Counts
## 1 0 0 0 0
## 2 0 0 1 2
## 3 0 1 0 4
## 4 0 1 1 0
## 5 1 0 0 259
## 6 1 0 1 35
## 7 1 1 0 0
## 8 1 1 1 0
## attr(,"class")
## [1] "VennCounts"
This package was made possible thanks to:
Code for creating the vignette
## Create the vignette
library("rmarkdown")
system.time(render("derfinderPlot.Rmd", "BiocStyle::html_document"))
## Extract the R code
library("knitr")
knit("derfinderPlot.Rmd", tangle = TRUE)
## Clean up
unlink("chr21", recursive = TRUE)
Date the vignette was generated.
## [1] "2024-12-13 16:53:51 EST"
Wallclock time spent generating the vignette.
## Time difference of 1.857 mins
R
session information.
## ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
## setting value
## version R Under development (unstable) (2024-10-21 r87258)
## os Ubuntu 24.04.1 LTS
## system x86_64, linux-gnu
## ui X11
## language (EN)
## collate C
## ctype en_US.UTF-8
## tz America/New_York
## date 2024-12-13
## pandoc 3.1.3 @ /usr/bin/ (via rmarkdown)
##
## ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
## package * version date (UTC) lib source
## abind 1.4-8 2024-09-12 [2] CRAN (R 4.5.0)
## AnnotationDbi * 1.69.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## AnnotationFilter 1.31.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## backports 1.5.0 2024-05-23 [2] CRAN (R 4.5.0)
## base64enc 0.1-3 2015-07-28 [2] CRAN (R 4.5.0)
## bibtex 0.5.1 2023-01-26 [2] CRAN (R 4.5.0)
## Biobase * 2.67.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## BiocFileCache 2.15.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## BiocGenerics * 0.53.3 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## BiocIO 1.17.1 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## BiocManager 1.30.25 2024-08-28 [2] CRAN (R 4.5.0)
## BiocParallel 1.41.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## BiocStyle * 2.35.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## biomaRt 2.63.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## Biostrings 2.75.2 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## biovizBase 1.55.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## bit 4.5.0.1 2024-12-03 [2] CRAN (R 4.5.0)
## bit64 4.5.2 2024-09-22 [2] CRAN (R 4.5.0)
## bitops 1.0-9 2024-10-03 [2] CRAN (R 4.5.0)
## blob 1.2.4 2023-03-17 [2] CRAN (R 4.5.0)
## bookdown 0.41 2024-10-16 [2] CRAN (R 4.5.0)
## BSgenome 1.75.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## bslib 0.8.0 2024-07-29 [2] CRAN (R 4.5.0)
## bumphunter * 1.49.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## cachem 1.1.0 2024-05-16 [2] CRAN (R 4.5.0)
## checkmate 2.3.2 2024-07-29 [2] CRAN (R 4.5.0)
## cli 3.6.3 2024-06-21 [2] CRAN (R 4.5.0)
## cluster 2.1.8 2024-12-11 [3] CRAN (R 4.5.0)
## codetools 0.2-20 2024-03-31 [3] CRAN (R 4.5.0)
## colorspace 2.1-1 2024-07-26 [2] CRAN (R 4.5.0)
## crayon 1.5.3 2024-06-20 [2] CRAN (R 4.5.0)
## curl 6.0.1 2024-11-14 [2] CRAN (R 4.5.0)
## data.table 1.16.4 2024-12-06 [2] CRAN (R 4.5.0)
## DBI 1.2.3 2024-06-02 [2] CRAN (R 4.5.0)
## dbplyr 2.5.0 2024-03-19 [2] CRAN (R 4.5.0)
## DelayedArray 0.33.3 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## derfinder * 1.41.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## derfinderData * 2.25.0 2024-12-12 [2] Bioconductor 3.21 (R 4.5.0)
## derfinderHelper 1.41.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## derfinderPlot * 1.41.0 2024-12-13 [1] Bioconductor 3.21 (R 4.5.0)
## dichromat 2.0-0.1 2022-05-02 [2] CRAN (R 4.5.0)
## digest 0.6.37 2024-08-19 [2] CRAN (R 4.5.0)
## doRNG 1.8.6 2023-01-16 [2] CRAN (R 4.5.0)
## dplyr 1.1.4 2023-11-17 [2] CRAN (R 4.5.0)
## ensembldb 2.31.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## evaluate 1.0.1 2024-10-10 [2] CRAN (R 4.5.0)
## fansi 1.0.6 2023-12-08 [2] CRAN (R 4.5.0)
## farver 2.1.2 2024-05-13 [2] CRAN (R 4.5.0)
## fastmap 1.2.0 2024-05-15 [2] CRAN (R 4.5.0)
## filelock 1.0.3 2023-12-11 [2] CRAN (R 4.5.0)
## foreach * 1.5.2 2022-02-02 [2] CRAN (R 4.5.0)
## foreign 0.8-87 2024-06-26 [3] CRAN (R 4.5.0)
## Formula 1.2-5 2023-02-24 [2] CRAN (R 4.5.0)
## generics * 0.1.3 2022-07-05 [2] CRAN (R 4.5.0)
## GenomeInfoDb * 1.43.2 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## GenomeInfoDbData 1.2.13 2024-10-23 [2] Bioconductor
## GenomicAlignments 1.43.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## GenomicFeatures * 1.59.1 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## GenomicFiles 1.43.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## GenomicRanges * 1.59.1 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## GGally 2.2.1 2024-02-14 [2] CRAN (R 4.5.0)
## ggbio 1.55.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## ggplot2 3.5.1 2024-04-23 [2] CRAN (R 4.5.0)
## ggstats 0.7.0 2024-09-22 [2] CRAN (R 4.5.0)
## glue 1.8.0 2024-09-30 [2] CRAN (R 4.5.0)
## graph 1.85.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## gridExtra 2.3 2017-09-09 [2] CRAN (R 4.5.0)
## gtable 0.3.6 2024-10-25 [2] CRAN (R 4.5.0)
## Hmisc 5.2-1 2024-12-02 [2] CRAN (R 4.5.0)
## hms 1.1.3 2023-03-21 [2] CRAN (R 4.5.0)
## htmlTable 2.4.3 2024-07-21 [2] CRAN (R 4.5.0)
## htmltools 0.5.8.1 2024-04-04 [2] CRAN (R 4.5.0)
## htmlwidgets 1.6.4 2023-12-06 [2] CRAN (R 4.5.0)
## httr 1.4.7 2023-08-15 [2] CRAN (R 4.5.0)
## httr2 1.0.7 2024-11-26 [2] CRAN (R 4.5.0)
## IRanges * 2.41.2 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## iterators * 1.0.14 2022-02-05 [2] CRAN (R 4.5.0)
## jquerylib 0.1.4 2021-04-26 [2] CRAN (R 4.5.0)
## jsonlite 1.8.9 2024-09-20 [2] CRAN (R 4.5.0)
## KEGGREST 1.47.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## knitr * 1.49 2024-11-08 [2] CRAN (R 4.5.0)
## labeling 0.4.3 2023-08-29 [2] CRAN (R 4.5.0)
## lattice 0.22-6 2024-03-20 [3] CRAN (R 4.5.0)
## lazyeval 0.2.2 2019-03-15 [2] CRAN (R 4.5.0)
## lifecycle 1.0.4 2023-11-07 [2] CRAN (R 4.5.0)
## limma 3.63.2 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## locfit * 1.5-9.10 2024-06-24 [2] CRAN (R 4.5.0)
## lubridate 1.9.4 2024-12-08 [2] CRAN (R 4.5.0)
## magick 2.8.5 2024-09-20 [2] CRAN (R 4.5.0)
## magrittr 2.0.3 2022-03-30 [2] CRAN (R 4.5.0)
## Matrix 1.7-1 2024-10-18 [3] CRAN (R 4.5.0)
## MatrixGenerics 1.19.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## matrixStats 1.4.1 2024-09-08 [2] CRAN (R 4.5.0)
## memoise 2.0.1 2021-11-26 [2] CRAN (R 4.5.0)
## munsell 0.5.1 2024-04-01 [2] CRAN (R 4.5.0)
## nnet 7.3-19 2023-05-03 [3] CRAN (R 4.5.0)
## org.Hs.eg.db * 3.20.0 2024-10-24 [2] Bioconductor
## OrganismDbi 1.49.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## pillar 1.9.0 2023-03-22 [2] CRAN (R 4.5.0)
## pkgconfig 2.0.3 2019-09-22 [2] CRAN (R 4.5.0)
## plyr 1.8.9 2023-10-02 [2] CRAN (R 4.5.0)
## png 0.1-8 2022-11-29 [2] CRAN (R 4.5.0)
## prettyunits 1.2.0 2023-09-24 [2] CRAN (R 4.5.0)
## progress 1.2.3 2023-12-06 [2] CRAN (R 4.5.0)
## ProtGenerics 1.39.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## purrr 1.0.2 2023-08-10 [2] CRAN (R 4.5.0)
## qvalue 2.39.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## R6 2.5.1 2021-08-19 [2] CRAN (R 4.5.0)
## rappdirs 0.3.3 2021-01-31 [2] CRAN (R 4.5.0)
## RBGL 1.83.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## RColorBrewer 1.1-3 2022-04-03 [2] CRAN (R 4.5.0)
## Rcpp 1.0.13-1 2024-11-02 [2] CRAN (R 4.5.0)
## RCurl 1.98-1.16 2024-07-11 [2] CRAN (R 4.5.0)
## RefManageR * 1.4.0 2022-09-30 [2] CRAN (R 4.5.0)
## reshape2 1.4.4 2020-04-09 [2] CRAN (R 4.5.0)
## restfulr 0.0.15 2022-06-16 [2] CRAN (R 4.5.0)
## rjson 0.2.23 2024-09-16 [2] CRAN (R 4.5.0)
## rlang 1.1.4 2024-06-04 [2] CRAN (R 4.5.0)
## rmarkdown 2.29 2024-11-04 [2] CRAN (R 4.5.0)
## rngtools 1.5.2 2021-09-20 [2] CRAN (R 4.5.0)
## rpart 4.1.23 2023-12-05 [3] CRAN (R 4.5.0)
## Rsamtools 2.23.1 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## RSQLite 2.3.9 2024-12-03 [2] CRAN (R 4.5.0)
## rstudioapi 0.17.1 2024-10-22 [2] CRAN (R 4.5.0)
## rtracklayer 1.67.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## S4Arrays 1.7.1 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## S4Vectors * 0.45.2 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## sass 0.4.9 2024-03-15 [2] CRAN (R 4.5.0)
## scales 1.3.0 2023-11-28 [2] CRAN (R 4.5.0)
## sessioninfo * 1.2.2 2021-12-06 [2] CRAN (R 4.5.0)
## SparseArray 1.7.2 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## statmod 1.5.0 2023-01-06 [2] CRAN (R 4.5.0)
## stringi 1.8.4 2024-05-06 [2] CRAN (R 4.5.0)
## stringr 1.5.1 2023-11-14 [2] CRAN (R 4.5.0)
## SummarizedExperiment 1.37.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## tibble 3.2.1 2023-03-20 [2] CRAN (R 4.5.0)
## tidyr 1.3.1 2024-01-24 [2] CRAN (R 4.5.0)
## tidyselect 1.2.1 2024-03-11 [2] CRAN (R 4.5.0)
## timechange 0.3.0 2024-01-18 [2] CRAN (R 4.5.0)
## tinytex 0.54 2024-11-01 [2] CRAN (R 4.5.0)
## TxDb.Hsapiens.UCSC.hg19.knownGene * 3.2.2 2024-10-25 [2] Bioconductor
## txdbmaker 1.3.1 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## UCSC.utils 1.3.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## utf8 1.2.4 2023-10-22 [2] CRAN (R 4.5.0)
## VariantAnnotation 1.53.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## vctrs 0.6.5 2023-12-01 [2] CRAN (R 4.5.0)
## withr 3.0.2 2024-10-28 [2] CRAN (R 4.5.0)
## xfun 0.49 2024-10-31 [2] CRAN (R 4.5.0)
## XML 3.99-0.17 2024-06-25 [2] CRAN (R 4.5.0)
## xml2 1.3.6 2023-12-04 [2] CRAN (R 4.5.0)
## XVector 0.47.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
## yaml 2.3.10 2024-07-26 [2] CRAN (R 4.5.0)
## zlibbioc 1.53.0 2024-12-13 [2] Bioconductor 3.21 (R 4.5.0)
##
## [1] /tmp/Rtmp18nEXT/Rinst6975f11e5d4bd
## [2] /home/biocbuild/bbs-3.21-bioc/R/site-library
## [3] /home/biocbuild/bbs-3.21-bioc/R/library
##
## ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
This vignette was generated using BiocStyle (Oleś, 2024) with knitr (Xie, 2014) and rmarkdown (Allaire, Xie, Dervieux et al., 2024) running behind the scenes.
Citations made with RefManageR (McLean, 2017).
[1] J. Allaire, Y. Xie, C. Dervieux, et al. rmarkdown: Dynamic Documents for R. R package version 2.29. 2024. URL: https://github.com/rstudio/rmarkdown.
[2] S. Arora, M. Morgan, M. Carlson, et al. GenomeInfoDb: Utilities for manipulating chromosome and other ‘seqname’ identifiers. 2017. DOI: 10.18129/B9.bioc.GenomeInfoDb.
[3] 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/.
[4] M. Carlson and B. P. Maintainer. TxDb.Hsapiens.UCSC.hg19.knownGene: Annotation package for TxDb object(s). R package version 3.2.2. 2015.
[5] L. Collado-Torres, A. E. Jaffe, and J. T. Leek. derfinderPlot: Plotting functions for derfinder. https://github.com/leekgroup/derfinderPlot - R package version 1.41.0. 2017. DOI: 10.18129/B9.bioc.derfinderPlot. URL: http://www.bioconductor.org/packages/derfinderPlot.
[6] L. Collado-Torres, A. Jaffe, and J. Leek. derfinderData: Processed BigWigs from BrainSpan for examples. R package version 2.25.0. 2024. DOI: 10.18129/B9.bioc.derfinderData. URL: https://bioconductor.org/packages/derfinderData.
[7] L. Collado-Torres, A. Nellore, A. C. Frazee, et al. “Flexible expressed region analysis for RNA-seq with derfinder”. In: Nucl. Acids Res. (2017). DOI: 10.1093/nar/gkw852. URL: http://nar.oxfordjournals.org/content/early/2016/09/29/nar.gkw852.
[8] A. E. Jaffe, P. Murakami, H. Lee, 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.
[9] A. E. Jaffe, P. Murakami, H. Lee, et al. “Bump hunting to identify differentially methylated regions in epigenetic epidemiology studies”. In: International Journal of Epidemiology (2012).
[10] M. Lawrence, W. Huber, H. Pagès, 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}.
[11] M. W. McLean. “RefManageR: Import and Manage BibTeX and BibLaTeX References in R”. In: The Journal of Open Source Software (2017). DOI: 10.21105/joss.00338.
[12] E. Neuwirth. RColorBrewer: ColorBrewer Palettes. R package version 1.1-3. 2022. DOI: 10.32614/CRAN.package.RColorBrewer. URL: https://CRAN.R-project.org/package=RColorBrewer.
[13] A. Oleś. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.35.0. 2024. DOI: 10.18129/B9.bioc.BiocStyle. URL: https://bioconductor.org/packages/BiocStyle.
[14] R Core Team. R: A Language and Environment for Statistical Computing. R Foundation for Statistical Computing. Vienna, Austria, 2024. URL: https://www.R-project.org/.
[15] H. Wickham. “Reshaping Data with the reshape Package”. In: Journal of Statistical Software 21.12 (2007), pp. 1–20. URL: http://www.jstatsoft.org/v21/i12/.
[16] H. Wickham. “The Split-Apply-Combine Strategy for Data Analysis”. In: Journal of Statistical Software 40.1 (2011), pp. 1–29. URL: https://www.jstatsoft.org/v40/i01/.
[17] H. Wickham. ggplot2: Elegant Graphics for Data Analysis. Springer-Verlag New York, 2016. ISBN: 978-3-319-24277-4. URL: https://ggplot2.tidyverse.org.
[18] 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.
[19] H. Wickham, W. Chang, R. Flight, et al. sessioninfo: R Session Information. R package version 1.2.2. 2021. DOI: 10.32614/CRAN.package.sessioninfo. URL: https://CRAN.R-project.org/package=sessioninfo.
[20] H. Wickham, T. Pedersen, and D. Seidel. scales: Scale Functions for Visualization. R package version 1.3.0. 2023. DOI: 10.32614/CRAN.package.scales. URL: https://CRAN.R-project.org/package=scales.
[21] 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.
[22] 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.
[23] T. Yin, M. Lawrence, and D. Cook. biovizBase: Basic graphic utilities for visualization of genomic data. R package version 1.55.0. 2024. DOI: 10.18129/B9.bioc.biovizBase. URL: https://bioconductor.org/packages/biovizBase.