Struct light_phylogeny::Options
source · pub struct Options {Show 33 fields
pub branch: bool,
pub gene_internal: bool,
pub species_internal: bool,
pub clado_flag: bool,
pub species_only_flag: bool,
pub real_length_flag: bool,
pub open_browser: bool,
pub verbose: bool,
pub disp_gene: usize,
pub scale: f32,
pub ratio: f32,
pub rotate: bool,
pub remove: bool,
pub thickness_flag: bool,
pub thickness_thresh: usize,
pub thickness_gene: usize,
pub thickness_disp_score: bool,
pub tidy: bool,
pub tidy_leaves_check: bool,
pub optimisation: bool,
pub height: f32,
pub width: f32,
pub support: bool,
pub free_living: bool,
pub free_living_shift: bool,
pub uniform: bool,
pub gthickness: usize,
pub sthickness: usize,
pub squaresize: f32,
pub trans_start: Option<String>,
pub trans_end: Option<String>,
pub mid_dist: bool,
pub gene_colors: Vec<String>,
}
Expand description
Structure Options: these are the drawing options for the svg. Concerning recPhyloXML file, we assume in the documentation that it describes a gene/species reconciliation.
Fields§
§branch: bool
Display branch length
gene_internal: bool
Display internal gene nodes.
species_internal: bool
Display internal species nodes (recPhyloXML).
clado_flag: bool
Display a cladogramme.
species_only_flag: bool
Only draw the species tree (recPhyloXML).
real_length_flag: bool
Use the real branch length.
open_browser: bool
Open the svg in the browser.
verbose: bool
Verbose mode.
disp_gene: usize
Only draw gene tree number # (recPhyloXML).
scale: f32
Scale to be applied to real branch length.
ratio: f32
Ratio between species pipe tree width and cumulated gene trees width (recPhyloXML).
rotate: bool
Rotate the svg 90 counter clockwise.
remove: bool
Not yet implemented.
thickness_flag: bool
Draw only one transfer to represent redundant transfers (recPhyloXML).
thickness_thresh: usize
Abundance threshold for displaying redundant transfers (recPhyloXML).
thickness_gene: usize
Number of the gene tree to display when displaying redundant transfers as one (recPhyloXML).
thickness_disp_score: bool
Display the abundance of the redundant transfers (recPhyloXML).
tidy: bool
Tidy option
tidy_leaves_check: bool
Tidy option
optimisation: bool
Optimise species branches left/right orientation in order to minimise transfer crossings (recPhyloXML, under development).
height: f32
Scale to be applied to the heigth of the tree.
width: f32
Width to be applied to the heigth of the tree.
support: bool
Display support.
free_living: bool
Support “free living” species (i.e. not associated to a species of the species tree) in the gene tree(s).
free_living_shift: bool
With free_living option : in case of multriple gene trees, shifting free living trees instead superposing.
uniform: bool
Uniformise the species tree nodes size.
gthickness: usize
thickness of the stroke for gene trees
sthickness: usize
thickness of the stroke for species trees
squaresize: f32
used for the size of square, circle, etc.
trans_start: Option<String>
start transfert
trans_end: Option<String>
end transfert
mid_dist: bool
place les duplication et les branchingout a mi distance de leur parent
gene_colors: Vec<String>
user-defined list of colors for genes
Implementations§
source§impl Options
impl Options
sourcepub fn new() -> Self
pub fn new() -> Self
Examples found in repository?
More examples
5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23
fn main() {
let transfers = vec![];
let mut options: Options = Options::new();
options.free_living = true;
let config: Config = Config::new();
// Version portrait
let mut sp_tree: ArenaTree<String> = ArenaTree::default();
let mut gene_trees:std::vec::Vec<ArenaTree<String>> = Vec::new();
let mut global_roots: std::vec::Vec<usize> = Vec::new();
read_recphyloxml_multi("examples/file_bug.recphylo".to_string(),
&mut sp_tree, &mut gene_trees, &mut global_roots);
recphyloxml_processing(&mut sp_tree, &mut gene_trees, &mut options, &config, true,
&transfers, "read_recphyloxml_portrait.svg".to_string());
println!("Please open output file 'read_recphyloxml_portrait.svg' with your browser");
}
5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24
fn main() {
let transfers = vec![];
let mut options: Options = Options::new();
let config: Config = Config::new();
// Version portrait
options.species_internal = true;
let mut sp_tree: ArenaTree<String> = ArenaTree::default();
let mut gene_trees:std::vec::Vec<ArenaTree<String>> = Vec::new();
let mut global_roots: std::vec::Vec<usize> = Vec::new();
read_recphyloxml_multi("examples/gene_parasite_page4.recphylo".to_string(),
&mut sp_tree, &mut gene_trees, &mut global_roots);
recphyloxml_processing(&mut sp_tree, &mut gene_trees, &mut options, &config, true,
&transfers, "read_recphyloxml_portrait.svg".to_string());
println!("Please open output file 'read_recphyloxml_portrait.svg' with your browser");
}
6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22
fn main() {
// ============================================================================================
// Version portrait
let transfers = vec![];
let mut options: Options = Options::new();
let config: Config = Config::new();
let mut sp_tree: ArenaTree<String> = ArenaTree::default();
let mut gene_trees:std::vec::Vec<ArenaTree<String>> = Vec::new();
let mut global_roots: std::vec::Vec<usize> = Vec::new();
read_recphyloxml_multi("examples/example_obsolete_specloss.xml".to_string(),
&mut sp_tree, &mut gene_trees, &mut global_roots);
recphyloxml_processing(&mut sp_tree, &mut gene_trees, &mut options, &config, true,
&transfers, "example_obsolete_specloss.svg".to_string());
println!("Please open output file 'example_obsolete_specloss.svg' with your browser");
}
5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25
fn main() {
let transfers = vec![];
let mut options: Options = Options::new();
let config: Config = Config::new();
options.free_living = true;
options.gene_internal = true;
options.species_internal = true;
let mut sp_tree: ArenaTree<String> = ArenaTree::default();
let mut gene_trees:std::vec::Vec<ArenaTree<String>> = Vec::new();
let mut global_roots: std::vec::Vec<usize> = Vec::new();
read_recphyloxml_multi("examples/free_living_reconciliated.recphylo".to_string(),
&mut sp_tree, &mut gene_trees, &mut global_roots);
recphyloxml_processing(&mut sp_tree, &mut gene_trees, &mut options, &config, true,
&transfers, "read_recphyloxml_free_living_reconciliated.svg".to_string());
println!("Please open output file 'read_recphyloxml_free_living_reconciliated.svg' with your browser");
}
- examples/read_recphyloxml_free_living_2.rs
- examples/read_recphyloxml_free_living_3.rs
- examples/read_recphyloxml_free_living_quadruple.rs
- examples/lca_2.rs
- examples/read_recphyloxml_copie.rs
- examples/read_recphyloxml_big_2.rs
- examples/read_recphyloxml_free_living_double.rs
- examples/read_recphyloxml_obsolete.rs
- examples/read_recphyloxml_switch.rs
- examples/lca.rs
- examples/test_colors.rs
- examples/read_phyloxml_stroke_thickness.rs
- examples/build_tree.rs
- examples/read_phyloxml_tidy.rs
- examples/read_recphyloxml_big.rs
- examples/read_phyloxml.rs
- examples/modify_tree.rs
- examples/read_recphyloxml_4.rs
- examples/read_recphyloxml_width_height.rs
- examples/read_recphyloxml_thickness_2.rs
- examples/read_recphyloxml_tidyspecies.rs
- examples/read_recphyloxml_2.rs
- examples/read_phyloxml_tidy_virus.rs
- examples/read_newick.rs
- examples/read_recphyloxml_thickness.rs
- examples/read_recphyloxml_stroke_thickness.rs
- examples/read_recphyloxml_3.rs
- examples/read_recphyloxml_multi.rs
- examples/read_recphyloxml_thirdlevel_5.rs
- examples/read_recphyloxml_thirdlevel_6.rs
- examples/read_recphyloxml_thirdlevel_bug.rs
- examples/read_recphyloxml_thirdlevel_1.rs
- examples/read_recphyloxml_thirdlevel_3.rs
- examples/read_recphyloxml_thirdlevel_1_midist.rs
- examples/read_recphyloxml_thirdlevel_2.rs
- examples/read_recphyloxml_thirdlevel_stroke.rs
- examples/read_recphyloxml_thirdlevel_fl.rs
- examples/read_recphyloxml_thirdlevel_alex.rs
- examples/read_recphyloxml_thirdlevel_4.rs
- examples/read_recphyloxml_thirdlevel_4_landscape.rs
- examples/read_recphyloxml.rs