Skip to contents

Reduce hit coverage for a seqname

Usage

reduceHitCoverage(x, min.match = 0.9)

Arguments

x

AlignmentPairs object

min.match

filter out hits with fraction of matching bases less than min.match

Value

reduced and filtered GRanges object

Details

Given an AlignmentPairs object, reduce overlapping query regions. For instance, if one hit spans coordinates 100-400 and another 200-500, we merge to 100-500. The function calculates how many overlaps that are merged, the rationale being that if there are >1 overlaps for a seqname, the sequence is associated with multiple distinct regions in the subject sequence. Thus, the unit of interest here is the query, and the function seeks to identify reduced disjoint regions of a query sequence that map to a subject.

The function returns a GRanges object with reduced regions,
i.e. overlaps are merged. Information about matches and
mismatches is currently dropped as it usually is not possible
to infer where mismatches occur. Instead, the data column
`coverage` simply holds the ratio of the width of the region
to the transcript length. Summing up coverages from disjoint
regions then gives total coverage of the transcript.

The cutoff is used to filter regions based on the ratio of
matches to the width of the region.

Since the association between query and subject regions is
removed, the return value is a GRanges object consisting of
the reduced query ranges with a revmap and coverage attribute.

Examples

ranges <- IRanges::IRanges(
          start=c(100, 200, 700),
          end=c(400, 500, 1000)
)
qry <- GenomicRanges::GRanges(
          ranges=ranges,
          seqnames=c("t1"),
          seqinfo=GenomeInfoDb::Seqinfo(seqnames=c("t1"),
                                        seqlengths=c(1050))
)
ranges <- IRanges::IRanges(
          start=c(1000, 1000, 1000),
          end=c(1300, 1300, 1300)
)
sbj <- GenomicRanges::GRanges(ranges=ranges, seqnames=c("c1", "c2", "c1"))
x <- AlignmentPairs(query=qry, subject=sbj, matches=c(300, 290, 280))
if (FALSE) gr <- reduceHitCoverage(x, 0.1)