interest | R Documentation |
A read summarization function that countsns all the reads mapping to the introns/exons based on the users detailed parameter settings. The process can be run in parallel on multiple computing cores to improve it performance.
interest( bamFileYieldSize=1000000, bamFile, isPaired,
isPairedDuplicate=FALSE, isSingleReadDuplicate= NA, reference,
referenceGeneNames, referenceIntronExon, repeatsTableToFilter=c(),
junctionReadsOnly=FALSE, outFile, logFile="",
returnObj= FALSE, method=c("ExEx", "IntRet", "IntSpan", "ExSkip"),
strandSpecific,
bpparam, appendLogFile=FALSE, sampleName="",
scaleLength= c(TRUE,FALSE), scaleFragment= c(TRUE,TRUE),
limitRanges=GRanges(),
excludeFusionReads=FALSE,
loadLimitRangesReads=FALSE, ...)
bamFileYieldSize |
Maximum number of pair reads in the temprorary files created as the result of dividing the input .bam file. |
bamFile |
Path of the input bam file. |
isPaired |
Whether the bam file is the result of a paired end sequencing read mapping (TRUE) or not (FALSE). |
isPairedDuplicate |
Whether extract only (if set TRUE), filter (FALSE) or include (if set NA) PCR
dupplicates for paired mapped reads. It uses the FLAG field in the bam file to
filter the duplicate read. If the mapping software does not support detection
and flaging the duplicate reads |
isSingleReadDuplicate |
Whether extract only (if set TRUE), filter (FALSE) or include (if set NA) PCR dupplicates for single mapped reads. |
reference |
Dataframe to be used as reference; It should at least contain three same-size
vectors with the tag names |
referenceGeneNames |
A vector with the same size as the row-size of the reference which includes the gene names of the reference. |
referenceIntronExon |
A vector with the same size as the row-size of the reference with values "intron" and "exon" describing which (intron or exon) each row of the reference represents. |
repeatsTableToFilter |
A data.frame table with similar stucture to the |
junctionReadsOnly |
The parameter is considered if the |
outFile |
The name or path of the result file. |
logFile |
The log file path; if defined log information are written to the log file. |
returnObj |
If set |
method |
A vector describing the summarization methods to use; i.e. whether count reads
mapping to the introns ( |
strandSpecific |
The description for strand specificity of the RNAseq data. The values are either "unstranded", "stranded", or "reverse".If the reads are not strand specific or directional use "unstranded". If the first read in paired-read sequencing or the reads single-read sequencing is in the same direction as the the transcript strand use "stranded". If the first read in paired-read sequencing or the reads in single-read sequencing is in the oposite direction to the transcript strand use "reverse". |
bpparam |
An optional |
appendLogFile |
Whether log information should be appended to the |
sampleName |
The name of the sample being analyzed. It will be included in the returned
object if |
scaleLength |
A vector constructed of TRUE/FALSE values, same size as the
|
scaleFragment |
A vector constructed of TRUE/FALSE values, same size as the
|
limitRanges |
A GRanges object. If defined it loads sequencing reads that
fall in the defined coordinates. It is similar to |
excludeFusionReads |
Only valid if limitRanges is defined. It filters the
defined by |
loadLimitRangesReads |
Boolean (TRUE/FALSE) variable. If set as
|
... |
Other parameter settings specific to |
If returnObj
is set TRUE
in addition to making result text files,
dependant on whether a single or two method
is defined, the results
would be returned as a single object of class SummarizedExperiment
or as
a list of size 2 which includes 2 objects of class SummarizedExperiment
one for IntRet and the other for ExEx.
Ali Oghabian
interest.sequential
.
# Creating temp directory to store the results
outDir<- file.path(tempdir(),"interestFolder")
dir.create(outDir)
outDir<- normalizePath(outDir)
# Loading suitable bam file
bamF <- system.file("extdata", "small_test_SRR1691637_ZRSR2Mut_RHBDD3.bam",
package="IntEREst", mustWork=TRUE)
# Choosing reference for the gene RHBDD3
ref= u12[u12[,"gene_name"]=="RHBDD3",]
test= interest(
bamFileYieldSize=10000,
bamFile=bamF,
isPaired=TRUE,
isPairedDuplicate=FALSE,
isSingleReadDuplicate=NA,
reference=ref,
referenceGeneNames=ref[,"ens_gene_id"],
referenceIntronExon=ref[,"int_ex"],
repeatsTableToFilter=c(),
outFile=paste(outDir,
"interestRes.tsv", sep="/"),
logFile=paste(outDir,
"log.txt", sep="/"),
method=c("IntRet", "IntSpan"),
strandSpecific="unstranded",
junctionReadsOnly=FALSE,
returnObj=TRUE,
scaleLength= c(TRUE,FALSE),
scaleFragment= c(TRUE,TRUE)
)
test
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.