findBreakpointOverlaps: Finding overlapping breakpoints between two breakpoint sets

Description Usage Arguments Details Value Examples

View source: R/BreakpointGRanges.R

Description

Finding overlapping breakpoints between two breakpoint sets

Usage

1
2
3
4
5
6
7
8
9
findBreakpointOverlaps(
  query,
  subject,
  maxgap = -1L,
  minoverlap = 0L,
  ignore.strand = FALSE,
  sizemargin = NULL,
  restrictMarginToSizeMultiple = NULL
)

Arguments

query, subject

Both of the input objects should be GRanges objects. Unlike findOverlaps(), subject cannot be ommitted. Each breakpoint must be accompanied with a partner breakend, which is also in the GRanges, with the partner's id recorded in the partner field. See GenomicRanges::findOverlaps-methods for details.

maxgap, minoverlap

Valid overlapping thresholds of a maximum gap and a minimum overlapping positions between breakend intervals. Both should be scalar integers. maxgap allows non-negative values, and minoverlap allows positive values. See GenomicRanges::findOverlaps-methods for details.

ignore.strand

Default value is FALSE. strand information is ignored when set to TRUE. See GenomicRanges::findOverlaps-methods for details.

sizemargin

Error margin in allowable size to prevent matching of events of different sizes, e.g. a 200bp event matching a 1bp event when maxgap is set to 200.

restrictMarginToSizeMultiple

Size restriction multiplier on event size. The default value of 0.5 requires that the breakpoint positions can be off by at maximum, half the event size. This ensures that small deletion do actually overlap at least one base pair.

Details

findBreakpointOverlaps() is an efficient adaptation of findOverlaps-methods() for breakend ranges. It searches for overlaps between breakpoint objects, and return a matrix including index of overlapping ranges as well as error stats. All breakends must have their partner breakend included in the partner field. A valid overlap requires that breakends on boths sides meets the overlapping requirements.

See GenomicRanges::findOverlaps-methods for details of overlap calculation.

Value

A dataframe containing index and error stats of overlapping breakpoints.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
#reading in VCF files
query.file <- system.file("extdata", "gridss-na12878.vcf", package = "StructuralVariantAnnotation")
subject.file <- system.file("extdata", "gridss.vcf", package = "StructuralVariantAnnotation")
query.vcf <- VariantAnnotation::readVcf(query.file, "hg19")
subject.vcf <- VariantAnnotation::readVcf(subject.file, "hg19")
#parsing vcfs to GRanges objects
query.gr <- breakpointRanges(query.vcf)
subject.gr <- breakpointRanges(subject.vcf)
#find overlapping breakpoint intervals
findBreakpointOverlaps(query.gr, subject.gr)
findBreakpointOverlaps(query.gr, subject.gr, ignore.strand=TRUE)
findBreakpointOverlaps(query.gr, subject.gr, maxgap=100, sizemargin=0.5)

StructuralVariantAnnotation documentation built on Nov. 8, 2020, 5:43 p.m.