find_overlaps | R Documentation |
Find overlap between two Ranges
find_overlaps(x, y, maxgap, minoverlap, suffix = c(".x", ".y"))
## S3 method for class 'IntegerRanges'
find_overlaps(x, y, maxgap = -1L, minoverlap = 0L, suffix = c(".x", ".y"))
## S3 method for class 'GenomicRanges'
find_overlaps(x, y, maxgap = -1L, minoverlap = 0L, suffix = c(".x", ".y"))
find_overlaps_within(x, y, maxgap, minoverlap, suffix = c(".x", ".y"))
## S3 method for class 'IntegerRanges'
find_overlaps_within(
x,
y,
maxgap = -1L,
minoverlap = 0L,
suffix = c(".x", ".y")
)
## S3 method for class 'GenomicRanges'
find_overlaps_within(
x,
y,
maxgap = -1L,
minoverlap = 0L,
suffix = c(".x", ".y")
)
find_overlaps_directed(x, y, maxgap, minoverlap, suffix = c(".x", ".y"))
## S3 method for class 'GenomicRanges'
find_overlaps_directed(
x,
y,
maxgap = -1L,
minoverlap = 0L,
suffix = c(".x", ".y")
)
find_overlaps_within_directed(x, y, maxgap, minoverlap, suffix = c(".x", ".y"))
## S3 method for class 'GenomicRanges'
find_overlaps_within_directed(x, y, maxgap, minoverlap, suffix = c(".x", ".y"))
group_by_overlaps(x, y, maxgap, minoverlap)
## S3 method for class 'IntegerRanges'
group_by_overlaps(x, y, maxgap = -1L, minoverlap = 0L)
## S3 method for class 'GenomicRanges'
group_by_overlaps(x, y, maxgap = -1L, minoverlap = 0L)
x , y |
Objects representing ranges |
maxgap , minoverlap |
The maximimum gap between intervals as an integer greater than or equal to negative one. The minimum amount of overlap between intervals as an integer greater than zero, accounting for the maximum gap. |
suffix |
A character vector of length two used to identify metadata columns coming from x and y. |
find_overlaps()
will search for any overlaps between ranges
x and y and return a Ranges object of length equal to the number of times x
overlaps y. This Ranges object will have additional metadata columns
corresponding to the metadata columns in y. find_overlaps_within()
is
the same but will only search for overlaps within y. For GRanges objects strand is
ignored, unless find_overlaps_directed()
is used. If the Ranges objects have no
metadata, one could use group_by_overlaps()
to be able to
identify the index of the input Range x that overlaps a Range in y.
Alternatively,
pair_overlaps()
could be used to place the x ranges next to the range
in y they overlap.
A Ranges object with rows corresponding to the
ranges in x that overlap y. In the case of group_by_overlaps()
, returns
a GroupedRanges object, grouped by the number of overlaps
of ranges in x that overlap y (stored in a column called query).
IRanges::findOverlaps()
,
GenomicRanges::findOverlaps()
query <- data.frame(start = c(5,10, 15,20), width = 5, gc = runif(4)) %>%
as_iranges()
subject <- data.frame(start = 2:6, width = 3:7, label = letters[1:5]) %>%
as_iranges()
find_overlaps(query, subject)
find_overlaps(query, subject, minoverlap = 5)
find_overlaps_within(query, subject) # same result as minoverlap
find_overlaps(query, subject, maxgap = 1)
# -- GRanges objects, strand is ignored by default
query <- data.frame(seqnames = "chr1",
start = c(11,101),
end = c(21, 200),
name = c("a1", "a2"),
strand = c("+", "-"),
score = c(1,2)) %>%
as_granges()
subject <- data.frame(seqnames = "chr1",
strand = c("+", "-", "+", "-"),
start = c(21,91,101,201),
end = c(30,101,110,210),
name = paste0("b", 1:4),
score = 1:4) %>%
as_granges()
# ignores strandedness
find_overlaps(query, subject, suffix = c(".query", ".subject"))
find_overlaps(query, subject, suffix = c(".query", ".subject"), minoverlap = 2)
# adding directed prefix includes strand
find_overlaps_directed(query, subject, suffix = c(".query", ".subject"))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.