GappedAlignments-class {GenomicRanges} | R Documentation |
The GappedAlignments class is a simple container which purpose is to store a set of alignments that will hold just enough information for supporting the operations described below.
WARNING! This is work-in-progress. Expect frequent changes in functionalities.
A GappedAlignments object is a vector-like object where each element describes an alignment i.e. how a given sequence (called "query" or "read", typically short) aligns to a reference sequence (typically long).
Most of the time, a GappedAlignments object will be created by loading records from a BAM (or SAM) file and each element in the resulting object will correspond to a record. BAM/SAM records generally contain a lot of information but only part of that information is loaded in the GappedAlignments object. In particular, we discard the query sequences (SEQ field), the query ids (QNAME field), the query qualities (QUAL), the mapping qualities (MAPQ) and any other information that is not needed in order to support the operations or methods described below.
This means that multi-reads (i.e. reads with multiple hits in the reference) won't receive any special treatment i.e. the various SAM/BAM records corresponding to a multi-read will show up in the GappedAlignments object as if they were coming from different/unrelated queries. Also paired-end reads will be treated as single-end reads and the pairing information will be lost.
Each element of a GappedAlignments object consists of:
Note that the last 2 items are not expicitly stored in the GappedAlignments object: they are inferred on-the-fly from the CIGAR and the "start".
The rest of this man page will focus on describing how to:
readGappedAlignments(file, format="BAM", ...)
:
Read a file as a GappedAlignments object.
The function is just a front-end that delegates to a format-specific
back-end function (any extra argument is passed to the back-end
function).
Only the BAM format is supported for now. Its back-end is the
readBamGappedAlignments
function defined
in the Rsamtools package.
See ?readBamGappedAlignments
for
more information (you might need to install and load the package
first).
In the code snippets below, x
is a GappedAlignments object.
length(x)
:
Returns the number of alignments in x
.
rname(x)
:
Returns a character factor of length length(x)
containing the name of the reference sequence for each alignment.
rname(x) <- value
:
Replace the name of the reference sequence for each alignment.
value
must be a character factor/vector, or a 'character' Rle,
or a 'factor' Rle, with the same length as x
.
strand(x)
:
Returns a character factor of length length(x)
(with levels +, - and *) containing the strand in the reference
sequence to which the query is aligned.
cigar(x)
:
Returns a character vector of length length(x)
containing the CIGAR string for each alignment.
qwidth(x)
:
Returns an integer vector of length length(x)
containing the length of the query *after* hard clipping
(i.e. the length of the query sequence that is stored in
the corresponding SAM/BAM record).
grglist(x)
, grg(x)
, rglist(x)
:
Returns a GRangesList, a GRanges or a
CompressedNormalIRangesList object of
length length(x)
where each element represents the regions
in the reference to which a query is aligned. See Details section
above for more information.
start(x)
, end(x)
:
Returns an integer vector of length length(x)
containing the "start" and "end" (respectively) of the query
for each alignment. See Details section above for the exact
definitions of the "start" and "end" of a query.
Note that start(x)
and end(x)
are equivalent
to start(grg(x))
and end(grg(x))
, respectively
(or, alternatively, to min(rglist(x))
and
max(rglist(x))
, respectively).
width(x)
:
Equivalent to width(grg(x))
(or, alternatively, to
end(x) - start(x) + 1L
).
Note that this is generally different from qwidth(x)
except for alignments with a trivial CIGAR string (i.e. a
string of the form "<n>M"
where <n> is a number).
ngap(x)
:
Returns an integer vector of length length(x)
containing the number of gaps for each alignment.
Equivalent to elementLengths(rglist(x)) - 1L
.
In the code snippets below, x
is a GappedAlignments object.
x[i]
:
Returns a new GappedAlignments object made of the selected
alignments. i
can be a numeric or logical vector.
qnarrow(x, start=NA, end=NA, width=NA)
:
x
is a GappedAlignments object.
Returns a new GappedAlignments object of the same length as x
describing how the narrowed query sequences align to the reference.
The start
/end
/width
arguments describe how
to narrow the query sequences. They must be vectors of integers.
NAs and negative values are accepted and "solved" according to the
rules of the SEW (Start/End/Width) interface (see
?solveUserSEW
for the details).
narrow(x, start=NA, end=NA, width=NA)
:
x
is a GappedAlignments object.
Returns a new GappedAlignments object of the same length as x
describing the narrowed alignments. Unlike with qnarrow
now the start
/end
/width
arguments describe
the narrowing on the reference side, not the query side.
Like with qnarrow
, they must be vectors of integers.
NAs and negative values are accepted and "solved" according to the
rules of the SEW (Start/End/Width) interface (see
?solveUserSEW
for the details).
pintersect(x, y)
:
Either x
is a GappedAlignments object and y
is a
GRanges object or x
is a GRanges object and y
is a
GappedAlignments object.
Returns a new GappedAlignments object of the same length as the
GappedAlignments input arguments. Like with narrow
, the
resulting "parallel" intersection is with respect to the reference.
coverage(x, shift=0L, width=NULL, weight=1L)
:
x
is a GappedAlignments object.
Returns a named RleList object with one element
(integer-Rle) per unique reference sequence. Each element represents
x
's coverage of the corresponding reference sequence, that is,
how many times each nucleotide position in the sequence is covered
by the alignments in x
.
Note that the semantic of the coverage
method for
GappedAlignments objects is different from the semantic of the method
for Ranges objects (the latter returns a single
integer-Rle object representing the coverage of all ranges
relatively to a unique imaginary reference sequence).
findOverlaps(query, subject, ...)
,
countOverlaps(query, subject, ...)
,
subsetByOverlaps(query, subject, ...)
,
match(x, table, nomatch=NA_integer_, incomparables=NULL)
,
x %in% table
:
query
or subject
or both are GappedAlignments objects.
findOverlaps(query, subject, ...)
is equivalent to
findOverlaps(grglist(query), subject, ...)
when
query
is a GappedAlignments object, or to
findOverlaps(query, grglist(subject), ...)
when
subject
is a GappedAlignments object, or to
findOverlaps(grglist(query), grglist(subject), ...)
when
both are GappedAlignments objects.
The same apply to countOverlaps(query, subject, ...)
and
subsetByOverlaps(query, subject, ...)
.
See ?`findOverlaps,GRangesList,GenomicRanges-method`
,
?`countOverlaps,GRangesList,GenomicRanges-method`
and
?`subsetByOverlaps,GRangesList,GenomicRanges-method`
for more information (in particular for descriptions of the
extra arguments and the returned object).
H. Pages and P. Aboyoun
http://samtools.sourceforge.net/
readBamGappedAlignments
,
GRangesList-class,
NormalIRanges-class,
CompressedNormalIRangesList-class,
coverage
,
RleList-class,
pintersect,GRanges,GRanges-method,
findOverlaps,GRangesList,GenomicRanges-method,
countOverlaps,GRangesList,GenomicRanges-method,
subsetByOverlaps,GRangesList,GenomicRanges-method,
library(Rsamtools) # the toy file below is there aln1_file <- system.file("extdata", "ex1.bam", package="Rsamtools") aln1 <- readGappedAlignments(aln1_file) aln1 ## --------------------------------------------------------------------- ## A. BASIC MANIPULATION ## --------------------------------------------------------------------- length(aln1) head(aln1) head(rname(aln1)) levels(rname(aln1)) ## Rename the reference sequences: rname(aln1) <- sub("seq", "chr", rname(aln1)) levels(rname(aln1)) head(strand(aln1)) head(cigar(aln1)) head(qwidth(aln1)) table(qwidth(aln1)) grglist(aln1) # a GRangesList object grg(aln1) # a GRanges object rglist(aln1) # a CompressedNormalIRangesList object stopifnot(identical(elementLengths(grglist(aln1)), elementLengths(rglist(aln1)))) head(start(aln1)) head(end(aln1)) head(width(aln1)) head(ngap(aln1)) ## --------------------------------------------------------------------- ## B. SUBSETTING ## --------------------------------------------------------------------- aln1[strand(aln1) == "-"] aln1[grep("I", cigar(aln1), fixed=TRUE)] aln1[grep("N", cigar(aln1), fixed=TRUE)] # no gaps ## A confirmation that all the queries map to the reference with no ## gaps: stopifnot(all(ngap(aln1) == 0)) ## Different ways to subset: aln1[6] # a GappedAlignments object of length 1 grglist(aln1)[[6]] # a GRanges object of length 1 rglist(aln1)[[6]] # a NormalIRanges object of length 1 ## Ds are NOT gaps: ii <- grep("D", cigar(aln1), fixed=TRUE) aln1[ii] ngap(aln1[ii]) grglist(aln1[ii]) ## qwidth() vs width(): aln1[qwidth(aln1) != width(aln1)] ## This MUST return an empty object: aln1[cigar(aln1) == "35M" & qwidth(aln1) != 35] ## but this doesn't have too: aln1[cigar(aln1) != "35M" & qwidth(aln1) == 35] ## --------------------------------------------------------------------- ## C. qnarrow()/narrow() ## --------------------------------------------------------------------- ## Note that there is no difference between qnarrow() and narrow() when ## all the alignments are simple and with no indels. ## This trims 3 nucleotides on the left and 5 nucleotides on the right ## of each alignment: qnarrow(aln1, start=4, end=-6) ## Note that the 'start' and 'end' arguments specify what part of each ## query sequence should be kept (negative values being relative to the ## right end of the query sequence), not what part should be trimmed. ## Trimming on the left doesn't change the "end" of the queries. qnarrow(aln1, start=21) stopifnot(identical(end(qnarrow(aln1, start=21)), end(aln1))) ## --------------------------------------------------------------------- ## D. coverage() ## --------------------------------------------------------------------- coverage(aln1) ## --------------------------------------------------------------------- ## E. findOverlaps()/countOverlaps() ## --------------------------------------------------------------------- findOverlaps(aln1, grg(aln1)[1]) sum(countOverlaps(aln1, grg(aln1)[1])) subsetByOverlaps(aln1, grg(aln1)[1]) table(match(aln1, grg(aln1)[1]), useNA = "ifany") table(aln1 %in% grg(aln1)[1])