Title: | Making 2D Hilbert Curve |
---|---|
Description: | Hilbert curve is a type of space-filling curves that fold one dimensional axis into a two dimensional space, but with still preserves the locality. This package aims to provide an easy and flexible way to visualize data through Hilbert curve. |
Authors: | Zuguang Gu [aut, cre] |
Maintainer: | Zuguang Gu <[email protected]> |
License: | MIT + file LICENSE |
Version: | 2.1.0 |
Built: | 2024-12-18 03:24:27 UTC |
Source: | https://github.com/bioc/HilbertCurve |
Default color overlay for adding new layers
default_overlay(r0, g0, b0, r, g, b, alpha = 1)
default_overlay(r0, g0, b0, r, g, b, alpha = 1)
r0 |
red channel for the layers that are already in the plot. |
g0 |
green channel for the layers that are already in the plot. |
b0 |
blue channel for the layers that are already in the plot. |
r |
red channel for the new layer |
g |
green channel for the new layer |
b |
blue channel for the new layer |
alpha |
alpha channel for the new layer |
The default overlay is (take red channel for example) r*alpha + r0*(1-alpha)
.
A list which contains overlayed RGB colors (values between 0 and 1).
Zuguang Gu <[email protected]>
Color overlay function is always used in hc_layer,HilbertCurve-method
or hc_layer,GenomicHilbertCurve-method
.
# red (1, 0, 0) overlay to the grey (0.5, 0.5, 0.5) with 0.5 transparency default_overlay(1, 0, 0, 0.5, 0.5, 0.5, 0.5)
# red (1, 0, 0) overlay to the grey (0.5, 0.5, 0.5) with 0.5 transparency default_overlay(1, 0, 0, 0.5, 0.5, 0.5, 0.5)
Initialize a Hilbert curve specifically for genomic data
GenomicHilbertCurve(chr = paste0("chr", c(1:22, "X", "Y")), species = "hg19", background = NULL, ...)
GenomicHilbertCurve(chr = paste0("chr", c(1:22, "X", "Y")), species = "hg19", background = NULL, ...)
chr |
a vector of chromosome names. Note it should have 'chr' prefix. This argument will be ignored when |
species |
abbreviation of species, e.g. 'hg19' or 'mm10'. |
background |
the background can be provided as a |
... |
common arguments in |
Multiple chromosomes can be visualized in a same Hilbert curve. All chromosomes are concatenated on after the other based on the order which is specified.
Since chromosomes will have irregular shapes on the curve, under 'pixel' mode,
users can set border
option in hc_map,GenomicHilbertCurve-method
to highlight
borders of chromosomes to identify their locations on the curve.
A GenomicHilbertCurve-class
object
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed() gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_points(hc, gr) hc = GenomicHilbertCurve(chr = c("chr1", "chr2")) hc_points(hc, gr) bg = GRanges(seqnames = c("chr1", "chr2"), ranges = IRanges(c(1,10000000), c(10000000,20000000))) hc = GenomicHilbertCurve(background = bg, level = 6) hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr)))) hc_map(hc, fill = NA, border = "grey", add = TRUE)
require(circlize) require(GenomicRanges) bed = generateRandomBed() gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_points(hc, gr) hc = GenomicHilbertCurve(chr = c("chr1", "chr2")) hc_points(hc, gr) bg = GRanges(seqnames = c("chr1", "chr2"), ranges = IRanges(c(1,10000000), c(10000000,20000000))) hc = GenomicHilbertCurve(background = bg, level = 6) hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr)))) hc_map(hc, fill = NA, border = "grey", add = TRUE)
The GenomicHilbertCurve class
The GenomicHilbertCurve-class
is inherited from the HilbertCurve-class
. Basically
the structure of this class is almost the same as the HilbertCurve-class
but with several
additional slots added to facilitate visualizing genomic data.
The GenomicHilbertCurve-class
provides following methods:
GenomicHilbertCurve
: constructor method;
hc_points,GenomicHilbertCurve-method
: add points;
hc_segments,GenomicHilbertCurve-method
: add lines;
hc_rect,GenomicHilbertCurve-method
: add rectangles;
hc_polygon,GenomicHilbertCurve-method
: add poygons;
hc_text,GenomicHilbertCurve-method
: add text;
hc_layer,GenomicHilbertCurve-method
: add layers undel "pixel" mode;
hc_map,GenomicHilbertCurve-method
: show the map of different categories on the curve. Works both for "normal" and "pixel" mode
The usage of above functions are almost same as those functions for the HilbertCurve-class
except that the second argument which specifies the intervals should be a GRanges
object.
Zuguang Gu <[email protected]>
NULL
NULL
Add text to the center of the block
## S4 method for signature 'HilbertCurve' hc_centered_text(object, ir = NULL, labels, x1 = NULL, x2 = NULL, gp = gpar(), ...)
## S4 method for signature 'HilbertCurve' hc_centered_text(object, ir = NULL, labels, x1 = NULL, x2 = NULL, gp = gpar(), ...)
object |
A |
ir |
an |
labels |
text corresponding to intervals in |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
gp |
graphic parameters for text. It should be specified by |
... |
pass to |
If the interval is long enough that it represents as a block in the 2D space, the corresponding label is put approximately at center (or at least inside) of the block.
Please use hc_text,HilbertCurve-method
directly.
NULL
Zuguang Gu <[email protected]>
It is used in hc_map,GenomicHilbertCurve-method
to put chromosome names in the center
of chromosomes.
hc = HilbertCurve(1, 10) hc_rect(hc, x1 = c(1, 3, 7), x2 = c(3, 7, 10), gp = gpar(fill = 2:5)) hc_centered_text(hc, x1 = 1, x2 = 3, labels = "A") hc_centered_text(hc, x1 = 3, x2 = 7, labels = "B") hc_centered_text(hc, x1 = 7, x2 = 10, labels = "C")
hc = HilbertCurve(1, 10) hc_rect(hc, x1 = c(1, 3, 7), x2 = c(3, 7, 10), gp = gpar(fill = 2:5)) hc_centered_text(hc, x1 = 1, x2 = 3, labels = "A") hc_centered_text(hc, x1 = 3, x2 = 7, labels = "B") hc_centered_text(hc, x1 = 7, x2 = 10, labels = "C")
Method dispatch page for hc_layer
.
hc_layer
can be dispatched on following classes:
hc_layer,GenomicHilbertCurve-method
, GenomicHilbertCurve-class
class method
hc_layer,HilbertCurve-method
, HilbertCurve-class
class method
# no example NULL
# no example NULL
Add a new layer to the Hilbert curve
## S4 method for signature 'GenomicHilbertCurve' hc_layer(object, gr, col = "red", border = NA, mean_mode = c("w0", "absolute", "weighted", "max_freq"), grid_line = 0, grid_line_col = "black", overlay = default_overlay)
## S4 method for signature 'GenomicHilbertCurve' hc_layer(object, gr, col = "red", border = NA, mean_mode = c("w0", "absolute", "weighted", "max_freq"), grid_line = 0, grid_line_col = "black", overlay = default_overlay)
object |
a |
gr |
a |
col |
a scalar or a vector of colors which correspond to regions in |
border |
a scalar or a vector of colors which correspond to the borders of regions. Set it to |
mean_mode |
Under 'pixel' mode, each pixel represents a small window. This argument provides methods to summarize value for the small window if the input genomic regions can not completely overlap with the window, pass to |
grid_line |
whether add grid lines to show blocks of the Hilber curve, pass to |
grid_line_col |
color for the grid lines, pass to |
overlay |
a self-defined function which defines how to overlay new layer to the plot, pass to |
It is basically a wrapper of hc_layer,HilbertCurve-method
.
Refer to hc_layer,HilbertCurve-method
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed() gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve(mode = "pixel", level = 9) hc_layer(hc, gr, col = rand_color(length(gr)))
require(circlize) require(GenomicRanges) bed = generateRandomBed() gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve(mode = "pixel", level = 9) hc_layer(hc, gr, col = rand_color(length(gr)))
Add a new layer to the Hilbert curve
## S4 method for signature 'HilbertCurve' hc_layer(object, ir = NULL, x1 = NULL, x2 = x1, col = "red", border = NA, mean_mode = c("w0", "absolute", "weighted", "max_freq"), grid_line = 0, grid_line_col = "black", overlay = default_overlay)
## S4 method for signature 'HilbertCurve' hc_layer(object, ir = NULL, x1 = NULL, x2 = x1, col = "red", border = NA, mean_mode = c("w0", "absolute", "weighted", "max_freq"), grid_line = 0, grid_line_col = "black", overlay = default_overlay)
object |
A |
ir |
an |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
col |
a scalar or a vector of colors which correspond to intervals in |
border |
a scalar or a vector of colors for the borders of intervals. Set it to |
mean_mode |
Under 'pixel' mode, each pixel represents a small window. This argument provides methods to summarize value for the small window if the input intervals can not completely overlap with the window. See explanation in |
grid_line |
whether add grid lines to show blocks of the Hilber curve. It should be an integer number and there will be |
grid_line_col |
color for the grid lines |
overlay |
a self-defined function which defines how to overlay new layer to the plot. By default it is |
This function only works under 'pixel' mode.
Under "pixel" mode, color is the only graphic representation of values in the input intervals.
To make a more precise and robust color mapping, users may consider colorRamp2
to create
a color mapping function.
If you want to add more than one layers to the curve, remember to set colors with transparency.
overlay
argument is useful for changing color themes for the overlapped areas, please refer to the vignette
to see examples of how to swith color themes in easy ways.
No value is returned.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 100, level = 9, mode = "pixel") x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_layer(hc, ir) hc = HilbertCurve(1, 100, level = 9, mode = "pixel") hc_layer(hc, ir, grid_line = 3) hc = HilbertCurve(1, 100, level = 9, mode = "pixel") hc_layer(hc, ir, border = "black")
hc = HilbertCurve(1, 100, level = 9, mode = "pixel") x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_layer(hc, ir) hc = HilbertCurve(1, 100, level = 9, mode = "pixel") hc_layer(hc, ir, grid_line = 3) hc = HilbertCurve(1, 100, level = 9, mode = "pixel") hc_layer(hc, ir, border = "black")
Level of the Hilbert curve
## S4 method for signature 'HilbertCurve' hc_level(object)
## S4 method for signature 'HilbertCurve' hc_level(object)
object |
A |
The level of the Hilbert curve.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 100) hc_level(hc) hc = HilbertCurve(1, 100, level = 5) hc_level(hc)
hc = HilbertCurve(1, 100) hc_level(hc) hc = HilbertCurve(1, 100, level = 5) hc_level(hc)
Draw a map which represents positions of different chromosomes on the curve
## S4 method for signature 'GenomicHilbertCurve' hc_map(object, level = 7, fill = rand_color(length(background), transparency = 0.5), border = NA, labels = names(object@background), show_labels = TRUE, labels_gp = gpar(), add = FALSE, ...)
## S4 method for signature 'GenomicHilbertCurve' hc_map(object, level = 7, fill = rand_color(length(background), transparency = 0.5), border = NA, labels = names(object@background), show_labels = TRUE, labels_gp = gpar(), add = FALSE, ...)
object |
a |
level |
Since a map does not need to have high resolution, a value of around 7 would be enough. If |
fill |
colors for different chromosomes, or more generally, for different 'seqnames'. |
border |
colors for the borders of chromosomes. Set it to |
labels |
label for each chromosome, or more generally, for different 'sequences' |
show_labels |
whether show text labels |
labels_gp |
graphic settings for labels |
add |
whether add the map to the current curve or draw it in a new graphic device. Notice if |
... |
pass to |
When multiple genomic categories (e.g. chromosomes) are drawn into one single Hilbert curve, a map which shows the positions of categories on the curve is necessary to distinguish different genomic categories.
Under "pixel" mode, if the map is directly added to the Hilbert curve, no chromosome name is drawn. The chromosome names are only drawn if the map is plotted in a new graphic device or added to the Hilbert curve under "normal" mode.
Just be careful if you directly overlay the map to the curve that the color of the map does not affect the original plot too much.
A GenomicHilbertCurve-class
object
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr)))) # add it in the same graphic device hc_map(hc, fill = rand_color(24, transparency = 0.5), add = TRUE) # add the map only with borders hc = GenomicHilbertCurve() hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr)))) hc_map(hc, fill = NA, border = "grey", add = TRUE) # or open a new graphic device hc_map(hc, fill = rand_color(24))
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr)))) # add it in the same graphic device hc_map(hc, fill = rand_color(24, transparency = 0.5), add = TRUE) # add the map only with borders hc = GenomicHilbertCurve() hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr)))) hc_map(hc, fill = NA, border = "grey", add = TRUE) # or open a new graphic device hc_map(hc, fill = rand_color(24))
Add points to the Hilbert curve
## S4 method for signature 'HilbertCurve' hc_normal_points(object, ir = NULL, x1 = NULL, x2 = x1, gp = gpar(), pch = 1, size = unit(1, "char"))
## S4 method for signature 'HilbertCurve' hc_normal_points(object, ir = NULL, x1 = NULL, x2 = x1, gp = gpar(), pch = 1, size = unit(1, "char"))
object |
A |
ir |
an |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
size |
size of the points. It should be a |
pch |
shape of points, pass to |
gp |
graphic parameters for points. It should be specified by |
Points are added at the middle of the intervals in ir
(or x1
and x2
),
so there is only one point for each interval.
This function is used internally. Please use hc_points,HilbertCurve-method
instead.
A data frame which contains coordinates (in the 2D space) of points.
Zuguang Gu <[email protected]>
# see documentation of hc_points NULL
# see documentation of hc_points NULL
Adjust positions
## S4 method for signature 'HilbertCurve' hc_offset(object, x)
## S4 method for signature 'HilbertCurve' hc_offset(object, x)
object |
A |
x |
positions. |
Since internally positions are transformed to positive integers, if input positions are specified as negative values when initializing the Hilbert curve, a shift will be recorded internally and positions are transformed to positive value automatically.
A positive numeric value.
Zuguang Gu <[email protected]>
hc = HilbertCurve(-100, 100) hc_offset(hc, c(-100, -50, 0, 50, 100))
hc = HilbertCurve(-100, 100) hc_offset(hc, c(-100, -50, 0, 50, 100))
Save Hilbert curve as a PNG figure
## S4 method for signature 'HilbertCurve' hc_png(object, file = "HilbertCurve.png")
## S4 method for signature 'HilbertCurve' hc_png(object, file = "HilbertCurve.png")
object |
A |
file |
file name. If the suffix of the file name is not |
A PNG figure with resolution of 2^level x 2^level
is generated.
Only the body of the Hilbert curve will be written to PNG file.
This function only works under 'pixel' mode.
No value is returned.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 100, level = 9, mode = "pixel") x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_layer(hc, ir) hc_png(hc, file = "test.png")
hc = HilbertCurve(1, 100, level = 9, mode = "pixel") x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_layer(hc, ir) hc_png(hc, file = "test.png")
Method dispatch page for hc_points
.
hc_points
can be dispatched on following classes:
hc_points,HilbertCurve-method
, HilbertCurve-class
class method
hc_points,GenomicHilbertCurve-method
, GenomicHilbertCurve-class
class method
# no example NULL
# no example NULL
Add points to the Hilbert curve
## S4 method for signature 'GenomicHilbertCurve' hc_points(object, gr, np = max(c(2, 10 - hc_level(object))), size = unit(1, "char"), pch = 1, gp = gpar(), mean_mode = c("w0", "absolute", "weighted", "max_freq"), shape = "circle")
## S4 method for signature 'GenomicHilbertCurve' hc_points(object, gr, np = max(c(2, 10 - hc_level(object))), size = unit(1, "char"), pch = 1, gp = gpar(), mean_mode = c("w0", "absolute", "weighted", "max_freq"), shape = "circle")
object |
a |
gr |
a |
np |
pass to |
size |
size of points when |
pch |
shape of the points when |
gp |
graphic parameters of the points when |
mean_mode |
pass to |
shape |
shape of the points when |
It is basically a wrapper of hc_points,HilbertCurve-method
.
Refer to hc_points,HilbertCurve-method
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr))))
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_points(hc, gr, gp = gpar(fill = rand_color(length(gr))))
Add points to the Hilbert curve
## S4 method for signature 'HilbertCurve' hc_points(object, ir = NULL, x1 = NULL, x2 = x1, np = max(c(2, 10 - hc_level(object))), size = unit(1, "char"), pch = 1, gp = gpar(), mean_mode = c("w0", "absolute", "weighted", "max_freq"), shape = "circle")
## S4 method for signature 'HilbertCurve' hc_points(object, ir = NULL, x1 = NULL, x2 = x1, np = max(c(2, 10 - hc_level(object))), size = unit(1, "char"), pch = 1, gp = gpar(), mean_mode = c("w0", "absolute", "weighted", "max_freq"), shape = "circle")
object |
A |
ir |
an |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
np |
number of points (a circle or a square, ...) that are put in a segment. |
size |
size of the points. It should be a |
pch |
shape of points, used for points if |
gp |
graphic parameters for points. It should be specified by |
mean_mode |
when |
shape |
shape of points, used for points if |
If np
is set to 1 or NULL
, points will be added in the middle for each interval in ir
(or x1
, x2
).
If np
is set to a value larger or equal to 2, every segment on the curve is split by np
points (e.g. circles).
In this case, each point actually represent a window on the curve and when the window is not fully covered by
the input intervals, there are three different metrics to average the values in the window.
Following illustrates different settings for mean_mode
:
100 80 60 values in ir ++++++ +++ +++++ ir ================ window (width = 16) 4 3 3 overlap absolute: (100 + 80 + 60)/3 weighted: (100*4 + 80*3 + 60*3)/(4 + 3 + 3) w0: (100*4 + 80*3 + 60*3 + 0*6)/16
So which mode to use depends on specific scenario. If the background is not of interest, absolute
and weighted
modes may be proper and if the value also needs to be averaged with background, w0
is the proper choice. Section "Averaging models"
in the vignette gives a more detailed explanation for this argument.
There is one more value for mean_mode
which is max_freq
. max_freq
is mainly for discrete signals and in a segment,
value with the highest frequency (or with the highest length) is selected for this segment
If np >= 2
, the value of np
also controls the size of points.
Graphic parameters is always represented as numeric values (e.g. colors can be converted into numeric RGB values) and they will be averaged according to above rules.
Internally, it will depatch to hc_normal_points,HilbertCurve-method
or hc_segmented_points,HilbertCurve-method
depending on the value of np
.
A data frame which contains coordinates (in the 2D space) of points.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_points(hc, ir) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, x1 = c(1.5, 50.5), x2 = c(10.5, 60.5)) require(circlize) value = runif(length(ir)) col_fun = colorRamp2(range(value), c("white", "red")) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, ir, np = 3, shape = "star", gp = gpar(fill = col_fun(value))) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, ir, np = 0) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, np = 0, x1 = c(1.5, 50.5), x2 = c(10.5, 60.5)) hc_points(hc, np = 0, x1 = 70.5, gp = gpar(col = "red"))
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_points(hc, ir) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, x1 = c(1.5, 50.5), x2 = c(10.5, 60.5)) require(circlize) value = runif(length(ir)) col_fun = colorRamp2(range(value), c("white", "red")) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, ir, np = 3, shape = "star", gp = gpar(fill = col_fun(value))) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, ir, np = 0) hc = HilbertCurve(1, 100, level = 4, reference = TRUE) hc_points(hc, np = 0, x1 = c(1.5, 50.5), x2 = c(10.5, 60.5)) hc_points(hc, np = 0, x1 = 70.5, gp = gpar(col = "red"))
Method dispatch page for hc_polygon
.
hc_polygon
can be dispatched on following classes:
hc_polygon,GenomicHilbertCurve-method
, GenomicHilbertCurve-class
class method
hc_polygon,HilbertCurve-method
, HilbertCurve-class
class method
# no example NULL
# no example NULL
Add text to Hilbert curve
## S4 method for signature 'GenomicHilbertCurve' hc_polygon(object, gr, gp = gpar(), end_type = c("average", "expanding", "shrinking"))
## S4 method for signature 'GenomicHilbertCurve' hc_polygon(object, gr, gp = gpar(), end_type = c("average", "expanding", "shrinking"))
object |
a |
gr |
a |
gp |
pass to |
end_type |
pass to |
It is basically a wrapper of hc_polygon,HilbertCurve-method
.
Refer to hc_polygon,HilbertCurve-method
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 20) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_polygon(hc, gr)
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 20) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_polygon(hc, gr)
Add polygons to Hilbert curve
## S4 method for signature 'HilbertCurve' hc_polygon(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(), end_type = c("expanding", "average", "shrinking"))
## S4 method for signature 'HilbertCurve' hc_polygon(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(), end_type = c("expanding", "average", "shrinking"))
object |
a |
ir |
an |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
gp |
graphic parameters. It should be specified by |
end_type |
since two ends of a continuous interval do not necessarily completely overlap with the Hilbert curve segments, this argument controls how to determine the ends of the interval which will be presented on the curve. |
Drawing polygons are quite visually similar as drawing rectangles. The major differences are: 1) for rectangles, colors for the ends of the interval can change if they are not completely covered by the Hilbert curve segments, and 2) polygons can have borders.
Normally polygons are used to mark areas in the Hilbert curve.
No value is returned.
Zuguang Gu <[email protected]>
require(IRanges) ir = IRanges(10, 40) hc = HilbertCurve(0, 100, level = 4, reference = TRUE) hc_segments(hc, ir) hc_text(hc, x1 = 10:40, labels = 10:40) hc_polygon(hc, ir, gp = gpar(fill = "#FF000080", col = 1))
require(IRanges) ir = IRanges(10, 40) hc = HilbertCurve(0, 100, level = 4, reference = TRUE) hc_segments(hc, ir) hc_text(hc, x1 = 10:40, labels = 10:40) hc_polygon(hc, ir, gp = gpar(fill = "#FF000080", col = 1))
Method dispatch page for hc_rect
.
hc_rect
can be dispatched on following classes:
hc_rect,GenomicHilbertCurve-method
, GenomicHilbertCurve-class
class method
hc_rect,HilbertCurve-method
, HilbertCurve-class
class method
# no example NULL
# no example NULL
Add rectangles on Hilbert curve
## S4 method for signature 'GenomicHilbertCurve' hc_rect(object, gr, gp = gpar(fill = "red", col = "red"), mean_mode = c("w0", "absolute", "weighted", "max_freq"))
## S4 method for signature 'GenomicHilbertCurve' hc_rect(object, gr, gp = gpar(fill = "red", col = "red"), mean_mode = c("w0", "absolute", "weighted", "max_freq"))
object |
a |
gr |
a |
gp |
pass to |
mean_mode |
pass to |
It is basically a wrapper of hc_rect,HilbertCurve-method
.
Refer to hc_rect,HilbertCurve-method
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_rect(hc, gr, gp = gpar(fill = rand_color(length(gr))))
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_rect(hc, gr, gp = gpar(fill = rand_color(length(gr))))
Add rectangles on Hilbert curve
## S4 method for signature 'HilbertCurve' hc_rect(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(fill = "red"), mean_mode = c("w0", "absolute", "weighted", "max_freq"))
## S4 method for signature 'HilbertCurve' hc_rect(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(fill = "red"), mean_mode = c("w0", "absolute", "weighted", "max_freq"))
object |
A |
ir |
an |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
gp |
graphic parameters for rectangles. It should be specified by |
mean_mode |
when a segment in the curve can not be overlapped with intervals in |
Rectangles are put if a segment in the Hilbert curve overlaps with the input intervals. You cannot set the width or height of the rectangles. It is always fixed (actually it is a square).
It can be thought as the low-resolution version of hc_layer,HilbertCurve-method
.
A data frame which contains coordinates (in the 2D space) of rectangles.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_rect(hc, ir)
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_rect(hc, ir)
Add points to the Hilbert curve
## S4 method for signature 'HilbertCurve' hc_segmented_points(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(), np = max(c(2, 10 - hc_level(object))), mean_mode = c("w0", "absolute", "weighted", "max_freq"), shape = "circle")
## S4 method for signature 'HilbertCurve' hc_segmented_points(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(), np = max(c(2, 10 - hc_level(object))), mean_mode = c("w0", "absolute", "weighted", "max_freq"), shape = "circle")
object |
A |
ir |
an |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
np |
number of points (a circle or a square, ...) that are put in a segment. |
gp |
graphic parameters for points. It should be specified by |
mean_mode |
when a segment in the curve overlaps with intervals in |
shape |
shape of points. Possible values are "circle", "square", "triangle", "hexagon", "star". |
Every segment on the curve is split by np
points.
This function is used internally, please use hc_points,HilbertCurve-method
directly.
A data frame which contains coordinates (in the 2D space) of points.
Zuguang Gu <[email protected]>
# see documentation of hc_points NULL
# see documentation of hc_points NULL
Method dispatch page for hc_segments
.
hc_segments
can be dispatched on following classes:
hc_segments,HilbertCurve-method
, HilbertCurve-class
class method
hc_segments,GenomicHilbertCurve-method
, GenomicHilbertCurve-class
class method
# no example NULL
# no example NULL
Add line segments to Hilbert curve
## S4 method for signature 'GenomicHilbertCurve' hc_segments(object, gr, gp = gpar(lty = 1, lwd = 1, col = 1))
## S4 method for signature 'GenomicHilbertCurve' hc_segments(object, gr, gp = gpar(lty = 1, lwd = 1, col = 1))
object |
a |
gr |
a |
gp |
pass to |
It is basically a wrapper of hc_segments,HilbertCurve-method
.
Refer to hc_segments,HilbertCurve-method
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_segments(hc, gr, gp = gpar(col = rand_color(length(gr))))
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 100) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_segments(hc, gr, gp = gpar(col = rand_color(length(gr))))
Add line segments to Hilbert curve
## S4 method for signature 'HilbertCurve' hc_segments(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(lty = 1, lwd = 1, col = 1))
## S4 method for signature 'HilbertCurve' hc_segments(object, ir = NULL, x1 = NULL, x2 = NULL, gp = gpar(lty = 1, lwd = 1, col = 1))
object |
A |
ir |
an |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
gp |
graphic parameters for lines. It should be specified by |
A data frame which contains coordinates (in the 2D space) of segments.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_segments(hc, ir)
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) hc_segments(hc, ir)
Method dispatch page for hc_text
.
hc_text
can be dispatched on following classes:
hc_text,HilbertCurve-method
, HilbertCurve-class
class method
hc_text,GenomicHilbertCurve-method
, GenomicHilbertCurve-class
class method
# no example NULL
# no example NULL
Add text to Hilbert curve
## S4 method for signature 'GenomicHilbertCurve' hc_text(object, gr, labels, gp = gpar(), centered_by = c("interval", "polygon"), ...)
## S4 method for signature 'GenomicHilbertCurve' hc_text(object, gr, labels, gp = gpar(), centered_by = c("interval", "polygon"), ...)
object |
a |
gr |
a |
labels |
pass to |
gp |
pass to |
centered_by |
how to define the "center" of the interval represented in Hilbert curve. Pass to |
... |
pass to |
It is basically a wrapper of hc_text,HilbertCurve-method
.
Refer to hc_text,HilbertCurve-method
Zuguang Gu <[email protected]>
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 20) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_text(hc, gr, labels = sample(letters, nrow(bed), replace = TRUE))
require(circlize) require(GenomicRanges) bed = generateRandomBed(nr = 20) gr = GRanges(seqnames = bed[[1]], ranges = IRanges(bed[[2]], bed[[3]])) hc = GenomicHilbertCurve() hc_text(hc, gr, labels = sample(letters, nrow(bed), replace = TRUE))
Add text to Hilbert curve
## S4 method for signature 'HilbertCurve' hc_text(object, ir = NULL, labels, x1 = NULL, x2 = x1, gp = gpar(), centered_by = c("interval", "polygon"), ...)
## S4 method for signature 'HilbertCurve' hc_text(object, ir = NULL, labels, x1 = NULL, x2 = x1, gp = gpar(), centered_by = c("interval", "polygon"), ...)
object |
A |
ir |
an |
labels |
text corresponding to intervals in |
x1 |
if start positions are not integers, they can be set by |
x2 |
if end positions are not integers, they can be set by |
gp |
graphic parameters for text. It should be specified by |
centered_by |
how to define the "center" of the interval represented in Hilbert curve. See Details section. |
... |
pass to |
If centered_by == "interval"
, the text is added correspoding to the middle of each interval in ir
,
while if centered_by == "polygon"
, the text is put in the visual center of the polygon of the interval
in the Hilbert curve.
A data frame which contains coordinates (in the 2D space) of text.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) labels = sample(letters, length(ir), replace = TRUE) hc_text(hc, ir, labels = labels)
hc = HilbertCurve(1, 100, level = 4, reference = TRUE) x = sort(sample(100, 20)) s = x[1:10*2 - 1] e = x[1:10*2] require(IRanges) ir = IRanges(s, e) labels = sample(letters, length(ir), replace = TRUE) hc_text(hc, ir, labels = labels)
Method dispatch page for hc_which
.
hc_which
can be dispatched on following classes:
hc_which,GenomicHilbertCurve-method
, GenomicHilbertCurve-class
class method
hc_which,HilbertCurve-method
, HilbertCurve-class
class method
# no example NULL
# no example NULL
Query regions
## S4 method for signature 'GenomicHilbertCurve' hc_which(object, ix, iy)
## S4 method for signature 'GenomicHilbertCurve' hc_which(object, ix, iy)
object |
a |
ix |
A single position on x-axis. |
iy |
A single position on y-axis. |
Values of ix
and iy
should be integers and take values in [1, 2^level].
A data frame with three columns chr
, start
and end
. The value corresponds to the genomic ranges.
# There is no example NULL
# There is no example NULL
Query regions
## S4 method for signature 'HilbertCurve' hc_which(object, ix, iy)
## S4 method for signature 'HilbertCurve' hc_which(object, ix, iy)
object |
A |
ix |
A single position on x-axis. |
iy |
A single position on y-axis. |
Values of ix
and iy
should be integers and take values in [1, 2^level].
A data frame with two columns start
and end
. The value corresponds to the range in data.
# There is no example NULL
# There is no example NULL
Initialize a Hilbert curve
HilbertCurve(s, e, level = 4, mode = c("normal", "pixel"), reference = FALSE, reference_gp = gpar(lty = 3, col = "#999999"), arrow = TRUE, zoom = NULL, newpage = TRUE, background_col = "transparent", background_border = NA, title = NULL, title_gp = gpar(fontsize = 16), start_from = c("bottomleft", "topleft", "bottomright", "topright"), first_seg = c("horizontal", "vertical"), legend = list(), padding = unit(2, "mm"))
HilbertCurve(s, e, level = 4, mode = c("normal", "pixel"), reference = FALSE, reference_gp = gpar(lty = 3, col = "#999999"), arrow = TRUE, zoom = NULL, newpage = TRUE, background_col = "transparent", background_border = NA, title = NULL, title_gp = gpar(fontsize = 16), start_from = c("bottomleft", "topleft", "bottomright", "topright"), first_seg = c("horizontal", "vertical"), legend = list(), padding = unit(2, "mm"))
s |
position that will be mapped as the start of the Hilbert curve. The value should a single numeric value. If it is a vector, the minimum is used. |
e |
position that will be mapped as the end of the Hilbert curve. The value should a single numeric value. If it is a vector, the maximum is used. |
level |
iteration level of the Hilbert curve. There will by |
mode |
"normal" mode is used for low |
reference |
whether add reference lines on the plot. Only works under 'normal' mode. The reference line is only used for illustrating how the curve folds. |
reference_gp |
graphic settings for the reference lines. It should be specified by |
arrow |
whether add arrows on the reference line. Only works under 'normal' mode. |
zoom |
Internally, position are stored as integer values. To better map the data to the Hilbert curve, the original positions are zoomed according to the range and the level of Hilbert curve. E.g. if the curve visualizes data ranging from 1 to 2 but level of the curve is set to 4, the positions will be zoomed by ~x2000 so that values like 1.5, 1.555 can be mapped to the curve with more accuracy. You don't need to care the zooming thing, proper zooming factor is calculated automatically. |
newpage |
whether call |
background_col |
background color. |
background_border |
background border border. |
title |
title of the plot. |
title_gp |
graphic parameters for the title. It should be specified by |
start_from |
which corner on the plot should the curve starts? |
first_seg |
the orientation of the first segment. |
legend |
a |
padding |
padding around the Hilbert curve. |
This funciton initializes a Hilbert curve with level level
which corresponds
to the range between s
and e
.
Under 'normal' mode, there is a visible Hilbert curve which plays like a folded axis and different low-level graphics can be added afterwards according to the coordinates. It works nice if the level of the Hilbert curve is small (say less than 6).
When the level is high (e.g. > 10), the whole 2D space will be almost completely filled by the curve and
it is impossible to add or visualize e.g. points on the curve. In this case, the 'pixel'
mode visualizes each tiny 'segment' as a pixel and maps values to colors. Internally, the whole plot
is represented as an RGB matrix and every time a new layer is added to the plot, the RGB matrix
will be updated according to the color overlay. When all the layers are added, normally a PNG figure is generated
directly from the RGB matrix. So the Hilbert
curve with level 11 will generate a PNG figure with 2048x2048 resolution. This is extremely
useful for visualize genomic data. E.g. If we make a Hilbert curve for human chromosome 1 with
level 11, then each pixel can represent 60bp (249250621/2048/2048
) which is of very high resolution.
Under 'pixel' mode, if the current device is an interactive deivce, every time a new layer is added,
the image will be add to the interactive device as a rastered image. But still you can use hc_png,HilbertCurve-method
to export the plot as PNG file.
To make it short and clear, under "normal" mode, you can use following low-level graphic functions:
And under "pixel" mode, you can use following functions:
Notice, s
and e
are not necessarily to be integers, it can be any values (e.g. numeric or even negative values).
A HilbertCurve-class
object.
Zuguang Gu <[email protected]>
HilbertCurve(1, 100, reference = TRUE) HilbertCurve(1, 100, level = 5, reference = TRUE) HilbertCurve(1, 100, title = "title", reference = TRUE) HilbertCurve(1, 100, start_from = "topleft", reference = TRUE) # plot with one legend require(ComplexHeatmap) legend = Legend(labels = c("a", "b"), title = "foo", legend_gp = gpar(fill = c("red", "blue"))) hc = HilbertCurve(1, 100, title = "title", legend = legend) hc_segments(hc, x1 = 20, x2 = 40) # plot with more than one legend require(circlize) legend1 = Legend(labels = c("a", "b"), title = "foo", legend_gp = gpar(fill = c("red", "blue"))) col_fun = colorRamp2(c(-1, 0, 1), c("green", "white", "red")) legend2 = Legend(col_fun = col_fun, title = "bar") hc = HilbertCurve(1, 100, title = "title", legend = list(legend1, legend2)) hc_segments(hc, x1 = 20, x2 = 40)
HilbertCurve(1, 100, reference = TRUE) HilbertCurve(1, 100, level = 5, reference = TRUE) HilbertCurve(1, 100, title = "title", reference = TRUE) HilbertCurve(1, 100, start_from = "topleft", reference = TRUE) # plot with one legend require(ComplexHeatmap) legend = Legend(labels = c("a", "b"), title = "foo", legend_gp = gpar(fill = c("red", "blue"))) hc = HilbertCurve(1, 100, title = "title", legend = legend) hc_segments(hc, x1 = 20, x2 = 40) # plot with more than one legend require(circlize) legend1 = Legend(labels = c("a", "b"), title = "foo", legend_gp = gpar(fill = c("red", "blue"))) col_fun = colorRamp2(c(-1, 0, 1), c("green", "white", "red")) legend2 = Legend(col_fun = col_fun, title = "bar") hc = HilbertCurve(1, 100, title = "title", legend = list(legend1, legend2)) hc_segments(hc, x1 = 20, x2 = 40)
The HilbertCurve class
Hilbert curve (https://en.wikipedia.org/wiki/Hilbert_curve ) is a type of space-filling curves that folds one-dimensional axis into a two-dimensional space, but still keeps the locality. It has advantages to visualize data with long axis with high resolution.
This package aims to provide an easy and flexible way to visualize data through Hilbert curve. The implementation and example figures are based on following sources:
The HilbertCurve-class
provides following methods:
HilbertCurve
: constructor method;
hc_points,HilbertCurve-method
: add points;
hc_segments,HilbertCurve-method
: add lines;
hc_rect,HilbertCurve-method
: add rectangles;
hc_polygon,HilbertCurve-method
: add polygons;
hc_text,HilbertCurve-method
: add text;
hc_layer,HilbertCurve-method
: add layers, works under "pixel" mode;
hc_png,HilbertCurve-method
: save plot as PNG format, works under "pixel" mode.
Zuguang Gu <[email protected]>
The GenomicHilbertCurve-class
inherits HilbertCurve-class
and is designed specifically for handling genomic data.
NULL
NULL
Whether the color is white
is_white(r, g, b, maxColorValue = 1)
is_white(r, g, b, maxColorValue = 1)
r |
Red channel. |
g |
Green channel. |
b |
Blue channel. |
maxColorValue |
1 or 255. |
# There is no example NULL
# There is no example NULL
Print the HilbertCurve object
## S4 method for signature 'HilbertCurve' show(object)
## S4 method for signature 'HilbertCurve' show(object)
object |
A |
No value is returned.
Zuguang Gu <[email protected]>
HilbertCurve(1, 100)
HilbertCurve(1, 100)
Transform zoomed positions to their original values
## S4 method for signature 'HilbertCurve' unzoom(object, x)
## S4 method for signature 'HilbertCurve' unzoom(object, x)
object |
A |
x |
positions. |
This is a reverse function of zoom,HilbertCurve-method
.
The function is used internally.
A numeric vector of original positions.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 2) z = zoom(hc, 1.5) unzoom(hc, z)
hc = HilbertCurve(1, 2) z = zoom(hc, 1.5) unzoom(hc, z)
Zoom original positions
## S4 method for signature 'HilbertCurve' zoom(object, x)
## S4 method for signature 'HilbertCurve' zoom(object, x)
object |
A |
x |
original positions. |
Internally, position are stored as integer values. To better map the data to the Hilbert curve, the original positions are zoomed according to the range and the level of Hilbert curve. E.g. if the curve visualizes data ranging from 1 to 2 but level of the curve is set to 4, the positions will be zoomed by ~x2000 so that values like 1.5, 1.555 can be mapped to the curve with more accuracy.
The function is used internally.
A numeric vector which is zoomed positions.
Zuguang Gu <[email protected]>
hc = HilbertCurve(1, 2) zoom(hc, 1.5)
hc = HilbertCurve(1, 2) zoom(hc, 1.5)