| BamFile {Rsamtools} | R Documentation |
Use BamFile() to create a reference to a BAM file (and
optionally its index). The reference remains open across calls to
methods, avoiding costly index re-loading.
BamFileList() provides a convenient way of managing a list of
BamFile instances.
## Constructors
BamFile(file, index=file)
BamFileList(...)
## Opening / closing
## S3 method for class 'BamFile'
open(con, ...)
## S3 method for class 'BamFile'
close(con, ...)
## accessors; also path(), index()
## S4 method for signature 'BamFile'
isOpen(con, rw="")
## actions
## S4 method for signature 'BamFile'
scanBamHeader(files, ...)
## S4 method for signature 'BamFile'
seqinfo(x)
## S4 method for signature 'BamFile'
scanBam(file, index=file, ..., param=ScanBamParam(what=scanBamWhat()))
## S4 method for signature 'BamFile'
countBam(file, index=file, ..., param=ScanBamParam())
## S4 method for signature 'BamFileList'
countBam(file, index=file, ..., param=ScanBamParam())
## S4 method for signature 'BamFile'
filterBam(file, destination, index=file, ...,
indexDestination=TRUE, param=ScanBamParam(what=scanBamWhat()))
## S4 method for signature 'BamFile'
indexBam(files, ...)
## S4 method for signature 'BamFile'
sortBam(file, destination, ..., byQname=FALSE, maxMemory=512)
## S4 method for signature 'BamFile'
readBamGappedAlignments(file, index=file, use.names=FALSE, param=NULL)
## S4 method for signature 'BamFile'
readBamGappedReads(file, index=file, use.names=FALSE, param=NULL)
## counting
## S4 method for signature 'GRanges,BamFileList'
summarizeOverlaps(
features, reads, mode, ignore.strand = FALSE, ..., param = ScanBamParam())
... |
Additional arguments. For |
con |
An instance of |
x, file, files |
A character vector of BAM file paths (for
|
index |
A character vector of indices (for |
destination |
character(1) file path to write filtered reads to. |
indexDestination |
logical(1) indicating whether the destination file should also be indexed. |
byQname, maxMemory |
See |
param |
An optional |
use.names |
Construct the names of the returned object from the query template names (QNAME field)? If not (the default), then the returned object has no names. |
rw |
Mode of file; ignored. |
reads |
A BamFileList that represents the data to be counted by summarizeOverlaps. |
features |
A GRanges or a GRangesList object of genomic regions of interest. When a GRanges is supplied, each row is considered a feature. When a GRangesList is supplied, each higher list-level is considered a feature. This distinction is important when defining an overlap between a read and a feature. See examples for details. |
mode |
A function that defines the method to be used when a read overlaps more than one feature. Pre-defined options are "Union", "IntersectionStrict", or "IntersectionNotEmpty" and are designed after the counting modes available in the HTSeq package by Simon Anders (see references).
|
ignore.strand |
A logical value indicating if strand should be considered when matching. |
Objects are created by calls of the form BamFile().
The BamFile class inherits fields from the
RsamtoolsFile class.
BamFileList inherits methods from
RsamtoolsFileList and SimpleList.
Opening / closing:
Opens the (local or remote) path and
index (if bamIndex is not character(0)),
files. Returns a BamFile instance.
Closes the BamFile con; returning
(invisibly) the updated BamFile. The instance may be
re-opened with open.BamFile.
Accessors:
Returns a character(1) vector of BAM path names.
Returns a character(1) vector of BAM index path names.
Methods:
Visit the path in path(file), returning
the information contained in the file header; see
scanBamHeader.
Visit the path in path(file), returning
a Seqinfo instance containing information on
the lengths of each sequence.
Visit the path in path(file), returning the
result of scanBam applied to the specified path.
Visit the path(s) in path(file), returning
the result of countBam applied to the specified
path.
Visit the path in path(file), returning
the result of filterBam applied to the specified
path.
Visit the path in path(file), returning
the result of indexBam applied to the specified
path.
Visit the path in path(file), returning the
result of sortBam applied to the specified path.
Visit
the path in path(file), returning the result of
readBamGappedAlignments or readBamGappedReads
applied to the specified path.
See readBamGappedAlignments.
Compactly display the object.
Martin Morgan and Marc Carlson
The GenomicRanges package is where the summarizeOverlaps
method originates.
fl <- system.file("extdata", "ex1.bam", package="Rsamtools")
bf <- open(BamFile(fl)) # implicit index
bf
identical(scanBam(bf), scanBam(fl))
rng <- GRanges(c("seq1", "seq2"), IRanges(1, c(1575, 1584)))
## repeatedly visit 'bf'
sapply(seq_len(length(rng)), function(i, bamFile, rng) {
param <- ScanBamParam(which=rng[i], what="seq")
bam <- scanBam(bamFile, param=param)[[1]]
alphabetFrequency(bam[["seq"]], baseOnly=TRUE, collapse=TRUE)
}, bf, rng)
##---------------------------------------------------------------------------##
## How to use summarizeOverlaps with a BamFileList object.
fls = list.files(system.file("extdata",package="GenomicRanges"),
recursive=TRUE, pattern="*bam$", full=TRUE)
bfs <- BamFileList(fls)
## "features" will be the argument for an "annotations" object (GRanges
## or GRangesList object.
group_id <- c("A", "B", "C", "C", "D", "D", "E", "F", "G", "H", "H")
features <- GRanges(
seqnames = Rle(c("chr2L", "chr2R", "chr2L", "chr2R", "chr2L", "chr2R",
"chr2L", "chr2R", "chr2R", "chr3L", "chr3L")),
strand = strand(rep("+", length(group_id))),
ranges = IRanges(
start=c(1000, 2000, 3000, 3600, 7000, 7500, 4000, 4000, 3000, 5000, 5400),
width=c(500, 900, 500, 300, 600, 300, 500, 900, 500, 500, 500)),
DataFrame(group_id)
)
## Then call the method:
summarizeOverlaps(features, bfs, mode = Union, ignore.strand=TRUE)