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,
...
)
A DESeqDataSet object with the results from running DESeq.
The title of the project.
interesting groups: a character vector of names in
colData(x)
to use for grouping. This parameter is passed to functions
such as plotPCA.
vector of colors used in heatmap. If NULL
, then a
a default set of colors will be used. This argument is passed to
pheatmap.
A DESeqResults object. If NULL
, then
results will be used on dds
with default parameters.
The number of features to include in the interactive table. Features are ordered by their adjusted p-values.
The number of best features to make plots of their counts. We recommend a small number, say 20.
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.
The name of output directory.
The name of output HTML file (without the html extension).
If TRUE
the HTML report is opened in your browser once
it's completed.
The graphical device used when knitting. See more at
http://yihui.name/knitr/options (dev
argument).
Template file to use for the report. If not provided, will use the default file found in DESeq2Exploration/DESeq2Exploration.Rmd within the package source.
A url used for searching the name of the features in
the web. By default http://www.ncbi.nlm.nih.gov/gene/?term=
is used
which is the recommended option when features are genes. It's only used
when the output is a HTML file.
A ggplot2 theme to use for the plots made with ggplot2.
The number of digits to round to in the interactive table of
the top nBestFeatures
. Note that p-values and adjusted p-values won't
be rounded.
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
The name of the package used for performing the
differential expression analysis. Either DESeq2
or edgeR
.
A DGEList object. NULL
by default and only
used by edgeReport.
The function call. NULL
by default and only used by
edgeReport.
Either html_document
, pdf_document
or
knitrBootstrap::bootstrap_document
unless you modify the YAML
template.
Logical, whether to clean the results or not. Passed to render.
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.
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.
## 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")
#> Loading required package: DEXSeq
#> Loading required package: BiocParallel
#> Loading required package: Biobase
#> 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, aperm, 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, saveRDS, setdiff, table,
#> tapply, union, unique, unsplit, which.max, which.min
#> Welcome to Bioconductor
#>
#> Vignettes contain introductory material; view with
#> 'browseVignettes()'. To cite Bioconductor, see
#> 'citation("Biobase")', and for packages 'citation("pkgname")'.
#> Loading required package: SummarizedExperiment
#> Loading required package: MatrixGenerics
#> Loading required package: matrixStats
#>
#> Attaching package: ‘matrixStats’
#> The following objects are masked from ‘package:Biobase’:
#>
#> anyMissing, rowMedians
#>
#> Attaching package: ‘MatrixGenerics’
#> The following objects are masked from ‘package:matrixStats’:
#>
#> colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
#> colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
#> colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
#> colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
#> colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
#> colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
#> colWeightedMeans, colWeightedMedians, colWeightedSds,
#> colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
#> rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
#> rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
#> rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
#> rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
#> rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
#> rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
#> rowWeightedSds, rowWeightedVars
#> The following object is masked from ‘package:Biobase’:
#>
#> rowMedians
#> Loading required package: GenomicRanges
#> Loading required package: stats4
#> Loading required package: S4Vectors
#>
#> Attaching package: ‘S4Vectors’
#> The following object is masked from ‘package:utils’:
#>
#> findMatches
#> The following objects are masked from ‘package:base’:
#>
#> I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: GenomeInfoDb
#> Loading required package: DESeq2
#> Loading required package: AnnotationDbi
#> Loading required package: RColorBrewer
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)
#> estimating size factors
#> estimating dispersions
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimates
#> fitting model and testing
## 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"
)
#> Writing 9 Bibtex entries ...
#> OK
#> Results written to file ‘DESeq2Report-example/DESeq2Exploration.bib’
#>
#>
#> processing file: DESeq2Exploration.Rmd
#> 1/47
#> 2/47 [docSetup]
#> 3/47
#> 4/47 [setup]
#> 5/47
#> 6/47 [PCA]
#> 7/47
#> 8/47 [sampleDist]
#> 9/47
#> 10/47 [MAplotalpha]
#> 11/47
#> 12/47 [MAplotalphaHalf]
#> 13/47
#> 14/47 [MAplotalpha-nBest]
#> 15/47
#> 16/47 [pvalueHistogram]
#> 17/47
#> 18/47 [pvalueSumm]
#> 19/47
#> 20/47 [pvalueTable]
#> 21/47
#> 22/47 [padjHistogram]
#> 23/47
#> 24/47 [padjSumm]
#> 25/47
#> 26/47 [padjTable]
#> 27/47
#> 28/47 [topFeatures]
#> 29/47
#> 30/47 [plotCounts]
#> 31/47
#> 32/47 [edgeR-BCV]
#> 33/47
#> 34/47 [edgeR-MDS]
#> 35/47
#> 36/47 [unnamed-chunk-1]
#> 37/47
#> 38/47 [thecall]
#> 39/47
#> 40/47 [reproducibility1]
#> 41/47
#> 42/47 [reproducibility2]
#> 43/47
#> 44/47 [reproducibility3]
#> 45/47
#> 46/47 [bibliography]
#> 47/47
#> output file: DESeq2Exploration.knit.md
#> /usr/bin/pandoc +RTS -K512m -RTS DESeq2Exploration.knit.md --to html4 --from markdown+autolink_bare_uris+tex_math_single_backslash --output DESeq2Exploration.html --lua-filter /__w/_temp/Library/bookdown/rmarkdown/lua/custom-environment.lua --lua-filter /usr/local/lib/R/site-library/rmarkdown/rmarkdown/lua/pagebreak.lua --lua-filter /usr/local/lib/R/site-library/rmarkdown/rmarkdown/lua/latex-div.lua --lua-filter /usr/local/lib/R/site-library/rmarkdown/rmarkdown/lua/table-classes.lua --embed-resources --standalone --wrap preserve --variable bs3=TRUE --section-divs --table-of-contents --toc-depth 3 --variable toc_float=1 --variable toc_selectors=h1,h2,h3 --variable toc_collapsed=1 --variable toc_smooth_scroll=1 --variable toc_print=1 --template /tmp/RtmplloSdM/BiocStyle/template.html --no-highlight --variable highlightjs=1 --number-sections --variable theme=bootstrap --css /__w/_temp/Library/BiocStyle/resources/html/bioconductor.css --mathjax --variable 'mathjax-url=https://mathjax.rstudio.com/latest/MathJax.js?config=TeX-AMS-MML_HTMLorMML' --include-in-header /tmp/RtmplloSdM/rmarkdown-str9b331c17d83.html --variable code_folding=hide --variable code_menu=1
#>
#> Output created: DESeq2Exploration.html
if (interactive()) {
## Browse the report
browseURL(report)
}
## See the example output at
## http://leekgroup.github.io/regionReport/reference/DESeq2Report-example/DESeq2Exploration.html
if (FALSE) { # \dontrun{
## Note that you can run the example using:
example("DESeq2Report", "regionReport", ask = FALSE)
} # }