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 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71
|
% Generated by roxygen2: do not edit by hand
% Please edit documentation in R/BreakpointGRanges.R
\name{findBreakpointOverlaps}
\alias{findBreakpointOverlaps}
\title{Finding overlapping breakpoints between two breakpoint sets}
\usage{
findBreakpointOverlaps(
query,
subject,
maxgap = -1L,
minoverlap = 0L,
ignore.strand = FALSE,
sizemargin = NULL,
restrictMarginToSizeMultiple = NULL
)
}
\arguments{
\item{query, subject}{Both of the input objects should be GRanges objects.
Unlike \code{findOverlaps()}, \code{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 \code{partner} field.
See GenomicRanges::findOverlaps-methods for details.}
\item{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.}
\item{ignore.strand}{Default value is FALSE. strand information is ignored when set to
TRUE.
See GenomicRanges::findOverlaps-methods for details.}
\item{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.}
\item{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.}
}
\value{
A dataframe containing index and error stats of overlapping breakpoints.
}
\description{
Finding overlapping breakpoints between two breakpoint sets
}
\details{
\code{findBreakpointOverlaps()} is an efficient adaptation of \code{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 \code{partner}
field. A valid overlap requires that breakends on boths sides meets the overlapping
requirements.
See GenomicRanges::findOverlaps-methods for details of overlap calculation.
}
\examples{
#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)
}
|