splitBam.Rd
shift the bam files by 5'ends and split the bam files.
splitBam( bamfile, tags, index = bamfile, outPath = NULL, txs, genome, conservation, positive = 4L, negative = 5L, breaks = c(0, 100, 180, 247, 315, 473, 558, 615, Inf), labels = c("NucleosomeFree", "inter1", "mononucleosome", "inter2", "dinucleosome", "inter3", "trinucleosome", "others"), seqlev = paste0("chr", c(1:22, "X", "Y")), cutoff = 0.8, flag = scanBamFlag(isSecondaryAlignment = FALSE, isUnmappedQuery = FALSE, isNotPassingQualityControls = FALSE, isSupplementaryAlignment = FALSE) )
bamfile | character(1). File name of bam. |
---|---|
tags | A vector of characters indicates the tags in bam file. |
index | The names of the index file of the 'BAM' file being processed; This is given without the '.bai' extension. |
outPath | Output file path. |
txs | GRanges of transcripts. |
genome | An object of BSgenome |
conservation | An object of GScores. |
positive | integer(1). the size to be shift for positive strand |
negative | integer(1). the size to be shift for negative strand |
breaks | A numeric vector for fragment size of nucleosome free, mononucleosome, dinucleosome and trinucleosome |
labels | A vector of characters indicates the labels for the levels of the resulting category. The length of labels = length of breaks - 1 |
seqlev | A vector of characters indicates the sequence levels. |
cutoff | numeric(1). Cutoff value for prediction by randomForest. |
flag | An integer(2) vector used to filter reads based on their 'flag' entry. |
an invisible list of GAlignments
Jianhong Ou
if(Sys.getenv("USER")=="jianhongou"){ bamfile <- system.file("extdata", "GL1.bam", package="ATACseqQC") tags <- c("AS", "XN", "XM", "XO", "XG", "NM", "MD", "YS", "YT") library(BSgenome.Hsapiens.UCSC.hg19) library(TxDb.Hsapiens.UCSC.hg19.knownGene) txs <- transcripts(TxDb.Hsapiens.UCSC.hg19.knownGene) library(phastCons100way.UCSC.hg19) objs <- splitBam(bamfile, tags, txs=txs, genome=Hsapiens, conservation=phastCons100way.UCSC.hg19, seqlev="chr1") }