source: trunk/GDE/PHYML20130708/phyml/src/main.c

Last change on this file was 10307, checked in by aboeckma, 11 years ago

added most recent version of phyml

File size: 13.2 KB
Line 
1/*
2
3PhyML:  a program that  computes maximum likelihood phyLOGenies from
4DNA or AA homoLOGous sequences.
5
6Copyright (C) Stephane Guindon. Oct 2003 onward.
7
8All parts of the source except where indicated are distributed under
9the GNU public licence. See http://www.opensource.org for details.
10
11*/
12
13#include "spr.h"
14#include "utilities.h"
15#include "lk.h"
16#include "optimiz.h"
17#include "bionj.h"
18#include "models.h"
19#include "free.h"
20#include "help.h"
21#include "simu.h"
22#include "eigen.h"
23#include "pars.h"
24#include "alrt.h"
25#include "mixt.h"
26#include "sergeii.h"
27//#include "geo.h"
28
29#ifdef MPI
30#include "mpi_boot.h"
31#endif
32
33
34
35
36#if (defined PHYML || EVOLVE)
37
38int main(int argc, char **argv)
39{ 
40  calign *cdata;
41  option *io;
42  t_tree *tree;
43  int num_data_set;
44  int num_tree,num_rand_tree;
45  t_mod *mod;
46  time_t t_beg,t_end;
47  phydbl best_lnL;
48  int r_seed;
49  char *most_likely_tree=NULL;
50
51 
52#ifdef MPI
53  int rc;
54  rc = MPI_Init(&argc,&argv);
55  if (rc != MPI_SUCCESS) {
56    PhyML_Printf("\n== Err starting MPI program. Terminating.\n");
57    MPI_Abort(MPI_COMM_WORLD, rc);
58  }
59  MPI_Comm_size(MPI_COMM_WORLD,&Global_numTask);
60  MPI_Comm_rank(MPI_COMM_WORLD,&Global_myRank);
61#endif
62
63#ifdef QUIET
64  setvbuf(stdout,NULL,_IOFBF,2048);
65#endif
66
67  tree             = NULL;
68  mod              = NULL;
69  best_lnL         = UNLIKELY;
70
71  io = (option *)Get_Input(argc,argv);
72  if(!io) return(-1);
73
74
75  #ifdef EVOLVE
76  io->colalias = NO;
77  #endif
78
79  r_seed = (io->r_seed < 0)?(time(NULL)):(io->r_seed);
80  srand(r_seed);
81  io->r_seed = r_seed;
82
83  if(io->in_tree == 2) Test_Multiple_Data_Set_Format(io);
84  else io->n_trees = 1;
85
86  if(io->n_trees == 0 && io->in_tree == 2)
87    {
88      PhyML_Printf("\n. The input tree file does not provide a tree in valid format.");
89      Exit("\n");
90    }
91
92  if((io->n_data_sets > 1) && (io->n_trees > 1))
93    {
94      io->n_data_sets = MIN(io->n_trees,io->n_data_sets);
95      io->n_trees     = MIN(io->n_trees,io->n_data_sets);
96    }
97
98  For(num_data_set,io->n_data_sets)
99    {
100      best_lnL = UNLIKELY;
101      Get_Seq(io);     
102      Make_Model_Complete(io->mod);
103      Set_Model_Name(io->mod);
104      Print_Settings(io);
105      mod = io->mod;
106       
107
108      if(io->data)
109        {
110          if(io->n_data_sets > 1) PhyML_Printf("\n. Data set [#%d]\n",num_data_set+1);
111          cdata = Compact_Data(io->data,io);
112
113          Free_Seq(io->data,cdata->n_otu);
114         
115          for(num_tree=(io->n_trees == 1)?(0):(num_data_set);num_tree < io->n_trees;num_tree++)
116            {
117              if(!io->mod->s_opt->random_input_tree) io->mod->s_opt->n_rand_starts = 1;
118
119              For(num_rand_tree,io->mod->s_opt->n_rand_starts)
120                {
121                  if((io->mod->s_opt->random_input_tree) && (io->mod->s_opt->topo_search != NNI_MOVE))
122                    if(!io->quiet) PhyML_Printf("\n. [Random start %3d/%3d]\n",num_rand_tree+1,io->mod->s_opt->n_rand_starts);
123
124                  Init_Model(cdata,mod,io);
125
126                  if(io->mod->use_m4mod) M4_Init_Model(mod->m4mod,cdata,mod);
127
128                  switch(io->in_tree)
129                    {
130                    case 0 : case 1 : { tree = Dist_And_BioNJ(cdata,mod,io); break; }
131                    case 2 :          { tree = Read_User_Tree(cdata,mod,io); break; }
132                    }
133
134                  if(io->fp_in_constraint_tree != NULL) 
135                    {
136                      char *s;
137                      io->cstr_tree = Read_Tree_File_Phylip(io->fp_in_constraint_tree);
138                      s = Add_Taxa_To_Constraint_Tree(io->fp_in_constraint_tree,cdata);
139                      fflush(NULL);
140                      Free_Tree(tree);
141                      tree = Read_Tree(&s);
142                      io->in_tree = 2;
143                      Free(s);
144                      Check_Constraint_Tree_Taxa_Names(io->cstr_tree,cdata);
145                      Alloc_Bip(io->cstr_tree); 
146                      Get_Bip(io->cstr_tree->a_nodes[0],
147                              io->cstr_tree->a_nodes[0]->v[0],
148                              io->cstr_tree);
149                      if(!tree->has_branch_lengths) Add_BioNJ_Branch_Lengths(tree,cdata,mod);
150                    }
151
152                  if(!tree) continue;
153
154                  time(&t_beg);
155                  time(&(tree->t_beg));
156                 
157                  tree->mod          = mod;
158                  tree->io           = io;
159                  tree->data         = cdata;
160                  tree->n_pattern    = tree->data->crunch_len;
161                  tree->n_root       = NULL;
162                  tree->e_root       = NULL;
163
164                 
165                  Set_Both_Sides(YES,tree);     
166                 
167                  if(mod->s_opt->random_input_tree) 
168                    {
169                      PhyML_Printf("\n");
170                      Random_Tree(tree);
171                    }
172
173                  if((!num_data_set) && (!num_tree) && (!num_rand_tree)) Check_Memory_Amount(tree);
174
175                  if(io->cstr_tree && !Check_Topo_Constraints(tree,io->cstr_tree))
176                    {
177                      PhyML_Printf("\n\n== The initial tree does not satisfy the topological constraint.");
178                      PhyML_Printf("\n== Please use the user input tree option with an adequate tree topology.");
179                      Exit("\n");
180                    }
181
182                  Prepare_Tree_For_Lk(tree);
183                  Br_Len_Not_Involving_Invar(tree);
184                  Unscale_Br_Len_Multiplier_Tree(tree);
185                 
186                 
187#ifdef PHYML
188                 
189                  if(io->in_tree == 1) Spr_Pars(tree);
190                 
191                  if(io->do_alias_subpatt)
192                    {
193                      MIXT_Set_Alias_Subpatt(YES,tree);
194                      Lk(NULL,tree);
195                      MIXT_Set_Alias_Subpatt(NO,tree);
196                    }
197
198                  if(tree->mod->s_opt->opt_topo)
199                    {
200                      if(tree->mod->s_opt->topo_search      == NNI_MOVE) Simu_Loop(tree);
201                      else if(tree->mod->s_opt->topo_search == SPR_MOVE) Speed_Spr_Loop(tree);
202                      else                                               Best_Of_NNI_And_SPR(tree);
203
204                      if(tree->n_root) Add_Root(tree->a_edges[0],tree);
205                    }
206                  else
207                    {
208                      if(tree->mod->s_opt->opt_subst_param || 
209                         tree->mod->s_opt->opt_bl)                       Round_Optimize(tree,tree->data,ROUND_MAX);
210                      else                                               Lk(NULL,tree);
211                    }
212
213                  if(tree->mod->gamma_mgf_bl) Optimum_Root_Position_IL_Model(tree);
214
215                  Set_Both_Sides(YES,tree);
216                  Lk(NULL,tree);
217                  Pars(NULL,tree);
218                  Get_Tree_Size(tree);
219                  PhyML_Printf("\n\n. Log likelihood of the current tree: %f.",tree->c_lnL);
220
221                  //
222                  /* ML_Ancestral_Sequences(tree); */
223                  //
224
225                  Check_Br_Lens(tree);
226                  Br_Len_Involving_Invar(tree);
227                  Rescale_Br_Len_Multiplier_Tree(tree);
228                 
229#elif defined EVOLVE
230                  Evolve(tree->data,tree->mod,tree);
231                  Exit("\n");
232#endif
233
234
235                  if(!tree->n_root) Get_Best_Root_Position(tree);
236
237                  /* Print the tree estimated using the current random (or BioNJ) starting tree */
238                  if(io->mod->s_opt->n_rand_starts > 1)
239                    {
240                      Print_Tree(io->fp_out_trees,tree);
241                      fflush(NULL);
242                    }
243
244                  /* Record the most likely tree in a string of characters */
245                  if(tree->c_lnL > best_lnL)
246                    {
247                      best_lnL = tree->c_lnL;
248                      if(most_likely_tree) Free(most_likely_tree);
249                      most_likely_tree = Write_Tree(tree,NO);
250                    }
251
252                  time(&t_end);
253
254
255                  Print_Fp_Out(io->fp_out_stats,t_beg,t_end,tree,
256                               io,num_data_set+1,
257                               (tree->mod->s_opt->n_rand_starts > 1)?
258                               (num_rand_tree):(num_tree),
259                               (num_rand_tree == io->mod->s_opt->n_rand_starts-1)?(YES):(NO));
260                 
261                  if(tree->io->print_site_lnl) Print_Site_Lk(tree,io->fp_out_lk);
262
263
264                  /* Start from BioNJ tree */
265                  if((num_rand_tree == io->mod->s_opt->n_rand_starts-1) && (tree->mod->s_opt->random_input_tree))
266                    {
267                      /* Do one more iteration in the loop, but don't randomize the tree */
268                      num_rand_tree--;
269                      tree->mod->s_opt->random_input_tree = NO;
270                    }
271                 
272                  if(io->fp_in_constraint_tree != NULL) Free_Tree(io->cstr_tree);
273                  Free_Spr_List(tree);
274                  Free_Triplet(tree->triplet_struct);
275                  Free_Tree_Pars(tree);
276                  Free_Tree_Lk(tree);
277                  Free_Tree(tree);
278                }
279             
280              /* Launch bootstrap analysis */
281              if(mod->bootstrap) 
282                {
283                  if(!io->quiet) PhyML_Printf("\n\n. Launch bootstrap analysis on the most likely tree...");
284
285                  #ifdef MPI
286                  MPI_Bcast (most_likely_tree, strlen(most_likely_tree)+1, MPI_CHAR, 0, MPI_COMM_WORLD);
287                  if(!io->quiet)  PhyML_Printf("\n\n. The bootstrap analysis will use %d CPUs.",Global_numTask);
288                  #endif
289
290                  most_likely_tree = Bootstrap_From_String(most_likely_tree,cdata,mod,io);
291                }
292              else if(io->ratio_test) 
293                {
294                  /* Launch aLRT */
295                  most_likely_tree = aLRT_From_String(most_likely_tree,cdata,mod,io);
296                }
297
298              /* Print the most likely tree in the output file */
299              if(!io->quiet) PhyML_Printf("\n\n. Printing the most likely tree in file '%s'.", Basename(io->out_tree_file));
300              if(io->n_data_sets == 1) rewind(io->fp_out_tree);
301              PhyML_Fprintf(io->fp_out_tree,"%s\n",most_likely_tree);
302             
303              if(io->n_trees > 1 && io->n_data_sets > 1) break;
304            }
305          Free_Cseq(cdata);
306        }
307      else
308        {
309          PhyML_Printf("\n== No data was found.\n");
310          PhyML_Printf("\n== Err in file %s at line %d\n",__FILE__,__LINE__);
311          Warn_And_Exit("");
312        }
313      Free_Model_Complete(mod);
314    }
315 
316  if(most_likely_tree) Free(most_likely_tree);
317
318  if(mod->s_opt->n_rand_starts > 1) PhyML_Printf("\n. Best log likelihood: %f\n",best_lnL);
319
320  Free_Optimiz(mod->s_opt);
321  M4_Free_M4_Model(mod->m4mod);
322  Free_Model_Basic(mod);
323 
324  if(io->fp_in_constraint_tree) fclose(io->fp_in_constraint_tree);
325  if(io->fp_in_align)           fclose(io->fp_in_align);
326  if(io->fp_in_tree)            fclose(io->fp_in_tree);
327  if(io->fp_out_lk)             fclose(io->fp_out_lk);
328  if(io->fp_out_tree)           fclose(io->fp_out_tree);
329  if(io->fp_out_trees)          fclose(io->fp_out_trees);
330  if(io->fp_out_stats)          fclose(io->fp_out_stats);
331  if(io->fp_out_trace)          fclose(io->fp_out_trace);
332
333  Free_Input(io);
334
335  time(&t_end);
336  Print_Time_Info(t_beg,t_end);
337
338#ifdef MPI
339  MPI_Finalize();
340#endif
341
342  return 0;
343}
344
345#elif(M4)
346#include "m4.h"
347int main(int argc, char **argv)
348{
349  M4_main(argc, argv);
350  return 1;
351}
352
353#elif(PART)
354#include "mg.h"
355int main(int argc, char **argv)
356{
357  PART_main(argc, argv);
358  return 1;
359}
360
361#elif(PHYTIME)
362#include "times.h"
363int main(int argc, char **argv)
364{
365  TIMES_main(argc, argv);
366  return 1;
367}
368
369#elif(PHYCONT)
370#include "continuous.h"
371int main(int argc, char **argv)
372{
373  CONT_main(argc, argv);
374  return 1;
375}
376
377#elif(RF)
378int main(int argc, char **argv)
379{
380  t_tree *tree1, *tree2;
381  FILE *fp_tree1, *fp_tree2;
382  int i,j;
383
384  fp_tree1 = (FILE *)fopen(argv[1],"r");
385  fp_tree2 = (FILE *)fopen(argv[2],"r");
386
387  tree1 = Read_Tree_File_Phylip(fp_tree1);
388  tree2 = Read_Tree_File_Phylip(fp_tree2);
389
390  Match_Nodes_In_Small_Tree(tree1,tree2);
391
392  For(i,2*tree1->n_otu-2)
393    {
394      printf("\n. Node %d in tree1 matches node %d in tree2",i,(tree1->noeud[i]->match_node)?(tree1->noeud[i]->match_node->num):(-1));
395    }
396
397
398
399
400/*   t_tree *tree1, *tree2; */
401/*   FILE *fp_tree1, *fp_tree2; */
402/*   int i,j,rf,n_edges,n_common,bip_size; */
403/*   phydbl thresh; */
404/*   t_edge *b; */
405
406
407/*   fp_tree1 = (FILE *)fopen(argv[1],"r"); */
408/*   fp_tree2 = (FILE *)fopen(argv[2],"r"); */
409/*   thresh = (phydbl)atof(argv[3]); */
410
411/*   tree1 = Read_Tree_File(fp_tree1); */
412/*   tree2 = Read_Tree_File(fp_tree2); */
413
414/*   Get_Rid_Of_Prefix('_',tree1); */
415
416/* /\*   Find_Common_Tips(tree1,tree2); *\/ */
417
418/*   Alloc_Bip(tree1); */
419/*   Alloc_Bip(tree2); */
420
421/*   Get_Bip(tree1->noeud[0],tree1->noeud[0]->v[0],tree1); */
422/*   Get_Bip(tree2->noeud[0],tree2->noeud[0]->v[0],tree2); */
423 
424/* /\*   PhyML_Printf("\n. rf=%f\n",Compare_Bip_On_Existing_Edges(thresh,tree1,tree2)); *\/ */
425/*   For(i,2*tree1->n_otu-3) tree1->a_edges[i]->bip_score = 0; */
426/*   For(i,2*tree2->n_otu-3) tree2->a_edges[i]->bip_score = 0; */
427 
428/*   rf = 0; */
429/*   n_edges = 0; */
430
431/*   /\* First tree *\/ */
432/*   For(i,2*tree1->n_otu-3)  */
433/*     { */
434/*       /\* Consider the branch only if the corresponding bipartition has size > 1 *\/ */
435/*       b = tree1->a_edges[i]; */
436/*       bip_size = MIN(b->left->bip_size[b->l_r],b->rght->bip_size[b->r_l]); */
437         
438/*       if(bip_size > 1) */
439/*      { */
440/*        /\* with non-zero length *\/ */
441/*        if(tree1->a_edges[i]->l > thresh)   */
442/*          { */
443/*            n_edges++; */
444/*            /\* This t_edge is not found in tree2 *\/ */
445/*            if(!tree1->a_edges[i]->bip_score) rf++; ; */
446/*          } */
447/*      } */
448/*     } */
449
450
451/*   /\* Second tree *\/ */
452/*   For(i,2*tree2->n_otu-3)  */
453/*     { */
454/*       b = tree2->a_edges[i]; */
455/*       bip_size = MIN(b->left->bip_size[b->l_r],b->rght->bip_size[b->r_l]); */
456
457/*       if(bip_size > 1) */
458/*      { */
459/*        if(tree2->a_edges[i]->l > thresh)   */
460/*          { */
461/*            n_edges++; */
462/*            /\* This t_edge is not found in tree1 *\/ */
463/*            if(!tree2->a_edges[i]->bip_score) rf++; ; */
464/*          } */
465/*      } */
466/*     } */
467
468/*   if(!n_edges) */
469/*     { */
470/*       Exit("\n. No comparable internal edges were found.\n"); */
471/*     } */
472/*   else */
473/*     { */
474/*       PhyML_Printf("\n. Robinson and Foulds distance: %f.",(double)rf/(n_edges)); */
475/* /\*       PhyML_Printf("\n. %d internal edges were processed (%d in the first tree, %d in the second).\n",n_edges,n_edges_t1,n_edges-n_edges_t1); *\/ */
476/*       PhyML_Printf("\n"); */
477/*     } */
478
479  return 1;
480}
481
482#elif(TIPORDER)
483#include "tiporder.h"
484int main(int argc, char **argv)
485{
486  TIPO_main(argc, argv);
487  return 1;
488}
489
490#elif(TEST)
491#include "xml.h"
492int main(int argc, char **argv)
493{
494}
495
496#elif(SERGEII)
497#include "sergeii.h"
498int main(int argc, char **argv)
499{
500  /*My_Function(argc, argv);*/
501  /* PhyTime_XML(argc, argv); */
502  Get_Input(argc,argv);
503  return 1;
504}
505
506#elif(GEO)
507#include "geo.h"
508int main(int argc, char **argv)
509{
510  GEO_Main(argc,argv);
511  return 1;
512}
513
514#elif(CHECKPOINT)
515#include "checkpoint.h"
516int main(int argc, char **argv)
517{
518  CHECK_Main(argc,argv);
519  return 1;
520}
521
522#endif
523
Note: See TracBrowser for help on using the repository browser.