Title: | Plotting functions for derfinder |
---|---|
Description: | This package provides plotting functions for results from the derfinder package. This helps separate the graphical dependencies required for making these plots from the core functionality of derfinder. |
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-12-14 03:16:16 UTC |
Source: | https://github.com/bioc/derfinderPlot |
For a given region found in calculatePvalues, plot the coverage for the cluster this region belongs to as well as some padding. The mean by group is shown to facilitate comparisons between groups. If annotation exists, you can plot the trancripts and exons (if any) overlapping in the vicinity of the region of interest.
plotCluster( idx, regions, annotation, coverageInfo, groupInfo, titleUse = "qval", txdb = NULL, p.ideogram = NULL, ... )
plotCluster( idx, regions, annotation, coverageInfo, groupInfo, titleUse = "qval", txdb = NULL, p.ideogram = NULL, ... )
idx |
A integer specifying the index number of the region of interest. This region is graphically highlighted by a red bar. |
regions |
The |
annotation |
The output from running matchGenes on the output from calculatePvalues. |
coverageInfo |
A DataFrame resulting from
loadCoverage using |
groupInfo |
A factor specifying the group membership of each sample. It will be used to color the samples by group. |
titleUse |
Whether to show the p-value ( |
txdb |
A transcript data base such as
|
p.ideogram |
If |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
|
See the parameter significantCut
in
calculatePvalues for how the significance cutoffs are
determined.
A ggplot2 plot that is ready to be printed out. Tecnically it is a ggbio object. The region with the red bar is the one whose information is shown in the title.
Leonardo Collado-Torres
loadCoverage, calculatePvalues, matchGenes, plotIdeogram
## Load data library("derfinder") ## Annotate the results with bumphunter::matchGenes() library("bumphunter") library("TxDb.Hsapiens.UCSC.hg19.knownGene") library("org.Hs.eg.db") genes <- annotateTranscripts( txdb = TxDb.Hsapiens.UCSC.hg19.knownGene, annotationPackage = "org.Hs.eg.db" ) annotation <- matchGenes(x = genomeRegions$regions, subject = genes) ## Make the plot plotCluster( idx = 1, regions = genomeRegions$regions, annotation = annotation, coverageInfo = genomeDataRaw$coverage, groupInfo = genomeInfo$pop, txdb = TxDb.Hsapiens.UCSC.hg19.knownGene ) ## Resize the plot window and the labels will look good. ## Not run: ## For a custom plot, check the ggbio and ggplot2 packages. ## Also feel free to look at the code for this function: plotCluster ## End(Not run)
## Load data library("derfinder") ## Annotate the results with bumphunter::matchGenes() library("bumphunter") library("TxDb.Hsapiens.UCSC.hg19.knownGene") library("org.Hs.eg.db") genes <- annotateTranscripts( txdb = TxDb.Hsapiens.UCSC.hg19.knownGene, annotationPackage = "org.Hs.eg.db" ) annotation <- matchGenes(x = genomeRegions$regions, subject = genes) ## Make the plot plotCluster( idx = 1, regions = genomeRegions$regions, annotation = annotation, coverageInfo = genomeDataRaw$coverage, groupInfo = genomeInfo$pop, txdb = TxDb.Hsapiens.UCSC.hg19.knownGene ) ## Resize the plot window and the labels will look good. ## Not run: ## For a custom plot, check the ggbio and ggplot2 packages. ## Also feel free to look at the code for this function: plotCluster ## End(Not run)
Plots an overview of the genomic locations of the identified regions (see calculatePvalues) in a karyotype view. The coloring can be done either by significant regions according to their p-values, significant by adjusted p-values, or by annotated region if using matchGenes.
plotOverview( regions, annotation = NULL, type = "pval", significantCut = c(0.05, 0.1), ... )
plotOverview( regions, annotation = NULL, type = "pval", significantCut = c(0.05, 0.1), ... )
regions |
The |
annotation |
The output from running matchGenes
on the output from calculatePvalues. It is only required
if |
type |
Must be either |
significantCut |
A vector of length two specifiying the cutoffs used to determine significance. The first element is used to determine significance for the p-values and the second element is used for the q-values. |
... |
Arguments passed to other methods and/or advanced arguments. Advanced arguments:
Passed to extendedMapSeqlevels. |
A ggplot2 plot that is ready to be printed out. Tecnically it is a ggbio object.
Leonardo Collado-Torres
## Construct toy data chrs <- paste0('chr', c(1:22, 'X', 'Y')) chrs <- factor(chrs, levels=chrs) library('GenomicRanges') regs <- GRanges(rep(chrs, 10), ranges=IRanges(runif(240, 1, 4e7), width=1e3), significant=sample(c(TRUE, FALSE), 240, TRUE, p=c(0.05, 0.95)), significantQval=sample(c(TRUE, FALSE), 240, TRUE, p=c(0.1, 0.9)), area=rnorm(240)) annotation <- data.frame(region=sample(c("upstream", "promoter", "overlaps 5'", "inside", "overlaps 3'", "close to 3'", "downstream"), 240, TRUE)) ## Type pval plotOverview(regs) ## Not run: ## Type qval plotOverview(regs, type='qval') ## Annotation plotOverview(regs, annotation, type='annotation') ## Resize the plots if needed. ## You might prefer to leave the legend at ggplot2's default option: right plotOverview(regs, legend.position='right') ## Although the legend looks better on the bottom plotOverview(regs, legend.position='bottom') ## Example knitr chunk for higher res plot using the CairoPNG device ```{r overview, message=FALSE, fig.width=7, fig.height=9, dev='CairoPNG', dpi=300} plotOverview(regs, base_size=30, areaRel=10, legend.position=c(0.95, 0.12)) ``` ## For more custom plots, take a look at the ggplot2 and ggbio packages ## and feel free to look at the code of this function: plotOverview ## End(Not run)
## Construct toy data chrs <- paste0('chr', c(1:22, 'X', 'Y')) chrs <- factor(chrs, levels=chrs) library('GenomicRanges') regs <- GRanges(rep(chrs, 10), ranges=IRanges(runif(240, 1, 4e7), width=1e3), significant=sample(c(TRUE, FALSE), 240, TRUE, p=c(0.05, 0.95)), significantQval=sample(c(TRUE, FALSE), 240, TRUE, p=c(0.1, 0.9)), area=rnorm(240)) annotation <- data.frame(region=sample(c("upstream", "promoter", "overlaps 5'", "inside", "overlaps 3'", "close to 3'", "downstream"), 240, TRUE)) ## Type pval plotOverview(regs) ## Not run: ## Type qval plotOverview(regs, type='qval') ## Annotation plotOverview(regs, annotation, type='annotation') ## Resize the plots if needed. ## You might prefer to leave the legend at ggplot2's default option: right plotOverview(regs, legend.position='right') ## Although the legend looks better on the bottom plotOverview(regs, legend.position='bottom') ## Example knitr chunk for higher res plot using the CairoPNG device ```{r overview, message=FALSE, fig.width=7, fig.height=9, dev='CairoPNG', dpi=300} plotOverview(regs, base_size=30, areaRel=10, legend.position=c(0.95, 0.12)) ``` ## For more custom plots, take a look at the ggplot2 and ggbio packages ## and feel free to look at the code of this function: plotOverview ## End(Not run)
This function takes the regions found in calculatePvalues and assigns them genomic states contructed with makeGenomicState. The main workhorse functions are countOverlaps and findOverlaps. For an alternative plot check plotCluster which is much slower and we recommend it's use only after quickly checking the results with this function.
plotRegionCoverage( regions, regionCoverage, groupInfo, nearestAnnotation, annotatedRegions, txdb = NULL, whichRegions = seq_len(min(100, length(regions))), colors = NULL, scalefac = 32, ask = interactive(), ylab = "Coverage", verbose = TRUE )
plotRegionCoverage( regions, regionCoverage, groupInfo, nearestAnnotation, annotatedRegions, txdb = NULL, whichRegions = seq_len(min(100, length(regions))), colors = NULL, scalefac = 32, ask = interactive(), ylab = "Coverage", verbose = TRUE )
regions |
The |
regionCoverage |
The output from getRegionCoverage
used on |
groupInfo |
A factor specifying the group membership of each sample. It will be used to color the samples by group. |
nearestAnnotation |
The output from matchGenes
used on |
annotatedRegions |
The output from annotateRegions
used on |
txdb |
A TxDb object. If specified, transcript annotation will be extracted from this object and used to plot the transcripts. |
whichRegions |
An integer vector with the index of the regions to plot. |
colors |
If |
scalefac |
The parameter used in preprocessCoverage. |
ask |
If |
ylab |
The name of the of the Y axis. |
verbose |
If |
A plot for every region showing the coverage of each sample at each base of the region as well as the summarized annotation information.
Andrew Jaffe, Leonardo Collado-Torres
calculatePvalues, getRegionCoverage, matchGenes, annotateRegions, plotCluster
## Load data library("derfinder") ## Annotate regions, first two regions only regions <- genomeRegions$regions[1:2] annotatedRegions <- annotateRegions( regions = regions, genomicState = genomicState$fullGenome, minoverlap = 1 ) ## Find nearest annotation with bumphunter::matchGenes() library("bumphunter") library("TxDb.Hsapiens.UCSC.hg19.knownGene") genes <- annotateTranscripts(txdb = TxDb.Hsapiens.UCSC.hg19.knownGene) nearestAnnotation <- matchGenes(x = regions, subject = genes) ## Obtain fullCov object fullCov <- list("21" = genomeDataRaw$coverage) ## Assign chr lengths using hg19 information library("GenomicRanges") seqlengths(regions) <- seqlengths(getChromInfoFromUCSC("hg19", as.Seqinfo = TRUE ))[ mapSeqlevels(names(seqlengths(regions)), "UCSC") ] ## Get the region coverage regionCov <- getRegionCoverage(fullCov = fullCov, regions = regions) # ## Make plots for the regions plotRegionCoverage( regions = regions, regionCoverage = regionCov, groupInfo = genomeInfo$pop, nearestAnnotation = nearestAnnotation, annotatedRegions = annotatedRegions, whichRegions = 1:2 ) ## Re-make plots with transcript information plotRegionCoverage( regions = regions, regionCoverage = regionCov, groupInfo = genomeInfo$pop, nearestAnnotation = nearestAnnotation, annotatedRegions = annotatedRegions, whichRegions = 1:2, txdb = TxDb.Hsapiens.UCSC.hg19.knownGene ) ## Not run: ## If you prefer, you can save the plots to a pdf file pdf("ders.pdf", h = 6, w = 9) plotRegionCoverage( regions = regions, regionCoverage = regionCov, groupInfo = genomeInfo$pop, nearestAnnotation = nearestAnnotation, annotatedRegions = annotatedRegions, whichRegions = 1:2, txdb = TxDb.Hsapiens.UCSC.hg19.knownGene, ask = FALSE ) dev.off() ## End(Not run)
## Load data library("derfinder") ## Annotate regions, first two regions only regions <- genomeRegions$regions[1:2] annotatedRegions <- annotateRegions( regions = regions, genomicState = genomicState$fullGenome, minoverlap = 1 ) ## Find nearest annotation with bumphunter::matchGenes() library("bumphunter") library("TxDb.Hsapiens.UCSC.hg19.knownGene") genes <- annotateTranscripts(txdb = TxDb.Hsapiens.UCSC.hg19.knownGene) nearestAnnotation <- matchGenes(x = regions, subject = genes) ## Obtain fullCov object fullCov <- list("21" = genomeDataRaw$coverage) ## Assign chr lengths using hg19 information library("GenomicRanges") seqlengths(regions) <- seqlengths(getChromInfoFromUCSC("hg19", as.Seqinfo = TRUE ))[ mapSeqlevels(names(seqlengths(regions)), "UCSC") ] ## Get the region coverage regionCov <- getRegionCoverage(fullCov = fullCov, regions = regions) # ## Make plots for the regions plotRegionCoverage( regions = regions, regionCoverage = regionCov, groupInfo = genomeInfo$pop, nearestAnnotation = nearestAnnotation, annotatedRegions = annotatedRegions, whichRegions = 1:2 ) ## Re-make plots with transcript information plotRegionCoverage( regions = regions, regionCoverage = regionCov, groupInfo = genomeInfo$pop, nearestAnnotation = nearestAnnotation, annotatedRegions = annotatedRegions, whichRegions = 1:2, txdb = TxDb.Hsapiens.UCSC.hg19.knownGene ) ## Not run: ## If you prefer, you can save the plots to a pdf file pdf("ders.pdf", h = 6, w = 9) plotRegionCoverage( regions = regions, regionCoverage = regionCov, groupInfo = genomeInfo$pop, nearestAnnotation = nearestAnnotation, annotatedRegions = annotatedRegions, whichRegions = 1:2, txdb = TxDb.Hsapiens.UCSC.hg19.knownGene, ask = FALSE ) dev.off() ## End(Not run)
Makes a venn diagram for the regions given the genomic state showing how many regions overlap introns, exons, intergenic regions, none or multiple groups.
vennRegions(annotatedRegions, subsetIndex = NULL, ...)
vennRegions(annotatedRegions, subsetIndex = NULL, ...)
annotatedRegions |
The output from annotateRegions
used on |
subsetIndex |
A vector of to use to subset the regions to use for the
venn diagram. It can be a logical vector of length equal to the number of
regions or an integer vector. If |
... |
Arguments passed to vennDiagram. |
Makes a venn diagram plot for the annotation given the genomic state and the actual venn counts used to make the plot.
Leonardo Collado-Torres
annotateRegions, vennCounts, vennDiagram
## Load data library("derfinder") ## Annotate regions annotatedRegions <- annotateRegions( regions = genomeRegions$regions, genomicState = genomicState$fullGenome, minoverlap = 1 ) ## Make venn diagram venn <- vennRegions(annotatedRegions) ## Add title and choose text color venn2 <- vennRegions(annotatedRegions, main = "Venn diagram", counts.col = "blue" ) ## Subset to only significant regions, so you don't have to annotate them ## again venn3 <- vennRegions(annotatedRegions, subsetIndex = genomeRegions$regions$significant == "TRUE", main = "Significant only" )
## Load data library("derfinder") ## Annotate regions annotatedRegions <- annotateRegions( regions = genomeRegions$regions, genomicState = genomicState$fullGenome, minoverlap = 1 ) ## Make venn diagram venn <- vennRegions(annotatedRegions) ## Add title and choose text color venn2 <- vennRegions(annotatedRegions, main = "Venn diagram", counts.col = "blue" ) ## Subset to only significant regions, so you don't have to annotate them ## again venn3 <- vennRegions(annotatedRegions, subsetIndex = genomeRegions$regions$significant == "TRUE", main = "Significant only" )