plotAVA | R Documentation |
This function takes PAF output file from minimap2 reporting all-versus-all alignments of multiple FASTA sequences and visualize the alignments in a miropeat style.
plotAVA(
paf.table,
seqnames.order = NULL,
min.deletion.size = NULL,
min.insertion.size = NULL,
highlight.sv = NULL,
binsize = NULL,
color.by = "direction",
perc.identity.breaks = c(90, 95, 99, 99.5, 99.6, 99.7, 99.8, 99.9),
color.palette = NULL,
outline.alignments = FALSE
)
paf.table |
A |
seqnames.order |
A user defined order sequence names to be plotted from top to the bottom. |
min.deletion.size |
A minimum size (in base pairs) of a deletion to be retained. |
min.insertion.size |
A minimum size (in base pairs) of an insertion to be retained. |
highlight.sv |
Visualize alignment embedded structural variation either as an outlined ('outline') or filled ('fill') miropeats. |
binsize |
A size of a bin in base pairs to split a PAF alignment into. |
color.by |
Color alignments either by directionality ('direction'), fraction of matched base pairs ('identity'), or a custom column name present in submitted 'paf.table'. |
perc.identity.breaks |
Set of percentage values to categorize alignment percent identity into a set of discrete colors. |
color.palette |
A discrete color palette defined as named character vector (elements = colors, names = discrete levels) to color alignment directionality, ‘[default: color.palette <- c(’-' = 'cornflowerblue', '+' = 'forestgreen')]'. |
outline.alignments |
Set to |
A list
of miropeat style plots.
David Porubsky
## Get PAF to plot
paf.file <- system.file("extdata", "test_ava.paf", package = "SVbyEye")
## Read in PAF
paf.table <- readPaf(paf.file = paf.file, include.paf.tags = TRUE, restrict.paf.tags = "cg")
## Make a plot
## Color by alignment directionality
plotAVA(paf.table = paf.table, color.by = "direction")
## Color by fraction of matched bases in each alignment
plotAVA(paf.table = paf.table, color.by = "identity")
## Define custom sample order
seqnames.order <- c("HG00438_2", "HG01358_2", "HG02630_2", "HG03453_2")
plotAVA(paf.table = paf.table, color.by = "direction", seqnames.order = seqnames.order)
## Only samples present in custom sample order are being plotted
seqnames.order <- c("HG00438_2", "HG01358_2", "HG03453_2")
plotAVA(paf.table = paf.table, color.by = "direction", seqnames.order = seqnames.order)
## Outline PAF alignments
plotAVA(paf.table = paf.table, outline.alignments = TRUE)
## Highlight structural variants
plotAVA(
paf.table = paf.table, min.deletion.size = 1000, min.insertion.size = 1000,
highlight.sv = "outline"
)
## Bin PAF alignments into user defined bin and color them by sequence identity (% of matched bases)
plotAVA(paf.table = paf.table, binsize = 10000)
## Add annotation to all-versus-all alignments ##
plt <- plotAVA(paf.table = paf.table, color.by = "direction")
annot.file <- system.file("extdata", "test_annot_ava.RData", package = "SVbyEye")
annot.gr <- get(load(annot.file))
addAnnotation(ggplot.obj = plt, annot.gr = annot.gr, coordinate.space = "self", y.label.id = "ID")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.