renderReport: Generate a HTML/PDF report exploring a set of genomic regions

Description Usage Arguments Format Details Value Author(s) Examples

View source: R/renderReport.R

Description

This function generates a HTML report with quality checks, genome location exploration, and an interactive table with the results. Other output formats are possible such as PDF but lose the interactivity. Users can easily append to the report by providing a R Markdown file to customCode, or can customize the entire template by providing an R Markdown file to template.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
renderReport(
  regions,
  project = "",
  pvalueVars = c(`P-values` = "pval"),
  densityVars = NULL,
  significantVar = mcols(regions)$pval <= 0.05,
  annotation = NULL,
  nBestRegions = 500,
  customCode = NULL,
  outdir = "regionExploration",
  output = "regionExploration",
  browse = interactive(),
  txdb = NULL,
  device = "png",
  densityTemplates = list(Pvalue = templatePvalueDensity, Common = templateDensity,
    Manhattan = templateManhattan),
  template = NULL,
  theme = NULL,
  digits = 2,
  ...
)

templatePvalueDensity

templateDensity

templateManhattan

templatePvalueHistogram

templateHistogram

Arguments

regions

The set of genomic regions of interest as a GRanges object. All sequence lengths must be provided.

project

The title of the project.

pvalueVars

The names of the variables with values between 0 and 1 to plot density values by chromosome and a table for commonly used cutoffs. Most commonly used to explore p-value distributions. If a named character vector is provided, the names are used in the plot titles.

densityVars

The names of variables to use for making density plots by chromosome. Commonly used to explore scores and other variables given by region. If a named character vector is provided, the names are used in the plot titles.

significantVar

A logical variable differentiating statistically significant regions from the rest. When provided, both types of regions are compared against each other to see differences in width, location, etc.

annotation

The output from matchGenes used on regions. Note that this can take time for a large set of regions so it's better to pre-compute this information and save it.

nBestRegions

The number of regions to include in the interactive table.

customCode

An absolute path to a child R Markdown file with code to be evaluated before the reproducibility section. Its useful for users who want to customize the report by adding conclusions derived from the data and/or further quality checks and plots.

outdir

The name of output directory.

output

The name of output HTML file (without the html extension).

browse

If TRUE the HTML report is opened in your browser once it's completed.

txdb

Specify the transcription database to use for identifying the closest genes via matchGenes. If NULL it will use TxDb.Hsapiens.UCSC.hg19.knownGene by default.

device

The graphical device used when knitting. See more at http://yihui.name/knitr/options (dev argument).

densityTemplates

A list of length 3 with templates for the p-value density plots (variables from pvalueVars), the continuous variables density plots (variables from densityVars), and Manhattan plots for the p-value variables (pvalueVars). These templates are processed by whisker.render. Check the default templates for more information. The densityTemplates argument is available for those users interested in customizing these plots. For example, to show histograms instead of density plots use templatePvalueHistogram and templateHistogram instead of templatePvalueDensity and templateDensity respectively.

template

Template file to use for the report. If not provided, will use the default file found in regionExploration/regionExploration.Rmd within the package source.

theme

A ggplot2 theme to use for the plots made with ggplot2.

digits

The number of digits to round to in the interactive table of the top nBestRegions. Note that p-values and adjusted p-values won't be rounded.

...

Arguments passed to other methods and/or advanced arguments. Advanced arguments:

overviewParams

A two element list with base_size and areaRel that control the text size for the genomic overview plots.

output_format

Either html_document, pdf_document or knitrBootstrap::bootstrap_document unless you modify the YAML template.

clean

Logical, whether to clean the results or not. Passed to render.

Format

An object of class character of length 1.

An object of class character of length 1.

An object of class character of length 1.

An object of class character of length 1.

An object of class character of length 1.

Details

Set output_format to 'knitrBootstrap::bootstrap_document' or 'pdf_document' if you want a HTML report styled by knitrBootstrap or a PDF report respectively. If using knitrBootstrap, we recommend the version available only via GitHub at https://github.com/jimhester/knitrBootstrap which has nicer features than the current version available via CRAN. You can also set the output_format to 'html_document' for a HTML report styled by rmarkdown. The default is set to 'BiocStyle::html_document'.

If you modify the YAML front matter of template, you can use other values for output_format.

The HTML report styled with knitrBootstrap can be smaller in size than the 'html_document' report.

Value

An HTML report with a basic exploration for the given set of genomic regions. See the example report at http://leekgroup.github.io/regionReport/reference/renderReport-example/regionExploration.html.

Author(s)

Leonardo Collado-Torres

Examples

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
## Load derfinder for an example set of regions
library("derfinder")
regions <- genomeRegions$regions

## Assign chr length
library("GenomicRanges")
seqlengths(regions) <- c("chr21" = 48129895)

## The output will be saved in the 'renderReport-example' directory
dir.create("renderReport-example", showWarnings = FALSE, recursive = TRUE)

## Generate the HTML report
report <- renderReport(regions, "Example run",
    pvalueVars = c(
        "Q-values" = "qvalues", "P-values" = "pvalues"
    ), densityVars = c(
        "Area" = "area", "Mean coverage" = "meanCoverage"
    ),
    significantVar = regions$qvalues <= 0.05, nBestRegions = 20,
    outdir = "renderReport-example"
)

if (interactive()) {
    ## Browse the report
    browseURL(report)
}

## See the example report at
## http://leekgroup.github.io/regionReport/reference/renderReport-example/regionExploration.html


## Check the default templates. For users interested in customizing these
## plots.
## For p-value variables:
cat(regionReport::templatePvalueDensity)

## For continous variables:
cat(regionReport::templateDensity)

## For Manhattan plots
cat(regionReport::templateManhattan)

##################################################
## bumphunter example mentioned in the vignette ##
##################################################

## Load bumphunter
library("bumphunter")

## 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)

## Explore data
lapply(tab, head)

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

## 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"
)

## See the example report at
## http://leekgroup.github.io/regionReport/reference/bumphunter-example/bumphunter-example.html

regionReport documentation built on Dec. 20, 2020, 2:01 a.m.