Last data update: 2014.03.03
R: Adding disjoint levels to a GenomicRanges object
addStepping-method R Documentation
Adding disjoint levels to a GenomicRanges object
Description
Adding disjoint levels to a GenomicRanges object
Usage
## S4 method for signature 'GenomicRanges'
addStepping(obj, group.name, extend.size = 0,
fix = "center",
group.selfish = TRUE)
Arguments
obj
A GenomicRanges object
group.name
Column name in the elementMetadata which specify
the grouping information of all the entries. If provided, this
will make sure all intervals belong to the same group will try to
be on the same level and nothing falls in between.
extend.size
Adding invisible buffered region to the
GenomicRanges object, if it's 10, then adding 5 at both end. This
make the close neighbors assigned to the different levels and make
your eyes easy to identify.
fix
"start", "end", or "center"(default) passed to resize
, denoting
what to use as an anchor for each element of obj, and add extend.size
to it's current width.
group.selfish
Passed to addStepping
, control whether to show each group as
unique level or not. If set to FALSE
, if two groups are not
overlapped with each other, they will probably be layout in the same
level to save space.
Details
This is a tricky question, for example, pair-end RNA-seq data
could be represented as two set of GenomicRanges object, one
indicates the read, one indicates the junction. At the same time,
we need to make sure pair-ended read are shown on the same level,
and nothing falls in between. For better visualization of the
data, we may hope to add invisible extended buffer to the reads, so
closely neighbored reads will be on the different levels.
Value
A modified GenmicRanges object with stepping
as one column.
Author(s)
Tengfei Yin
Examples
library(GenomicRanges)
set.seed(1)
N <- 500
## sample GRanges
gr <- GRanges(seqnames =
sample(c("chr1", "chr2", "chr3", "chrX", "chrY"),
size = N, replace = TRUE),
IRanges(
start = sample(1:300, size = N, replace = TRUE),
width = sample(70:75, size = N,replace = TRUE)),
strand = sample(c("+", "-", "*"), size = N,
replace = TRUE),
value = rnorm(N, 10, 3), score = rnorm(N, 100, 30),
group = sample(c("Normal", "Tumor"),
size = N, replace = TRUE),
pair = sample(letters, size = N,
replace = TRUE))
## grouping and extending
head(addStepping(gr))
head(addStepping(gr, group.name = "pair"))
gr.close <- GRanges(c("chr1", "chr1"), IRanges(c(10, 20), width = 9))
addStepping(gr.close)
addStepping(gr.close, extend.size = 5)
Results
R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)
R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.
R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.
Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.
> library(biovizBase)
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/biovizBase/addStepping-method.Rd_%03d_medium.png", width=480, height=480)
> ### Name: addStepping-method
> ### Title: Adding disjoint levels to a GenomicRanges object
> ### Aliases: addStepping addStepping,GenomicRanges-method
>
> ### ** Examples
>
> library(GenomicRanges)
Loading required package: BiocGenerics
Loading required package: parallel
Attaching package: 'BiocGenerics'
The following objects are masked from 'package:parallel':
clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
clusterExport, clusterMap, parApply, parCapply, parLapply,
parLapplyLB, parRapply, parSapply, parSapplyLB
The following objects are masked from 'package:stats':
IQR, mad, xtabs
The following objects are masked from 'package:base':
Filter, Find, Map, Position, Reduce, anyDuplicated, append,
as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
rbind, rownames, sapply, setdiff, sort, table, tapply, union,
unique, unsplit
Loading required package: S4Vectors
Loading required package: stats4
Attaching package: 'S4Vectors'
The following objects are masked from 'package:base':
colMeans, colSums, expand.grid, rowMeans, rowSums
Loading required package: IRanges
Loading required package: GenomeInfoDb
> set.seed(1)
> N <- 500
> ## sample GRanges
> gr <- GRanges(seqnames =
+ sample(c("chr1", "chr2", "chr3", "chrX", "chrY"),
+ size = N, replace = TRUE),
+ IRanges(
+ start = sample(1:300, size = N, replace = TRUE),
+ width = sample(70:75, size = N,replace = TRUE)),
+ strand = sample(c("+", "-", "*"), size = N,
+ replace = TRUE),
+ value = rnorm(N, 10, 3), score = rnorm(N, 100, 30),
+ group = sample(c("Normal", "Tumor"),
+ size = N, replace = TRUE),
+ pair = sample(letters, size = N,
+ replace = TRUE))
>
> ## grouping and extending
> head(addStepping(gr))
GRanges object with 6 ranges and 5 metadata columns:
seqnames ranges strand | value score group pair
<Rle> <IRanges> <Rle> | <numeric> <numeric> <character> <character>
chr1 chr1 [241, 313] * | 8.120639 102.31909 Tumor n
chr1 chr1 [ 54, 126] + | 7.493114 64.50273 Tumor w
chr1 chr1 [273, 343] - | 13.374793 87.88790 Normal m
chr1 chr1 [138, 207] * | 9.951429 87.32883 Normal p
chr1 chr1 [ 25, 99] + | 12.831509 111.22355 Normal z
chr1 chr1 [ 33, 103] * | 12.346409 108.72000 Tumor m
stepping
<numeric>
chr1 2
chr1 29
chr1 9
chr1 18
chr1 20
chr1 21
-------
seqinfo: 5 sequences from an unspecified genome; no seqlengths
> head(addStepping(gr, group.name = "pair"))
GRanges object with 6 ranges and 5 metadata columns:
seqnames ranges strand | value score group pair
<Rle> <IRanges> <Rle> | <numeric> <numeric> <character> <character>
chr1 chr1 [241, 313] * | 8.120639 102.31909 Tumor n
chr1 chr1 [ 54, 126] + | 7.493114 64.50273 Tumor w
chr1 chr1 [273, 343] - | 13.374793 87.88790 Normal m
chr1 chr1 [138, 207] * | 9.951429 87.32883 Normal p
chr1 chr1 [ 25, 99] + | 12.831509 111.22355 Normal z
chr1 chr1 [ 33, 103] * | 12.346409 108.72000 Tumor m
stepping
<numeric>
chr1 14
chr1 22
chr1 13
chr1 16
chr1 25
chr1 13
-------
seqinfo: 5 sequences from an unspecified genome; no seqlengths
> gr.close <- GRanges(c("chr1", "chr1"), IRanges(c(10, 20), width = 9))
> addStepping(gr.close)
GRanges object with 2 ranges and 1 metadata column:
seqnames ranges strand | stepping
<Rle> <IRanges> <Rle> | <numeric>
chr1 chr1 [10, 18] * | 1
chr1 chr1 [20, 28] * | 1
-------
seqinfo: 1 sequence from an unspecified genome; no seqlengths
> addStepping(gr.close, extend.size = 5)
GRanges object with 2 ranges and 1 metadata column:
seqnames ranges strand | stepping
<Rle> <IRanges> <Rle> | <numeric>
chr1 chr1 [10, 18] * | 1
chr1 chr1 [20, 28] * | 2
-------
seqinfo: 1 sequence from an unspecified genome; no seqlengths
>
>
>
>
>
> dev.off()
null device
1
>