View source: R/phylo.rcandy.plot.R
RCandyVis | R Documentation |
This function reads a reference genome in GFF format and then plots the genetic features (coding sequences) on both forward and reverse strands.
RCandyVis( tree.file.name, taxon.metadata.file = NULL, taxon.metadata.columns = NULL, gubbins.gff.file = NULL, recom.input.type = "Gubbins", ref.genome.name = NULL, metadata.column.label.angle = 90, show.gene.label = FALSE, gene.label.angle = 45, show.metadata.columns = TRUE, subtree.tips = NULL, color.pallette = "inferno", taxon.metadata.columns.colors = NULL, taxon.id.column = NULL, show.genome.ticks = TRUE, show.genome.axis = TRUE, rec.heatmap.color = c("red", "blue"), tree.scale.length = NULL, show.rec.events = TRUE, show.metadata.label = TRUE, taxon.metadata.label.cex = 0.95, taxon.metadata.delimiter = "\t", ref.genome.length = NULL, show.rec.freq.per.base = FALSE, show.rec.freq.per.genome = TRUE, show.rec.per.genome.scale = FALSE, rec.events.per.base.as.heatmap = TRUE, ladderize.tree.right = NULL, midpoint.root = FALSE, rec.plot.bg.transparency = 0.1, show.genome.annot = TRUE, show.rec.plot.tracks = FALSE, show.rec.plot.border = FALSE, ace.model.name = "ER", trait.for.ancestral.reconstr = NULL, save.to.this.file = NULL, plot.width = 12, plot.height = 9.5, show.tip.label = FALSE, align.tip.label = FALSE, show.fig.legend = TRUE, genome.start = NULL, genome.end = NULL, color.tree.tips.by.column = NULL, tree.tip.node.cex = 0.35, tree.tip.label.cex = 0.35, tree.node.cex = 0.6 )
tree.file.name |
File name or ape phylo object containing the phylogenetic tree in Newick format. |
taxon.metadata.file |
File name or data frame containing taxon metadata. |
taxon.metadata.columns |
Vector containing metadata columns to plotted. |
gubbins.gff.file |
Gubbins output recombination file in GFF format. |
recom.input.type |
Type of input recombination data, either "Gubbins" GFF or "BRATNextGen" tabular data. |
ref.genome.name |
Reference genome file name in GFF format. |
metadata.column.label.angle |
Angle for the metadata column labels. |
show.gene.label |
A Boolean value indicating whether or not to show gene labels in the reference genome. |
gene.label.angle |
Angle for the gene labels. |
show.metadata.columns |
A Boolean indicating whether or not to show metadata in the figure. |
subtree.tips |
A vector containing a subset of taxons/taxa used to generate a subtree from the main phylogenetic tree. |
color.pallette |
A vector containing names of the viridis colour palletes for visualisation. Choose from "plasma", "cividis", "viridis", "magma" and "inferno". |
taxon.metadata.columns.colors |
A vector containing column names with custom colours (should be equal in size and same order as the vector specified for taxon.metadata.columns option) |
taxon.id.column |
Character or string for the column name containing the strain/taxon name in the metadata file. |
show.genome.ticks |
A Boolean indicating whether to show the xticks for the recombination events diagram/heatmap. |
show.genome.axis |
A Boolean indicating whether to show the axis for the recombination events diagram/heatmap. |
rec.heatmap.color |
A two-value vector containing colour names to use for the recombination event diagram/heatmap. |
tree.scale.length |
A positive number showing the length of the phylogenetic tree branches. |
show.rec.events |
A Boolean indicating whether to show the recombination event diagram/diagram. |
show.metadata.label |
A Boolean indicating whether to show labels for the selected metadata columns. |
taxon.metadata.label.cex |
A number for the size of the labels for the selected matadata columns |
taxon.metadata.delimiter |
A delimeter separating metadata columns. |
ref.genome.length |
An optional reference genome length, otherwise it's read from the reference genome GFF file or data frame. |
show.rec.freq.per.base |
A Boolean indicating whether to show the frequency of recombination per genomic position/base. |
show.rec.freq.per.genome |
A Boolean indicating whether to show the frequency of recombination events per genome/taxon. |
show.rec.per.genome.scale |
A boolean indicating whether to show the barplot scale for the number of recombination events per genome. |
rec.events.per.base.as.heatmap |
A Boolean indicating whether to show the frequency of recombination events per genome/taxon as a barchart or colour scale (heatmap). |
ladderize.tree.right |
A Boolean indicating whether to ladderize the phylogenetic tree to the right. |
midpoint.root |
A Boolean indicating whether to root the phylogenetic tree at midpoint. |
rec.plot.bg.transparency |
A value between 0 and 1 indicating the transparency level of the background for the recombination events diagram/heatmap. |
show.genome.annot |
A Boolean indicating whether to show genome annotation above the recombination events diagram/heatmap. |
show.rec.plot.tracks |
A Boolean indicating whether to plot genome tracks for each taxa. |
show.rec.plot.border |
A Boolean indicating whether to show the border for the recombination events diagram/heatmap. |
ace.model.name |
A character or string for the model name used for the discrete ancestral character reconstruction. Choose from "ARD", "ER" and "SYM". |
trait.for.ancestral.reconstr |
A character or string for the column in the metadata file or data frame used for discrete ancestral character reconstruction. |
save.to.this.file |
If speficified save the plot to this filename, otherwise show the plot in R. |
plot.width |
Width of the figure |
plot.height |
Height of the figure |
show.tip.label |
A Boolean indicating whether to show the phylogenetic tip labels. |
align.tip.label |
A Boolean indicating whether to align the phylogenetic tip labels. |
show.fig.legend |
A Boolean indicating whether to show the legend for the selected metadata columns. |
genome.start |
A positive number indicating start position in the genome to zoom in. |
genome.end |
A positive number indicating end position in the genome to zoom in. |
color.tree.tips.by.column |
Character or string for the column name in the metadata file for colouring the phylogenetic tree tips or terminal nodes. |
tree.tip.node.cex |
A number for the terminal node or tip size in the phylogenetic tree. |
tree.tip.label.cex |
A number for the tip label size in the phylogenetic tree. |
tree.node.cex |
A number for the size of the nodes phylogenetic tree. |
None
Chrispin Chaguza, Chrispin.Chaguza@gmail.com
https://github.com/ChrispinChaguza/RCandy
## Not run: Read phylogenetic tree in Newick format, reference genome in GFF formatted file, generated usign readseq, and Gubbins GFF file to plot the genomic features metadata.file<-system.file("extdata", "ST320.tsv", package = "RCandy", mustWork = TRUE) tree.file<-system.file("extdata", "ST320.final_tree.tre", package = "RCandy", mustWork = TRUE) gubbins.gff<-system.file("extdata", "ST320.recombination_predictions.gff", package = "RCandy",mustWork = TRUE) ref.genome.gff<-system.file("extdata", "Hungary19A-6.gff", package = "RCandy", mustWork = TRUE) RCandyVis(tree.file.name = tree.file, taxon.metadata.file = metadata.file, taxon.metadata.columns = c("Source","Country"),ref.genome.name = ref.genome.gff, gubbins.gff.file = gubbins.gff,color.tree.tips.by.column = "Country", show.rec.freq.per.base = FALSE,show.gene.label = FALSE,ladderize.tree.right = TRUE, midpoint.root = TRUE) RCandyVis(tree.file.name = tree.file, taxon.metadata.file = metadata.file, taxon.metadata.columns = c("Source","Country"),ref.genome.name = ref.genome.gff, gubbins.gff.file = gubbins.gff,color.tree.tips.by.column = "Country", show.rec.freq.per.base = FALSE,show.gene.label = FALSE,ladderize.tree.right = TRUE, midpoint.root = TRUE,genome.start = 30000, genome.end = 60000,show.gene.label=TRUE, save.to.this.file = "RCandy.output.pdf",) ## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.