Last data update: 2014.03.03

R: Import bedgraph to GRanges
bedgraph2normR Documentation

Import bedgraph to GRanges

Description

Function importing data from bedgraph format compatible with UCSC Genome Browser to norm_GR data frame. Warning: Compatible only with bedgraph files generated by norm2bedgraph function (bedgraph needs to have 2 tracks, first for plus strand, second for minus strand). May be used for transforming normalized data to another different annotation sets.

Usage

bedgraph2norm(bedgraph_file, fasta_file, txDb, bed_file,
  column_name = "bedgraph_score", add_to, track_strand = "+")

Arguments

bedgraph_file

path to compatible bedgraph file

fasta_file

path to fasta file which is used for a) choosing which transcripts to use (transcripts absent from fasta are not reported), b) providing sequence for to display in GRanges metadata

txDb

TranscriptDb object with transcript definitions. Names must match those in fasta_file

bed_file

character containing file path to BED file with transcript definitions. Supply txDb XOR bedfile

column_name

How to name imported metadata in GRanges

add_to

GRanges object made by other normalization function (dtcr(), slograt(), swinsor(), compdata()) to which values from bedgraph should be added.

track_strand

specifies which genomic strand the supplied bedgraph describes ("+" or "-"). Used only if the bedgraph file is composed of only one track.

Value

Function creates GRanges object or (if add_to specified) adds metadata to already existing object

Author(s)

Lukasz Jan Kielpinski, Nikos Sidiropoulos

See Also

norm2bedgraph, GR2norm_df, plotRNA, BED2txDb, dtcr, slograt, swinsor, compdata

Examples

dummy_euc_GR_control <- GRanges(seqnames="DummyRNA",
IRanges(start=round(runif(100)*100), width=round(runif(100)*100+1)),
        strand="+", EUC=round(runif(100)*100))
dummy_euc_GR_treated <- GRanges(seqnames="DummyRNA",
                                IRanges(start=round(runif(100)*100),
                                        width=round(runif(100)*100+1)),
                                strand="+", EUC=round(runif(100)*100))
dummy_comp_GR_control <- comp(dummy_euc_GR_control)
dummy_comp_GR_treated <- comp(dummy_euc_GR_treated)
dummy_norm <- dtcr(control_GR=dummy_comp_GR_control,
                   treated_GR=dummy_comp_GR_treated)

write(strwrap("chr1\t134212702\t134229870\tDummyRNA\t0\t+
              \t134212806\t134228958\t0\t8\t347,121,24,152,66,120,133,1973,
              \t0,8827,10080,11571,12005,13832,14433,15195,", width = 300),
      file="dummy.bed")
norm2bedgraph(norm_GR = dummy_norm, bed_file = "dummy.bed")

write(c(">DummyRNA", paste(sample(c("A","C","G","T"), 100, replace=TRUE),
                           collapse="")), file="dummy.fa")
bedgraph2norm(bedgraph_file = "out_file.bedgraph", fasta_file = "dummy.fa",
              bed_file = "dummy.bed")

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(RNAprobR)
Loading required package: GenomicFeatures
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
Loading required package: GenomicRanges
Loading required package: AnnotationDbi
Loading required package: Biobase
Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

Loading required package: plyr

Attaching package: 'plyr'

The following object is masked from 'package:IRanges':

    desc

The following object is masked from 'package:S4Vectors':

    rename

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/RNAprobR/bedgraph2norm.Rd_%03d_medium.png", width=480, height=480)
> ### Name: bedgraph2norm
> ### Title: Import bedgraph to GRanges
> ### Aliases: bedgraph2norm
> 
> ### ** Examples
> 
> dummy_euc_GR_control <- GRanges(seqnames="DummyRNA",
+ IRanges(start=round(runif(100)*100), width=round(runif(100)*100+1)),
+         strand="+", EUC=round(runif(100)*100))
> dummy_euc_GR_treated <- GRanges(seqnames="DummyRNA",
+                                 IRanges(start=round(runif(100)*100),
+                                         width=round(runif(100)*100+1)),
+                                 strand="+", EUC=round(runif(100)*100))
> dummy_comp_GR_control <- comp(dummy_euc_GR_control)
Fasta file not specified.
0 % of EUCs removed due to cutoff
> dummy_comp_GR_treated <- comp(dummy_euc_GR_treated)
Fasta file not specified.
0 % of EUCs removed due to cutoff
> dummy_norm <- dtcr(control_GR=dummy_comp_GR_control,
+                    treated_GR=dummy_comp_GR_treated)
> 
> write(strwrap("chr1\t134212702\t134229870\tDummyRNA\t0\t+
+               \t134212806\t134228958\t0\t8\t347,121,24,152,66,120,133,1973,
+               \t0,8827,10080,11571,12005,13832,14433,15195,", width = 300),
+       file="dummy.bed")
> norm2bedgraph(norm_GR = dummy_norm, bed_file = "dummy.bed")
Warning: normalization method to convert not specified. dtcr chosen.
Warning message:
Using togroup() on a CompressedIRangesList object is deprecated. Please
  use togroup(PartitioningByWidth(...)) instead. 
> 
> write(c(">DummyRNA", paste(sample(c("A","C","G","T"), 100, replace=TRUE),
+                            collapse="")), file="dummy.fa")
> bedgraph2norm(bedgraph_file = "out_file.bedgraph", fasta_file = "dummy.fa",
+               bed_file = "dummy.bed")
For RNA  DummyRNA positions outside FASTA annotation exist. N's added
GRanges object with 189 ranges and 2 metadata columns:
               seqnames     ranges strand |          nt bedgraph_score
                  <Rle>  <IRanges>  <Rle> | <character>      <numeric>
    DummyRNA.1 DummyRNA     [1, 1]      + |           C      0.3333333
    DummyRNA.2 DummyRNA     [2, 2]      + |           A      0.5460317
    DummyRNA.3 DummyRNA     [3, 3]      + |           A      0.3342521
    DummyRNA.4 DummyRNA     [4, 4]      + |           A      0.3757211
    DummyRNA.5 DummyRNA     [5, 5]      + |           T      0.2265148
           ...      ...        ...    ... .         ...            ...
  DummyRNA.185 DummyRNA [185, 185]      + |           N              0
  DummyRNA.186 DummyRNA [186, 186]      + |           N              0
  DummyRNA.187 DummyRNA [187, 187]      + |           N              0
  DummyRNA.188 DummyRNA [188, 188]      + |           N              0
  DummyRNA.189 DummyRNA [189, 189]      + |           N              0
  -------
  seqinfo: 1 sequence from an unspecified genome; no seqlengths
Warning messages:
1: In asMethod(object) : NAs introduced by coercion
2: Using togroup() on a CompressedIRangesList object is deprecated. Please
  use togroup(PartitioningByWidth(...)) instead. 
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>