kpPlotRegions: kpPlotRegions

Description Usage Arguments Details Value See Also Examples

View source: R/kpPlotRegions.R

Description

Plots rectangles along the genome representing the regions (or intervals) specified by a GRanges object

Usage

1
kpPlotRegions(karyoplot, data, data.panel=1, r0=NULL, r1=NULL, col="black", border=NULL, avoid.overlapping=TRUE, num.layers=NULL, layer.margin=0.05, clipping=TRUE, ...)

Arguments

karyoplot

(a KaryoPlot object) This is the first argument to all data plotting functions of karyoploteR. A KaryoPlot object referring to the currently active plot.

data

(a GRanges or equivalent) It can be any of the formats accepted by the toGRanges function from the package regioneR.

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 plotKaryotype. (defaults to 1)

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 background color of the regions. (defaults to black)

border

(color) The color used to draw the border of the regions. If NULL, no border is drawn. (defaults to NULL)

avoid.overlapping

(boolean) Whether overlapping regions should be drawn as stacks (TRUE) on drawing one occluding the other in a single layer (FALSE). (defaults to TRUE)

num.layers

(numeric) The number of layers the plotting space should be divided into to allow for plotting overlapping regions. The lotting region will be divided into this many pieces regardless if any overlapping regions actually exist. If NULL, the maximum number of regions overlapping a single point in the genome. (defaults to NULL)

layer.margin

(numeric) The blank space left between layers of regions. (defaults to 0.05)

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.

Details

This is one of the high-level, or specialized, plotting functions of karyoploteR. It takes a GRanges object and plots its content. Overlapping regions can be stacked and the number of layers for overlapping regions can be set. In contrast with the low-level functions such as kpRect, it is not possible to specify the data using independent numeric vectors and the function only takes in GRanges.

Value

Returns the original karyoplot object, unchanged.

See Also

plotKaryotype, kpRect, kpSegments

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
 
 
 set.seed(1000)
 
 #Example 1: create 20 sets of non-overlapping random regions and plot them all. Add a coverage plot on top.
 kp <- plotKaryotype("hg19", plot.type=1, chromosomes=c("chr1", "chr2"))
 
 all.regs <- GRanges()
 
 nreps <- 20
 for(i in 1:nreps) {
   regs <- createRandomRegions(nregions = 100, length.mean = 10000000, length.sd = 1000000,
                               non.overlapping = TRUE, genome = "hg19", mask=NA)
   all.regs <- c(all.regs, regs)
   kpPlotRegions(kp, regs, r0 = (i-1)*(0.8/nreps), r1 = (i)*(0.8/nreps), col="#AAAAAA")
 }
 
 kpPlotCoverage(kp, all.regs, ymax = 20, r0=0.8,  r1=1, col="#CCCCFF")
 kpAxis(kp, ymin = 0, ymax= 20, numticks = 2, r0 = 0.8, r1=1)
 
 
 #Example 2: Do the same with a single bigger set of possibly overlapping regions
 
 kp <- plotKaryotype("hg19", plot.type=1, chromosomes=c("chr1", "chr2"))
 
 regs <- createRandomRegions(nregions = 1000, length.mean = 10000000, length.sd = 1000000,
                             non.overlapping = FALSE, genome = "hg19", mask=NA)
                             
 kpPlotRegions(kp, regs, r0 = 0, r1 = 0.8, col="#AAAAAA")
 
 kpPlotCoverage(kp, regs, ymax = 20, r0=0.8,  r1=1, col="#CCCCFF")
 kpAxis(kp, ymin = 0, ymax= 20, numticks = 2, r0 = 0.8, r1=1)
 
 
 

karyoploteR documentation built on Nov. 8, 2020, 5:52 p.m.