readsEndPlot | R Documentation |
Plot the reads shifted from start/stop position of CDS.
readsEndPlot(
bamfile,
CDS,
toStartCodon = TRUE,
fiveEnd = TRUE,
shift = 0,
window = c(-29, 30),
readLen = 25:30,
ignore.seqlevelsStyle = FALSE
)
bamfile |
A BamFile object. |
CDS |
Output of prepareCDS |
toStartCodon |
What to search: start or end codon |
fiveEnd |
Search from five or three ends of the reads. |
shift |
number(1). Search from 5' end or 3' end of given number. if fiveEnd set to false, please set the shift as a negative number. |
window |
The window of CDS region to plot |
readLen |
The reads length used to plot |
ignore.seqlevelsStyle |
Ignore the sequence name style detection or not. |
The invisible list with counts numbers and reads in GRanges.
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"))
re <- readsEndPlot(bamfile, CDS, toStartCodon=TRUE)
readsEndPlot(re$reads, CDS, toStartCodon=TRUE, fiveEnd=FALSE)
#re <- readsEndPlot(bamfile, CDS, toStartCodon=FALSE)
#readsEndPlot(re$reads, CDS, toStartCodon=FALSE, fiveEnd=FALSE)
readsEndPlot(bamfile, CDS, shift=13)
#readsEndPlot(bamfile, CDS, fiveEnd=FALSE, shift=-16)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.