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:
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!
## 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-10-30 05:45:47.290979 fullCoverage: processing chromosome chr21
## 2024-10-30 05:45:47.302928 loadCoverage: finding chromosome lengths
## 2024-10-30 05:45:47.323242 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB103.bw
## 2024-10-30 05:45:47.46993 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB114.bw
## 2024-10-30 05:45:47.606987 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB123.bw
## 2024-10-30 05:45:47.702125 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB126.bw
## 2024-10-30 05:45:47.768818 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB130.bw
## 2024-10-30 05:45:47.858552 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB135.bw
## 2024-10-30 05:45:47.922148 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB136.bw
## 2024-10-30 05:45:47.989344 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB145.bw
## 2024-10-30 05:45:48.076263 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB149.bw
## 2024-10-30 05:45:48.156914 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB150.bw
## 2024-10-30 05:45:48.232102 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB154.bw
## 2024-10-30 05:45:48.329043 loadCoverage: loading BigWig file /github/workspace/pkglib/derfinderData/extdata/A1C/HSB178.bw
## 2024-10-30 05:45:49.223264 loadCoverage: applying the cutoff to the merged data
## 2024-10-30 05:45:49.243991 filterData: originally there were 48129895 rows, now there are 48129895 rows. Meaning that 0 percent was filtered.
## user system elapsed
## 1.970 0.028 1.998
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.
## 2024-10-30 05:45:49.321298 sampleDepth: Calculating sample quantiles
## 2024-10-30 05:45:49.401524 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.
## 2024-10-30 05:45:49.641668 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-10-30 05:45:50.571039 analyzeChr: Pre-processing the coverage data
## 2024-10-30 05:45:51.866829 analyzeChr: Calculating statistics
## 2024-10-30 05:45:51.869407 calculateStats: calculating the F-statistics
## 2024-10-30 05:45:51.998141 analyzeChr: Calculating pvalues
## 2024-10-30 05:45:51.998847 analyzeChr: Using the following theoretical cutoff for the F-statistics 5.31765507157871
## 2024-10-30 05:45:52.000106 calculatePvalues: identifying data segments
## 2024-10-30 05:45:52.007249 findRegions: segmenting information
## 2024-10-30 05:45:52.038273 findRegions: identifying candidate regions
## 2024-10-30 05:45:52.089237 findRegions: identifying region clusters
## 2024-10-30 05:45:52.186609 calculatePvalues: calculating F-statistics for permutation 1 and seed 20140924
## 2024-10-30 05:45:52.299255 findRegions: segmenting information
## 2024-10-30 05:45:52.326331 findRegions: identifying candidate regions
## 2024-10-30 05:45:52.384212 calculatePvalues: calculating F-statistics for permutation 2 and seed 20140925
## 2024-10-30 05:45:52.502136 findRegions: segmenting information
## 2024-10-30 05:45:52.528964 findRegions: identifying candidate regions
## 2024-10-30 05:45:52.570482 calculatePvalues: calculating F-statistics for permutation 3 and seed 20140926
## 2024-10-30 05:45:52.697095 findRegions: segmenting information
## 2024-10-30 05:45:52.724653 findRegions: identifying candidate regions
## 2024-10-30 05:45:52.76721 calculatePvalues: calculating F-statistics for permutation 4 and seed 20140927
## 2024-10-30 05:45:52.891298 findRegions: segmenting information
## 2024-10-30 05:45:52.918767 findRegions: identifying candidate regions
## 2024-10-30 05:45:52.961642 calculatePvalues: calculating F-statistics for permutation 5 and seed 20140928
## 2024-10-30 05:45:53.082798 findRegions: segmenting information
## 2024-10-30 05:45:53.109737 findRegions: identifying candidate regions
## 2024-10-30 05:45:53.151353 calculatePvalues: calculating F-statistics for permutation 6 and seed 20140929
## 2024-10-30 05:45:53.278686 findRegions: segmenting information
## 2024-10-30 05:45:53.305925 findRegions: identifying candidate regions
## 2024-10-30 05:45:53.348425 calculatePvalues: calculating F-statistics for permutation 7 and seed 20140930
## 2024-10-30 05:45:53.472527 findRegions: segmenting information
## 2024-10-30 05:45:53.499543 findRegions: identifying candidate regions
## 2024-10-30 05:45:53.541557 calculatePvalues: calculating F-statistics for permutation 8 and seed 20140931
## 2024-10-30 05:45:53.670149 findRegions: segmenting information
## 2024-10-30 05:45:53.697205 findRegions: identifying candidate regions
## 2024-10-30 05:45:53.739056 calculatePvalues: calculating F-statistics for permutation 9 and seed 20140932
## 2024-10-30 05:45:53.866207 findRegions: segmenting information
## 2024-10-30 05:45:53.894046 findRegions: identifying candidate regions
## 2024-10-30 05:45:53.937964 calculatePvalues: calculating F-statistics for permutation 10 and seed 20140933
## 2024-10-30 05:45:54.061508 findRegions: segmenting information
## 2024-10-30 05:45:54.088749 findRegions: identifying candidate regions
## 2024-10-30 05:45:54.130578 calculatePvalues: calculating F-statistics for permutation 11 and seed 20140934
## 2024-10-30 05:45:54.252267 findRegions: segmenting information
## 2024-10-30 05:45:54.279226 findRegions: identifying candidate regions
## 2024-10-30 05:45:54.321364 calculatePvalues: calculating F-statistics for permutation 12 and seed 20140935
## 2024-10-30 05:45:54.448265 findRegions: segmenting information
## 2024-10-30 05:45:54.475616 findRegions: identifying candidate regions
## 2024-10-30 05:45:54.519212 calculatePvalues: calculating F-statistics for permutation 13 and seed 20140936
## 2024-10-30 05:45:54.64929 findRegions: segmenting information
## 2024-10-30 05:45:54.676754 findRegions: identifying candidate regions
## 2024-10-30 05:45:54.719918 calculatePvalues: calculating F-statistics for permutation 14 and seed 20140937
## 2024-10-30 05:45:55.647879 findRegions: segmenting information
## 2024-10-30 05:45:55.674853 findRegions: identifying candidate regions
## 2024-10-30 05:45:55.71657 calculatePvalues: calculating F-statistics for permutation 15 and seed 20140938
## 2024-10-30 05:45:55.837188 findRegions: segmenting information
## 2024-10-30 05:45:55.866058 findRegions: identifying candidate regions
## 2024-10-30 05:45:55.907748 calculatePvalues: calculating F-statistics for permutation 16 and seed 20140939
## 2024-10-30 05:45:56.029791 findRegions: segmenting information
## 2024-10-30 05:45:56.057525 findRegions: identifying candidate regions
## 2024-10-30 05:45:56.098849 calculatePvalues: calculating F-statistics for permutation 17 and seed 20140940
## 2024-10-30 05:45:56.21142 findRegions: segmenting information
## 2024-10-30 05:45:56.2381 findRegions: identifying candidate regions
## 2024-10-30 05:45:56.278898 calculatePvalues: calculating F-statistics for permutation 18 and seed 20140941
## 2024-10-30 05:45:56.39836 findRegions: segmenting information
## 2024-10-30 05:45:56.42607 findRegions: identifying candidate regions
## 2024-10-30 05:45:56.467728 calculatePvalues: calculating F-statistics for permutation 19 and seed 20140942
## 2024-10-30 05:45:56.591923 findRegions: segmenting information
## 2024-10-30 05:45:56.61913 findRegions: identifying candidate regions
## 2024-10-30 05:45:56.661023 calculatePvalues: calculating F-statistics for permutation 20 and seed 20140943
## 2024-10-30 05:45:56.77408 findRegions: segmenting information
## 2024-10-30 05:45:56.80163 findRegions: identifying candidate regions
## 2024-10-30 05:45:56.865969 calculatePvalues: calculating the p-values
## 2024-10-30 05:45:56.912657 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
## 42.932 7.161 42.032
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 @ref(fig:plotOverview) 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-10-30 05:46:32.721407 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 @ref(fig:plotOverview2) 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-10-30 05:46:34.133483 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-10-30 05:46:34.763665 annotateRegions: counting
## 2024-10-30 05:46:34.816212 annotateRegions: annotating
## 2024-10-30 05:46:34.895705 getRegionCoverage: processing chr21
## 2024-10-30 05:46:34.935381 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 @ref(fig:plotRegCov1),
@ref(fig:plotRegCov2), @ref(fig:plotRegCov3), @ref(fig:plotRegCov4),
@ref(fig:plotRegCov5), @ref(fig:plotRegCov6), @ref(fig:plotRegCov7),
@ref(fig:plotRegCov8), @ref(fig:plotRegCov9),
@ref(fig:plotRegCov10)).
## 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 @ref(fig:plotCluster) 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 @ref(fig:plotCluster2)) 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
@ref(fig:vennRegions).
## 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)
Date the vignette was generated.
## [1] "2024-10-30 05:46:49 UTC"
Wallclock time spent generating the vignette.
## Time difference of 1.247 mins
R
session information.
## ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
## setting value
## version R version 4.4.1 (2024-06-14)
## os Ubuntu 24.04.1 LTS
## system x86_64, linux-gnu
## ui X11
## language (EN)
## collate C
## ctype en_US.UTF-8
## tz Etc/UTC
## date 2024-10-30
## pandoc 3.2.1 @ /usr/local/bin/ (via rmarkdown)
##
## ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
## package * version date (UTC) lib source
## abind 1.4-8 2024-09-12 [2] RSPM (R 4.4.0)
## AnnotationDbi * 1.69.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## AnnotationFilter 1.31.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## backports 1.5.0 2024-05-23 [2] RSPM (R 4.4.0)
## base64enc 0.1-3 2015-07-28 [2] RSPM (R 4.4.0)
## bibtex 0.5.1 2023-01-26 [2] RSPM (R 4.4.0)
## Biobase * 2.67.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## BiocFileCache 2.15.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## BiocGenerics * 0.53.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## BiocIO 1.17.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## BiocManager 1.30.25 2024-08-28 [2] RSPM (R 4.4.0)
## BiocParallel 1.39.0 2024-10-23 [2] https://bioc.r-universe.dev (R 4.4.1)
## BiocStyle * 2.35.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## biomaRt 2.63.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## Biostrings 2.75.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## biovizBase 1.55.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## bit 4.5.0 2024-09-20 [2] RSPM (R 4.4.0)
## bit64 4.5.2 2024-09-22 [2] RSPM (R 4.4.0)
## bitops 1.0-9 2024-10-03 [2] RSPM (R 4.4.0)
## blob 1.2.4 2023-03-17 [2] RSPM (R 4.4.0)
## BSgenome 1.73.1 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## bslib 0.8.0 2024-07-29 [2] RSPM (R 4.4.0)
## buildtools 1.0.0 2024-10-28 [3] local (/pkg)
## bumphunter * 1.49.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
## cachem 1.1.0 2024-05-16 [2] RSPM (R 4.4.0)
## checkmate 2.3.2 2024-07-29 [2] RSPM (R 4.4.0)
## cli 3.6.3 2024-06-21 [2] RSPM (R 4.4.0)
## cluster 2.1.6 2023-12-01 [2] RSPM (R 4.4.0)
## codetools 0.2-20 2024-03-31 [2] RSPM (R 4.4.0)
## colorspace 2.1-1 2024-07-26 [2] RSPM (R 4.4.0)
## crayon 1.5.3 2024-06-20 [2] RSPM (R 4.4.0)
## curl 5.2.3 2024-09-20 [2] RSPM (R 4.4.0)
## data.table 1.16.2 2024-10-10 [2] RSPM (R 4.4.0)
## DBI 1.2.3 2024-06-02 [2] RSPM (R 4.4.0)
## dbplyr 2.5.0 2024-03-19 [2] RSPM (R 4.4.0)
## DelayedArray 0.31.14 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## derfinder * 1.39.0 2024-10-16 [2] https://bioc.r-universe.dev (R 4.4.1)
## derfinderData * 2.23.0 2024-05-02 [2] Bioconductor 3.20 (R 4.4.1)
## derfinderHelper 1.39.0 2024-10-17 [2] https://bioc.r-universe.dev (R 4.4.1)
## derfinderPlot * 1.41.0 2024-10-30 [1] https://bioc.r-universe.dev (R 4.4.1)
## dichromat 2.0-0.1 2022-05-02 [2] RSPM (R 4.4.0)
## digest 0.6.37 2024-08-19 [2] RSPM (R 4.4.0)
## doRNG 1.8.6 2023-01-16 [2] RSPM (R 4.4.0)
## dplyr 1.1.4 2023-11-17 [2] RSPM (R 4.4.0)
## ensembldb 2.29.1 2024-10-21 [2] https://bioc.r-universe.dev (R 4.4.1)
## evaluate 1.0.1 2024-10-10 [2] RSPM (R 4.4.0)
## fansi 1.0.6 2023-12-08 [2] RSPM (R 4.4.0)
## farver 2.1.2 2024-05-13 [2] RSPM (R 4.4.0)
## fastmap 1.2.0 2024-05-15 [2] RSPM (R 4.4.0)
## filelock 1.0.3 2023-12-11 [2] RSPM (R 4.4.0)
## foreach * 1.5.2 2022-02-02 [2] RSPM (R 4.4.0)
## foreign 0.8-87 2024-06-26 [2] RSPM (R 4.4.0)
## Formula 1.2-5 2023-02-24 [2] RSPM (R 4.4.0)
## generics 0.1.3 2022-07-05 [2] RSPM (R 4.4.0)
## GenomeInfoDb * 1.41.2 2024-10-02 [2] https://bioc.r-universe.dev (R 4.4.1)
## GenomeInfoDbData 1.2.13 2024-10-30 [2] Bioconductor
## GenomicAlignments 1.41.0 2024-10-28 [2] https://bioc.r-universe.dev (R 4.4.1)
## GenomicFeatures * 1.57.1 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## GenomicFiles 1.41.0 2024-10-28 [2] https://bioc.r-universe.dev (R 4.4.1)
## GenomicRanges * 1.57.2 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## GGally 2.2.1 2024-02-14 [2] RSPM (R 4.4.0)
## ggbio 1.53.0 2024-10-28 [2] https://bioc.r-universe.dev (R 4.4.1)
## ggplot2 3.5.1 2024-04-23 [2] RSPM (R 4.4.0)
## ggstats 0.7.0 2024-09-22 [2] RSPM (R 4.4.0)
## glue 1.8.0 2024-09-30 [2] RSPM (R 4.4.0)
## graph 1.83.0 2024-10-01 [2] https://bioc.r-universe.dev (R 4.4.1)
## gridExtra 2.3 2017-09-09 [2] RSPM (R 4.4.0)
## gtable 0.3.6 2024-10-25 [2] RSPM (R 4.4.0)
## highr 0.11 2024-05-26 [2] RSPM (R 4.4.0)
## Hmisc 5.2-0 2024-10-28 [2] RSPM (R 4.4.0)
## hms 1.1.3 2023-03-21 [2] RSPM (R 4.4.0)
## htmlTable 2.4.3 2024-07-21 [2] RSPM (R 4.4.0)
## htmltools 0.5.8.1 2024-04-04 [2] RSPM (R 4.4.0)
## htmlwidgets 1.6.4 2023-12-06 [2] RSPM (R 4.4.0)
## httr 1.4.7 2023-08-15 [2] RSPM (R 4.4.0)
## httr2 1.0.5 2024-09-26 [2] RSPM (R 4.4.0)
## IRanges * 2.39.2 2024-10-25 [2] https://bioc.r-universe.dev (R 4.4.1)
## iterators * 1.0.14 2022-02-05 [2] RSPM (R 4.4.0)
## jquerylib 0.1.4 2021-04-26 [2] RSPM (R 4.4.0)
## jsonlite 1.8.9 2024-09-20 [2] RSPM (R 4.4.0)
## KEGGREST 1.45.1 2024-10-16 [2] https://bioc.r-universe.dev (R 4.4.1)
## knitr * 1.48 2024-07-07 [2] RSPM (R 4.4.0)
## labeling 0.4.3 2023-08-29 [2] RSPM (R 4.4.0)
## lattice 0.22-6 2024-03-20 [2] RSPM (R 4.4.0)
## lazyeval 0.2.2 2019-03-15 [2] RSPM (R 4.4.0)
## lifecycle 1.0.4 2023-11-07 [2] RSPM (R 4.4.0)
## limma 3.61.12 2024-10-01 [2] https://bioc.r-universe.dev (R 4.4.1)
## locfit * 1.5-9.10 2024-06-24 [2] RSPM (R 4.4.0)
## lubridate 1.9.3 2023-09-27 [2] RSPM (R 4.4.0)
## magrittr 2.0.3 2022-03-30 [2] RSPM (R 4.4.0)
## maketools 1.3.1 2024-10-28 [3] Github (jeroen/maketools@d46f92c)
## Matrix 1.7-1 2024-10-18 [2] RSPM (R 4.4.0)
## MatrixGenerics 1.17.1 2024-10-23 [2] https://bioc.r-universe.dev (R 4.4.1)
## matrixStats 1.4.1 2024-09-08 [2] RSPM (R 4.4.0)
## memoise 2.0.1 2021-11-26 [2] RSPM (R 4.4.0)
## munsell 0.5.1 2024-04-01 [2] RSPM (R 4.4.0)
## nnet 7.3-19 2023-05-03 [2] RSPM (R 4.4.0)
## org.Hs.eg.db * 3.20.0 2024-10-30 [2] Bioconductor
## OrganismDbi 1.47.0 2024-10-11 [2] https://bioc.r-universe.dev (R 4.4.1)
## pillar 1.9.0 2023-03-22 [2] RSPM (R 4.4.0)
## pkgconfig 2.0.3 2019-09-22 [2] RSPM (R 4.4.0)
## plyr 1.8.9 2023-10-02 [2] RSPM (R 4.4.0)
## png 0.1-8 2022-11-29 [2] RSPM (R 4.4.0)
## prettyunits 1.2.0 2023-09-24 [2] RSPM (R 4.4.0)
## progress 1.2.3 2023-12-06 [2] RSPM (R 4.4.0)
## ProtGenerics 1.37.1 2024-10-22 [2] https://bioc.r-universe.dev (R 4.4.1)
## purrr 1.0.2 2023-08-10 [2] RSPM (R 4.4.0)
## qvalue 2.37.0 2024-10-01 [2] https://bioc.r-universe.dev (R 4.4.1)
## R6 2.5.1 2021-08-19 [2] RSPM (R 4.4.0)
## rappdirs 0.3.3 2021-01-31 [2] RSPM (R 4.4.0)
## RBGL 1.81.0 2024-10-18 [2] https://bioc.r-universe.dev (R 4.4.1)
## RColorBrewer 1.1-3 2022-04-03 [2] RSPM (R 4.4.0)
## Rcpp 1.0.13 2024-07-17 [2] RSPM (R 4.4.0)
## RCurl 1.98-1.16 2024-07-11 [2] RSPM (R 4.4.0)
## RefManageR * 1.4.0 2022-09-30 [2] RSPM (R 4.4.0)
## reshape2 1.4.4 2020-04-09 [2] RSPM (R 4.4.0)
## restfulr 0.0.15 2022-06-16 [2] RSPM (R 4.4.1)
## rjson 0.2.23 2024-09-16 [2] RSPM (R 4.4.0)
## rlang 1.1.4 2024-06-04 [2] RSPM (R 4.4.0)
## rmarkdown 2.28 2024-08-17 [2] RSPM (R 4.4.0)
## rngtools 1.5.2 2021-09-20 [2] RSPM (R 4.4.0)
## rpart 4.1.23 2023-12-05 [2] RSPM (R 4.4.0)
## Rsamtools 2.21.2 2024-10-26 [2] https://bioc.r-universe.dev (R 4.4.1)
## RSQLite 2.3.7 2024-05-27 [2] RSPM (R 4.4.0)
## rstudioapi 0.17.1 2024-10-22 [2] RSPM (R 4.4.0)
## rtracklayer 1.65.0 2024-10-23 [2] https://bioc.r-universe.dev (R 4.4.1)
## S4Arrays 1.5.11 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## S4Vectors * 0.43.2 2024-10-17 [2] https://bioc.r-universe.dev (R 4.4.1)
## sass 0.4.9 2024-03-15 [2] RSPM (R 4.4.0)
## scales 1.3.0 2023-11-28 [2] RSPM (R 4.4.0)
## sessioninfo * 1.2.2 2021-12-06 [2] RSPM (R 4.4.0)
## SparseArray 1.5.45 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## statmod 1.5.0 2023-01-06 [2] RSPM (R 4.4.0)
## stringi 1.8.4 2024-05-06 [2] RSPM (R 4.4.0)
## stringr 1.5.1 2023-11-14 [2] RSPM (R 4.4.0)
## SummarizedExperiment 1.35.5 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## sys 3.4.3 2024-10-04 [2] RSPM (R 4.4.0)
## tibble 3.2.1 2023-03-20 [2] RSPM (R 4.4.0)
## tidyr 1.3.1 2024-01-24 [2] RSPM (R 4.4.0)
## tidyselect 1.2.1 2024-03-11 [2] RSPM (R 4.4.0)
## timechange 0.3.0 2024-01-18 [2] RSPM (R 4.4.0)
## TxDb.Hsapiens.UCSC.hg19.knownGene * 3.2.2 2024-10-30 [2] Bioconductor
## txdbmaker 1.1.2 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## UCSC.utils 1.1.0 2024-10-29 [2] https://bioc.r-universe.dev (R 4.4.1)
## utf8 1.2.4 2023-10-22 [2] RSPM (R 4.4.0)
## VariantAnnotation 1.51.2 2024-10-24 [2] https://bioc.r-universe.dev (R 4.4.1)
## vctrs 0.6.5 2023-12-01 [2] RSPM (R 4.4.0)
## withr 3.0.2 2024-10-28 [2] RSPM (R 4.4.0)
## xfun 0.48 2024-10-03 [2] RSPM (R 4.4.0)
## XML 3.99-0.17 2024-06-25 [2] RSPM (R 4.4.0)
## xml2 1.3.6 2023-12-04 [2] RSPM (R 4.4.0)
## XVector 0.45.0 2024-10-02 [2] https://bioc.r-universe.dev (R 4.4.1)
## yaml 2.3.10 2024-07-26 [2] RSPM (R 4.4.0)
## zlibbioc 1.51.2 2024-10-21 [2] Bioconductor 3.20 (R 4.4.1)
##
## [1] /tmp/RtmpzCrfqY/Rinst22ec68edb3bc
## [2] /github/workspace/pkglib
## [3] /usr/local/lib/R/site-library
## [4] /usr/lib/R/site-library
## [5] /usr/lib/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.28. 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.23.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.
[13] A. Oleś. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.35.0. 2024. URL: https://github.com/Bioconductor/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, https://r-lib.github.io/sessioninfo/. 2021. URL: https://github.com/r-lib/sessioninfo#readme.
[20] H. Wickham, T. Pedersen, and D. Seidel. scales: Scale Functions for Visualization. R package version 1.3.0, https://github.com/r-lib/scales. 2023. URL: https://scales.r-lib.org.
[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.