Title: | Generate HTML or PDF reports for a set of genomic regions or DESeq2/edgeR results |
---|---|
Description: | Generate HTML or PDF reports to explore a set of regions such as the results from annotation-agnostic expression analysis of RNA-seq data at base-pair resolution performed by derfinder. You can also create reports for DESeq2 or edgeR results. |
Authors: | Leonardo Collado-Torres [aut, cre] , Andrew E. Jaffe [aut] , Jeffrey T. Leek [aut, ths] |
Maintainer: | Leonardo Collado-Torres <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.41.0 |
Built: | 2024-10-31 04:50:56 UTC |
Source: | https://github.com/bioc/regionReport |
This function generates a HTML report exploring the basic results from single base-level approach derfinder analysis results (<www.bioconductor.org/packages/derfinder>). The HTML report itself is generated using rmarkdown (http://rmarkdown.rstudio.com/). It works best after using mergeResults.
derfinderReport( prefix, outdir = "basicExploration", output = "basicExploration", project = prefix, browse = interactive(), nBestRegions = 100, makeBestClusters = TRUE, nBestClusters = 2, fullCov = NULL, hg19 = TRUE, p.ideos = NULL, txdb = NULL, device = "png", significantVar = "qvalue", customCode = NULL, template = NULL, theme = NULL, digits = 2, ... )
derfinderReport( prefix, outdir = "basicExploration", output = "basicExploration", project = prefix, browse = interactive(), nBestRegions = 100, makeBestClusters = TRUE, nBestClusters = 2, fullCov = NULL, hg19 = TRUE, p.ideos = NULL, txdb = NULL, device = "png", significantVar = "qvalue", customCode = NULL, template = NULL, theme = NULL, digits = 2, ... )
prefix |
The main data directory path where
mergeResults was run. It should be the same as
|
outdir |
The name of output directory relative to |
output |
The name of output HTML file (without the html extension). |
project |
The title of the project. |
browse |
If |
nBestRegions |
The number of region plots to make, ordered by area. |
makeBestClusters |
If |
nBestClusters |
The number of region cluster plots to make by taking
the |
fullCov |
A list where each element is the result from
loadCoverage used with |
hg19 |
If |
p.ideos |
A list where each element is the result of
plotIdeogram. If it's |
txdb |
Specify the transcription database to use for making the plots
for the top regions by area. If |
device |
The graphical device used when knitting. See more at
http://yihui.name/knitr/options ( |
significantVar |
A character variable specifying whether to use the
p-values, the FDR adjusted p-values or the FWER adjusted p-values to
determine significance. Has to be either |
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. |
template |
Template file to use for the report. If not provided, will use the default file found in basicExploration/basicExploration.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 |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
Passed to extendedMapSeqlevels. |
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.
An HTML report with a basic exploration of the derfinder results. See the example output at http://leekgroup.github.io/regionReport/reference/derfinderReport-example/basicExploration/basicExploration.html.
Leonardo Collado-Torres
mergeResults, analyzeChr, fullCoverage
## Load derfinder library("derfinder") ## The output will be saved in the 'derfinderReport-example' directory dir.create("derfinderReport-example", showWarnings = FALSE, recursive = TRUE) ## For convenience, the derfinder output has been pre-computed file.copy(system.file(file.path("extdata", "chr21"), package = "derfinder", mustWork = TRUE ), "derfinderReport-example", recursive = TRUE) ## Not run: ## If you prefer, you can generate the output from derfinder initialPath <- getwd() setwd(file.path(initialPath, "derfinderReport-example")) ## Collapse the coverage information collapsedFull <- collapseFullCoverage(list(genomeData$coverage), verbose = TRUE ) ## Calculate library size adjustments sampleDepths <- sampleDepth(collapsedFull, probs = c(0.5), nonzero = TRUE, verbose = TRUE ) ## Build the models group <- genomeInfo$pop adjustvars <- data.frame(genomeInfo$gender) models <- makeModels(sampleDepths, testvars = group, adjustvars = adjustvars) ## Analyze chromosome 21 analyzeChr( chr = "21", coverageInfo = genomeData, models = models, cutoffFstat = 1, cutoffType = "manual", seeds = 20140330, groupInfo = group, mc.cores = 1, writeOutput = TRUE, returnOutput = FALSE ) ## Change the directory back to the original one setwd(initialPath) ## End(Not run) ## Merge the results from the different chromosomes. In this case, there's ## only one: chr21 mergeResults( chrs = "21", prefix = "derfinderReport-example", genomicState = genomicState$fullGenome ) ## Load the options used for calculating the statistics load(file.path("derfinderReport-example", "chr21", "optionsStats.Rdata")) ## Generate the HTML report report <- derfinderReport( prefix = "derfinderReport-example", browse = FALSE, nBestRegions = 15, makeBestClusters = TRUE, fullCov = list("21" = genomeDataRaw$coverage), optionsStats = optionsStats ) if (interactive()) { ## Browse the report browseURL(report) } ## See the example output at ## http://leekgroup.github.io/regionReport/reference/derfinderReport-example/basicExploration/basicExploration.html ## Not run: ## Note that you can run the example using: example("derfinderReport", "regionReport", ask = FALSE) ## End(Not run)
## Load derfinder library("derfinder") ## The output will be saved in the 'derfinderReport-example' directory dir.create("derfinderReport-example", showWarnings = FALSE, recursive = TRUE) ## For convenience, the derfinder output has been pre-computed file.copy(system.file(file.path("extdata", "chr21"), package = "derfinder", mustWork = TRUE ), "derfinderReport-example", recursive = TRUE) ## Not run: ## If you prefer, you can generate the output from derfinder initialPath <- getwd() setwd(file.path(initialPath, "derfinderReport-example")) ## Collapse the coverage information collapsedFull <- collapseFullCoverage(list(genomeData$coverage), verbose = TRUE ) ## Calculate library size adjustments sampleDepths <- sampleDepth(collapsedFull, probs = c(0.5), nonzero = TRUE, verbose = TRUE ) ## Build the models group <- genomeInfo$pop adjustvars <- data.frame(genomeInfo$gender) models <- makeModels(sampleDepths, testvars = group, adjustvars = adjustvars) ## Analyze chromosome 21 analyzeChr( chr = "21", coverageInfo = genomeData, models = models, cutoffFstat = 1, cutoffType = "manual", seeds = 20140330, groupInfo = group, mc.cores = 1, writeOutput = TRUE, returnOutput = FALSE ) ## Change the directory back to the original one setwd(initialPath) ## End(Not run) ## Merge the results from the different chromosomes. In this case, there's ## only one: chr21 mergeResults( chrs = "21", prefix = "derfinderReport-example", genomicState = genomicState$fullGenome ) ## Load the options used for calculating the statistics load(file.path("derfinderReport-example", "chr21", "optionsStats.Rdata")) ## Generate the HTML report report <- derfinderReport( prefix = "derfinderReport-example", browse = FALSE, nBestRegions = 15, makeBestClusters = TRUE, fullCov = list("21" = genomeDataRaw$coverage), optionsStats = optionsStats ) if (interactive()) { ## Browse the report browseURL(report) } ## See the example output at ## http://leekgroup.github.io/regionReport/reference/derfinderReport-example/basicExploration/basicExploration.html ## Not run: ## Note that you can run the example using: example("derfinderReport", "regionReport", ask = FALSE) ## End(Not run)
This function generates a HTML report with exploratory data analysis plots
for DESeq2 results created with DESeq. 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
.
DESeq2Report( dds, project = "", intgroup, colors = NULL, res = NULL, nBest = 500, nBestFeatures = 20, customCode = NULL, outdir = "DESeq2Exploration", output = "DESeq2Exploration", browse = interactive(), device = "png", template = NULL, searchURL = "http://www.ncbi.nlm.nih.gov/gene/?term=", theme = NULL, digits = 2, ... )
DESeq2Report( dds, project = "", intgroup, colors = NULL, res = NULL, nBest = 500, nBestFeatures = 20, customCode = NULL, outdir = "DESeq2Exploration", output = "DESeq2Exploration", browse = interactive(), device = "png", template = NULL, searchURL = "http://www.ncbi.nlm.nih.gov/gene/?term=", theme = NULL, digits = 2, ... )
dds |
A DESeqDataSet object with the results from running DESeq. |
project |
The title of the project. |
intgroup |
interesting groups: a character vector of names in
|
colors |
vector of colors used in heatmap. If |
res |
A DESeqResults object. If |
nBest |
The number of features to include in the interactive table. Features are ordered by their adjusted p-values. |
nBestFeatures |
The number of best features to make plots of their counts. We recommend a small number, say 20. |
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 |
device |
The graphical device used when knitting. See more at
http://yihui.name/knitr/options ( |
template |
Template file to use for the report. If not provided, will use the default file found in DESeq2Exploration/DESeq2Exploration.Rmd within the package source. |
searchURL |
A url used for searching the name of the features in
the web. By default |
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 |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
|
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.
An HTML report with a basic exploration for the given set of DESeq2 results. See an example at http://leekgroup.github.io/regionReport/reference/DESeq2Report-example/DESeq2Exploration.html.
Leonardo Collado-Torres
## Load example data from the pasilla package as done in the DESeq2 vignette ## at ## <https://bioconductor.org/packages/release/bioc/vignettes/DESeq2/inst/doc/DESeq2.html#count-matrix-input>. library("pasilla") pasCts <- system.file("extdata", "pasilla_gene_counts.tsv", package = "pasilla", mustWork = TRUE ) pasAnno <- system.file("extdata", "pasilla_sample_annotation.csv", package = "pasilla", mustWork = TRUE ) cts <- as.matrix(read.csv(pasCts, sep = "\t", row.names = "gene_id")) coldata <- read.csv(pasAnno, row.names = 1) coldata <- coldata[, c("condition", "type")] coldata$condition <- factor(coldata$condition) coldata$type <- factor(coldata$type) rownames(coldata) <- sub("fb", "", rownames(coldata)) cts <- cts[, rownames(coldata)] ## Create DESeqDataSet object from the pasilla package library("DESeq2") dds <- DESeqDataSetFromMatrix( countData = cts, colData = coldata, design = ~condition ) dds <- DESeq(dds) ## The output will be saved in the 'DESeq2Report-example' directory dir.create("DESeq2Report-example", showWarnings = FALSE, recursive = TRUE) ## Generate the HTML report report <- DESeq2Report(dds, "DESeq2-example", c("condition", "type"), outdir = "DESeq2Report-example" ) if (interactive()) { ## Browse the report browseURL(report) } ## See the example output at ## http://leekgroup.github.io/regionReport/reference/DESeq2Report-example/DESeq2Exploration.html ## Not run: ## Note that you can run the example using: example("DESeq2Report", "regionReport", ask = FALSE) ## End(Not run)
## Load example data from the pasilla package as done in the DESeq2 vignette ## at ## <https://bioconductor.org/packages/release/bioc/vignettes/DESeq2/inst/doc/DESeq2.html#count-matrix-input>. library("pasilla") pasCts <- system.file("extdata", "pasilla_gene_counts.tsv", package = "pasilla", mustWork = TRUE ) pasAnno <- system.file("extdata", "pasilla_sample_annotation.csv", package = "pasilla", mustWork = TRUE ) cts <- as.matrix(read.csv(pasCts, sep = "\t", row.names = "gene_id")) coldata <- read.csv(pasAnno, row.names = 1) coldata <- coldata[, c("condition", "type")] coldata$condition <- factor(coldata$condition) coldata$type <- factor(coldata$type) rownames(coldata) <- sub("fb", "", rownames(coldata)) cts <- cts[, rownames(coldata)] ## Create DESeqDataSet object from the pasilla package library("DESeq2") dds <- DESeqDataSetFromMatrix( countData = cts, colData = coldata, design = ~condition ) dds <- DESeq(dds) ## The output will be saved in the 'DESeq2Report-example' directory dir.create("DESeq2Report-example", showWarnings = FALSE, recursive = TRUE) ## Generate the HTML report report <- DESeq2Report(dds, "DESeq2-example", c("condition", "type"), outdir = "DESeq2Report-example" ) if (interactive()) { ## Browse the report browseURL(report) } ## See the example output at ## http://leekgroup.github.io/regionReport/reference/DESeq2Report-example/DESeq2Exploration.html ## Not run: ## Note that you can run the example using: example("DESeq2Report", "regionReport", ask = FALSE) ## End(Not run)
This function generates a HTML report with exploratory data analysis plots
for edgeR results created. Other output formats are possible such as PDF
reports but they 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
.
edgeReport( dge, object, project = "", intgroup, colors = NULL, pAdjustMethod = "BH", alpha = 0.1, independentFiltering = FALSE, filter, theta, filterFun, nBest = 500, nBestFeatures = 20, customCode = NULL, outdir = "edgeRexploration", output = "edgeRexploration", browse = interactive(), device = "png", template = NULL, searchURL = "http://www.ncbi.nlm.nih.gov/gene/?term=", theme = NULL, digits = 2, ... )
edgeReport( dge, object, project = "", intgroup, colors = NULL, pAdjustMethod = "BH", alpha = 0.1, independentFiltering = FALSE, filter, theta, filterFun, nBest = 500, nBestFeatures = 20, customCode = NULL, outdir = "edgeRexploration", output = "edgeRexploration", browse = interactive(), device = "png", template = NULL, searchURL = "http://www.ncbi.nlm.nih.gov/gene/?term=", theme = NULL, digits = 2, ... )
dge |
A DGEList object. |
object |
A DGEExact or
DGELRT object that contains p-values stored in
|
project |
The title of the project. |
intgroup |
interesting groups: a character vector of names in
|
colors |
vector of colors used in heatmap. If |
pAdjustMethod |
the method to use for adjusting p-values, see p.adjust. This argument will be passed to results. |
alpha |
the significance cutoff used for optimizing the independent filtering (by default 0.1). If the adjusted p-value cutoff (FDR) will be a value other than 0.1, alpha should be set to that value. This argument will be passed to results. |
independentFiltering |
logical, whether independent filtering should be
applied automatically. By default it's set to |
filter |
the vector of filter statistics over which the independent filtering
will be optimized. By default the logCPM will be used if
|
theta |
the quantiles at which to assess the number of rejections from independent filtering. This argument is passed results. |
filterFun |
an optional custom function as described in results. |
nBest |
The number of features to include in the interactive table. Features are ordered by their adjusted p-values. |
nBestFeatures |
The number of best features to make plots of their counts. We recommend a small number, say 20. |
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 |
device |
The graphical device used when knitting. See more at
http://yihui.name/knitr/options ( |
template |
Template file to use for the report. If not provided, will use the default file found in DESeq2Exploration/DESeq2Exploration.Rmd within the package source. |
searchURL |
A url used for searching the name of the features in
the web. By default |
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 |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
|
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.
If you modify the YAML front matter of template
, you can use other
values for output_format
.
This report is similar to the one created by DESeq2Report with two additional plots exclusive for edgeR results. We designed the reports to be very similar intentionally and use the Bioconductor package DEFormats to achieve this goal.
An HTML report with a basic exploration for the given set of edgeR results. See the example report at http://leekgroup.github.io/regionReport/reference/edgeReport-example/edgeRexploration.html.
Leonardo Collado-Torres
## Create example data using DEFormats library("DEFormats") set.seed(20160407) counts <- simulateRnaSeqData() group <- rep(c("A", "B"), each = 3) ## Create DGEList object library("edgeR") dge <- DGEList(counts, group = group) ## Perform DE analysis with edgeR design <- model.matrix(~group) dge <- estimateDisp(dge, design) fit <- glmFit(dge, design) lrt <- glmLRT(fit, coef = 2) ## The output will be saved in the 'edgeReport-example' directory dir.create("edgeReport-example", showWarnings = FALSE, recursive = TRUE) ## Generate the HTML report report <- edgeReport(dge, lrt, project = "edgeR-example", intgroup = "group", outdir = "edgeReport-example" ) if (interactive()) { ## Browse the report browseURL(report) } ## See the example report at ## http://leekgroup.github.io/regionReport/reference/edgeReport-example/edgeRexploration.html ## Not run: ## Note that you can run the example using: example("edgeReport", "regionReport", ask = FALSE) ## End(Not run)
## Create example data using DEFormats library("DEFormats") set.seed(20160407) counts <- simulateRnaSeqData() group <- rep(c("A", "B"), each = 3) ## Create DGEList object library("edgeR") dge <- DGEList(counts, group = group) ## Perform DE analysis with edgeR design <- model.matrix(~group) dge <- estimateDisp(dge, design) fit <- glmFit(dge, design) lrt <- glmLRT(fit, coef = 2) ## The output will be saved in the 'edgeReport-example' directory dir.create("edgeReport-example", showWarnings = FALSE, recursive = TRUE) ## Generate the HTML report report <- edgeReport(dge, lrt, project = "edgeR-example", intgroup = "group", outdir = "edgeReport-example" ) if (interactive()) { ## Browse the report browseURL(report) } ## See the example report at ## http://leekgroup.github.io/regionReport/reference/edgeReport-example/edgeRexploration.html ## Not run: ## Note that you can run the example using: example("edgeReport", "regionReport", ask = FALSE) ## End(Not run)
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
.
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
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
regions |
The set of genomic regions of interest as a |
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 |
annotation |
The output from matchGenes used on
|
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 |
txdb |
Specify the transcription database to use for identifying the
closest genes via matchGenes. If |
device |
The graphical device used when knitting. See more at
http://yihui.name/knitr/options ( |
densityTemplates |
A list of length 3 with templates for the p-value
density plots (variables from |
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 |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
|
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.
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.
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.
Leonardo Collado-Torres
## 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
## 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