kpPlotBigWig | R Documentation |
Plots the wiggle values in a BigWig file. This function does not work on windows.
kpPlotBigWig(karyoplot, data, ymin=NULL, ymax="global", data.panel=1, r0=NULL, r1=NULL, col=NULL, border=NULL, clipping=TRUE, ...)
karyoplot |
(a |
data |
(a |
ymin |
(numeric) The minimum value to be plotted on the data panel. If NULL, the minimum between 0 and the minimum value in the WHOLE GENOME will be used. (deafults to NULL) |
ymax |
(numeric or c("global", "per.chr", "visible.region")) The maximum value to be plotted on the data.panel. It can be either a numeric value or one of c("global", "per.chr", "per.region"). "Global" will set ymax to the maximum value in the whole data file. "per.chr" will set ymax to the maximum value in each chromosome. "visible.region" will set ymax to the maximum value in the visible region. (defaults to "global") |
data.panel |
(numeric) The identifier of the data panel where the data is to be plotted. The available data panels depend on the plot type selected in the call to |
r0 |
(numeric) r0 and r1 define the vertical range of the data panel to be used to draw this plot. They can be used to split the data panel in different vertical ranges (similar to tracks in a genome browser) to plot differents data. If NULL, they are set to the min and max of the data panel, it is, to use all the available space. (defaults to NULL) |
r1 |
(numeric) r0 and r1 define the vertical range of the data panel to be used to draw this plot. They can be used to split the data panel in different vertical ranges (similar to tracks in a genome browser) to plot differents data. If NULL, they are set to the min and max of the data panel, it is, to use all the available space. (defaults to NULL) |
col |
(color) The fill color of the area. A single color. If NULL the color will be assigned automatically, either a lighter version of the color used for the outer line or gray if the line color is not defined. If NA no area will be drawn. (defaults to NULL) |
border |
(color) The color of the line enclosing the area. A single color. If NULL the color will be assigned automatically, either a darker version of the color used for the area or black if col=NA. If NA no border will be drawn. (Defaults to NULL) |
clipping |
(boolean) Only used if zooming is active. If TRUE, the data representation will be not drawn out of the drawing area (i.e. in margins, etc) even if the data overflows the drawing area. If FALSE, the data representation may overflow into the margins of the plot. (defaults to TRUE) |
... |
The ellipsis operator can be used to specify any additional graphical parameters. Any additional parameter will be passed to the internal calls to the R base plotting functions. |
kpPlotBigWig
plots the data contained in a binary file format called
BigWig. BigWig are used to efficiently store numeric values computed for
windows covering the whole genome, ususally the coverage from an NGS
experiment such as ChIP-seq. Only data required for the plotted region
is loaded, and when more than one chromosome is visible, it will load the
data for one crhomosome at a time.
The function accepts either a BigWigFile
oject or a
character
with the path to a valid big wig file. The character can
also be a URL to a remote server. In this case data will be loaded
transparently using the import
function from
rtracklayer
.
The data is plotted using kpArea
and therefore it is possible
to plot as a single line, a line with shaded area below or as a shaded area
only adjusting the col
and border
parameters.
Returns the original karyoplot object with the data computed (ymax and ymin values used) stored at karyoplot$latest.plot
Since this functions uses rtracklayer
BigWig infrastructure and it does not work on windows, this function won't work on windows either.
plotKaryotype
, kpArea
, kpPlotBAMDensity
if (.Platform$OS.type != "windows") {
bigwig.file <- system.file("extdata", "BRCA.genes.hg19.bw", package = "karyoploteR")
brca.genes.file <- system.file("extdata", "BRCA.genes.hg19.txt", package = "karyoploteR")
brca.genes <- toGRanges(brca.genes.file)
seqlevelsStyle(brca.genes) <- "UCSC"
kp <- plotKaryotype(zoom = brca.genes[1])
kp <- kpPlotBigWig(kp, data=bigwig.file, r0=0, r1=0.2)
kp <- kpPlotBigWig(kp, data=bigwig.file, r0=0.25, r1=0.45, border="red", lwd=2)
kp <- kpPlotBigWig(kp, data=bigwig.file, r0=0.5, r1=0.7, ymin=0, ymax=1000, border="gold", col=NA)
kpAxis(kp, r0=0.5, r1=0.7, ymin=0, ymax=1000)
kp <- kpPlotBigWig(kp, data=bigwig.file, r0=0.75, r1=0.95, ymin=0, ymax="visible.region", border="orchid", col=NA)
kpAxis(kp, r0=0.75, r1=0.95, ymin=0, ymax=kp$latest.plot$computed.values$ymax)
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.