estimatePsite: Estimate P site position

View source: R/estimatePsite.R

estimatePsiteR Documentation

Estimate P site position

Description

Estimate P site position from a subset reads.

Usage

estimatePsite(
  bamfile,
  CDS,
  genome,
  anchor = "5end",
  readLen = c(25:30),
  ignore.seqlevelsStyle = FALSE
)

Arguments

bamfile

A BamFile object.

CDS

Output of prepareCDS

genome

A BSgenome object.

anchor

5end or 3end. Default is 5end.

readLen

The reads length used to estimate.

ignore.seqlevelsStyle

Ignore the sequence name style detection or not.

Value

A best P site position.

References

1: Bazzini AA, Johnstone TG, Christiano R, Mackowiak SD, Obermayer B, Fleming ES, Vejnar CE, Lee MT, Rajewsky N, Walther TC, Giraldez AJ. Identification of small ORFs in vertebrates using ribosome footprinting and evolutionary conservation. EMBO J. 2014 May 2;33(9):981-93. doi: 10.1002/embj.201488411. Epub 2014 Apr 4. PubMed PMID: 24705786; PubMed Central PMCID: PMC4193932.

Examples

library(Rsamtools)
bamfilename <- system.file("extdata", "RPF.WT.1.bam",
                           package="ribosomeProfilingQC")
yieldSize <- 10000000
bamfile <- BamFile(bamfilename, yieldSize = yieldSize)
#library(GenomicFeatures)
library(BSgenome.Drerio.UCSC.danRer10)
#txdb <- makeTxDbFromGFF(system.file("extdata",
 #         "Danio_rerio.GRCz10.91.chr1.gtf.gz",
 #         package="ribosomeProfilingQC"),
 #         organism = "Danio rerio",
 #         chrominfo = seqinfo(Drerio)["chr1"],
 #         taxonomyId = 7955)
#CDS <- prepareCDS(txdb)
CDS <- readRDS(system.file("extdata", "CDS.rds",
                           package="ribosomeProfilingQC"))
estimatePsite(bamfile, CDS, Drerio)


jianhong/ribosomeProfilingQC documentation built on April 15, 2024, 7:10 p.m.