biovizBase (version 1.20.0)

transformToGenome: Transform GRanges to different coordinates and layout

Description

Used for coordiante genome transformation, other transformation in circular view.

Usage

"transformToGenome"(data, space.skip = 0.1, chr.weight = NULL) "transformToGenome"(data, space.skip = 0.1, chr.weight = NULL)
"transformToArch"(data, width = 1) transformToCircle(data, x = NULL, y = NULL, ylim = NULL, radius = 10, trackWidth =10, direction = c("clockwise", "anticlockwise"), mul = 0.05)
transformToRectInCircle(data, y = NULL, space.skip = 0.1, trackWidth = 10, radius = 10, direction = c("clockwise", "anticlockwise"), n = 100, mul = 0.05, chr.weight = NULL)
transformToBarInCircle(data, y = NULL, space.skip = 0.1, trackWidth = 10, radius = 10, direction = c("clockwise", "anticlockwise"), n = 100, mul = 0.05, chr.weight = NULL)
transformToSegInCircle(data, y = NULL, space.skip = 0.1, trackWidth = 10, radius = 10, direction = c("clockwise", "anticlockwise"), n = 100, chr.weight = NULL)
transformToLinkInCircle(data, linked.to, space.skip = 0.1, trackWidth = 10, radius = 10, link.fun = function(x, y, n = 100) bezier(x, y, evaluation = n), direction = c("clockwise", "anticlockwise"), chr.weight = NULL)
transformDfToGr(data, seqnames = NULL, start = NULL, end = NULL, width = NULL, strand = NULL, to.seqnames = NULL, to.start = NULL, to.end = NULL, to.width = NULL, to.strand = NULL, linked.to  = to.gr)
"transformToDf"(data)
is_coord_genome(data)

Arguments

data
a GRanges object.

for function transformDfToGr it's data.frame.

x
character for variable as x axis used for transformation.
y
character for variable as y axis used for transformation.
ylim
numeric range to control the ylimits on tracks when 'y' information is involved.
space.skip
numeric values indicates skipped ratio of whole space, not skipped space is identical between each space.
radius
numeric value, indicates radius when transform to a circle.
trackWidth
numeric value, for track width.
direction
"clockwise" or "counterclockwise", for layout or transform direction to circle.
mul
numeric value, passed to expand_range function, to control margin of y in the track.
n
integer value, control interpolated points numbers.
linked.to
a column name of GRanges object, indicate the linked line's end point which represented as a GRanges too..
link.fun
function used to generate linking lines.
seqnames
character or integer values for column name or index indicate variable mapped to seqnames, default NULL use "seqnames".
start
character or integer values for column name or index indicate variable mapped to start, default NULL use "start".
end
character or integer values for column name or index indicate variable mapped to end, default NULL use "end".
width
character or integer values for column name or index indicate variable mapped to width, default NULL use "width".
strand
character or integer values for column name or index indicate variable mapped to strand, default NULL use "strand".
to.seqnames
character or integer values for column name or index indicate variable mapped to linked seqnames, default NULL, create GRanges without new GRanges attached as column. If this varialbe is not NULL, this mean you try to parse linked GRanges object.
to.start
character or integer values for column name or index indicate variable mapped to start of linked GRanges, default NULL use "to.start".
to.end
character or integer values for column name or index indicate variable mapped to end of linked GRanges, default NULL use "to.end".
to.width
character or integer values for column name or index indicate variable mapped to width of linked GRanges, default NULL use "to.width".
to.strand
character or integer values for column name or index indicate variable mapped to strand, default NULL use "to.strand" or just use *.
chr.weight
numeric vectors which sum to

Value

A GRanges object, with calculated new variables, including ".circle.x" for transformed x position, ".circle.y" for transformed y position, ".circle.angle" for transformed angle.

Examples

Run this code
library(biovizBase)
library(GenomicRanges)
set.seed(1)
gr1 <- GRanges("chr1", IRanges(start = as.integer(runif(20, 1, 2e9)),
width = 5))
gr2 <- GRanges("chr2", IRanges(start = as.integer(runif(20, 1, 2e9)),
width = 5))
gr <- c(gr1, gr2)
gr.t <- transformToGenome(gr, space.skip = 0.1)
is_coord_genome(gr.t)
transformToCircle(gr.t)
transformToRectInCircle(gr)
transformToSegInCircle(gr)
values(gr1)$to.gr <- gr2
transformToLinkInCircle(gr1, linked.to = "to.gr")

Run the code above in your browser using DataCamp Workspace