findOverlaps-methods | R Documentation |
Various methods for finding/counting overlaps between objects containing genomic ranges. This man page describes the methods that operate on GenomicRanges and GRangesList objects.
NOTE: The findOverlaps
generic function
and methods for IntegerRanges and
IntegerRangesList objects are defined and
documented in the IRanges package.
The methods for GAlignments,
GAlignmentPairs, and
GAlignmentsList objects are defined and
documented in the GenomicAlignments package.
GenomicRanges and GRangesList objects also support
countOverlaps
, overlapsAny
, and subsetByOverlaps
thanks to the default methods defined in the IRanges package and
to the findOverlaps
and countOverlaps
methods defined in
this package and documented below.
## S4 method for signature 'GenomicRanges,GenomicRanges'
findOverlaps(query, subject,
maxgap=-1L, minoverlap=0L,
type=c("any", "start", "end", "within", "equal"),
select=c("all", "first", "last", "arbitrary"),
ignore.strand=FALSE)
## S4 method for signature 'GRangesList,GenomicRanges'
findOverlaps(query, subject,
maxgap=-1L, minoverlap=0L,
type=c("any", "start", "end", "within", "equal"),
select=c("all", "first", "last", "arbitrary"),
ignore.strand=FALSE)
## S4 method for signature 'GenomicRanges,GRangesList'
findOverlaps(query, subject,
maxgap=-1L, minoverlap=0L,
type=c("any", "start", "end", "within", "equal"),
select=c("all", "first", "last", "arbitrary"),
ignore.strand=FALSE)
## S4 method for signature 'GRangesList,GRangesList'
findOverlaps(query, subject,
maxgap=-1L, minoverlap=0L,
type=c("any", "start", "end", "within", "equal"),
select=c("all", "first", "last", "arbitrary"),
ignore.strand=FALSE)
## S4 method for signature 'GenomicRanges,GenomicRanges'
countOverlaps(query, subject,
maxgap=-1L, minoverlap=0L,
type=c("any", "start", "end", "within", "equal"),
ignore.strand=FALSE)
query , subject |
A GRanges or GRangesList object. |
maxgap , minoverlap , type |
See IMPORTANT NOTE about how |
select |
When |
ignore.strand |
When set to |
When the query and the subject are GRanges or
GRangesList objects, findOverlaps
uses the triplet
(sequence name, range, strand) to determine which features (see
paragraph below for the definition of feature) from the query
overlap which features in the subject
, where a strand value
of "*"
is treated as occurring on both the "+"
and
"-"
strand.
An overlap is recorded when a feature in the query
and a feature
in the subject
have the same sequence name, have a compatible
pairing of strands (e.g. "+"
/"+"
, "-"
/"-"
,
"*"
/"+"
, "*"
/"-"
, etc.), and satisfy the
interval overlap requirements.
In the context of findOverlaps
, a feature is a collection of
ranges that are treated as a single entity. For GRanges objects,
a feature is a single range; while for GRangesList objects,
a feature is a list element containing a set of ranges. In the results,
the features are referred to by number, which run from 1 to
length(query)
/length(subject)
.
For type="equal"
with GRangesList objects, query[[i]]
matches subject[[j]]
iff for each range in query[[i]]
there is an identical range in subject[[j]]
, and vice versa.
For findOverlaps
: either a Hits object when
select="all"
or an integer vector otherwise.
For countOverlaps
: an integer vector containing the tabulated
query overlap hits.
P. Aboyoun, S. Falcon, M. Lawrence, and H. Pagès
The Hits class in the S4Vectors package for representing a set of hits between 2 vector-like objects.
The findOverlaps
generic function defined
in the IRanges package.
The GNCList constructor and class for preprocessing and representing a GenomicRanges or object as a data structure based on Nested Containment Lists.
The GRanges and GRangesList classes.
## ---------------------------------------------------------------------
## BASIC EXAMPLES
## ---------------------------------------------------------------------
## GRanges object:
gr <- GRanges(
seqnames=Rle(c("chr1", "chr2", "chr1", "chr3"), c(1, 3, 2, 4)),
ranges=IRanges(1:10, width=10:1, names=head(letters,10)),
strand=Rle(strand(c("-", "+", "*", "+", "-")), c(1, 2, 2, 3, 2)),
score=1:10,
GC=seq(1, 0, length=10)
)
gr
## GRangesList object:
gr1 <- GRanges(seqnames="chr2", ranges=IRanges(4:3, 6),
strand="+", score=5:4, GC=0.45)
gr2 <- GRanges(seqnames=c("chr1", "chr1"),
ranges=IRanges(c(7,13), width=3),
strand=c("+", "-"), score=3:4, GC=c(0.3, 0.5))
gr3 <- GRanges(seqnames=c("chr1", "chr2"),
ranges=IRanges(c(1, 4), c(3, 9)),
strand=c("-", "-"), score=c(6L, 2L), GC=c(0.4, 0.1))
grl <- GRangesList("gr1"=gr1, "gr2"=gr2, "gr3"=gr3)
## Overlapping two GRanges objects:
table(!is.na(findOverlaps(gr, gr1, select="arbitrary")))
countOverlaps(gr, gr1)
findOverlaps(gr, gr1)
subsetByOverlaps(gr, gr1)
countOverlaps(gr, gr1, type="start")
findOverlaps(gr, gr1, type="start")
subsetByOverlaps(gr, gr1, type="start")
findOverlaps(gr, gr1, select="first")
findOverlaps(gr, gr1, select="last")
findOverlaps(gr1, gr)
findOverlaps(gr1, gr, type="start")
findOverlaps(gr1, gr, type="within")
findOverlaps(gr1, gr, type="equal")
## ---------------------------------------------------------------------
## MORE EXAMPLES
## ---------------------------------------------------------------------
table(!is.na(findOverlaps(gr, gr1, select="arbitrary")))
countOverlaps(gr, gr1)
findOverlaps(gr, gr1)
subsetByOverlaps(gr, gr1)
## Overlaps between a GRanges and a GRangesList object:
table(!is.na(findOverlaps(grl, gr, select="first")))
countOverlaps(grl, gr)
findOverlaps(grl, gr)
subsetByOverlaps(grl, gr)
countOverlaps(grl, gr, type="start")
findOverlaps(grl, gr, type="start")
subsetByOverlaps(grl, gr, type="start")
findOverlaps(grl, gr, select="first")
table(!is.na(findOverlaps(grl, gr1, select="first")))
countOverlaps(grl, gr1)
findOverlaps(grl, gr1)
subsetByOverlaps(grl, gr1)
countOverlaps(grl, gr1, type="start")
findOverlaps(grl, gr1, type="start")
subsetByOverlaps(grl, gr1, type="start")
findOverlaps(grl, gr1, select="first")
## Overlaps between two GRangesList objects:
countOverlaps(grl, rev(grl))
findOverlaps(grl, rev(grl))
subsetByOverlaps(grl, rev(grl))
## ---------------------------------------------------------------------
## INTERPRETATION OF 'minoverlap' WHEN 'query' OR 'subject' IS A
## GRangesList OBJECT
## ---------------------------------------------------------------------
gr1 <- GRanges("chr5:1-26")
gr2 <- GRanges("chr5:31-40")
gr3 <- c(GRanges("chr5:11-20"), gr2)
grl123 <- GRangesList(gr1, gr2, gr3)
grl123
query <- GRanges("chr5:17-35")
findOverlaps(query, grl123[[1]], minoverlap=8) # 1 hit
findOverlaps(query, grl123[[2]], minoverlap=8) # no hit
findOverlaps(query, grl123[[3]], minoverlap=8) # no hit
## Using GRangesList object 'grl123' as the subject:
findOverlaps(query, grl123, minoverlap=8)
## As we can see, a hit is reported with the 3rd element in the subject.
## That's because the total number of positions in this overlap is 9:
## - positions 17 to 20 in the first range of grl123[[3]], so 4 positions
## - positions 31 to 35 in its second range, so 5 positions
## Sanity checks:
hits <- findOverlaps(query, grl123[[1]], minoverlap=8)
stopifnot(length(hits) == 1)
hits <- findOverlaps(query, grl123[[2]], minoverlap=8)
stopifnot(length(hits) == 0)
hits <- findOverlaps(query, grl123[[3]], minoverlap=8)
stopifnot(length(hits) == 0)
hits <- findOverlaps(query, grl123, minoverlap=8)
stopifnot(identical(subjectHits(hits), c(1L, 3L)))
hits <- findOverlaps(query, grl123, minoverlap=9)
stopifnot(identical(subjectHits(hits), c(1L, 3L)))
hits <- findOverlaps(query, grl123, minoverlap=10)
stopifnot(identical(subjectHits(hits), 1L))
hits <- findOverlaps(query, grl123, minoverlap=11)
stopifnot(length(hits) == 0)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.