findOverlaps: Find overlaps between GenomicInteractions

Description Usage Arguments Value One-dimensional overlaps Two-dimensional overlaps Author(s) Examples

Description

Find one-dimensional overlaps between GenomicInteractions and other genomic intervals, or two-dimensional overlaps between two GenomicInteractions objects.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
## S4 method for signature 'GenomicInteractions,GenomicRanges'
findOverlaps(query, subject,
  maxgap = -1L, minoverlap = 0L, type = c("any", "start", "end",
  "within", "equal"), select = c("all", "first", "last", "arbitrary"),
  ignore.strand = TRUE, ..., use.region = c("any", "first", "second",
  "both"))

## S4 method for signature 'GenomicRanges,GenomicInteractions'
findOverlaps(query, subject,
  maxgap = -1L, minoverlap = 0L, type = c("any", "start", "end",
  "within", "equal"), select = c("all", "first", "last", "arbitrary"),
  ignore.strand = TRUE, ..., use.region = c("any", "first", "second",
  "both"))

## S4 method for signature 'GenomicInteractions,GenomicInteractions'
findOverlaps(query,
  subject, maxgap = -1L, minoverlap = 0L, type = c("any", "start",
  "end", "within", "equal"), select = c("all", "first", "last",
  "arbitrary"), ignore.strand = TRUE, ..., use.region = c("any",
  "match", "reverse"))

Arguments

query

A GenomicInteractions instance or any Vector-like object that can be overlapped with a GenomicRanges via findOverlaps.

subject

Same as query.

maxgap

Integer scalar specifying the maximum gap between two genomic intervals to consider them as overlapping, see ?findOverlaps for details.

minoverlap

Integer scalar specifying the minimum overlap between two genomic intervals to consider them as overlapping, see ?findOverlaps for details.

type

String specifying the type of overlap, see ?findOverlaps for details.

select

String specifying what kind of output to return, see ?findOverlaps for details.

ignore.strand

Logical scalar indicating whether strand information should be ignored. If FALSE, stranded intervals must be on the same strand to be considered as overlapping.

...

Further arguments to pass to findOverlaps.

use.region

String specifying which regions should be used to perform the overlap, see below.

Value

If select="all", a Hits object is returned specifying the overlaps between query and subject.

Otherwise, an integer vector is returned of length equal to length(query), containing the selected index of subject that is overlapped by each entry of query (or NA, if there are no overlaps).

One-dimensional overlaps

If only one of query or subject is an GenomicInteractions object, a one-dimensional overlap is performed. This involves identifying overlaps between the individual anchor regions without consideration of the pairing of anchor regions in each interaction. One-dimensional overlaps are useful for identifying any kind of overlap between the interactions' anchor regions and the query/subject of interest.

Let's say that query is the GenomoicInteractions object, in which case:

The same principles apply when subject is the GenomicInteractions object.

Overlaps between genomic regions are defined based on the various parameters passed to findOverlaps, e.g., maxgap, minoverlap. If query is the GenomicInteractions object, its anchor regions will also be the query in the findOverlaps call. Conversely, if subject is the GenomicInteractions, the anchor regions will be the subject. This has implications for overlap settings that are not symmetric, e.g., type="within".

By default, all overlaps with GenomicInteractions objects are performed with ignore.strand=TRUE. This reflects the fact that interactions generally occur between unstranded genomic loci rather than stranded transcripts.

Two-dimensional overlaps

If both query and subject are GenomicInteractions objects, a two-dimensional overlap can be performed. An interation in query only overlaps an interaction in subject if both of the query's anchor regions overlaps both of the subject's anchor regions. This is useful for identifying interactions that span the same part of the two-dimensional interaction space.

The exact nature of the overlap can be controlled with use.region:

Again, overlaps between genomic regions are defined based on the various parameters passed to findOverlaps, e.g., maxgap, minoverlap. The query's anchor regions will be used as the query in the findOverlaps call, which ensures that asymmetric modes like type="within" are correctly handled.

Author(s)

Aaron Lun

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
######################
## One-dimensional ###
######################

anchor1 <- GRanges(c("chr1", "chr1", "chr1", "chr1"), 
    IRanges(c(10, 20, 30, 20), width=5))
anchor2 <- GRanges(c("chr1", "chr1", "chr1", "chr2"), 
    IRanges(c(100, 200, 300, 50), width=5))
test <- GenomicInteractions(anchor1, anchor2)
test

findOverlaps(test, GRanges("chr1:25-100"))
findOverlaps(test, GRanges("chr1:25-100"), use.region="first")
findOverlaps(test, GRanges("chr1:25-100"), use.region="second")

findOverlaps(GRanges("chr1:25-100"), test)
findOverlaps(GRanges("chr1:25-100"), test, use.region="first")
findOverlaps(GRanges("chr1:25-100"), test, use.region="second")

######################
## Two-dimensional ###
######################

alt <- GenomicInteractions(
   GRanges("chr1:1-20"), GRanges("chr1:100-150")
)
findOverlaps(test, alt)
findOverlaps(test, alt, use.region="match")
findOverlaps(test, alt, use.region="reverse")

alt2 <- GenomicInteractions(
   GRanges("chr2:1-60"), GRanges("chr1:1-30")
)
suppressWarnings(findOverlaps(alt2, test))
suppressWarnings(findOverlaps(alt2, test, use.region="match"))
findOverlaps(alt2, test, use.region="reverse")

LTLA/GenomicInteractions documentation built on June 24, 2019, 2:09 p.m.