readBamFile {ATACseqQC}R Documentation

read in bam files

Description

wraper for readGAlignments/readGAlignmentsList to read in bam files.

Usage

readBamFile(bamFile, which, tag = character(0), what = c("qname",
  "flag", "mapq", "isize", "seq", "qual", "mrnm"),
  flag = scanBamFlag(isSecondaryAlignment = FALSE, isUnmappedQuery =
  FALSE, isNotPassingQualityControls = FALSE), asMates = FALSE,
  bigFile = FALSE, ...)

Arguments

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. This is most easily created with the

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

Value

A GAlignmentsList object when asMats=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.

Author(s)

Jianhong Ou

Examples

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)

[Package ATACseqQC version 1.8.1 Index]