vignettes/derfinderPlot.Rmd
derfinderPlot.Rmd
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.35.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, 2023) 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"))
## 2023-05-07 05:32:04.630076 fullCoverage: processing chromosome chr21
## 2023-05-07 05:32:04.649001 loadCoverage: finding chromosome lengths
## 2023-05-07 05:32:04.683889 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB103.bw
## 2023-05-07 05:32:04.931039 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB114.bw
## 2023-05-07 05:32:05.151544 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB123.bw
## 2023-05-07 05:32:05.319494 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB126.bw
## 2023-05-07 05:32:05.435475 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB130.bw
## 2023-05-07 05:32:05.583077 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB135.bw
## 2023-05-07 05:32:05.699146 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB136.bw
## 2023-05-07 05:32:05.827028 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB145.bw
## 2023-05-07 05:32:05.963995 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB149.bw
## 2023-05-07 05:32:06.120738 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB150.bw
## 2023-05-07 05:32:06.239938 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB154.bw
## 2023-05-07 05:32:06.393255 loadCoverage: loading BigWig file /__w/_temp/Library/derfinderData/extdata/A1C/HSB178.bw
## 2023-05-07 05:32:06.555567 loadCoverage: applying the cutoff to the merged data
## 2023-05-07 05:32:06.593495 filterData: originally there were 48129895 rows, now there are 48129895 rows. Meaning that 0 percent was filtered.
## user system elapsed
## 2.009 0.028 2.037
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)
## 2023-05-07 05:32:06.898357 sampleDepth: Calculating sample quantiles
## 2023-05-07 05:32:07.046805 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)
## 2023-05-07 05:32:07.43633 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)
))
## 2023-05-07 05:32:08.493556 analyzeChr: Pre-processing the coverage data
## 2023-05-07 05:32:10.323049 analyzeChr: Calculating statistics
## 2023-05-07 05:32:10.326123 calculateStats: calculating the F-statistics
## 2023-05-07 05:32:10.524869 analyzeChr: Calculating pvalues
## 2023-05-07 05:32:10.525218 analyzeChr: Using the following theoretical cutoff for the F-statistics 5.31765507157871
## 2023-05-07 05:32:10.526452 calculatePvalues: identifying data segments
## 2023-05-07 05:32:10.534732 findRegions: segmenting information
## 2023-05-07 05:32:11.61813 findRegions: identifying candidate regions
## 2023-05-07 05:32:11.681973 findRegions: identifying region clusters
## 2023-05-07 05:32:11.849668 calculatePvalues: calculating F-statistics for permutation 1 and seed 20140924
## 2023-05-07 05:32:12.005713 findRegions: segmenting information
## 2023-05-07 05:32:12.030672 findRegions: identifying candidate regions
## 2023-05-07 05:32:12.097174 calculatePvalues: calculating F-statistics for permutation 2 and seed 20140925
## 2023-05-07 05:32:12.26795 findRegions: segmenting information
## 2023-05-07 05:32:12.293012 findRegions: identifying candidate regions
## 2023-05-07 05:32:12.34187 calculatePvalues: calculating F-statistics for permutation 3 and seed 20140926
## 2023-05-07 05:32:12.49784 findRegions: segmenting information
## 2023-05-07 05:32:12.522959 findRegions: identifying candidate regions
## 2023-05-07 05:32:12.580178 calculatePvalues: calculating F-statistics for permutation 4 and seed 20140927
## 2023-05-07 05:32:12.73575 findRegions: segmenting information
## 2023-05-07 05:32:12.760945 findRegions: identifying candidate regions
## 2023-05-07 05:32:12.810367 calculatePvalues: calculating F-statistics for permutation 5 and seed 20140928
## 2023-05-07 05:32:12.979759 findRegions: segmenting information
## 2023-05-07 05:32:13.005017 findRegions: identifying candidate regions
## 2023-05-07 05:32:13.053943 calculatePvalues: calculating F-statistics for permutation 6 and seed 20140929
## 2023-05-07 05:32:13.223524 findRegions: segmenting information
## 2023-05-07 05:32:13.260629 findRegions: identifying candidate regions
## 2023-05-07 05:32:13.310692 calculatePvalues: calculating F-statistics for permutation 7 and seed 20140930
## 2023-05-07 05:32:13.466949 findRegions: segmenting information
## 2023-05-07 05:32:13.491972 findRegions: identifying candidate regions
## 2023-05-07 05:32:13.540413 calculatePvalues: calculating F-statistics for permutation 8 and seed 20140931
## 2023-05-07 05:32:13.71283 findRegions: segmenting information
## 2023-05-07 05:32:13.737769 findRegions: identifying candidate regions
## 2023-05-07 05:32:13.786609 calculatePvalues: calculating F-statistics for permutation 9 and seed 20140932
## 2023-05-07 05:32:13.942104 findRegions: segmenting information
## 2023-05-07 05:32:13.967305 findRegions: identifying candidate regions
## 2023-05-07 05:32:14.025037 calculatePvalues: calculating F-statistics for permutation 10 and seed 20140933
## 2023-05-07 05:32:14.178943 findRegions: segmenting information
## 2023-05-07 05:32:14.204043 findRegions: identifying candidate regions
## 2023-05-07 05:32:14.25239 calculatePvalues: calculating F-statistics for permutation 11 and seed 20140934
## 2023-05-07 05:32:14.417238 findRegions: segmenting information
## 2023-05-07 05:32:14.442406 findRegions: identifying candidate regions
## 2023-05-07 05:32:14.491484 calculatePvalues: calculating F-statistics for permutation 12 and seed 20140935
## 2023-05-07 05:32:14.655616 findRegions: segmenting information
## 2023-05-07 05:32:14.680612 findRegions: identifying candidate regions
## 2023-05-07 05:32:14.728853 calculatePvalues: calculating F-statistics for permutation 13 and seed 20140936
## 2023-05-07 05:32:14.89555 findRegions: segmenting information
## 2023-05-07 05:32:14.920547 findRegions: identifying candidate regions
## 2023-05-07 05:32:14.968822 calculatePvalues: calculating F-statistics for permutation 14 and seed 20140937
## 2023-05-07 05:32:15.145993 findRegions: segmenting information
## 2023-05-07 05:32:15.171136 findRegions: identifying candidate regions
## 2023-05-07 05:32:15.220175 calculatePvalues: calculating F-statistics for permutation 15 and seed 20140938
## 2023-05-07 05:32:15.37665 findRegions: segmenting information
## 2023-05-07 05:32:15.401826 findRegions: identifying candidate regions
## 2023-05-07 05:32:15.450926 calculatePvalues: calculating F-statistics for permutation 16 and seed 20140939
## 2023-05-07 05:32:15.621174 findRegions: segmenting information
## 2023-05-07 05:32:15.646229 findRegions: identifying candidate regions
## 2023-05-07 05:32:15.695088 calculatePvalues: calculating F-statistics for permutation 17 and seed 20140940
## 2023-05-07 05:32:15.867494 findRegions: segmenting information
## 2023-05-07 05:32:15.892528 findRegions: identifying candidate regions
## 2023-05-07 05:32:15.941264 calculatePvalues: calculating F-statistics for permutation 18 and seed 20140941
## 2023-05-07 05:32:16.113336 findRegions: segmenting information
## 2023-05-07 05:32:16.138349 findRegions: identifying candidate regions
## 2023-05-07 05:32:16.187688 calculatePvalues: calculating F-statistics for permutation 19 and seed 20140942
## 2023-05-07 05:32:16.361938 findRegions: segmenting information
## 2023-05-07 05:32:16.387105 findRegions: identifying candidate regions
## 2023-05-07 05:32:16.435997 calculatePvalues: calculating F-statistics for permutation 20 and seed 20140943
## 2023-05-07 05:32:16.605025 findRegions: segmenting information
## 2023-05-07 05:32:16.630204 findRegions: identifying candidate regions
## 2023-05-07 05:32:16.703005 calculatePvalues: calculating the p-values
## 2023-05-07 05:32:16.776075 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
## 58.483 2.352 57.359
## 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 @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")
## 2023-05-07 05:33:06.041294 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"
)
## 2023-05-07 05:33:07.773298 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)
## 2023-05-07 05:33:08.577075 annotateRegions: counting
## 2023-05-07 05:33:08.662564 annotateRegions: annotating
regionCov <- getRegionCoverage(fullCov, regions)
## 2023-05-07 05:33:08.789818 getRegionCoverage: processing chr21
## 2023-05-07 05:33:08.847713 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: Transformation introduced infinite values in continuous y-axis
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).
## 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] "2023-05-07 05:33:34 UTC"
Wallclock time spent generating the vignette.
## Time difference of 1.794 mins
R
session information.
## ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
## setting value
## version R version 4.3.0 (2023-04-21)
## os Ubuntu 22.04.2 LTS
## system x86_64, linux-gnu
## ui X11
## language en
## collate en_US.UTF-8
## ctype en_US.UTF-8
## tz UTC
## date 2023-05-07
## pandoc 2.19.2 @ /usr/local/bin/ (via rmarkdown)
##
## ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
## package * version date (UTC) lib source
## AnnotationDbi * 1.62.1 2023-05-02 [1] Bioconductor
## AnnotationFilter 1.24.0 2023-04-25 [1] Bioconductor
## backports 1.4.1 2021-12-13 [1] CRAN (R 4.3.0)
## base64enc 0.1-3 2015-07-28 [2] RSPM (R 4.3.0)
## bibtex 0.5.1 2023-01-26 [1] RSPM (R 4.3.0)
## Biobase * 2.60.0 2023-04-25 [1] Bioconductor
## BiocFileCache 2.8.0 2023-04-25 [1] Bioconductor
## BiocGenerics * 0.46.0 2023-04-25 [1] Bioconductor
## BiocIO 1.10.0 2023-04-25 [1] Bioconductor
## BiocManager 1.30.20 2023-02-24 [2] CRAN (R 4.3.0)
## BiocParallel 1.34.1 2023-05-05 [1] Bioconductor
## BiocStyle * 2.28.0 2023-04-25 [1] Bioconductor
## biomaRt 2.56.0 2023-04-25 [1] Bioconductor
## Biostrings 2.68.0 2023-04-25 [1] Bioconductor
## biovizBase 1.48.0 2023-04-25 [1] Bioconductor
## bit 4.0.5 2022-11-15 [1] CRAN (R 4.3.0)
## bit64 4.0.5 2020-08-30 [1] CRAN (R 4.3.0)
## bitops 1.0-7 2021-04-24 [1] CRAN (R 4.3.0)
## blob 1.2.4 2023-03-17 [1] RSPM (R 4.3.0)
## bookdown 0.33 2023-03-06 [1] RSPM (R 4.3.0)
## BSgenome 1.68.0 2023-04-25 [1] Bioconductor
## bslib 0.4.2 2022-12-16 [2] RSPM (R 4.3.0)
## bumphunter * 1.42.0 2023-04-25 [1] Bioconductor
## cachem 1.0.8 2023-05-01 [2] RSPM (R 4.3.0)
## checkmate 2.2.0 2023-04-27 [1] RSPM (R 4.3.0)
## cli 3.6.1 2023-03-23 [2] RSPM (R 4.3.0)
## cluster 2.1.4 2022-08-22 [3] CRAN (R 4.3.0)
## codetools 0.2-19 2023-02-01 [3] CRAN (R 4.3.0)
## colorspace 2.1-0 2023-01-23 [1] RSPM (R 4.3.0)
## crayon 1.5.2 2022-09-29 [2] RSPM (R 4.3.0)
## curl 5.0.0 2023-01-12 [2] RSPM (R 4.3.0)
## data.table 1.14.8 2023-02-17 [1] RSPM (R 4.3.0)
## DBI 1.1.3 2022-06-18 [1] CRAN (R 4.3.0)
## dbplyr 2.3.2 2023-03-21 [1] RSPM (R 4.3.0)
## DelayedArray 0.26.2 2023-05-05 [1] Bioconductor
## derfinder * 1.34.0 2023-04-25 [1] Bioconductor
## derfinderData * 2.18.0 2023-04-27 [1] Bioconductor
## derfinderHelper 1.34.0 2023-04-25 [1] Bioconductor
## derfinderPlot * 1.35.0 2023-05-07 [1] Bioconductor
## desc 1.4.2 2022-09-08 [2] RSPM (R 4.3.0)
## dichromat 2.0-0.1 2022-05-02 [1] CRAN (R 4.3.0)
## digest 0.6.31 2022-12-11 [2] RSPM (R 4.3.0)
## doRNG 1.8.6 2023-01-16 [1] RSPM (R 4.3.0)
## dplyr 1.1.2 2023-04-20 [1] RSPM (R 4.3.0)
## ensembldb 2.24.0 2023-04-25 [1] Bioconductor
## evaluate 0.20 2023-01-17 [2] RSPM (R 4.3.0)
## fansi 1.0.4 2023-01-22 [2] RSPM (R 4.3.0)
## farver 2.1.1 2022-07-06 [1] CRAN (R 4.3.0)
## fastmap 1.1.1 2023-02-24 [2] RSPM (R 4.3.0)
## filelock 1.0.2 2018-10-05 [1] CRAN (R 4.3.0)
## foreach * 1.5.2 2022-02-02 [1] CRAN (R 4.3.0)
## foreign 0.8-84 2022-12-06 [3] CRAN (R 4.3.0)
## Formula 1.2-5 2023-02-24 [1] RSPM (R 4.3.0)
## fs 1.6.2 2023-04-25 [2] RSPM (R 4.3.0)
## generics 0.1.3 2022-07-05 [1] CRAN (R 4.3.0)
## GenomeInfoDb * 1.36.0 2023-04-25 [1] Bioconductor
## GenomeInfoDbData 1.2.10 2023-05-07 [1] Bioconductor
## GenomicAlignments 1.36.0 2023-04-25 [1] Bioconductor
## GenomicFeatures * 1.52.0 2023-04-25 [1] Bioconductor
## GenomicFiles 1.36.0 2023-04-25 [1] Bioconductor
## GenomicRanges * 1.52.0 2023-04-25 [1] Bioconductor
## GGally 2.1.2 2021-06-21 [1] CRAN (R 4.3.0)
## ggbio 1.48.0 2023-04-25 [1] Bioconductor
## ggplot2 3.4.2 2023-04-03 [1] RSPM (R 4.3.0)
## glue 1.6.2 2022-02-24 [2] RSPM (R 4.3.0)
## graph 1.78.0 2023-04-25 [1] Bioconductor
## gridExtra 2.3 2017-09-09 [1] CRAN (R 4.3.0)
## gtable 0.3.3 2023-03-21 [1] RSPM (R 4.3.0)
## highr 0.10 2022-12-22 [2] RSPM (R 4.3.0)
## Hmisc 5.0-1 2023-03-08 [1] RSPM (R 4.3.0)
## hms 1.1.3 2023-03-21 [1] RSPM (R 4.3.0)
## htmlTable 2.4.1 2022-07-07 [1] CRAN (R 4.3.0)
## htmltools 0.5.5 2023-03-23 [2] RSPM (R 4.3.0)
## htmlwidgets 1.6.2 2023-03-17 [2] RSPM (R 4.3.0)
## httr 1.4.5 2023-02-24 [2] RSPM (R 4.3.0)
## IRanges * 2.34.0 2023-04-25 [1] Bioconductor
## iterators * 1.0.14 2022-02-05 [1] CRAN (R 4.3.0)
## jquerylib 0.1.4 2021-04-26 [2] RSPM (R 4.3.0)
## jsonlite 1.8.4 2022-12-06 [2] RSPM (R 4.3.0)
## KEGGREST 1.40.0 2023-04-25 [1] Bioconductor
## knitr * 1.42 2023-01-25 [2] RSPM (R 4.3.0)
## labeling 0.4.2 2020-10-20 [1] CRAN (R 4.3.0)
## lattice 0.21-8 2023-04-05 [3] CRAN (R 4.3.0)
## lazyeval 0.2.2 2019-03-15 [1] CRAN (R 4.3.0)
## lifecycle 1.0.3 2022-10-07 [2] RSPM (R 4.3.0)
## limma 3.56.0 2023-04-25 [1] Bioconductor
## locfit * 1.5-9.7 2023-01-02 [1] RSPM (R 4.3.0)
## lubridate 1.9.2 2023-02-10 [1] RSPM (R 4.3.0)
## magrittr 2.0.3 2022-03-30 [2] RSPM (R 4.3.0)
## Matrix 1.5-4 2023-04-04 [3] CRAN (R 4.3.0)
## MatrixGenerics 1.12.0 2023-04-25 [1] Bioconductor
## matrixStats 0.63.0 2022-11-18 [1] CRAN (R 4.3.0)
## memoise 2.0.1 2021-11-26 [2] RSPM (R 4.3.0)
## munsell 0.5.0 2018-06-12 [1] CRAN (R 4.3.0)
## nnet 7.3-19 2023-05-03 [3] RSPM (R 4.3.0)
## org.Hs.eg.db * 3.17.0 2023-05-07 [1] Bioconductor
## OrganismDbi 1.42.0 2023-04-25 [1] Bioconductor
## pillar 1.9.0 2023-03-22 [2] RSPM (R 4.3.0)
## pkgconfig 2.0.3 2019-09-22 [2] RSPM (R 4.3.0)
## pkgdown 2.0.7 2022-12-14 [2] RSPM (R 4.3.0)
## plyr 1.8.8 2022-11-11 [1] CRAN (R 4.3.0)
## png 0.1-8 2022-11-29 [1] CRAN (R 4.3.0)
## prettyunits 1.1.1 2020-01-24 [2] RSPM (R 4.3.0)
## progress 1.2.2 2019-05-16 [1] CRAN (R 4.3.0)
## ProtGenerics 1.32.0 2023-04-25 [1] Bioconductor
## purrr 1.0.1 2023-01-10 [2] RSPM (R 4.3.0)
## qvalue 2.32.0 2023-04-25 [1] Bioconductor
## R6 2.5.1 2021-08-19 [2] RSPM (R 4.3.0)
## ragg 1.2.5 2023-01-12 [2] RSPM (R 4.3.0)
## rappdirs 0.3.3 2021-01-31 [2] RSPM (R 4.3.0)
## RBGL 1.76.0 2023-04-25 [1] Bioconductor
## RColorBrewer 1.1-3 2022-04-03 [1] CRAN (R 4.3.0)
## Rcpp 1.0.10 2023-01-22 [2] RSPM (R 4.3.0)
## RCurl 1.98-1.12 2023-03-27 [1] RSPM (R 4.3.0)
## RefManageR * 1.4.0 2022-09-30 [1] CRAN (R 4.3.0)
## reshape 0.8.9 2022-04-12 [1] CRAN (R 4.3.0)
## reshape2 1.4.4 2020-04-09 [1] CRAN (R 4.3.0)
## restfulr 0.0.15 2022-06-16 [1] CRAN (R 4.3.0)
## rjson 0.2.21 2022-01-09 [1] CRAN (R 4.3.0)
## rlang 1.1.1 2023-04-28 [2] RSPM (R 4.3.0)
## rmarkdown 2.21 2023-03-26 [2] RSPM (R 4.3.0)
## rngtools 1.5.2 2021-09-20 [1] CRAN (R 4.3.0)
## rpart 4.1.19 2022-10-21 [3] CRAN (R 4.3.0)
## rprojroot 2.0.3 2022-04-02 [2] RSPM (R 4.3.0)
## Rsamtools 2.16.0 2023-04-25 [1] Bioconductor
## RSQLite 2.3.1 2023-04-03 [1] RSPM (R 4.3.0)
## rstudioapi 0.14 2022-08-22 [2] RSPM (R 4.3.0)
## rtracklayer 1.60.0 2023-04-25 [1] Bioconductor
## S4Arrays 1.0.1 2023-05-01 [1] Bioconductor
## S4Vectors * 0.38.1 2023-05-02 [1] Bioconductor
## sass 0.4.6 2023-05-03 [2] RSPM (R 4.3.0)
## scales 1.2.1 2022-08-20 [1] CRAN (R 4.3.0)
## sessioninfo * 1.2.2 2021-12-06 [2] RSPM (R 4.3.0)
## stringi 1.7.12 2023-01-11 [2] RSPM (R 4.3.0)
## stringr 1.5.0 2022-12-02 [2] RSPM (R 4.3.0)
## SummarizedExperiment 1.30.1 2023-05-01 [1] Bioconductor
## systemfonts 1.0.4 2022-02-11 [2] RSPM (R 4.3.0)
## textshaping 0.3.6 2021-10-13 [2] RSPM (R 4.3.0)
## tibble 3.2.1 2023-03-20 [2] RSPM (R 4.3.0)
## tidyselect 1.2.0 2022-10-10 [1] CRAN (R 4.3.0)
## timechange 0.2.0 2023-01-11 [1] RSPM (R 4.3.0)
## TxDb.Hsapiens.UCSC.hg19.knownGene * 3.2.2 2022-12-06 [1] Bioconductor
## utf8 1.2.3 2023-01-31 [2] RSPM (R 4.3.0)
## VariantAnnotation 1.46.0 2023-04-25 [1] Bioconductor
## vctrs 0.6.2 2023-04-19 [2] RSPM (R 4.3.0)
## withr 2.5.0 2022-03-03 [2] RSPM (R 4.3.0)
## xfun 0.39 2023-04-20 [2] RSPM (R 4.3.0)
## XML 3.99-0.14 2023-03-19 [1] RSPM (R 4.3.0)
## xml2 1.3.4 2023-04-27 [2] RSPM (R 4.3.0)
## XVector 0.40.0 2023-04-25 [1] Bioconductor
## yaml 2.3.7 2023-01-23 [2] RSPM (R 4.3.0)
## zlibbioc 1.46.0 2023-04-25 [1] Bioconductor
##
## [1] /__w/_temp/Library
## [2] /usr/local/lib/R/site-library
## [3] /usr/local/lib/R/library
##
## ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
This vignette was generated using BiocStyle (Oleś, 2023) with knitr (Xie, 2014) and rmarkdown (Allaire, Xie, Dervieux et al., 2023) 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.21. 2023. 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.35.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.18.0. 2023. 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. URL: https://CRAN.R-project.org/package=RColorBrewer.
[13] A. Oleś. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.28.0. 2023. 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, 2023. URL: https://www.R-project.org/.
[15] H. Wickham. ggplot2: Elegant Graphics for Data Analysis. Springer-Verlag New York, 2016. ISBN: 978-3-319-24277-4. URL: https://ggplot2.tidyverse.org.
[16] 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/.
[17] 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.
[18] 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/.
[19] H. Wickham, W. Chang, R. Flight, et al. sessioninfo: R Session Information. https://github.com/r-lib/sessioninfo#readme, https://r-lib.github.io/sessioninfo/. 2021.
[20] H. Wickham and D. Seidel. scales: Scale Functions for Visualization. R package version 1.2.1. 2022. 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.48.0. 2023. DOI: 10.18129/B9.bioc.biovizBase. URL: https://bioconductor.org/packages/biovizBase.