regionReport 1.30.0
The bumphunter package can be used for methylation analyses where you are interested in identifying differentially methylated regions. The vignette explains in greater detail the data set we are using in this example.
## Load bumphunter
library("bumphunter")
## Loading required package: S4Vectors
## Loading required package: stats4
## Loading required package: BiocGenerics
##
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:stats':
##
## IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
##
## Filter, Find, Map, Position, Reduce, anyDuplicated, append,
## as.data.frame, basename, cbind, colnames, dirname, do.call,
## duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
## lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
## pmin.int, rank, rbind, rownames, sapply, setdiff, sort, table,
## tapply, union, unique, unsplit, which.max, which.min
##
## Attaching package: 'S4Vectors'
## The following objects are masked from 'package:base':
##
## I, expand.grid, unname
## Loading required package: IRanges
## Loading required package: GenomeInfoDb
## Loading required package: GenomicRanges
## Loading required package: foreach
## Loading required package: iterators
## Loading required package: parallel
## Loading required package: locfit
## locfit 1.5-9.5 2022-03-01
## Create data from the vignette
pos <- list(
pos1 = seq(1, 1000, 35),
pos2 = seq(2001, 3000, 35),
pos3 = seq(1, 1000, 50)
)
chr <- rep(paste0("chr", c(1, 1, 2)), times = sapply(pos, length))
pos <- unlist(pos, use.names = FALSE)
## Find clusters
cl <- clusterMaker(chr, pos, maxGap = 300)
## Build simulated bumps
Indexes <- split(seq_along(cl), cl)
beta1 <- rep(0, length(pos))
for (i in seq(along = Indexes)) {
ind <- Indexes[[i]]
x <- pos[ind]
z <- scale(x, median(x), max(x) / 12)
beta1[ind] <- i * (-1)^(i + 1) * pmax(1 - abs(z)^3, 0)^3 ## multiply by i to vary size
}
## Build data
beta0 <- 3 * sin(2 * pi * pos / 720)
X <- cbind(rep(1, 20), rep(c(0, 1), each = 10))
set.seed(23852577)
error <- matrix(rnorm(20 * length(beta1), 0, 1), ncol = 20)
y <- t(X[, 1]) %x% beta0 + t(X[, 2]) %x% beta1 + error
## Perform bumphunting
tab <- bumphunter(y, X, chr, pos, cl, cutoff = .5)
## [bumphunterEngine] Using a single core (backend: doSEQ, version: 1.5.2).
## [bumphunterEngine] Computing coefficients.
## [bumphunterEngine] Finding regions.
## [bumphunterEngine] Found 15 bumps.
## Explore data
lapply(tab, head)
## $table
## chr start end value area cluster indexStart indexEnd L
## 10 chr1 2316 2631 -1.5814747 15.8147473 2 39 48 10
## 7 chr2 451 551 1.5891293 4.7673878 3 68 70 3
## 2 chr1 456 526 1.0678828 3.2036485 1 14 16 3
## 5 chr1 2176 2211 0.7841794 1.5683589 2 35 36 2
## 6 chr1 2841 2841 1.2010184 1.2010184 2 54 54 1
## 4 chr1 771 771 0.7780902 0.7780902 1 23 23 1
## clusterL
## 10 29
## 7 20
## 2 29
## 5 29
## 6 29
## 4 29
##
## $coef
## [,1]
## [1,] 0.60960932
## [2,] -0.09052769
## [3,] -0.21482638
## [4,] 0.13053755
## [5,] -0.21723642
## [6,] 0.39934961
##
## $fitted
## [,1]
## [1,] 0.60960932
## [2,] -0.09052769
## [3,] -0.21482638
## [4,] 0.13053755
## [5,] -0.21723642
## [6,] 0.39934961
##
## $pvaluesMarginal
## [1] NA
Once we have the regions we can proceed to build the required GRanges
object.
library("GenomicRanges")
## Build GRanges with sequence lengths
regions <- GRanges(
seqnames = tab$table$chr,
IRanges(start = tab$table$start, end = tab$table$end),
strand = "*", value = tab$table$value, area = tab$table$area,
cluster = tab$table$cluster, L = tab$table$L, clusterL = tab$table$clusterL
)
## Assign chr lengths
seqlengths(regions) <- seqlengths(
getChromInfoFromUCSC("hg19", as.Seqinfo = TRUE)
)[
names(seqlengths(regions))
]
## Explore the regions
regions
## GRanges object with 15 ranges and 5 metadata columns:
## seqnames ranges strand | value area cluster L
## <Rle> <IRanges> <Rle> | <numeric> <numeric> <numeric> <numeric>
## [1] chr1 2316-2631 * | -1.581475 15.81475 2 10
## [2] chr2 451-551 * | 1.589129 4.76739 3 3
## [3] chr1 456-526 * | 1.067883 3.20365 1 3
## [4] chr1 2176-2211 * | 0.784179 1.56836 2 2
## [5] chr1 2841 * | 1.201018 1.20102 2 1
## ... ... ... ... . ... ... ... ...
## [11] chr1 631 * | 0.618603 0.618603 1 1
## [12] chr1 1 * | 0.609609 0.609609 1 1
## [13] chr1 2911 * | -0.576423 0.576423 2 1
## [14] chr2 251 * | -0.556160 0.556160 3 1
## [15] chr1 2806 * | -0.521606 0.521606 2 1
## clusterL
## <integer>
## [1] 29
## [2] 20
## [3] 29
## [4] 29
## [5] 29
## ... ...
## [11] 29
## [12] 29
## [13] 29
## [14] 20
## [15] 29
## -------
## seqinfo: 2 sequences from an unspecified genome
Now that we have identified a set of differentially methylated regions we can proceed to creating the HTML report. Note that this report has less information than the DiffBind example because we don’t have a p-value variable.
## Load regionReport
library("regionReport")
## Now create the report
report <- renderReport(regions, "Example bumphunter",
pvalueVars = NULL,
densityVars = c(
"Area" = "area", "Value" = "value",
"Cluster Length" = "clusterL"
), significantVar = NULL,
output = "bumphunter-example", outdir = "bumphunter-example",
device = "png"
)
You can view the final report at bumphunter-example/bumphunter-example.html
or here.
In case the link does not work, a pre-compiled version of this document and its corresponding report are available at leekgroup.github.io/regionReportSupp/.
## Date generated:
Sys.time()
## [1] "2022-04-26 17:49:17 EDT"
## Time spent making this page:
proc.time()
## user system elapsed
## 13.567 1.180 15.118
## R and packages info:
options(width = 120)
sessioninfo::session_info()
## ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
## setting value
## version R version 4.2.0 RC (2022-04-19 r82224)
## os Ubuntu 20.04.4 LTS
## system x86_64, linux-gnu
## ui X11
## language (EN)
## collate C
## ctype en_US.UTF-8
## tz America/New_York
## date 2022-04-26
## pandoc 2.5 @ /usr/bin/ (via rmarkdown)
##
## ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
## package * version date (UTC) lib source
## annotate 1.74.0 2022-04-26 [2] Bioconductor
## AnnotationDbi 1.58.0 2022-04-26 [2] Bioconductor
## assertthat 0.2.1 2019-03-21 [2] CRAN (R 4.2.0)
## backports 1.4.1 2021-12-13 [2] CRAN (R 4.2.0)
## base64enc 0.1-3 2015-07-28 [2] CRAN (R 4.2.0)
## Biobase 2.56.0 2022-04-26 [2] Bioconductor
## BiocFileCache 2.4.0 2022-04-26 [2] Bioconductor
## BiocGenerics * 0.42.0 2022-04-26 [2] Bioconductor
## BiocIO 1.6.0 2022-04-26 [2] Bioconductor
## BiocManager 1.30.17 2022-04-22 [2] CRAN (R 4.2.0)
## BiocParallel 1.30.0 2022-04-26 [2] Bioconductor
## BiocStyle * 2.24.0 2022-04-26 [2] Bioconductor
## biomaRt 2.52.0 2022-04-26 [2] Bioconductor
## Biostrings 2.64.0 2022-04-26 [2] Bioconductor
## bit 4.0.4 2020-08-04 [2] CRAN (R 4.2.0)
## bit64 4.0.5 2020-08-30 [2] CRAN (R 4.2.0)
## bitops 1.0-7 2021-04-24 [2] CRAN (R 4.2.0)
## blob 1.2.3 2022-04-10 [2] CRAN (R 4.2.0)
## bookdown 0.26 2022-04-15 [2] CRAN (R 4.2.0)
## BSgenome 1.64.0 2022-04-26 [2] Bioconductor
## bslib 0.3.1 2021-10-06 [2] CRAN (R 4.2.0)
## bumphunter * 1.38.0 2022-04-26 [2] Bioconductor
## cachem 1.0.6 2021-08-19 [2] CRAN (R 4.2.0)
## checkmate 2.1.0 2022-04-21 [2] CRAN (R 4.2.0)
## cli 3.3.0 2022-04-25 [2] CRAN (R 4.2.0)
## cluster 2.1.3 2022-03-28 [2] CRAN (R 4.2.0)
## codetools 0.2-18 2020-11-04 [2] CRAN (R 4.2.0)
## colorspace 2.0-3 2022-02-21 [2] CRAN (R 4.2.0)
## crayon 1.5.1 2022-03-26 [2] CRAN (R 4.2.0)
## curl 4.3.2 2021-06-23 [2] CRAN (R 4.2.0)
## data.table 1.14.2 2021-09-27 [2] CRAN (R 4.2.0)
## DBI 1.1.2 2021-12-20 [2] CRAN (R 4.2.0)
## dbplyr 2.1.1 2021-04-06 [2] CRAN (R 4.2.0)
## DEFormats 1.24.0 2022-04-26 [2] Bioconductor
## DelayedArray 0.22.0 2022-04-26 [2] Bioconductor
## derfinder 1.30.0 2022-04-26 [2] Bioconductor
## derfinderHelper 1.30.0 2022-04-26 [2] Bioconductor
## DESeq2 1.36.0 2022-04-26 [2] Bioconductor
## digest 0.6.29 2021-12-01 [2] CRAN (R 4.2.0)
## doRNG 1.8.2 2020-01-27 [2] CRAN (R 4.2.0)
## dplyr 1.0.8 2022-02-08 [2] CRAN (R 4.2.0)
## edgeR 3.38.0 2022-04-26 [2] Bioconductor
## ellipsis 0.3.2 2021-04-29 [2] CRAN (R 4.2.0)
## evaluate 0.15 2022-02-18 [2] CRAN (R 4.2.0)
## fansi 1.0.3 2022-03-24 [2] CRAN (R 4.2.0)
## fastmap 1.1.0 2021-01-25 [2] CRAN (R 4.2.0)
## filelock 1.0.2 2018-10-05 [2] CRAN (R 4.2.0)
## foreach * 1.5.2 2022-02-02 [2] CRAN (R 4.2.0)
## foreign 0.8-82 2022-01-16 [2] CRAN (R 4.2.0)
## Formula 1.2-4 2020-10-16 [2] CRAN (R 4.2.0)
## genefilter 1.78.0 2022-04-26 [2] Bioconductor
## geneplotter 1.74.0 2022-04-26 [2] Bioconductor
## generics 0.1.2 2022-01-31 [2] CRAN (R 4.2.0)
## GenomeInfoDb * 1.32.0 2022-04-26 [2] Bioconductor
## GenomeInfoDbData 1.2.8 2022-04-21 [2] Bioconductor
## GenomicAlignments 1.32.0 2022-04-26 [2] Bioconductor
## GenomicFeatures 1.48.0 2022-04-26 [2] Bioconductor
## GenomicFiles 1.32.0 2022-04-26 [2] Bioconductor
## GenomicRanges * 1.48.0 2022-04-26 [2] Bioconductor
## ggplot2 3.3.5 2021-06-25 [2] CRAN (R 4.2.0)
## glue 1.6.2 2022-02-24 [2] CRAN (R 4.2.0)
## gridExtra 2.3 2017-09-09 [2] CRAN (R 4.2.0)
## gtable 0.3.0 2019-03-25 [2] CRAN (R 4.2.0)
## Hmisc 4.7-0 2022-04-19 [2] CRAN (R 4.2.0)
## hms 1.1.1 2021-09-26 [2] CRAN (R 4.2.0)
## htmlTable 2.4.0 2022-01-04 [2] CRAN (R 4.2.0)
## htmltools 0.5.2 2021-08-25 [2] CRAN (R 4.2.0)
## htmlwidgets 1.5.4 2021-09-08 [2] CRAN (R 4.2.0)
## httr 1.4.2 2020-07-20 [2] CRAN (R 4.2.0)
## IRanges * 2.30.0 2022-04-26 [2] Bioconductor
## iterators * 1.0.14 2022-02-05 [2] CRAN (R 4.2.0)
## jpeg 0.1-9 2021-07-24 [2] CRAN (R 4.2.0)
## jquerylib 0.1.4 2021-04-26 [2] CRAN (R 4.2.0)
## jsonlite 1.8.0 2022-02-22 [2] CRAN (R 4.2.0)
## KEGGREST 1.36.0 2022-04-26 [2] Bioconductor
## knitr 1.38 2022-03-25 [2] CRAN (R 4.2.0)
## knitrBootstrap 1.0.2 2018-05-24 [2] CRAN (R 4.2.0)
## lattice 0.20-45 2021-09-22 [2] CRAN (R 4.2.0)
## latticeExtra 0.6-29 2019-12-19 [2] CRAN (R 4.2.0)
## lifecycle 1.0.1 2021-09-24 [2] CRAN (R 4.2.0)
## limma 3.52.0 2022-04-26 [2] Bioconductor
## locfit * 1.5-9.5 2022-03-03 [2] CRAN (R 4.2.0)
## lubridate 1.8.0 2021-10-07 [2] CRAN (R 4.2.0)
## magrittr 2.0.3 2022-03-30 [2] CRAN (R 4.2.0)
## markdown 1.1 2019-08-07 [2] CRAN (R 4.2.0)
## Matrix 1.4-1 2022-03-23 [2] CRAN (R 4.2.0)
## MatrixGenerics 1.8.0 2022-04-26 [2] Bioconductor
## matrixStats 0.62.0 2022-04-19 [2] CRAN (R 4.2.0)
## memoise 2.0.1 2021-11-26 [2] CRAN (R 4.2.0)
## munsell 0.5.0 2018-06-12 [2] CRAN (R 4.2.0)
## nnet 7.3-17 2022-01-16 [2] CRAN (R 4.2.0)
## pillar 1.7.0 2022-02-01 [2] CRAN (R 4.2.0)
## pkgconfig 2.0.3 2019-09-22 [2] CRAN (R 4.2.0)
## plyr 1.8.7 2022-03-24 [2] CRAN (R 4.2.0)
## png 0.1-7 2013-12-03 [2] CRAN (R 4.2.0)
## prettyunits 1.1.1 2020-01-24 [2] CRAN (R 4.2.0)
## progress 1.2.2 2019-05-16 [2] CRAN (R 4.2.0)
## purrr 0.3.4 2020-04-17 [2] CRAN (R 4.2.0)
## qvalue 2.28.0 2022-04-26 [2] Bioconductor
## R6 2.5.1 2021-08-19 [2] CRAN (R 4.2.0)
## rappdirs 0.3.3 2021-01-31 [2] CRAN (R 4.2.0)
## RColorBrewer 1.1-3 2022-04-03 [2] CRAN (R 4.2.0)
## Rcpp 1.0.8.3 2022-03-17 [2] CRAN (R 4.2.0)
## RCurl 1.98-1.6 2022-02-08 [2] CRAN (R 4.2.0)
## RefManageR 1.3.0 2020-11-13 [2] CRAN (R 4.2.0)
## regionReport * 1.30.0 2022-04-26 [1] Bioconductor
## reshape2 1.4.4 2020-04-09 [2] CRAN (R 4.2.0)
## restfulr 0.0.13 2017-08-06 [2] CRAN (R 4.2.0)
## rjson 0.2.21 2022-01-09 [2] CRAN (R 4.2.0)
## rlang 1.0.2 2022-03-04 [2] CRAN (R 4.2.0)
## rmarkdown 2.14 2022-04-25 [2] CRAN (R 4.2.0)
## rngtools 1.5.2 2021-09-20 [2] CRAN (R 4.2.0)
## rpart 4.1.16 2022-01-24 [2] CRAN (R 4.2.0)
## Rsamtools 2.12.0 2022-04-26 [2] Bioconductor
## RSQLite 2.2.12 2022-04-02 [2] CRAN (R 4.2.0)
## rstudioapi 0.13 2020-11-12 [2] CRAN (R 4.2.0)
## rtracklayer 1.56.0 2022-04-26 [2] Bioconductor
## S4Vectors * 0.34.0 2022-04-26 [2] Bioconductor
## sass 0.4.1 2022-03-23 [2] CRAN (R 4.2.0)
## scales 1.2.0 2022-04-13 [2] CRAN (R 4.2.0)
## sessioninfo 1.2.2 2021-12-06 [2] CRAN (R 4.2.0)
## stringi 1.7.6 2021-11-29 [2] CRAN (R 4.2.0)
## stringr 1.4.0 2019-02-10 [2] CRAN (R 4.2.0)
## SummarizedExperiment 1.26.0 2022-04-26 [2] Bioconductor
## survival 3.3-1 2022-03-03 [2] CRAN (R 4.2.0)
## tibble 3.1.6 2021-11-07 [2] CRAN (R 4.2.0)
## tidyselect 1.1.2 2022-02-21 [2] CRAN (R 4.2.0)
## utf8 1.2.2 2021-07-24 [2] CRAN (R 4.2.0)
## VariantAnnotation 1.42.0 2022-04-26 [2] Bioconductor
## vctrs 0.4.1 2022-04-13 [2] CRAN (R 4.2.0)
## xfun 0.30 2022-03-02 [2] CRAN (R 4.2.0)
## XML 3.99-0.9 2022-02-24 [2] CRAN (R 4.2.0)
## xml2 1.3.3 2021-11-30 [2] CRAN (R 4.2.0)
## xtable 1.8-4 2019-04-21 [2] CRAN (R 4.2.0)
## XVector 0.36.0 2022-04-26 [2] Bioconductor
## yaml 2.3.5 2022-02-21 [2] CRAN (R 4.2.0)
## zlibbioc 1.42.0 2022-04-26 [2] Bioconductor
##
## [1] /tmp/RtmpNAXEQ8/Rinstceb9a24cb4a3d
## [2] /home/biocbuild/bbs-3.15-bioc/R/library
##
## ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────