Skip to content

Commit b8a3cd5

Browse files
committed
fix poverlaps()
1 parent 8b97cb2 commit b8a3cd5

File tree

2 files changed

+15
-1
lines changed

2 files changed

+15
-1
lines changed

R/findOverlaps-methods.R

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -258,5 +258,5 @@ setMethod("poverlaps", c("GenomicRanges", "GenomicRanges"),
258258
seqnames(query) == seqnames(subject) &
259259
(if (ignore.strand) TRUE
260260
else compatibleStrand(strand(query), strand(subject))) &
261-
poverlaps(ranges(query), ranges(subject), maxgap, minoverlaps, type)
261+
poverlaps(ranges(query), ranges(subject), maxgap, minoverlap, type)
262262
})

inst/unitTests/test_findOverlaps-methods.R

Lines changed: 14 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -318,3 +318,17 @@ test_findOverlaps_with_circular_sequences <- function()
318318
.checkHits(1:4, 1:4, 4, 4, current5, select="all")
319319
}
320320

321+
test_poverlaps <- function() {
322+
ans <- poverlaps(GRanges(), GRanges())
323+
checkIdentical(ans, Rle())
324+
325+
ans <- poverlaps(GRanges("chr1:11-15"), GRanges("chr1:16-20"))
326+
checkIdentical(ans, Rle(FALSE))
327+
328+
ans <- poverlaps(GRanges("chr1:11-16"), GRanges("chr1:16-20"))
329+
checkIdentical(ans, Rle(TRUE))
330+
331+
ans <- poverlaps(GRanges(c("chr1:11-15", "chr1:11-16")),
332+
GRanges("chr1:16-20"))
333+
checkIdentical(ans, Rle(c(FALSE, TRUE)))
334+
}

0 commit comments

Comments
 (0)