Last data update: 2014.03.03

R: Returns the flank size around the TSS for the x % CDSs
aroundPromoterR Documentation

Returns the flank size around the TSS for the x % CDSs

Description

Returns the flank size around the TSS for the x % CDSs

Usage

aroundPromoter(txdb, alnGRanges, percBestExpressed, flankSize)

Arguments

txdb

a TxDb object containing the annotations info to intersect with the alignment files.

alnGRanges

A GRanges object containing the alignment information. In order to improve the performance of this function the GAlignments BAM object should be transformed into a GRanges object containing the cigar match size information as metadata.

percBestExpressed

a numeric [between 0 and 1]. The percentage of the best expressed CDSs on which to plot the coverage around the TSS. Necessary if the shiftValue parameter must be estimated. Default value 0.03 (3%).

flankSize

a numeric positive integer. How many bp left and right of the TSS should the coverage be performed? Ex. flankSize <- 20

Value

A GRanges object containing the 1 bp ranges for the selected CDSs in the TSS defined flanking region.

Examples

#read the BAM into a GAlignments object using
#GenomicAlignments::readGAlignments
#the GAlignments object should be similar to ctrlGAlignments
data(ctrlGAlignments)
aln <- ctrlGAlignments
#transform the GAlignments object into a GRanges object (faster processing)
alnGRanges <- readsToStartOrEnd(aln, what="start")
#make a txdb object containing the annotations for the specified species.
#In this case hg19.
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene::TxDb.Hsapiens.UCSC.hg19.knownGene
#Please make sure that seqnames of txdb correspond to
#the seqnames of the alignment files ("chr" particle)
#if not rename the txdb seqlevels
#renameSeqlevels(txdb, sub("chr", "",seqlevels(txdb)))
#get the flanking region around the promoter of the best expressed CDSs
oneBinRanges <- aroundPromoter(txdb, alnGRanges)

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(RiboProfiling)
Loading required package: Biostrings
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: XVector
Warning messages:
1: replacing previous import 'BiocGenerics::Position' by 'ggplot2::Position' when loading 'RiboProfiling' 
2: replacing previous import 'ggplot2::Position' by 'BiocGenerics::Position' when loading 'ggbio' 
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/RiboProfiling/aroundPromoter.Rd_%03d_medium.png", width=480, height=480)
> ### Name: aroundPromoter
> ### Title: Returns the flank size around the TSS for the x % CDSs
> ### Aliases: aroundPromoter
> 
> ### ** Examples
> 
> #read the BAM into a GAlignments object using
> #GenomicAlignments::readGAlignments
> #the GAlignments object should be similar to ctrlGAlignments
> data(ctrlGAlignments)
> aln <- ctrlGAlignments
> #transform the GAlignments object into a GRanges object (faster processing)
> alnGRanges <- readsToStartOrEnd(aln, what="start")
> #make a txdb object containing the annotations for the specified species.
> #In this case hg19.
> txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene::TxDb.Hsapiens.UCSC.hg19.knownGene
> #Please make sure that seqnames of txdb correspond to
> #the seqnames of the alignment files ("chr" particle)
> #if not rename the txdb seqlevels
> #renameSeqlevels(txdb, sub("chr", "",seqlevels(txdb)))
> #get the flanking region around the promoter of the best expressed CDSs
> oneBinRanges <- aroundPromoter(txdb, alnGRanges)
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>