1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541
/// name = "light_phylogeny"
/// version = "1.4.1"
/// authors = ["Simon Penel <simon.penel@univ-lyon1.fr>"]
/// license = "CECILL-2.1"
// Convention: "pipe" trees are equivalent to "upper" trees, "path" trees are equivalenet to "lower" trees
use log::{info};
use std::fs;
use crate::arena::Options;
use crate::arena::ArenaTree;
use crate::arena::Config;
use crate::arena::PIPEBLOCK;
use crate::arena::{newick2tree,xml2tree};
use crate::arena::{knuth_layout,cladogramme,check_contour_postorder,
check_contour_postorder_tidy_tree,shift_mod_xy,set_middle_postorder,real_length,
set_leaves_y_values,shift_nodes_y_values};
use crate::arena::{map_species_trees,set_species_width,check_vertical_contour_postorder,
bilan_mappings,center_gene_nodes,move_dupli_mappings,move_species_mappings,
species_uniformisation,process_fl,uniformise_gene_leaves_y_values};
use crate::arena::{find_sptrees,find_rgtrees,check_for_obsolete,scale_heigth,scale_width};
use crate::thirdlevel::{get_gtransfer,optimisation,check_optimisation,classify_transfer,reorder_transfers};
use crate::drawing::{draw_tree,draw_sptree_gntrees};
// Functions
// =========
/// Read a newick file and store the tree into ArenaTree structure.
pub fn read_newick(filename:String, tree: &mut ArenaTree<String>) {
let contents = fs::read_to_string(filename);
let contents = match contents {
Ok(contents) => contents,
Err(err) => {
eprintln!("\nERROR:\nSomething went wrong when reading the newick file.");
eprintln!("{}",err);
eprintln!("Please check file name and path.");
std::process::exit(1);
},
};
let root = tree.new_node("Root".to_string());
newick2tree(contents, tree, root, &mut 0);
}
/// Read a phyloxml file and store the tree into a ArenaTree structure.
pub fn read_phyloxml(filename: String, tree: &mut ArenaTree<String>) {
let contents = fs::read_to_string(filename);
let contents = match contents {
Ok(contents) => contents,
Err(err) => {
eprintln!("\nERROR:\nSomething went wrong when reading the phyloxml file.");
eprintln!("{}",err);
eprintln!("Please check file name and path.");
std::process::exit(1);
},
};
let doc = roxmltree::Document::parse(&contents);
let doc = match doc {
Ok(xml) => xml,
Err(_err) => {
eprintln!("\nERROR:\nSomething went wrong when parsing the input file.");
eprintln!("It seems this is not a correct xml file:");
eprintln!("{}",_err);
std::process::exit(1);
},
};
let descendants = doc.root().descendants();
// Search for the first occurence of the "clade" tag
for node in descendants {
if node.has_tag_name("clade"){
// The first clade is the root
// Initialize the index used for defining the value
let mut index = &mut 0;
// Val of the root
let root = "N".to_owned()+&index.to_string();
// Create the node, get its associated index and store it in root again
let root = tree.new_node(root.to_string());
// Call xlm2tree on the root
xml2tree(node, root, &mut index, tree);
break;
}
}
info!("Tree : {:?}",tree);
}
/// Adding invisible parent nodes in order to display multiple species trees.
pub fn add_virtual_roots(
global_roots: &mut std::vec::Vec<String>,
pere: &mut usize,
level: &mut usize,
global_pipe: &mut ArenaTree<String>,
) {
if global_roots.len() == 1 {
let nom_droite = "VIRTUAL_".to_string().to_owned() + &level.to_string() + "_"+&(*pere-1).to_string();
let nom_pere = "ROOT_".to_string().to_owned() + &level.to_string() + "_"+&(*pere).to_string();
let node_pere = global_pipe.new_node(nom_pere.clone());
global_pipe.arena[node_pere].name = nom_pere;
global_pipe.arena[node_pere].visible = false;
let node_droite = global_pipe.node(nom_droite.clone());
global_pipe.arena[node_droite].name = nom_droite;
global_pipe.arena[node_droite].visible = false;
let node_left = global_pipe.node(global_roots[0].to_string());
global_pipe.arena[node_pere].children.push(node_left);
global_pipe.arena[node_pere].children.push(node_droite);
global_pipe.arena[node_left].parent = Some(node_pere);
global_pipe.arena[node_droite].parent = Some(node_pere);
if *pere > 0 {
let mut my_vect : Vec<String> = (0..(*pere+1)).map(|n| "ROOT_".to_string().to_owned() + &level.to_string() + "_" + &n.to_string()).collect();
*level = *level + 1;
*pere = 0;
add_virtual_roots(&mut my_vect, pere, level, global_pipe);
}
}
if global_roots.len() == 2 {
let nom_pere = "ROOT_".to_string().to_owned() + &level.to_string() + "_"+&(*pere).to_string();
let node_pere = global_pipe.node(nom_pere.clone());
global_pipe.arena[node_pere].name = nom_pere;
global_pipe.arena[node_pere].visible = false;
let node_droite = global_pipe.node(global_roots[1].to_string());
let node_left = global_pipe.node(global_roots[0].to_string());
global_pipe.arena[node_pere].children.push(node_left);
global_pipe.arena[node_pere].children.push(node_droite);
global_pipe.arena[node_left].parent = Some(node_pere);
global_pipe.arena[node_droite].parent = Some(node_pere);
if *pere > 0 {
let mut my_vect : Vec<String> = (0..(*pere+1)).map(|n| "ROOT_".to_string().to_owned() + &level.to_string() + "_" + &n.to_string()).collect();
*level = *level + 1;
*pere = 0;
add_virtual_roots(&mut my_vect, pere, level, global_pipe);
}
}
if global_roots.len() > 2 {
let (tete, queue) = global_roots.split_at(2);
let nom_pere = "ROOT_".to_string().to_owned() + &level.to_string() + "_"+&(*pere).to_string();
let node_pere = global_pipe.node(nom_pere.clone());
global_pipe.arena[node_pere].name = nom_pere;
global_pipe.arena[node_pere].visible = false;
let node_droite = global_pipe.node(tete[1].to_string());
let node_left = global_pipe.node(tete[0].to_string());
global_pipe.arena[node_pere].children.push(node_left);
global_pipe.arena[node_pere].children.push(node_droite);
global_pipe.arena[node_left].parent = Some(node_pere);
global_pipe.arena[node_droite].parent = Some(node_pere);
let mut queue_clone = queue.to_vec();
*pere = *pere + 1;
add_virtual_roots(&mut queue_clone, pere,level, global_pipe);
}
}
/// Read a recphyloxml file and store the species and gene trees into several ArenaTree structures.
pub fn read_recphyloxml_multi(
filename: String,
global_pipe: &mut ArenaTree<String>,
gene_trees: &mut std::vec::Vec<ArenaTree<String>>,
global_roots: &mut std::vec::Vec<usize>
) {
let contents = fs::read_to_string(filename);
let contents = match contents {
Ok(contents) => contents,
Err(err) => {
eprintln!("\nERROR:\nSomething went wrong when reading the recPhyloXML file.");
eprintln!("{}",err);
eprintln!("Please check file name and path.");
std::process::exit(1);
},
};
// let doc = &mut roxmltree::Document::parse(&contents).unwrap();
let doc = roxmltree::Document::parse(&contents);
let doc = &mut match doc {
Ok(contents) => contents,
Err(_err) => {
eprintln!("\nERROR:\nSomething went wrong when parsing the input file.");
eprintln!("It seems this is not a correct xml file:");
eprintln!("{}",_err);
std::process::exit(1);
},
};
// Get the species trees:
// Get the list of nodes associated to the "spTree" tag
// let spnodes = find_sptrees(doc).expect("No clade spTree has been found in xml");
let spnodes = find_sptrees(doc);
let spnodes = match spnodes {
Ok(indexes) => indexes,
Err(_err) => {
eprintln!("\nERROR:\nNo clade spTree has been found in xml.\
\nIt seems that the input file is not a recPhyloXML file.\
\nUse option -F to force to use phyloXML or newick format.");
std::process::exit(1);
},
};
let mut index = &mut 0;
let mut global_root_names: std::vec::Vec<String> = Vec::new();
for spnode in spnodes {
info!("Search spTree node {:?}",spnode);
let spnode = doc.get_node(spnode).expect("Unable to get the Node associated to this nodeId");
info!("Associated spTree : {:?}",spnode);
// Analyse the gene tree
let descendants = spnode.descendants();
// Search for the first occurence of the "clade" tag
// index = &(index + 1) ;
*index += 1;
for node in descendants {
if node.has_tag_name("clade"){
let globalroot = "G".to_owned() + &index.to_string();
global_root_names.push(globalroot.clone());
// Create the node, get its associated index and store it in root again
info!("Create {}","G".to_owned() + &index.to_string());
let globalroot = global_pipe.new_node(globalroot.to_string());
global_roots.push(globalroot);
// Call xlm2tree on the root
xml2tree(node, globalroot, &mut index, global_pipe);
break;
}
}
}
if global_roots.len() > 1 {
add_virtual_roots(&mut global_root_names, &mut 0, &mut 0, global_pipe);
}
let nb_sptree = global_roots.len().clone();
println!("Number of species trees : {}",nb_sptree);
info!("List of species roots : {:?}",global_roots);
// Get the gene trees:
// Get the list of nodes associated to the "recGeneTree" tag
// let rgnodes = find_rgtrees(doc).expect("No clade recGeneTree has been found in xml");
let rgnodes = find_rgtrees(doc);
let rgnodes = match rgnodes {
Ok(indexes) => indexes,
Err(_err) => {
eprintln!("\nERROR:\nNo clade recGeneTree has been found in xml.\
\nIt seems that the input file is not a recPhyloXML file.\
\nUse option -F to force to use phyloXML or newick format.");
std::process::exit(1);
},
};
for rgnode in rgnodes {
let mut gene_tree: ArenaTree<String> = ArenaTree::default();
info!("Search recGeneTree node {:?}",rgnode);
let rgnode = doc.get_node(rgnode).expect("Unable to get the Node associated to this nodeId");
info!("Associated recGeneTree : {:?}",rgnode);
// Analyse the gene tree
let descendants = rgnode.descendants();
// Search for the first occurence of the "clade" tag
for node in descendants {
if node.has_tag_name("clade"){
// The first clade is the root
// Initialize the index used for defining the value
let mut index = &mut 0;
// Val of the root
let root = "N".to_owned()+&index.to_string();
// Create the node, get its associated index and store it in root again
let root = gene_tree.new_node(root.to_string());
// Call xlm2tree on the root
xml2tree(node, root, &mut index, &mut gene_tree);
break;
}
}
// Traitement des balises obsoletes potentielles (ancien format recPhyloXML)
check_for_obsolete(&mut gene_tree, global_pipe);
// Ajoute l'arbre de gene
gene_trees.push(gene_tree);
}
let nb_gntree = gene_trees.len().clone();
println!("Number of gene trees : {}",nb_gntree);
info!("List of gene trees : {:?}",gene_trees);
}
/// Create a svg of the tree in phyloxml context.
pub fn phyloxml_processing(
mut tree: &mut ArenaTree<String>, // tree
options: &Options, // display options
config: &Config, // svg configuration
outfile: String // output file
) {
info!("Tree : {:?}",tree);
// -----------------------
// Traitement en 4 étapes
// -----------------------
// Au départ l'arbre est orienté du haut vers le bas (i.e. selon Y)
// Le svg sera tourné de -90 a la fin.
//
//----------------------------------------------------------
// 1ère étape :initialisation des x,y de l'arbre :
// profondeur => Y, left => X= 0, right X=1
// ---------------------------------------------------------
let root = tree.get_root();
knuth_layout(&mut tree, root, &mut 1);
// ---------------------------------------------------------
// Option : Cladogramme
// ---------------------------------------------------------
if options.clado_flag {
cladogramme(&mut tree);
}
// ---------------------------------------------------------
// 2ème étape : Vérifie les contours
// ---------------------------------------------------------
check_contour_postorder(&mut tree, root);
// ---------------------------------------------------------
// 3eme etape : Decale toutes les valeurs de x en fonction
// de xmod
// ---------------------------------------------------------
shift_mod_xy(&mut tree, root, &mut 0.0, &mut 0.0);
// ---------------------------------------------------------
// 4ème étape : Place le parent entre les enfants
// ---------------------------------------------------------
set_middle_postorder(&mut tree, root);
// ---------------------------------------------------------
// OPTIONAL Scale the heigt if needed
// ---------------------------------------------------------
if options.height != 1.0 {scale_heigth(&mut tree,options.height)};
// ---------------------------------------------------------
// OPTIONAL Scale the width if needed
// ---------------------------------------------------------
if options.width != 1.0 {scale_width(&mut tree,options.width)};
// ---------------------------------------------------------
// Option : real_length
// ---------------------------------------------------------
if options.real_length_flag {
real_length(&mut tree, root, &mut 0.0, & options);
if options.tidy{
check_contour_postorder_tidy_tree(&mut tree, root, & options, & config);
set_middle_postorder(&mut tree, root);
}
}
// ---------------------------------------------------------
// Fin: Ecriture du fichier svg
// ---------------------------------------------------------
println!("Output filename is {}",outfile);
draw_tree(&mut tree, outfile, & options, & config);
}
/// Create a svg of the tree in recphyloxml context.
pub fn recphyloxml_processing(
mut sp_tree: &mut ArenaTree<String>, // species tree
mut gene_trees: &mut std::vec::Vec<ArenaTree<String>>, // gene trees
options: &mut Options, // display options
config: &Config, // svg configuration
mapping: bool, // map gene and species
transfers: & std::vec::Vec<(String,String)>, // optional additional transfers
outfile: String, // output file
) {
// -----------------------
// Traitement en 12 etapes
// -----------------------
// Au depart l'arbre est orienté du haut vers le bas (i.e. selon Y)
// Le svg sera tourné de -90 a la fin.
//
// Option : ajout d'une branche free_living
let initial_root = sp_tree.get_root();
let mut free_root = 0;
if options.free_living {
// Check if a free_living is already in the data
let test_root = sp_tree.get_index("FREE_LIVING_ROOT".to_string());
match test_root {
Ok(root) => {
info!("[recphyloxml_processing] FREE_LIVING_ROOT already in the tree.");
let free_children = & sp_tree.arena[root].children;
free_root = free_children[1];
info!("FREE_LIVING_ROOT : {}",free_root);
},
Err(_err) => {
info!("[recphyloxml_processing] No FREE_LIVING_ROOT in the tree, I create one.");
let left = sp_tree.get_root();
let right = sp_tree.new_node("free_living".to_string());
sp_tree.arena[right].name = "FREE_LIVING".to_string();
free_root = right;
let fl_root = sp_tree.new_node("free_living_root".to_string());
sp_tree.arena[fl_root].name = "FREE_LIVING_ROOT".to_string();
sp_tree.arena[fl_root].visible = false;
sp_tree.arena[fl_root].children.push(left);
sp_tree.arena[fl_root].children.push(right);
sp_tree.arena[right].parent = Some(fl_root);
sp_tree.arena[left].parent = Some(fl_root);
info!("[recphyloxml_processing] FREE_LIVING_ROOT : {}",free_root);
},
};
}
//----------------------------------------------------------
// 1ere étape :initialisation des x,y de l'arbre d'espèces :
// profondeur => Y, left => X= 0, right X=1
// ---------------------------------------------------------
let root = sp_tree.get_root();
knuth_layout(&mut sp_tree,root, &mut 1);
// --------------------
// OPTIONAL Cladogramme
// --------------------
if options.clado_flag {
cladogramme(&mut sp_tree);
}
// ---------------------------------------------------------
// 2eme étape : mapping des genes sur l'espèce pour
// connaître le nombre de noeuds d'arbre de gènes associés à
// chaque noeud de l'arbre d'espèces
// ---------------------------------------------------------
if mapping {
map_species_trees(&mut sp_tree, &mut gene_trees);
info!("Species tree after mapping : {:?}",sp_tree);
}
// Option: uniformise les noeuds de l'arbre d'espece
if options.uniform {
species_uniformisation(&mut sp_tree);
// Il fait decaler le tout (mais je ne comprend pas pourquoi?)
let smallest_y = sp_tree.get_smallest_y();
let root_width = sp_tree.arena[root].nbg as f32 * PIPEBLOCK ;
shift_nodes_y_values(sp_tree, root, -smallest_y + root_width);
}
// ---------------------------------------------------------
// Option : utilise les longueurs de branches
// ---------------------------------------------------------
if options.real_length_flag {
// Tout d'abord il faut uniformiser
species_uniformisation(&mut sp_tree);
let min_dist = sp_tree.get_smallest_l();
if min_dist == 0.0 {
eprintln!("\nERROR:\nFound a branch with a distance equal to 0.0 in the 'pipe' tree.");
eprintln!("It is not possible create an 'upper' tree presenting branches of zero length.");
std::process::exit(1);
}
// On veut que la longueur minimum soit un peu superieure a la moitie de l'epaisser des noods
// On utilise le nombre de gene , qui a etet uniformise precedemment
// options.scale = (sp_tree.arena[root].nbg as f32 / 2.0 + 0.25) / min_dist ;
let optimised_factor = (sp_tree.arena[root].nbg as f32 / 2.0 + 0.25) / min_dist ;
options.scale = optimised_factor * options.scale ;
real_length(&mut sp_tree, root, &mut 0.0, & options);
// Il faut decaler le tout (mais je ne comprend pas pourquoi?)
let smallest_y = sp_tree.get_smallest_y();
let root_width = sp_tree.arena[root].nbg as f32 * PIPEBLOCK ;
shift_nodes_y_values(sp_tree, root, -smallest_y + root_width);
}
// ---------------------------------------------------------
// OPTIONAL Optimisation if needed
// ---------------------------------------------------------
if options.optimisation {
println!("Optimisation of orientation according to transfers");
if gene_trees.len() > 1 {
eprintln!("Optimisation is working only with single gene tree reconciliation.\nExit");
std::process::exit(1);
}
// Analyse des transfers
// ---------------------
let gene_transfers = get_gtransfer(&mut gene_trees[0]);
let selected_transfers = gene_transfers;
println!("Optimisation: Gene transfers = {:?}",selected_transfers);
let nbt = selected_transfers.len();
println!("Optimisation: Number of transfers = {:?}",nbt);
let mut numt = 0 ;
while numt < nbt {
classify_transfer(&selected_transfers[numt], &mut sp_tree, & numt);
numt = numt + 1;
}
let mut ordered = Vec::<usize>::new();
reorder_transfers(&mut sp_tree, root, &mut ordered);
ordered.reverse();
println!("Optimisation: Reordered transfers : {:?}",ordered);
for numt in ordered {
let tr_root = optimisation(&selected_transfers[numt], &mut sp_tree);
check_optimisation(&selected_transfers[numt], &mut sp_tree, tr_root);
}
}
// ---------------------------------------------------------
// 3eme étape : Vérifie les conflits dans l'arbre d'espèces
// au niveau horizontal -> valeurs xmod
// ---------------------------------------------------------
check_contour_postorder(&mut sp_tree, root);
// ---------------------------------------------------------
// 4eme étape : Décale toutes les valeurs de x en fonction
// de xmod dans l'arbre d'espèces
// ---------------------------------------------------------
shift_mod_xy(&mut sp_tree, root, &mut 0.0, &mut 0.0);
// ---------------------------------------------------------
// 5eme étape : Place le parent entre les enfants dans
// l'arbre d'espèces
// ---------------------------------------------------------
set_middle_postorder(&mut sp_tree, root);
// ---------------------------------------------------------
// 6ème etape : Fixe l'épaisseur de l'arbre d'espèces
// ---------------------------------------------------------
set_species_width(&mut sp_tree, &mut gene_trees);
// ---------------------------------------------------------
// 7ème étape : Vérifie les conflits verticaux dans
// l'arbre d'espèces
// ---------------------------------------------------------
if ! options.real_length_flag {
check_vertical_contour_postorder(&mut sp_tree, root, 0.0);
}
// ---------------------------------------------------------
// Egalise les feuilles
// ---------------------------------------------------------
if ! options.real_length_flag {
// Egalise les feuilles
let largest_y = sp_tree.get_largest_y();
set_leaves_y_values(sp_tree, root, largest_y);
} else {
if options.tidy{
check_contour_postorder_tidy_tree(&mut sp_tree, root, & options, & config);
}
}
// ---------------------------------------------------------
// OPTIONAL Scale the heigt if needed
// ---------------------------------------------------------
if options.height != 1.0 {scale_heigth(&mut sp_tree, options.height)};
// ---------------------------------------------------------
// OPTIONAL Scale the width if needed
// ---------------------------------------------------------
if options.width != 1.0 { scale_width(&mut sp_tree, options.width)};
// ---------------------------------------------------------
// 8ème etape : décale les noeuds de gene associés à un
// noeud d'especes pour éviter qu'ils soit superposés
// ---------------------------------------------------------
bilan_mappings(&mut sp_tree, &mut gene_trees, initial_root, & options);
// ---------------------------------------------------------
// 9eme etape : centre les noeuds de genes dans le noeud de l'espece
// ---------------------------------------------------------
center_gene_nodes(&mut sp_tree,&mut gene_trees, initial_root);
// ---------------------------------------------------------
// 10eme etape traite spécifiquement les duplications et les feuilles
// ---------------------------------------------------------
move_dupli_mappings(&mut sp_tree, &mut gene_trees, initial_root, & options);
// 11eme etape : Ca des speciations dont les fils sont
// dans le meme noeud "tuyeau" que le pere
// Cela n'arrice que quand on mappe des genes sur des hotes
// via les parasites (thirdlevel)
// ---------------------------------------------------------
move_species_mappings(&mut sp_tree, &mut gene_trees, initial_root);
// Option : traitement specifique des gene "free living"
if options.free_living {
process_fl(&mut sp_tree, &mut gene_trees, free_root, &options);
}
if ! options.real_length_flag {
// Egalise les feuilles 2
uniformise_gene_leaves_y_values(sp_tree, gene_trees);
}
// ---------------------------------------------------------
// Fin: Ecriture du fichier svg
// ---------------------------------------------------------
println!("Output filename is {}",outfile);
match options.species_only_flag {
true => {
if options.species_internal {
options.gene_internal = true;
}
draw_tree(&mut sp_tree, outfile, &options, &config);
},
false => draw_sptree_gntrees(&mut sp_tree, &mut gene_trees, outfile,
&options, &config, &transfers),
};
}