View source: R/disjoinPafAlignments.R
disjoinPafAlignments | R Documentation |
This function takes loaded PAF alignments using readPaf
function and identify overlapping
genomic ranges either in target or query coordinates and then split PAF alignments at the positions of
these overlaps into a disjoined set of genomic ranges. Alternatively, user defined set of genomic ranges
in parameter 'disjoin.gr' can be used to disjoin PAF alignments.
disjoinPafAlignments(
paf.table,
min.overlap = 1000,
coordinates = "target",
disjoin.gr = NULL
)
paf.table |
A |
min.overlap |
A minimum number of overlapping base pairs to disjoin PAF alignments. |
coordinates |
PAF coordinates to disjoin, either 'target' or 'query'. |
disjoin.gr |
A |
A tibble
of disjoined PAF alignments.
David Porubsky
## Get PAF to plot
paf.file <- system.file("extdata", "test1.paf", package = "SVbyEye")
## Read in PAF
paf.table <- readPaf(paf.file = paf.file, include.paf.tags = TRUE, restrict.paf.tags = "cg")
## Disjoin PAF at target coordinates
disj.paf.table <- disjoinPafAlignments(paf.table = paf.table, coordinates = "target")
## Disjoin PAF at user defined coordinates
disj.gr <- GenomicRanges::GRanges(
seqnames = "query.region",
ranges = IRanges::IRanges(start = 16300000, end = 16350000)
)
disj.paf.table <- disjoinPafAlignments(
paf.table = paf.table,
coordinates = "query", disjoin.gr = disj.gr
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.