-
Notifications
You must be signed in to change notification settings - Fork 0
Home
Read, build, modify and displays as svg phylogenetic trees.
light_phylogeny allows to draw one or several reconciled gene trees mapped to a species tree (recPhyloXML files). It allows as well to draw simple gene trees from phyloXML and newick files, however its main purpose is to deal with recPhyloXML reconciled trees. Ultimately you may just want to use it to build and modify phylogenic trees.
thirdkind (see https://github.com/simonpenel/thirdkind/wiki) is using light_phylogeny.
Note: the library is intended for binary trees only.
Some examples of use:
Build a svg representation of a phylogenetic reconciled (or not) tree with events (loss, duplication, speciation, transfer).
Read a recphyloxml file: create a svg representation of -for example- the gene tree(s) and the associated species tree.
Read 2 nested recphyloxml files: create svg representations of -for example- the gene tree(s), the associated symbiot tree(s) and the associated host tree ("3-levels visualisation").
Read a phylogenetic tree and modify it by adding children, removing nodes, etc.
Build a phylogenetic tree from your data.
Keywords: phylogeny, reconciled trees, phylogenetic trees
phyloXML, recPhyloXML, rooted newick ( NHX balises will not be considered ).
https://crates.io/crates/light_phylogeny
You may find many code examples in the "examples" directory (https://github.com/simonpenel/light_phylogeny/tree/master/examples).
use light_phylogeny::{ArenaTree,Options,Config,read_newick,phyloxml_processing};
fn main() {
let mut tree: ArenaTree<String> = ArenaTree::default();
let options: Options = Options::new();
let config: Config = Config::new();
read_newick("examples/newick.txt".to_string(), &mut tree);
phyloxml_processing(&mut tree, &options, &config,"read_newick-clado.svg".to_string());
println!("Please open output file 'read_newick-clado.svg' with your browser");
let mut tree: ArenaTree<String> = ArenaTree::default();
let mut options: Options = Options::new();
options.real_length_flag = true;
let config: Config = Config::new();
read_newick("examples/newick.txt".to_string(), &mut tree);
phyloxml_processing(&mut tree, &options, &config,"read_newick-real-clado.svg".to_string());
println!("Please open output file 'read_newick-real-clado.svg' with your browser");
}
Some newick examples are available here : https://github.com/simonpenel/light_phylogeny/tree/master/newick_examples
use light_phylogeny::{ArenaTree,Options,Config,Event,add_child,move_child,phyloxml_processing,
summary,reset_pos};
fn main() {
let mut tree: ArenaTree<String> = ArenaTree::default();
let mut options: Options = Options::new();
let config: Config = Config::new();
// Create a new node root
let root = tree.new_node("root".to_string());
// Create new nodes
let a1 = tree.new_node("a1".to_string());
let a2 = tree.new_node("a2".to_string());
let a = tree.new_node("a".to_string());
let b1 = tree.new_node("b1".to_string());
let b2 = tree.new_node("b2".to_string());
let b = tree.new_node("b".to_string());
let c = tree.new_node("c".to_string());
let d = tree.new_node("d".to_string());
println!("Initial tree :");
summary(&mut tree);
// Set names
tree.arena[root].name = "MyRoot".to_string();
tree.arena[a].name = "Gene A".to_string();
tree.arena[a1].name = "Gene A1".to_string();
tree.arena[a2].name = "Gene A2".to_string();
tree.arena[b].name = "Gene B".to_string();
tree.arena[b1].name = "Gene B1".to_string();
tree.arena[b2].name = "Gene B2".to_string();
tree.arena[c].name = "Gene C".to_string();
tree.arena[d].name = "Gene D".to_string();
println!("");
println!("Tree after name attribution:");
summary(&mut tree);
// Set hierarchy
// a1 and a2 are children of a
add_child(&mut tree,a,a1);
add_child(&mut tree,a,a2);
// a1 and a2 are children of a
add_child(&mut tree,b,b1);
add_child(&mut tree,b,b2);
// a and b are children of c
add_child(&mut tree,c,a);
add_child(&mut tree,c,b);
// c and d are children of root
add_child(&mut tree,root,c);
add_child(&mut tree,root,d);
println!("");
println!("Tree after hierarchy attribution:");
summary(&mut tree);
// Display internal nodes
options.gene_internal = true ;
phyloxml_processing(&mut tree, &options, &config,"modify_tree_ini.svg".to_string());
println!("Add a loss to C");
let loss = tree.new_node("loss".to_string());
tree.arena[loss].name = "Loss".to_string();
tree.arena[loss].e = Event::Loss;
add_child(&mut tree,c,loss);
println!("Add a node up to B");
let add = tree.new_node("add".to_string());
tree.arena[add].name = "Added up to B".to_string();
println!("Move A to new node ");
move_child(&mut tree, a, add);
println!("Move B to new node ");
move_child(&mut tree, b, add);
println!("Move new node to C ");
add_child(&mut tree, c, add);
println!("Tree after hierarchy modification:");
summary(&mut tree);
reset_pos(&mut tree);
phyloxml_processing(&mut tree, &options, &config,"modify_tree_mod.svg".to_string());
println!("Please open output files 'modify_tree_ini.svg' and 'modify_tree_mod.svg' with your browser");
println!("OK.");
}
You may try the codes in the 'examples' directory:
cargo run --example read_newick
cargo run --example build_tree
cargo run --example lca
cargo run --example modify_tree
cargo run --example read_recphyloxml
Read and display a reconciled tree from a recPhyloXML file:
https://github.com/simonpenel/light_phylogeny/blob/master/examples/read_recphyloxml.rs
See Rust documentation : https://docs.rs/light_phylogeny/
For output examples, please see the thirdkind home page https://github.com/simonpenel/thirdkind/wiki
gene_internal:false [display internal node in gene trees]
species_internal:false [display internal node in gene trees (recphyloxml)]
clado_flag:true [display a cladogramme]
species_only_flag:false [display only the species tree (recphyloxml)]
real_length_flag:false [use real length branch]
open_browser:false [open the svg in the browser]
verbose:false [set log verbosity]
disp_gene:n=0 [if n > 0 display only the n th gene tree (recphyloxml)]
scale:1.0 [real length branch scale (if real_length_flag is true) ]
ratio:1.0 [ratio between width of species tree and gene trees (recphyloxml)]
rotate:true [rotate 90 degree, i.e. view as potrait]
thickness_flag: false [display redundant transfers as one, with an opacity according to transfer abundance. Only 1 gene tree is displayed ]
thickness_thresh:n=1 [with thickness_flag=true, abundance threshold for displaying a redundant transfer]
thickness_gene:n=1 [with thickness_flag=true, number of the gene tree to display]
thickness_disp_score:false [with thickness_flag=true, display the redundant transfer abundance]
height:1.0 [scaling svg height]
width:1.0 [scaling svg width]
support:false [display node support]
free_living:false [set free_living mode for species node named "FREE_LIVING"]
uniform:false [species nodes uniformisation]
species_color:"pink"
species_opacity:"0.9"
single_gene_color:"blue"
gene_opacity:"0.9"
species_police_color:"orange"
species_police_size:"12"
gene_police_size:"10"
bezier:"1"
thirdkind uses the library to display reconciled phylogenic trees, in particular nested reconciliations :
https://github.com/simonpenel/thirdkind/wiki
See http://phylariane.univ-lyon1.fr/recphyloxml/
recPhyloXML paper: https://www.ncbi.nlm.nih.gov/pmc/articles/PMC6198865/
phyloXML paper: https://www.ncbi.nlm.nih.gov/pmc/articles/PMC2774328/
- Possible problem with the obsolete version of recPhyloXML format (speciationLoss is supported, speciationOutLoss and speciationOut are not supported yet)
"Arena" Tree structure is inspired by the code proposed here