readBamFile.Rd
wraper for readGAlignments/readGAlignmentsList to read in bam files.
readBamFile( bamFile, which, tag = character(0), what = c("qname", "flag", "mapq", "isize", "seq", "qual", "mrnm"), flag = scanBamFlag(isSecondaryAlignment = FALSE, isUnmappedQuery = FALSE, isNotPassingQualityControls = FALSE, isSupplementaryAlignment = FALSE), asMates = FALSE, bigFile = FALSE, ... )
bamFile | character(1). Bam file name. |
---|---|
which | A GRanges, IntegerRangesList, or any object that can be coerced to a RangesList, or missing object, from which a IRangesList instance will be constructed. See ScanBamParam. |
tag | A vector of characters indicates the tag names to be read. See ScanBamParam. |
what | A character vector naming the fields to return. Fields are described on the Rsamtools[scanBam] help page. |
flag | An integer(2) vector used to filter reads based on their 'flag' entry. |
asMates | logical(1). Paired ends or not |
bigFile | If the file take too much memory, set it to true to avoid read the reads into memory. scanBamFlag helper function. |
... | parameters used by readGAlignmentsList or readGAlignments |
A GAlignmentsList object when asMates=TRUE, otherwise A GAlignments object. If bigFile is set to TRUE, no reads will be read into memory at this step and empty GAlignments/GAlignmentsList will be returned.
Jianhong Ou
library(BSgenome.Hsapiens.UCSC.hg19) which <- as(seqinfo(Hsapiens)["chr1"], "GRanges") bamfile <- system.file("extdata", "GL1.bam", package="ATACseqQC", mustWork=TRUE) readBamFile(bamfile, which=which, asMates=TRUE)#> GAlignmentsList object of length 44357: #> [[1]] #> GAlignments object with 2 alignments and 7 metadata columns: #> seqnames strand cigar qwidth start end width #> <Rle> <Rle> <character> <integer> <integer> <integer> <integer> #> [1] chr1 + 47M 47 567414 567460 47 #> [2] chr1 - 50M 50 567544 567593 50 #> njunc | qname flag mapq isize #> <integer> | <character> <integer> <integer> <integer> #> [1] 0 | 100 99 30 180 #> [2] 0 | 100 147 30 -180 #> seq qual mrnm #> <DNAStringSet> <PhredQuality> <factor> #> [1] ATTTAGCTGA...AAATGATCTG CCCFFFFFHH...JJJJJJJJJJ chr1 #> [2] GACATCGTAC...TCCACTATGT JJJJJJJHHJ...HHFFFFFCCC chr1 #> ------- #> seqinfo: 25 sequences from an unspecified genome #> #> [[2]] #> GAlignments object with 2 alignments and 7 metadata columns: #> seqnames strand cigar qwidth start end width #> <Rle> <Rle> <character> <integer> <integer> <integer> <integer> #> [1] chr1 - 50M 50 101701613 101701662 50 #> [2] chr1 + 50M 50 101701398 101701447 50 #> njunc | qname flag mapq isize #> <integer> | <character> <integer> <integer> <integer> #> [1] 0 | 100003 83 42 -265 #> [2] 0 | 100003 163 42 265 #> seq qual mrnm #> <DNAStringSet> <PhredQuality> <factor> #> [1] TGAATCTCTC...ACCTCTCCAT FFFFFFGHHH...FHFDEDFCB@ chr1 #> [2] GACCAAGTCC...CTAAGCAAAA @CCFFFFEHG...JJIGJIJJIJ chr1 #> ------- #> seqinfo: 25 sequences from an unspecified genome #> #> [[3]] #> GAlignments object with 2 alignments and 7 metadata columns: #> seqnames strand cigar qwidth start end width #> <Rle> <Rle> <character> <integer> <integer> <integer> <integer> #> [1] chr1 + 50M 50 101701631 101701680 50 #> [2] chr1 - 50M 50 101701766 101701815 50 #> njunc | qname flag mapq isize #> <integer> | <character> <integer> <integer> <integer> #> [1] 0 | 100005 99 42 185 #> [2] 0 | 100005 147 42 -185 #> seq qual mrnm #> <DNAStringSet> <PhredQuality> <factor> #> [1] GCAGCAGCCC...ATCTGGGACA CC@FFFFFHH...JIJJIJJIJE chr1 #> [2] AAAGCCACTG...GGATCTGGGG DHHIIIIJIH...HHFFFFFCCC chr1 #> ------- #> seqinfo: 25 sequences from an unspecified genome #> #> ... #> <44354 more elements>