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

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

added most recent version of phyml

File size: 73.1 KB
Line 
1#include "spr.h"
2#include "utilities.h"
3#include "lk.h"
4#include "optimiz.h"
5#include "bionj.h"
6#include "models.h"
7#include "free.h"
8#include "help.h"
9#include "simu.h"
10#include "eigen.h"
11#include "pars.h"
12#include "alrt.h"
13#include "mixt.h"
14#include "sergeii.h"
15#ifdef MPI
16#include "mpi_boot.h"
17#endif
18#include <stdlib.h>
19
20
21
22//////////////////////////////////////////////////////////////
23//////////////////////////////////////////////////////////////
24void PhyTime_XML(char *xml_file)
25{
26
27  FILE *f;
28  char **clade, *clade_name, **mon_list;
29  phydbl low, up;
30  int i, j, n_taxa, clade_size, node_num, n_mon; //rnd_num
31  xml_node *n_r, *n_t, *n_m, *n_cur;
32  t_cal *last_calib; 
33  /* t_cal *cur; */
34  align **data, **c_seq;
35  option *io;
36  calign *cdata;
37  t_opt *s_opt;
38  t_mod *mod;
39  time_t t_beg,t_end;
40  int r_seed;
41  char *most_likely_tree;
42  int user_lk_approx;
43  t_tree *tree;
44  t_node **a_nodes; //*node;
45  m4 *m4mod;
46 
47  srand(time(NULL));
48
49  i = 0;
50  j = 0;
51  last_calib       = NULL;
52  mod              = NULL;
53  most_likely_tree = NULL;
54  n_taxa = 0; 
55  node_num = -1;
56
57  //file can/cannot be open:
58  if ((f =(FILE *)fopen(xml_file, "r")) == NULL)
59    {
60      PhyML_Printf("\n== File '%s' can not be opened...\n",xml_file);
61      Exit("\n");
62    }
63
64  n_r = XML_Load_File(f);
65
66  //XML_Search_Node_Attribute_Value_Clade("id", "clade1", NO, n_r -> child);
67  //Exit("\n");
68
69  //memory allocation for model parameters:
70  io = (option *)Make_Input();
71  mod   = (t_mod *)Make_Model_Basic();
72  s_opt = (t_opt *)Make_Optimiz();
73  m4mod = (m4 *)M4_Make_Light();
74  Set_Defaults_Input(io);                                                                         
75  Set_Defaults_Model(mod);
76  Set_Defaults_Optimiz(s_opt); 
77  io -> mod = mod;
78  mod = io -> mod;
79  mod -> s_opt = s_opt;
80  clade_size = -1;
81
82  ////////////////////////////////////////////////////////////////////////////
83  //////////////////////reading tree topology:////////////////////////////////
84
85  //looking for a node <topology>
86  n_t = XML_Search_Node_Name("topology", YES, n_r);
87
88  //setting tree:
89  tree        = (t_tree *)mCalloc(1,sizeof(t_tree));
90  n_cur = XML_Search_Node_Name("instance", YES, n_t);
91  if(n_cur != NULL)
92    {
93      if(XML_Search_Attribute(n_cur, "user.tree") != NULL) 
94        {
95          strcpy(io -> out_tree_file, XML_Search_Attribute(n_cur, "user.tree") -> value); 
96          io -> fp_out_tree  = Openfile(io -> out_tree_file, 1);
97          io -> tree = Read_Tree_File_Phylip(io -> fp_in_tree);
98        }
99      else
100        {
101          PhyML_Printf("\n==Tree was not found. \n");
102          PhyML_Printf("\n==Either specify tree file name or enter the whole tree. \n");
103          Exit("\n");
104        }
105    }
106  else io -> tree  = Read_Tree(&n_t -> value);
107  io -> n_otu = io -> tree -> n_otu;
108  tree = io -> tree;
109
110
111  //setting initial values to n_calib:
112  For(i, 2 * tree -> n_otu - 2) 
113    {
114      tree -> a_nodes[i] -> n_calib = 0;
115      //PhyML_Printf("\n. '%d' \n", tree -> a_nodes[i] -> n_calib);
116    }
117
118
119  ////////////////////////////////////////////////////////////////////////////
120  ////////////////////////////////////////////////////////////////////////////
121
122  //memory for nodes:
123  a_nodes   = (t_node **)mCalloc(2 * io -> n_otu - 1,sizeof(t_node *));
124  For(i, 2 * io -> n_otu - 2) a_nodes[i] = (t_node *)mCalloc(1,sizeof(t_node));
125 
126  //setting a model:
127  tree -> rates = RATES_Make_Rate_Struct(io -> n_otu);                                   
128  RATES_Init_Rate_Struct(tree -> rates, io -> rates, tree -> n_otu);
129
130  //reading seed:
131  if(XML_Search_Attribute(n_r, "seed")) io -> r_seed = String_To_Dbl(XML_Search_Attribute(n_r, "seed") -> value);
132   
133  //TO DO: check that the tree has a root...
134  Update_Ancestors(io -> tree -> n_root, io -> tree -> n_root -> v[2], io -> tree);
135  Update_Ancestors(io -> tree -> n_root, io -> tree -> n_root -> v[1], io -> tree);
136                 
137
138  ////////////////////////////////////////////////////////////////////////////
139  //////////////////////memory allocation for temp parameters/////////////////
140
141  //memory for monitor flag:
142  io -> mcmc -> monitor = (int *)mCalloc(2 * io -> n_otu - 1,sizeof(int));
143
144 
145  //memory for sequences:
146  n_cur = XML_Search_Node_Name("alignment", YES, n_r);
147
148  data   = (align **)mCalloc(io -> n_otu,sizeof(align *));
149  For(i, io -> n_otu)
150    {
151      data[i]          = (align *)mCalloc(1,sizeof(align));
152      data[i] -> name  = (char *)mCalloc(T_MAX_NAME,sizeof(char));
153      if(n_cur -> child -> value != NULL) data[i] -> state = (char *)mCalloc(strlen(n_cur -> child -> value) + 1,sizeof(char));
154      else data[i] -> state = (char *)mCalloc(T_MAX_SEQ,sizeof(char));
155    }
156  io -> data = data;
157  //tree -> data = data;
158
159 
160  //memory for clade:
161  clade_name = (char *)mCalloc(T_MAX_NAME,sizeof(char));
162  clade = (char **)mCalloc(tree -> n_otu, sizeof(char *));
163  For(i, tree -> n_otu)
164    {
165      clade[i] = (char *)mCalloc(T_MAX_NAME,sizeof(char));
166    }
167
168  //memory for list of clades to be monitored
169  mon_list = (char **)mCalloc(T_MAX_FILE,sizeof(char *));
170  For(i, T_MAX_FILE)
171    {
172      mon_list[i] = (char *)mCalloc(T_MAX_NAME,sizeof(char));
173    }
174
175 ////////////////////////////////////////////////////////////////////////////
176 ////////////////////////////////////////////////////////////////////////////
177 
178  //reading monitor node:
179  i = 0;
180  n_m = XML_Search_Node_Name("monitor", YES, n_r);
181  if(n_m != NULL)
182    {
183      do
184        {
185          strcpy(mon_list[i], n_m -> child -> attr -> value);
186          i++;
187          if(n_m -> child) n_m -> child = n_m -> child -> next;
188          else break;
189        }
190      while(n_m -> child);
191      n_mon = i;
192    }
193  else
194    {
195      n_mon = 0;
196      PhyML_Printf("\n. There is no clade to be monitored. \n");
197    }
198
199 ////////////////////////////////////////////////////////////////////////////
200 ////////////////////////////////////////////////////////////////////////////
201
202  //chekcing for calibration node (upper or lower bound) to exist:
203  n_cur = XML_Search_Node_Name("calibration", YES, n_r);
204
205  if(n_cur == NULL)
206    {
207      PhyML_Printf("\n==There is no calibration information provided. \n");
208      PhyML_Printf("\n==Please check your data. \n");
209      Exit("\n");
210    }
211  else
212    {
213      if(XML_Search_Node_Name("upper", NO, n_cur -> child) == NULL && XML_Search_Node_Name("lower", NO, n_cur -> child) == NULL)
214        {
215          PhyML_Printf("\n==There is no calibration information provided. \n");
216          PhyML_Printf("\n==Please check your data. \n");
217          Exit("\n");
218        }
219    }
220
221 ////////////////////////////////////////////////////////////////////////////
222 ////////////////////////////////////////////////////////////////////////////
223  if(XML_Search_Attribute(n_r, "use_data") != NULL) 
224    {
225      if(XML_Search_Attribute(n_r, "use_data") -> value != NULL && (!strcmp(XML_Search_Attribute(n_r, "use_data") -> value, "YES")))
226        io -> mcmc -> use_data = YES;
227      if(XML_Search_Attribute(n_r, "use_data") -> value != NULL && (!strcmp(XML_Search_Attribute(n_r, "use_data") -> value, "NO")))
228        io -> mcmc -> use_data = NO;
229    }
230
231  n_r = n_r -> child; 
232  tree -> rates -> tot_num_cal = 0;
233  do
234    {
235      if(!strcmp(n_r -> name, "alignment"))//looking for a node <alignment>.
236        {
237          if(!n_r -> attr -> value)
238          {
239                PhyML_Printf("\n==Not found sequence type (nt / aa). \n");
240                PhyML_Printf("\n==Please, include data to node <%s> attribute value. \n", n_r -> name);
241                Exit("\n");
242          }     
243          if(!strcmp(To_Upper_String(n_r -> attr -> value), "NT")) 
244            {
245              io -> datatype = 0;
246              io -> mod -> ns = 4;
247            }
248          if(!strcmp(To_Upper_String(n_r -> attr -> value), "AA")) 
249            {
250              io -> datatype = 1;
251              io -> mod -> ns = 20;
252            }
253
254          n_cur = XML_Search_Node_Name("instance", YES, n_r);
255          if(n_cur != NULL)
256            {
257              if(XML_Search_Attribute(n_cur, "user.alignment") != NULL) 
258                {
259                  strcpy(io -> in_align_file, XML_Search_Attribute(n_cur, "user.alignment") -> value); 
260                  io -> fp_in_align  = Openfile(io -> in_align_file, 1);
261                  Detect_Align_File_Format(io);
262                  io -> data = Get_Seq(io);
263                }
264            }
265          else
266            {   
267              i = 0;
268              do
269                {
270                  strcpy(io -> in_align_file, "sergeii"); 
271                  strcpy(io -> data[i] -> name, n_r -> child -> attr -> value);
272                  strcpy(io -> data[i] -> state, To_Upper_String(n_r -> child -> value));
273                  i++;
274                  if(n_r -> child -> next) n_r -> child = n_r -> child -> next;
275                  else break;
276                }
277              while(n_r -> child);
278              n_taxa = i;
279            }
280
281          //checking if a sequences of the same lengths:
282
283          i = 1;
284          For(i, n_taxa) if(strlen(io -> data[0] -> state) != strlen(io -> data[i] -> state))
285            {
286              printf("\n. Sequences are of different length. Please check your data...\n");
287              Exit("\n");
288              break;
289            } 
290
291          //checking sequence names:
292          i = 0;
293          For(i, n_taxa) Check_Sequence_Name(io -> data[i] -> name);
294
295          //check if a number of tips is equal to a number of taxa:
296          if(n_taxa != io -> n_otu)
297            {
298              PhyML_Printf("\n==Number of taxa is not the same as a number of tips. Check your data...\n");
299              Exit("\n");
300            }
301
302          //deleting '-', etc. from sequences:
303          io -> data[0] -> len = strlen(io -> data[0] -> state);
304          Post_Process_Data(io);   
305          n_r = n_r -> next;
306        }
307      else if(!strcmp(n_r -> name, "calibration"))//looking for a node <calibration>.
308        {
309          tree -> rates -> tot_num_cal++;
310          if (tree -> rates -> calib == NULL) tree -> rates -> calib = Make_Calib(tree -> n_otu);
311          if(last_calib)
312            {
313              last_calib -> next = tree -> rates -> calib;
314              tree -> rates -> calib -> prev = last_calib;
315            }
316          last_calib = tree -> rates -> calib;
317
318          low = -BIG;
319          up  = BIG;
320          n_cur = XML_Search_Node_Name("lower", YES, n_r);
321          if(n_cur != NULL) low = String_To_Dbl(n_cur -> value); 
322          n_cur = XML_Search_Node_Name("upper", YES, n_r);
323          if(n_cur != NULL) up = String_To_Dbl(n_cur -> value);
324          do
325            {
326              if(!strcmp("appliesto", n_r -> child -> name)) 
327                {
328                  //case of internal node:
329                  strcpy(clade_name, n_r -> child -> attr -> value);//reached clade names n_r -> child -> attr -> value
330                  if(!strcmp("root", clade_name))
331                    {
332                      node_num = io -> tree -> n_root -> num;
333                      //printf("\n. Node number [%d] \n", node_num);
334                    }
335                  else if(strcmp("NO_CLADE", clade_name) && strcmp("root", clade_name))
336                    {
337                      xml_node *n_clade, *nd2;
338                      nd2 = n_r -> parent;
339                      /* n_clade = XML_Search_Node_Attribute_Value_Clade("id", clade_name, NO, n_r -> parent -> child); */
340                      n_clade = XML_Search_Node_Generic("clade", "id", clade_name, YES, nd2);
341                      if(n_clade) //found clade with a given name
342                        {
343                          i = 0;
344                          xml_node *nd;
345                          nd = n_clade -> child;
346                          /* do */
347                          /*   { */
348                          /*     strcpy(clade[i], nd -> attr -> value);  */
349                          /*     i++; */
350                          /*     nd = nd -> next; */
351                          /*     if(!nd) break; */
352                          /*   } */
353                          /* while(n_clade -> child); */
354                          /* clade_size = i; */
355                          clade = XML_Reading_Clade(nd, tree);
356                          clade_size = XML_Number_Of_Taxa_In_Clade(nd);
357                          node_num = Find_Clade(clade, clade_size, io -> tree);
358                          /* printf("\n. Clade size [%d] \n", clade_size); */
359                          /* printf("\n. Node number [%d] \n", node_num); */
360                         }
361                      else
362                        {
363                          PhyML_Printf("==Calibration information on the clade [%s] was not found. \n", clade_name);
364                          PhyML_Printf("\n. Err in file %s at line %d\n",__FILE__,__LINE__);
365                          Exit("\n");
366                        }
367                    }
368                  if(strcmp("NO_CLADE", clade_name))
369                    {
370                      For(j, n_mon)
371                        {
372                          if(!strcmp(clade_name, mon_list[j])) io -> mcmc -> monitor[node_num] = YES;
373                        }
374                    }
375                  /* For(i, clade_size) PhyML_Printf("\n. Clade name [%s] Taxon name: [%s]", clade_name, clade[i]); */
376                  if(strcmp("NO_CLADE", clade_name)) 
377                    {
378                      tree -> rates -> calib -> proba[node_num] = String_To_Dbl(n_r -> child -> attr -> next -> value);
379                     
380                      if(!n_r -> child -> attr -> next && n_r -> child -> next == NULL) tree -> rates -> calib -> proba[node_num] = 1.;
381                      if(!n_r -> child -> attr -> next && n_r -> child -> next)
382                        {
383                          PhyML_Printf("==You either need to provide information about probability with which calibration \n");
384                          PhyML_Printf("==applies to a node or you need to apply calibartion only to one node. \n");
385                          PhyML_Printf("\n. Err in file %s at line %d\n",__FILE__,__LINE__);
386                          Exit("\n");
387                        }
388                     
389                      tree -> rates -> calib -> all_applies_to[tree -> rates -> calib -> n_all_applies_to] -> num = node_num;
390                    }
391                  else tree -> rates -> calib -> proba[2 * tree -> n_otu - 1] = String_To_Dbl(n_r -> child -> attr -> next -> value);
392                  tree -> rates -> calib -> n_all_applies_to++;
393                  tree -> rates -> calib -> lower = low;
394                  tree -> rates -> calib -> upper = up;
395                  /* printf("\n. Porbability [%f] \n", String_To_Dbl(n_r -> child -> attr -> next -> value)); */
396                 
397                  /////////////////////////////////////////////////////////////////////////////////////////////////////
398                  PhyML_Printf("\n. .......................................................................");
399                  PhyML_Printf("\n");
400                  if(strcmp(clade_name, "NO_CLADE")) PhyML_Printf("\n. Clade name: [%s]", clade_name);
401                  else 
402                    {
403                      PhyML_Printf("\n. Calibration with:");
404                      PhyML_Printf("\n. Lower bound set to: %15f time units.", low);
405                      PhyML_Printf("\n. Upper bound set to: %15f time units.", up);
406                      PhyML_Printf("\n. DOES NOT apply with probability [%f]", String_To_Dbl(n_r -> child -> attr -> next -> value));
407                      PhyML_Printf("\n. .......................................................................");
408                    }
409       
410                  if(strcmp(clade_name, "root") && strcmp(clade_name, "NO_CLADE"))
411                    {
412                      For(i, clade_size) PhyML_Printf("\n. Taxon name: [%s]", clade[i]);
413                    }
414                  if(strcmp(clade_name, "NO_CLADE")) 
415                    {
416                      PhyML_Printf("\n. Node number to which calibration applies to is [%d] with probability [%f]", node_num, String_To_Dbl(n_r -> child -> attr -> next -> value));
417                      PhyML_Printf("\n. Lower bound set to: %15f time units.", low);
418                      PhyML_Printf("\n. Upper bound set to: %15f time units.", up);
419                      PhyML_Printf("\n. .......................................................................");
420                    }
421                  /////////////////////////////////////////////////////////////////////////////////////////////////////
422                  if(n_r -> child -> next) n_r -> child = n_r -> child -> next;
423                  else break;   
424                }
425              else if(n_r -> child -> next) n_r -> child = n_r -> child -> next;
426              else break; 
427            }
428          while(n_r -> child);     
429          //PhyML_Printf("\n. '%d'\n", tree -> rates -> calib -> n_all_applies_to);
430          tree -> rates -> calib = tree -> rates -> calib -> next;         
431          n_r = n_r -> next;
432        }     
433      else if(!strcmp(n_r -> name, "ratematrices"))//initializing rate matrix:
434        {
435          if(n_r -> child) 
436            {
437              Make_Ratematrice_From_XML_Node(n_r -> child, io, mod);
438              n_r = n_r -> next;
439            } 
440          else n_r = n_r -> next;
441        }
442      else if(!strcmp(n_r -> name, "equfreqs"))//initializing frequencies:
443        {
444          if(n_r -> child) 
445            {
446               Make_Efrq_From_XML_Node(n_r -> child , io, mod);
447               n_r = n_r -> next;
448            }
449          else n_r = n_r -> next;
450        }
451      else if(!strcmp(n_r -> name, "siterates"))//initializing site rates:
452        {
453          if(n_r -> child) 
454            {
455              Make_RAS_From_XML_Node(n_r, io -> mod);
456              n_r = n_r -> next;
457            }
458          else n_r = n_r -> next;
459        }
460      else if (n_r -> next) n_r = n_r -> next;
461      else break;
462    }
463  while(1);
464 
465  tree -> rates -> calib = last_calib;
466  while(tree -> rates -> calib -> prev) tree -> rates -> calib = tree -> rates -> calib -> prev;
467  ////////////////////////////////////////////////////////////////////////////////////////////////
468  //Check for the sum of probabilities for one calibration add up to one
469  do
470    {
471      phydbl p = 0.0;
472      for(i = tree -> n_otu; i < 2 * tree -> n_otu; i++) 
473        {
474          p = p + tree -> rates -> calib -> proba[i];
475          /* PhyML_Printf("\n. # applies to %d \n", tree -> rates -> calib -> n_all_applies_to); */
476          /* PhyML_Printf("\n. %f \n", tree -> rates -> calib -> proba[i]); */
477        }
478      if(!Are_Equal(p, 1.0, 1.E-10))
479        {
480          PhyML_Printf("\n. .......................................................................");
481          PhyML_Printf("\n. WARNING! The sum of the probabilities for the calibration with:");
482          PhyML_Printf("\n. Lower bound set to: %15f time units.", tree -> rates -> calib -> lower);
483          PhyML_Printf("\n. Upper bound set to: %15f time units.", tree -> rates -> calib -> upper);
484          PhyML_Printf("\n. IS NOT equal to 1.");
485          PhyML_Printf("\n. The probability of NOT applying this calibration will be set to [%f].", 1.0 - p);
486          PhyML_Printf("\n. .......................................................................");
487          tree -> rates -> calib -> proba[2 * tree -> n_otu - 1] = 1.0 - p;
488          tree -> rates -> calib -> n_all_applies_to++;
489          /* PhyML_Printf("\n==You need to provide proper probabilities for the calibration. \n"); */
490          /* PhyML_Printf("\n. Err in file %s at line %d\n",__FILE__,__LINE__); */
491          /* Exit("\n"); */
492        } 
493      if(tree -> rates -> calib -> next) tree -> rates -> calib = tree -> rates -> calib -> next;
494      else break;
495    }
496  while(tree -> rates -> calib);
497  while(tree -> rates -> calib -> prev) tree -> rates -> calib = tree -> rates -> calib -> prev;
498  //////////////////////////////////////////////////////////////////////////////////////////////// 
499  //Set_Current_Calibration(0, tree);
500  //TIMES_Set_All_Node_Priors(tree);
501  //Slicing_Calibrations(tree);
502  ////////////////////////////////////////////////////////////////////////////
503  ////////////////////////////////////////////////////////////////////////////   
504  //clear memory:
505  free(clade_name);
506  For(i, tree -> n_otu)
507    {
508      free(clade[i]);
509    }
510  free(clade);
511  For(i, T_MAX_FILE)
512    {
513      free(mon_list[i]);
514    }
515  free(mon_list);
516
517  //Exit("\n");
518  ////////////////////////////////////////////////////////////////////////////
519  ////////////////////////////////////////////////////////////////////////////   
520  //START analysis:
521  r_seed = (io -> r_seed < 0)?(time(NULL)):(io -> r_seed);
522  srand(r_seed); 
523  rand();       
524  PhyML_Printf("\n. Seed: %d\n", r_seed);
525  PhyML_Printf("\n. Pid: %d\n",getpid()); 
526  PhyML_Printf("\n. Compressing sequences...\n");
527  data = io -> data;
528  data[0] -> len = strlen(data[0] -> state); 
529  ////////////////////////////////////////////////////////////////////////////
530  //memory for compressed sequences:
531  cdata         = (calign *)mCalloc(1,sizeof(calign));
532  c_seq   = (align **)mCalloc(io -> n_otu,sizeof(align *));
533  For(i, io -> n_otu)
534    {
535      c_seq[i]          = (align *)mCalloc(1,sizeof(align));
536      c_seq[i] -> name  = (char *)mCalloc(T_MAX_NAME,sizeof(char));
537      c_seq[i] -> state = (char *)mCalloc(data[0] -> len + 1,sizeof(char));
538    }
539  cdata -> c_seq = c_seq;
540  ////////////////////////////////////////////////////////////////////////////
541  cdata = Compact_Data(data, io);
542  Free_Seq(io -> data, cdata -> n_otu);
543  io -> mod -> io = io;
544  Check_Ambiguities(cdata, io -> mod -> io -> datatype, io -> state_len);           
545  Make_Model_Complete(mod);                           
546  Init_Model(cdata, mod, io);
547  if(io -> mod -> use_m4mod) M4_Init_Model(mod -> m4mod, cdata, mod);
548  time(&(t_beg));
549
550  tree -> mod         = mod;                                                                   
551  tree -> io          = io;
552  tree -> data        = cdata;
553  tree -> n_pattern   = tree -> data -> crunch_len / tree -> io -> state_len;
554 
555  Set_Both_Sides(YES, tree);
556  Prepare_Tree_For_Lk(tree);
557
558  //calculate the probabilities of each combination of calibrations:
559  TIMES_Calib_Partial_Proba(tree);
560  int cal_numb = 0;
561  do
562    {
563      if(!Are_Equal(tree -> rates -> times_partial_proba[cal_numb], 0.0, 1.E-10)) break;
564      else cal_numb += 1;
565    }
566  while(1); 
567  /* printf("\n. Calib number [%d] \n", cal_numb); */
568  Set_Current_Calibration(cal_numb, tree);
569  int tot_num_comb;
570  tot_num_comb = Number_Of_Comb(tree -> rates -> calib);                                                                                                               
571  PhyML_Printf("\n. The total number of calibration combinations is going to be considered is %d.\n", tot_num_comb);
572  TIMES_Set_All_Node_Priors(tree);
573
574
575  //set initial value for Hastings ratio for conditional jump:
576  tree -> rates -> c_lnL_Hastings_ratio = 0.0;
577 
578  TIMES_Get_Number_Of_Time_Slices(tree);
579  TIMES_Label_Edges_With_Calibration_Intervals(tree);
580
581
582  tree -> write_br_lens = NO;   
583
584                                                                       
585  PhyML_Printf("\n. Input tree with calibration information ('almost' compatible with MCMCtree).\n");
586  PhyML_Printf("\n. %s \n", Write_Tree(tree, YES));
587
588  tree -> write_br_lens = YES;
589
590
591  // Work with log of branch lengths?
592  if(tree -> mod -> log_l == YES) Log_Br_Len(tree);
593
594  if(io -> mcmc -> use_data == YES)                                                                                                                             
595    { 
596      // Force the exact likelihood score
597      user_lk_approx = tree -> io -> lk_approx;                                                                                                 
598      tree -> io -> lk_approx = EXACT;
599                     
600      // MLE for branch lengths                                                                                                                                                       
601      /* printf("\n. %s",Write_Tree(tree,NO)); */
602      /* printf("\n. alpha %f",tree->mod->ras->alpha->v); */
603      /* printf("\n.  %f %f %f %f",tree->mod->e_frq->pi->v[0],tree->mod->e_frq->pi->v[1],tree->mod->e_frq->pi->v[2],tree->mod->e_frq->pi->v[3]); */
604      /* Lk(NULL,tree); */
605      /* printf("\n. %f",tree->c_lnL); */
606      /* Exit("\n"); */
607
608      PhyML_Printf("\n");
609      Round_Optimize(tree, tree -> data, ROUND_MAX);
610                     
611      // Set vector of mean branch lengths for the Normal approximation of the likelihood
612      RATES_Set_Mean_L(tree);
613     
614      // Estimate the matrix of covariance for the Normal approximation of the likelihood
615      PhyML_Printf("\n");
616      PhyML_Printf("\n. Computing Hessian...");                                                                                             
617      tree -> rates -> bl_from_rt = 0;                                                                                                                                         
618      Free(tree -> rates -> cov_l);                                                                                                                                                     
619      tree -> rates -> cov_l = Hessian_Seo(tree);
620                                                                                                                         
621      // tree->rates->cov_l = Hessian_Log(tree);
622      For(i, (2 * tree -> n_otu - 3) * (2 * tree -> n_otu - 3)) tree -> rates -> cov_l[i] *= -1.0;
623      if(!Iter_Matinv(tree -> rates -> cov_l, 2 * tree -> n_otu - 3, 2 * tree -> n_otu - 3, YES)) Exit("\n");                   
624      tree -> rates -> covdet = Matrix_Det(tree -> rates -> cov_l, 2 * tree -> n_otu - 3, YES);                                                 
625      For(i,(2 * tree -> n_otu - 3) * (2 * tree -> n_otu - 3)) tree -> rates -> invcov[i] = tree -> rates -> cov_l[i]; 
626      if(!Iter_Matinv(tree -> rates -> invcov, 2 * tree -> n_otu - 3, 2 * tree -> n_otu - 3, YES)) Exit("\n");
627      tree -> rates -> grad_l = Gradient(tree);                                 
628
629      // Pre-calculation of conditional variances to speed up calculations                                     
630      RATES_Bl_To_Ml(tree);
631      RATES_Get_Conditional_Variances(tree);
632      RATES_Get_All_Reg_Coeff(tree);
633      RATES_Get_Trip_Conditional_Variances(tree);
634      RATES_Get_All_Trip_Reg_Coeff(tree);
635                     
636      Lk(NULL, tree);
637      PhyML_Printf("\n");
638      PhyML_Printf("\n. p(data|model) [exact ] ~ %.2f",tree -> c_lnL);                                                                                                                                                 
639      tree -> io -> lk_approx = NORMAL;                                                                                                                 
640      For(i,2 * tree -> n_otu - 3) tree -> rates -> u_cur_l[i] = tree -> rates -> mean_l[i] ;                           
641      tree -> c_lnL = Lk(NULL,tree);
642      PhyML_Printf("\n. p(data|model) [approx] ~ %.2f",tree -> c_lnL);
643
644      tree -> io -> lk_approx = user_lk_approx; 
645                                                                               
646    }
647
648 
649  tree -> rates -> model = io -> rates -> model;
650                                                                                                         
651  PhyML_Printf("\n. Selected model '%s' \n", RATES_Get_Model_Name(io -> rates -> model));
652
653  if(tree -> rates -> model == GUINDON) tree -> mod -> gamma_mgf_bl = YES;
654                                                                                                               
655  tree -> rates -> bl_from_rt = YES;
656
657  if(tree -> io -> cstr_tree) Find_Surviving_Edges_In_Small_Tree(tree, tree -> io -> cstr_tree);                                                                                                                                                                 
658  time(&t_beg);
659
660  tree -> mcmc = MCMC_Make_MCMC_Struct();
661 
662  MCMC_Copy_MCMC_Struct(tree -> io -> mcmc, tree -> mcmc, "phytime"); 
663
664  tree -> mod -> m4mod = m4mod;
665               
666  MCMC_Complete_MCMC(tree -> mcmc, tree);
667
668  tree -> mcmc -> is_burnin = NO;
669
670  //PhyML_Printf("\n");
671  //PhyML_Printf("\n. Computing Normalizing Constant(s) for the Node Times Prior Density...\n");
672  //tree -> K = Norm_Constant_Prior_Times(tree);
673  //Exit("\n");
674
675  MCMC(tree);   
676                                                                                                                                                               
677  MCMC_Close_MCMC(tree -> mcmc);                                                                                                                                       
678  MCMC_Free_MCMC(tree -> mcmc);                                                                                                                                         
679  PhyML_Printf("\n");   
680 
681  Free_Tree_Pars(tree);
682  Free_Tree_Lk(tree);
683  Free_Tree(tree);
684  Free_Cseq(cdata);
685  Free_Model(mod);
686  if(io -> fp_in_align)   fclose(io -> fp_in_align);
687  if(io -> fp_in_tree)    fclose(io -> fp_in_tree);
688  if(io -> fp_out_lk)     fclose(io -> fp_out_lk);
689  if(io -> fp_out_tree)   fclose(io -> fp_out_tree);
690  if(io -> fp_out_trees)  fclose(io -> fp_out_trees);
691  if(io -> fp_out_stats)  fclose(io -> fp_out_stats);
692  fclose(f);
693  Free(most_likely_tree);
694  Free_Input(io);
695  Free_Calib(tree -> rates -> calib);
696  time(&t_end);
697  Print_Time_Info(t_beg,t_end); 
698       
699  /* return 1;    */
700}
701
702
703//////////////////////////////////////////////////////////////
704//////////////////////////////////////////////////////////////
705//Calculate the prior probability for node times taking into account the
706//probailitis with which each calibration applies to the particular node.
707phydbl TIMES_Calib_Cond_Prob(t_tree *tree)
708{
709
710  phydbl times_lk, *Yule_val, *times_partial_proba, times_tot_proba, *t_prior_min, *t_prior_max, c, constant, ln_t; 
711  short int *t_has_prior;
712  int i, j, k, tot_num_comb;
713  t_cal *calib;
714 
715
716  times_tot_proba = 0.0;
717  calib = tree -> rates -> calib;
718  t_prior_min = tree -> rates -> t_prior_min;
719  t_prior_max = tree -> rates -> t_prior_max;
720  t_has_prior = tree -> rates -> t_has_prior;
721  times_partial_proba = tree -> rates -> times_partial_proba;
722  /* constant = tree -> K; */
723
724  tot_num_comb = Number_Of_Comb(calib);
725
726 
727  Yule_val = (phydbl *)mCalloc(tot_num_comb, sizeof(phydbl));   
728  /* times_partial_proba = (phydbl *)mCalloc(tot_num_comb, sizeof(phydbl));    */
729
730  For(i, tot_num_comb)
731    {
732      for(j = tree -> n_otu; j < 2 * tree -> n_otu - 1; j++) 
733        {
734          t_prior_min[j] = -BIG;
735          t_prior_max[j] = BIG;
736          t_has_prior[j] = NO; 
737        }
738      do
739        {
740          k = (i % Number_Of_Comb(calib)) / Number_Of_Comb(calib -> next);
741          if(calib -> all_applies_to[k] -> num)
742            {
743              t_prior_min[calib -> all_applies_to[k] -> num] = MAX(t_prior_min[calib -> all_applies_to[k] -> num], calib -> lower);
744              t_prior_max[calib -> all_applies_to[k] -> num] = MIN(t_prior_max[calib -> all_applies_to[k] -> num], calib -> upper);
745              t_has_prior[calib -> all_applies_to[k] -> num] = YES;
746              /* if((t_prior_min[calib -> all_applies_to[k] -> num] > t_prior_max[calib -> all_applies_to[k] -> num])) times_partial_proba[i] = 0.0;  */
747              /* else times_partial_proba[i] *= calib -> proba[calib -> all_applies_to[k] -> num];  */
748            }
749          else
750            {
751              if(calib -> next) calib = calib -> next;
752              else break;
753            }
754          if(calib -> next) calib = calib -> next;
755          else break;
756        }
757      while(calib);
758
759      int result;
760      result = TRUE;
761      TIMES_Set_All_Node_Priors_S(&result, tree);
762      /* printf("\n\n"); */
763      /* For(j, 2 * tree -> n_otu - 1) printf("\n. [1] Node [%d] min [%f] max [%f] node time [%f]\n", j, tree -> rates -> t_prior_min[j], tree -> rates -> t_prior_max[j], tree -> rates -> nd_t[j]); */
764      /* printf("\n. p[%i] = %f \n", i + 1, times_partial_proba[i]); */
765      /* printf("\n\n"); */
766
767      //tree -> rates -> birth_rate = 4.0;
768      times_lk = TIMES_Lk_Yule_Order(tree);
769      /* if(result != FALSE) times_lk = TIMES_Lk_Yule_Order(tree); */
770      /* else times_lk = 1.0; */
771
772      constant = 1.0; 
773      if(times_lk > -INFINITY && result != FALSE) constant = Slicing_Calibrations(tree);     
774      /* else */
775      /*   { */
776      /*     times_lk = 0.0; */
777      /*     times_partial_proba[i] = 0.0; */
778      /*   } */
779
780      /* printf("\n. K = [%f] \n", constant); */
781      /* K = Norm_Constant_Prior_Times(tree); */
782      /* Yule_val[i] = K[i] * TIMES_Lk_Yule_Order(tree); */
783 
784      /* For(j, 2 * tree -> n_otu - 1) printf("\n. [2] Node [%d] time [%f]\n", j, tree -> rates -> nd_t[j]); */
785 
786      Yule_val[i] = LOG(constant) + times_lk;
787
788      /* printf("\n. Yule = %f \n", Yule_val[i]);   */
789 
790      while(calib -> prev) calib = calib -> prev;
791    }
792 
793  /* min_value = 0.0; */
794  /* For(i, tot_num_comb) if(Yule_val[i] < min_value && Yule_val[i] > -INFINITY) min_value = Yule_val[i]; */
795  /* c = -600. - min_value; */ 
796
797 
798  c = .0;
799  times_tot_proba = 0.0;
800  For(i, tot_num_comb)
801    {
802      times_tot_proba += times_partial_proba[i] * EXP(Yule_val[i] + c);
803      /* printf("\n. Proba = [%f] \n", times_partial_proba[i]); */
804    } 
805
806  For(i, 2 * tree -> n_otu - 1) t_has_prior[i] = NO;
807
808  ln_t = -c + LOG(times_tot_proba);
809  /* printf("\n. Prior for node times = [%f] \n", ln_t); */
810  /* Set_Current_Calibration(1, tree); */
811  /* printf("\n\n"); */
812  free(Yule_val); 
813  /* free(times_partial_proba); */
814  /* Exit("\n"); */
815  return(ln_t);
816}
817
818
819
820////////////////////////////////////////////////////////////////////////////
821////////////////////////////////////////////////////////////////////////////
822//Function calculates the normalizing constant K of the joint distribution Yule_Order.
823//Use the fact that density Yule_Order can be used streight forward only in case of the complete
824//overlap of the calibration intervals for all of the nodes or in case of no overlap. 
825phydbl Slicing_Calibrations(t_tree *tree)
826{
827  int i, j, k, f, n_otu, *indic, *n_slice, *slice_numbers;
828  phydbl K, buf, chop_bound, *t_prior_min, *t_prior_max, *t_slice, *t_slice_min, *t_slice_max;
829
830
831  t_prior_min = tree -> rates -> t_prior_min;
832  t_prior_max = tree -> rates -> t_prior_max;
833  n_otu = tree -> n_otu;
834
835  t_slice        = (phydbl *)mCalloc(2 * (n_otu - 1), sizeof(phydbl)); //the vector of the union of lower and upper bounds, lined up in incresing order.
836  t_slice_min    = (phydbl *)mCalloc(2 * n_otu - 3, sizeof(phydbl));   //vector of the lower bounds of the sliced intervals.
837  t_slice_max    = (phydbl *)mCalloc(2 * n_otu - 3, sizeof(phydbl));   //vector of the upper bounds of the sliced intervals.
838  indic          = (int *)mCalloc((n_otu - 1) * (2 * n_otu - 3), sizeof(int));  //vector of the indicators, columns - node numbers (i + n_otu), rows - the number of the sliced interval.
839  slice_numbers  = (int *)mCalloc((n_otu - 1) * (2 * n_otu - 3), sizeof(int )); //vecor of the slice intervals numbers, columns node numbers (i + n_otu), rows - the number of the sliced interval.
840  n_slice        = (int *)mCalloc(n_otu - 1, sizeof(int));                      //vector of the numbers of sliced intervals that apply to one node with number (i + n_otu).
841 
842  i = 0;
843  K = 0;
844  j = n_otu;
845  ////////////////////////////////////////////////////////////////////////////
846  //Put prior bounds in one vector t_slice. Excluding tips.
847  For(i, n_otu - 1) 
848    {
849      t_slice[i] = t_prior_min[j];
850      j++;
851    }
852
853  j = n_otu; 
854  for(i = n_otu - 1; i < 2 * n_otu - 3; i++) 
855    {
856      t_slice[i] = t_prior_max[j];
857      j++;
858    }
859  if(tree -> rates -> nd_t[tree -> n_root -> num] > t_prior_min[tree -> n_root -> num]) chop_bound =  MIN(tree -> rates -> nd_t[tree -> n_root -> num], t_prior_max[tree -> n_root -> num]);
860  else chop_bound = t_prior_min[tree -> n_root -> num];
861  t_slice[2 * n_otu - 3] = chop_bound;
862  //printf("\n. Chop bound [%f] \n", chop_bound);
863  //t_slice[2 * n_otu - 3] = -1.1;
864  //For(j, 2 * n_otu - 2) printf("\n. Slice bound [%f] \n", t_slice[j]);
865  ////////////////////////////////////////////////////////////////////////////
866  //Get slices in increasing order. Excluding tips.
867  do
868    {
869      f = NO;
870      For(j, 2 * n_otu - 3)
871        {
872          if(t_slice[j] > t_slice[j + 1])
873            {
874              buf = t_slice[j];
875              t_slice[j] = t_slice[j + 1];
876              t_slice[j + 1] = buf;
877              f = YES;
878            }
879        }
880    }
881  while(f);
882  //For(j, 2 * n_otu - 2) printf("\n. [1] Slice bound [%f] \n", t_slice[j]);
883  for(j = 1; j < 2 * n_otu - 2; j++) t_slice[j] = MAX(chop_bound, t_slice[j]);
884  //for(j = 1; j < 2 * n_otu - 2; j++) t_slice[j] = MAX(-1.1, t_slice[j]);
885  //For(j, 2 * n_otu - 2) printf("\n. [2] Slice bound [%f] \n", t_slice[j]);
886  ////////////////////////////////////////////////////////////////////////////
887  //Get the intervals with respect to slices. Total number of t_slice_min(max) - 2 * n_otu - 3. Excluding tips.
888  i = 0;
889  For(j, 2 * n_otu - 3)
890    {
891      t_slice_min[j] = t_slice[i];
892      t_slice_max[j] = t_slice[i + 1];
893      i++;
894    }
895
896  //For(j, 2 * n_otu - 3) printf("\n. The interval number [%d] min [%f] max[%f] \n", j, t_slice_min[j], t_slice_max[j]);
897
898  ////////////////////////////////////////////////////////////////////////////
899  //Getting indicators for the node number [i + n_otu] to have slice. i = i + n_otu is the node number on the tree and j is the slice number, total
900  //number of intervals is 2 * n_otu - 3. Excluding tips.
901  For(i, n_otu - 1)
902    { 
903      For(j, 2 * n_otu - 3)
904        {
905
906          if(Are_Equal(t_prior_min[i + n_otu], t_slice_min[j], 1.E-10) && t_prior_max[i + n_otu] > t_slice_max[j] && t_prior_min[i + n_otu] < t_slice_max[j] && !Are_Equal(t_slice_max[j], t_slice_min[j], 1.E-10)) indic[i * (2 * n_otu - 3) + j] = 1;
907          else if(Are_Equal(t_prior_max[i + n_otu], t_slice_max[j], 1.E-10) && t_prior_min[i + n_otu] < t_slice_min[j] && t_prior_max[i + n_otu] > t_slice_min[j] && !Are_Equal(t_slice_max[j], t_slice_min[j], 1.E-10)) indic[i * (2 * n_otu - 3) + j] = 1;
908          else if(t_prior_min[i + n_otu] < t_slice_min[j] && t_prior_max[i + n_otu] > t_slice_max[j] && !Are_Equal(t_slice_max[j], t_slice_min[j], 1.E-10)) indic[i * (2 * n_otu - 3) + j] = 1;
909          else if(Are_Equal(t_prior_min[i + n_otu], t_slice_min[j], 1.E-10) && Are_Equal(t_prior_max[i + n_otu], t_slice_max[j], 1.E-10)) indic[i * (2 * n_otu - 3) + j] = 1;
910        }
911    }
912
913
914  For(i, n_otu - 2)
915    {
916      indic[i * (2 * n_otu - 3)] = 0;
917    }
918 
919  for(j = 1; j <  2 * n_otu - 3; j++)
920    {
921      indic[(n_otu - 2) * (2 * n_otu - 3) + j] = 0;
922    }
923
924  /* For(i, n_otu - 2) */
925  /*   { */
926  /*     indic[i * (2 * n_otu - 3)] = 0; */
927  /*   } */
928
929  /* For(i, n_otu - 1) */
930  /*   { */
931  /*     indic[i * (2 * n_otu - 3) + 1] = 0; */
932  /*   } */
933 
934 
935  /* printf("\n"); */
936  /* For(i, n_otu - 1) */
937  /*   { */
938  /*     printf(" ['%d]' ", i + n_otu); */
939  /*     For(j, 2 * n_otu - 3) */
940  /*       { */
941  /*         printf(". '%d' ", indic[i * (2 * n_otu - 3) + j]); */
942  /*       } */
943  /*     printf("\n"); */
944  /*   } */
945
946
947  ////////////////////////////////////////////////////////////////////////////
948  //Get the number of slices that can be applied for each node and the vectors of slice numbers for each node.
949  For(i, n_otu - 1)
950    {
951      k = 0;
952      For(j, 2 * n_otu - 3)
953        {
954          if(indic[i * (2 * n_otu - 3) + j] == 1) 
955            {
956              slice_numbers[i * (2 * n_otu - 3) + k] = j; //printf("\n. Node [%d] slice'%d' ", i + n_otu, slice_numbers[i * (2 * n_otu - 3) + j]);
957              n_slice[i]++; 
958              k++;
959            }
960        }
961      //printf(" Number of slices'%d' \n", n_slice[i]);
962    }
963  /*
964  printf("\n");
965  For(i, n_otu - 1)
966    { 
967      printf(" ['%d]' ", i + n_otu);
968       For(j, n_slice[i]) 
969        {         
970          printf(". '%d' ", slice_numbers[i * (2 * n_otu - 3) + j]);         
971        }
972      printf("\n");
973    }
974  */
975
976  ////////////////////////////////////////////////////////////////////////////
977  //Running through all of the combinations of slices
978  int l, tot_num_comb, *cur_slices, *cur_slices_shr, shr_num_slices;
979  phydbl P, *t_cur_slice_min, *t_cur_slice_max;
980  phydbl k_part;
981
982  tot_num_comb = 1;
983  P = 0.0; 
984 
985
986  t_cur_slice_min    = (phydbl *)mCalloc(n_otu - 1, sizeof(phydbl));
987  t_cur_slice_max    = (phydbl *)mCalloc(n_otu - 1, sizeof(phydbl));
988  cur_slices         = (int *)mCalloc(n_otu - 1, sizeof(int)); //the vector of the current slices with repetition.
989  cur_slices_shr     = (int *)mCalloc(n_otu - 1, sizeof(int)); //the vector of the current slices without repetition.
990
991  For(i, n_otu - 1) tot_num_comb = tot_num_comb * n_slice[i]; //printf("\n. Total number of combinations of slices [%d] \n", tot_num_comb);
992 
993  For(k, tot_num_comb)
994    {
995      shr_num_slices = 0;
996      //printf("\n");
997      For(i, n_otu - 1) //node number i + n_otu
998        {
999          //printf(" ['%d]' ", i + n_otu);
1000          l = (k % Number_Of_Comb_Slices(i, n_otu - 1, n_slice)) / Number_Of_Comb_Slices(i+1, n_otu - 1, n_slice); //printf(" Slice number'%d' ", slice_numbers[i * (2 * n_otu - 3) + l]);
1001          t_cur_slice_min[i] = t_slice_min[slice_numbers[i * (2 * n_otu - 3) + l]]; //printf(" '%f' ", t_cur_slice_min[i]);
1002          t_cur_slice_max[i] = t_slice_max[slice_numbers[i * (2 * n_otu - 3) + l]]; //printf(" '%f' ", t_cur_slice_max[i]);
1003          cur_slices[i] = slice_numbers[i * (2 * n_otu - 3) + l];
1004          //printf("\n");
1005        }
1006      //printf("\n");
1007      //For(i, n_otu - 1) printf(" Slice number'%d' ", cur_slices[i]);
1008      //printf("\n");
1009
1010      ///////////////////////////////////////////////////////////////////////////
1011      //Taking away duplicated slices
1012      For(i, n_otu - 1)
1013        {
1014          for(j = i + 1; j < n_otu - 1; j++)
1015            {
1016              if(cur_slices[i] == cur_slices[j]) cur_slices[j] = -1;
1017            } 
1018        }
1019      //For(i, n_otu - 1) printf(" Slice number'%d' \n", cur_slices[i]);
1020
1021      ///////////////////////////////////////////////////////////////////////////
1022      //Getting a vector of all of the slices without duplicates.
1023      For(i, n_otu -1)
1024        {
1025          if(cur_slices[i] >= 0) 
1026            {
1027              cur_slices_shr[shr_num_slices] = cur_slices[i];
1028              shr_num_slices++; 
1029            }
1030        }
1031      //printf("\n");
1032      //For(i, shr_num_slices) printf("\n. Slice number'%d' \n", cur_slices_shr[i]);
1033      //printf("\n");
1034
1035      ////////////////////////////////////////////////////////////////////////////
1036      //Check for the time slices to be set properly
1037      int result;
1038
1039      result = TRUE;
1040
1041      Check_Time_Slices(tree -> n_root, tree -> n_root -> v[1], &result, t_cur_slice_min, t_cur_slice_max, tree);
1042      Check_Time_Slices(tree -> n_root, tree -> n_root -> v[2], &result, t_cur_slice_min, t_cur_slice_max, tree);
1043      //printf("\n. '%d' \n", result);
1044
1045      //For(i, n_otu - 1) printf("\n. Node [%d] min [%f] max [%f] \n", i + n_otu, t_cur_slice_min[i], t_cur_slice_max[i]);
1046
1047      ////////////////////////////////////////////////////////////////////////////
1048      //Calculating k_part
1049
1050      k_part = 1.0;
1051
1052      if(result != TRUE) k_part = 0.0; 
1053      else
1054        {
1055          int n_1, n_2;
1056
1057          ////////////////////////////////////////////////////////////////////////////
1058          //Getting the root node in a given slice
1059          int *root_nodes;
1060          int num_elem;
1061           
1062          num_elem = 0;
1063           
1064          root_nodes = (int *)mCalloc(n_otu - 1, sizeof(int)); 
1065         
1066          //printf("\n. Number of slices shrinked [%d] \n", shr_num_slices);
1067         
1068          For(i, shr_num_slices)
1069            {
1070              //printf("\n. Hello [%d] \n", i);
1071              //printf("\n. The number of the shrinked interval [%d] min [%f] max [%f] \n", cur_slices_shr[i], t_slice_min[cur_slices_shr[i]], t_slice_max[cur_slices_shr[i]]);
1072              Search_Root_Node_In_Slice(tree -> n_root, tree -> n_root -> v[1], root_nodes, &num_elem, t_slice_min[cur_slices_shr[i]], t_slice_max[cur_slices_shr[i]], t_cur_slice_min, t_cur_slice_max, tree);
1073              Search_Root_Node_In_Slice(tree -> n_root, tree -> n_root -> v[2], root_nodes, &num_elem, t_slice_min[cur_slices_shr[i]], t_slice_max[cur_slices_shr[i]], t_cur_slice_min, t_cur_slice_max, tree); 
1074            }
1075          //printf("\n. Number of elements in a vector of the nodes [%d] \n", num_elem);
1076          //For(j, num_elem) printf("\n. Root node number [%d] \n", root_nodes[j] + n_otu);
1077          For(j, num_elem) 
1078            {
1079              //printf("\n. Root node number [%d] \n", tree -> a_nodes[root_nodes[j] + n_otu] -> num);
1080             
1081              n_1 = 0;
1082              n_2 = 0;
1083              Number_Of_Nodes_In_Slice(tree -> a_nodes[root_nodes[j] + n_otu], tree -> a_nodes[root_nodes[j] + n_otu] -> v[1], &n_1, t_cur_slice_min, t_cur_slice_max, tree);
1084              Number_Of_Nodes_In_Slice(tree -> a_nodes[root_nodes[j] + n_otu], tree -> a_nodes[root_nodes[j] + n_otu] -> v[2], &n_2, t_cur_slice_min, t_cur_slice_max, tree);
1085              //printf("\n. n_1 [%d] n_2 [%d]\n", n_1, n_2); 
1086              k_part = k_part * Factorial(n_1 + n_2) / ((phydbl)Factorial(n_1 + n_2 + 1) * Factorial(n_1) * Factorial(n_2));
1087            }
1088          //printf("\n. k_part [%f] \n", k_part);       
1089          ////////////////////////////////////////////////////////////////////////////
1090          //Calculating PRODUCT over all of the time slices k_part * (exp(-lmbd*l) - exp(-lmbd*u))/(exp(-lmbd*l) - exp(-lmbd*u))
1091          phydbl num, denom, lmbd;
1092         
1093          lmbd = tree -> rates -> birth_rate;
1094          num = 1;
1095          denom = 1;
1096          //lmbd = 4.0;
1097          /* For(j, n_otu - 1) num = num * (EXP(-lmbd * t_cur_slice_min[j]) - EXP(-lmbd * t_cur_slice_max[j]));  */
1098          /* for(j = n_otu; j < 2 * n_otu - 1; j++) denom = denom * (EXP(-lmbd * t_prior_min[j]) - EXP(-lmbd * t_prior_max[j]));  */
1099          For(j, n_otu - 2) num = num * (EXP(lmbd * t_cur_slice_max[j]) - EXP(lmbd * t_cur_slice_min[j])); 
1100          for(j = n_otu; j < 2 * n_otu - 2; j++) denom = denom * (EXP(lmbd * t_prior_max[j]) - EXP(lmbd * t_prior_min[j])); 
1101          k_part = (k_part * num) / denom; 
1102          //printf("\n. [2] k_part of the tree for one combination of slices [%f] \n", k_part);
1103        } 
1104      P = P + k_part;     
1105    }
1106  //printf("\n. [P] of the tree for one combination of slices [%f] \n", P);
1107  K = 1 / P;
1108  //printf("\n. [K] of the tree for one combination of slices [%f] \n", K);
1109  ////////////////////////////////////////////////////////////////////////////
1110  free(t_cur_slice_min);
1111  free(t_cur_slice_max);
1112  free(cur_slices);
1113  free(cur_slices_shr); 
1114  free(t_slice);       
1115  free(t_slice_min);   
1116  free(t_slice_max);   
1117  free(indic);         
1118  free(slice_numbers); 
1119  free(n_slice);
1120  //free(tree -> K);
1121 
1122  return(K); 
1123}
1124
1125////////////////////////////////////////////////////////////////////////////
1126////////////////////////////////////////////////////////////////////////////
1127int Number_Of_Comb_Slices(int m, int num_elem, int *n_slice)
1128{
1129  int i, num_comb;
1130
1131  i = 0;
1132  num_comb = 1;       
1133
1134  for(i = m; i < num_elem; i++) num_comb = num_comb * n_slice[i]; 
1135 
1136  return(num_comb);
1137}
1138
1139
1140
1141//////////////////////////////////////////////////////////////
1142//////////////////////////////////////////////////////////////
1143//Check the combination of the time slices to be set correctly.
1144void Check_Time_Slices(t_node *a, t_node *d, int *result, phydbl *t_cur_slice_min, phydbl *t_cur_slice_max, t_tree *tree)
1145{
1146  int n_otu;
1147
1148  n_otu = tree -> n_otu;
1149 
1150
1151  d -> anc = a;
1152  if(d -> tax) return;
1153  else
1154    { 
1155      if(t_cur_slice_max[d -> num - n_otu] < t_cur_slice_max[a -> num - n_otu])
1156        {
1157          *result = FALSE; 
1158        }
1159
1160      int i;
1161      For(i,3) 
1162        if((d -> v[i] != d -> anc) && (d -> b[i] != tree -> e_root))
1163             Check_Time_Slices(d, d -> v[i], result, t_cur_slice_min, t_cur_slice_max, tree);
1164    }
1165}
1166
1167//////////////////////////////////////////////////////////////
1168/////////////////////////////////////////////////////////////
1169//Getting the number of nodes on both sides from the node d_start, that are in the slice of that node.
1170void Number_Of_Nodes_In_Slice(t_node *d_start, t_node *d, int *n, phydbl *t_cur_slice_min, phydbl *t_cur_slice_max, t_tree *tree)
1171{
1172  int n_otu;
1173 
1174  n_otu = tree -> n_otu;
1175
1176
1177  if(d -> tax) return;
1178  else
1179    { 
1180      if(Are_Equal(t_cur_slice_max[d_start -> num - n_otu], t_cur_slice_max[d -> num - n_otu], 1.E-10) && Are_Equal(t_cur_slice_min[d_start -> num - n_otu], t_cur_slice_min[d -> num - n_otu], 1.E-10))
1181        {
1182          (*n)++; 
1183          int i;
1184          For(i,3) 
1185            if((d -> v[i] != d -> anc) && (d -> b[i] != tree -> e_root))
1186              Number_Of_Nodes_In_Slice(d_start, d -> v[i], n, t_cur_slice_min, t_cur_slice_max, tree);
1187        }     
1188    }
1189}
1190
1191
1192//////////////////////////////////////////////////////////////
1193/////////////////////////////////////////////////////////////
1194//Returnig the root node in the given slice.
1195void Search_Root_Node_In_Slice(t_node *d_start, t_node *d, int *root_nodes, int *num_elem, phydbl t_slice_min, phydbl t_slice_max, phydbl *t_cur_slice_min, phydbl *t_cur_slice_max, t_tree *tree)
1196{
1197  int j, n_otu, f;
1198  j = 0;
1199  f = FALSE;
1200  n_otu = tree -> n_otu; //printf("\n. Node number 1  [%d] \n", d_start -> num);printf("\n. Node number 1  [%d] \n", d -> num);
1201
1202  if(Are_Equal(t_cur_slice_max[d_start -> num - n_otu], t_slice_max, 1.E-10) && Are_Equal(t_cur_slice_min[d_start -> num - n_otu], t_slice_min, 1.E-10))
1203    {
1204      For(j, *num_elem) if(d_start -> num - n_otu == (root_nodes[j])) f = TRUE;
1205      if(f != TRUE)
1206        {
1207          (root_nodes[(*num_elem)]) = d_start -> num - n_otu; //printf("\n. Node number 2_2  [%d] \n", root_nodes[(*num_elem)] + n_otu);
1208          (*num_elem)++;  //printf("\n. Number of elements 2_2  [%d] \n", *num_elem);
1209          return;
1210        }
1211     
1212    }
1213  else
1214    {
1215      d -> anc = d_start;
1216      if(d -> tax) return;
1217      else
1218        { 
1219          if(Are_Equal(t_cur_slice_max[d -> num - n_otu], t_slice_max, 1.E-10) && Are_Equal(t_cur_slice_min[d -> num - n_otu], t_slice_min, 1.E-10))
1220            {
1221              (root_nodes[*num_elem]) = d -> num - n_otu; //printf("\n. Node number 3_2  [%d] \n",  root_nodes[*num_elem] + n_otu);
1222              (*num_elem)++; //printf("\n. Number of elements 3_2  [%d] \n", *num_elem); 
1223              return; 
1224            }
1225         
1226          int i;
1227          For(i,3) 
1228            if((d -> v[i] != d -> anc) && (d -> b[i] != tree -> e_root))
1229              Search_Root_Node_In_Slice(d, d -> v[i], root_nodes, num_elem, t_slice_min, t_slice_max, t_cur_slice_min, t_cur_slice_max, tree);
1230        }
1231    }
1232}
1233
1234
1235//////////////////////////////////////////////////////////////
1236//////////////////////////////////////////////////////////////
1237int Factorial(int base)
1238{
1239  if(base == 0) return(1); 
1240  if(base == 1) return(1);
1241  return(base * Factorial(base-1)); 
1242}
1243
1244
1245//////////////////////////////////////////////////////////////
1246//////////////////////////////////////////////////////////////
1247//Calculate the vector of the norm.constants for prior on node times.
1248//The length of the vector is the total number of combinations of calibrations.
1249phydbl *Norm_Constant_Prior_Times(t_tree *tree)
1250{
1251
1252  phydbl *t_prior_min, *t_prior_max, *K; 
1253  short int *t_has_prior;
1254  int i, j, k, tot_num_comb;
1255  t_cal *calib;
1256 
1257
1258
1259
1260  calib = tree -> rates -> calib;
1261  t_prior_min = tree -> rates -> t_prior_min;
1262  t_prior_max = tree -> rates -> t_prior_max;
1263  t_has_prior = tree -> rates -> t_has_prior;
1264
1265  tot_num_comb = Number_Of_Comb(calib);
1266 
1267  //PhyML_Printf("\n. The total number of calibration combinations: [%d]\n", tot_num_comb);
1268
1269  K = (phydbl *)mCalloc(tot_num_comb, sizeof(phydbl));   
1270
1271  For(i, tot_num_comb)
1272    {
1273      for(j = tree -> n_otu; j < 2 * tree -> n_otu - 1; j++) 
1274        {
1275          t_prior_min[j] = -BIG;
1276          t_prior_max[j] = BIG;
1277          t_has_prior[j] = NO; 
1278        }
1279      do
1280        {
1281          k = (i % Number_Of_Comb(calib)) / Number_Of_Comb(calib -> next);
1282          if(calib -> all_applies_to[k] -> num)
1283            {
1284              t_prior_min[calib -> all_applies_to[k] -> num] = MAX(t_prior_min[calib -> all_applies_to[k] -> num], calib -> lower);
1285              t_prior_max[calib -> all_applies_to[k] -> num] = MIN(t_prior_max[calib -> all_applies_to[k] -> num], calib -> upper);
1286              t_has_prior[calib -> all_applies_to[k] -> num] = YES;
1287            }
1288          else
1289            {
1290              if(calib -> next) calib = calib -> next;
1291              else break;
1292            }
1293          if(calib -> next) calib = calib -> next;
1294          else break;
1295        }
1296      while(calib);
1297      TIMES_Set_All_Node_Priors(tree);
1298      //for(j = tree -> n_otu; j < 2 * tree -> n_otu - 1; j++) printf("\n. [1] Node [%d] min [%f] max[%f]\n", j, tree -> rates -> t_prior_min[j], tree -> rates -> t_prior_max[j]);
1299      //tree -> rates -> birth_rate = 4.0;
1300      K[i] = Slicing_Calibrations(tree);     
1301      //PhyML_Printf("\n. Number [%d] normolizing constant [%f] \n", i+1, K[i]);
1302      while(calib -> prev) calib = calib -> prev;
1303    }
1304  return(K);
1305}
1306
1307
1308//////////////////////////////////////////////////////////////
1309//////////////////////////////////////////////////////////////
1310//Sets a vector of the partial probabilities for each combination of calibrations
1311void TIMES_Calib_Partial_Proba(t_tree *tree)
1312{
1313
1314  phydbl *times_partial_proba, proba, *t_prior_min, *t_prior_max; 
1315  int i, j, k, tot_num_comb;
1316  t_cal *calib;
1317  short int *t_has_prior;
1318
1319  proba = 0.0;
1320 
1321  times_partial_proba = tree -> rates -> times_partial_proba; 
1322  calib = tree -> rates -> calib;
1323  t_prior_min = tree -> rates -> t_prior_min;
1324  t_prior_max = tree -> rates -> t_prior_max;
1325  t_has_prior = tree -> rates -> t_has_prior;
1326
1327  tot_num_comb = Number_Of_Comb(calib);
1328
1329  For(i, tot_num_comb)
1330    {
1331      times_partial_proba[i] = 1.0;
1332      for(j = tree -> n_otu; j < 2 * tree -> n_otu - 1; j++) 
1333        {
1334          t_prior_min[j] = -BIG;
1335          t_prior_max[j] = BIG;
1336          t_has_prior[j] = NO; 
1337        }
1338      do
1339        {
1340          k = (i % Number_Of_Comb(calib)) / Number_Of_Comb(calib -> next);
1341          if(calib -> all_applies_to[k] -> num)
1342            {
1343              t_prior_min[calib -> all_applies_to[k] -> num] = MAX(t_prior_min[calib -> all_applies_to[k] -> num], calib -> lower);
1344              t_prior_max[calib -> all_applies_to[k] -> num] = MIN(t_prior_max[calib -> all_applies_to[k] -> num], calib -> upper);
1345              t_has_prior[calib -> all_applies_to[k] -> num] = YES;
1346              proba = calib -> proba[calib -> all_applies_to[k] -> num]; 
1347              times_partial_proba[i] *= proba;
1348              /* printf("\n. [1] Proba [%f] \n", proba); */
1349              /* if((t_prior_min[calib -> all_applies_to[k] -> num] > t_prior_max[calib -> all_applies_to[k] -> num])) times_partial_proba[i] = 0.0;  */
1350              /* else times_partial_proba[i] *= calib -> proba[calib -> all_applies_to[k] -> num];  */
1351            }
1352          else
1353            {
1354              proba = calib -> proba[2 * tree -> n_otu - 1]; 
1355              times_partial_proba[i] *= proba;
1356              /* printf("\n. [2] Proba [%f] \n", proba); */
1357              if(calib -> next) calib = calib -> next;
1358              else break;
1359            }
1360         
1361          if(calib -> next) calib = calib -> next;
1362          else break;
1363        }
1364      while(calib);
1365
1366      int result;
1367
1368      result = TRUE;
1369
1370      /* printf("\n. [3] Partial Proba [%f] \n", times_partial_proba[i]); */
1371
1372      TIMES_Set_All_Node_Priors_S(&result, tree);
1373
1374      if(result != TRUE) times_partial_proba[i] = 0; /* printf("\n. [4] Partial Proba [%f] \n", times_partial_proba[i]); */
1375      /* times_partial_proba[i] = 1.0; */
1376      /* do */
1377      /*   { */
1378      /*     k = (i % Number_Of_Comb(calib)) / Number_Of_Comb(calib -> next); */
1379      /*     if(calib -> all_applies_to[k] -> num) proba = calib -> proba[calib -> all_applies_to[k] -> num];  */
1380      /*     else proba = calib -> proba[2 * tree -> n_otu - 1]; */
1381      /*     times_partial_proba[i] *= proba; */
1382      /*     if(calib -> next) calib = calib -> next; */
1383      /*     else break; */
1384      /*   } */
1385      /* while(calib); */
1386      while(calib -> prev) calib = calib -> prev;
1387    }
1388
1389  phydbl sum_proba;
1390  sum_proba = 0.0;
1391  /* For(i, tot_num_comb) printf("\n. [1] Partial Proba [%f] \n", times_partial_proba[i]);  */
1392  For(i, tot_num_comb) sum_proba += times_partial_proba[i];
1393  if(!Are_Equal(sum_proba, 1.0, 1.E-10)) 
1394    {
1395      For(i, tot_num_comb) times_partial_proba[i] = times_partial_proba[i] / sum_proba;
1396    } 
1397  /* For(i, tot_num_comb) printf("\n. [2] Partial Proba [%f] \n", times_partial_proba[i]);  */
1398  /* Exit("\n"); */
1399}
1400
1401
1402
1403//////////////////////////////////////////////////////////////
1404//////////////////////////////////////////////////////////////
1405//Function checks if the randomized node times are within the
1406//upper and lower time limits, taken into account the times of
1407//the ancestor and descendent.
1408void Check_Node_Time(t_node *a, t_node *d, int *result, t_tree *tree)
1409{
1410  phydbl t_low, t_up;
1411  phydbl *t_prior_min, *t_prior_max, *nd_t;
1412
1413  t_prior_min = tree -> rates -> t_prior_min;
1414  t_prior_max = tree -> rates -> t_prior_max;
1415  nd_t = tree -> rates -> nd_t;
1416
1417  if(a == tree -> n_root && (nd_t[a -> num] > MIN(t_prior_max[a -> num], MIN(nd_t[a -> v[1] -> num], nd_t[a -> v[2] -> num])))) 
1418    {
1419      *result = FALSE;
1420      return;
1421    }
1422  if(d -> tax) return;
1423  else
1424    { 
1425      t_low = MAX(t_prior_min[d -> num], nd_t[d -> anc -> num]);
1426      t_up = MIN(t_prior_max[d -> num], MIN(nd_t[d -> v[1] -> num], nd_t[d -> v[2] -> num]));
1427      if(nd_t[d -> num] < t_low || nd_t[d -> num] > t_up)
1428        {
1429          *result = FALSE; 
1430        }
1431
1432      int i;
1433      For(i,3) 
1434        if((d -> v[i] != d -> anc) && (d -> b[i] != tree -> e_root))
1435             Check_Node_Time(d, d -> v[i], result, tree);
1436    }
1437}
1438
1439//////////////////////////////////////////////////////////////
1440//////////////////////////////////////////////////////////////
1441//Function calculates the TOTAL number of calibration combinations,
1442//given the number of nodes to which each calibartion applies to.
1443int Number_Of_Comb(t_cal *calib)
1444{
1445
1446  int num_comb;
1447
1448  if(!calib) return(1);
1449  num_comb = 1;
1450  do
1451    {
1452      num_comb *= calib -> n_all_applies_to;
1453      if(calib -> next) calib = calib -> next;
1454      else break;
1455    }
1456  while(calib);
1457  return(num_comb);
1458}
1459
1460
1461//////////////////////////////////////////////////////////////
1462//////////////////////////////////////////////////////////////
1463//Function calculates the TOTAL number of calibration combinations,
1464//given the number of nodes to which each calibartion applies to.
1465int Number_Of_Calib(t_cal *calib)
1466{
1467
1468  int num_calib;
1469
1470
1471  num_calib = 0;
1472  do
1473    {
1474      num_calib++;
1475      if(calib -> next) calib = calib -> next;
1476      else break;
1477    }
1478  while(calib);
1479  return(num_calib);
1480}
1481
1482
1483//////////////////////////////////////////////////////////////
1484//////////////////////////////////////////////////////////////
1485//Function sets current calibartion in the following way:
1486//Suppose we have a vector of calibrations C=(C1, C2, C3), each calibration 
1487//applies to a set of nodes. we can reach each node number through the indeces (corresponds
1488//to the number the information was read). C1={0,1,2}, C2={0,1}, C3={0};
1489//The total number of combinations is 3*2*1=6. The first combination with row number 0
1490//will be {0,0,0}, the second row will be {0,1,0} and so on. Calling the node numbers with
1491//the above indeces will return current calibration. Also sets the vector of the probabilities
1492//for current calibration combination.   
1493void Set_Current_Calibration(int row, t_tree *tree)
1494{
1495
1496  t_cal *calib;
1497  phydbl *t_prior_min, *t_prior_max; 
1498  short int *t_has_prior;
1499  int k, i, j, *curr_nd_for_cal;
1500
1501  calib = tree -> rates -> calib;
1502  t_prior_min = tree -> rates -> t_prior_min;
1503  t_prior_max = tree -> rates -> t_prior_max;
1504  t_has_prior = tree -> rates -> t_has_prior;
1505  curr_nd_for_cal = tree -> rates -> curr_nd_for_cal;
1506
1507  for(j = tree -> n_otu; j < 2 * tree -> n_otu - 1; j++) 
1508    {
1509      t_prior_min[j] = -BIG;
1510      t_prior_max[j] = BIG;
1511      t_has_prior[j] = NO; 
1512    }
1513 
1514  k = -1;
1515  i = 0;
1516  do
1517    {
1518      k = (row % Number_Of_Comb(calib)) / Number_Of_Comb(calib -> next); 
1519      if(calib -> all_applies_to[k] -> num) 
1520        {
1521          t_prior_min[calib -> all_applies_to[k] -> num] = MAX(t_prior_min[calib -> all_applies_to[k] -> num], calib -> lower);
1522          t_prior_max[calib -> all_applies_to[k] -> num] = MIN(t_prior_max[calib -> all_applies_to[k] -> num], calib -> upper);     
1523          t_has_prior[calib -> all_applies_to[k] -> num] = YES; 
1524          curr_nd_for_cal[i] = calib -> all_applies_to[k] -> num;
1525          i++;
1526        }
1527      else 
1528        {
1529          if(calib->next) calib = calib->next;
1530          else break;   
1531        }     
1532      if(calib->next) calib = calib->next;
1533      else break;   
1534    }
1535  while(calib);
1536  //while(calib -> prev) calib = calib -> prev;
1537}
1538
1539
1540
1541
1542
1543
1544//////////////////////////////////////////////////////////////
1545//////////////////////////////////////////////////////////////
1546//Randomly choose a combination of calibrations drawing an index of calibration combination,
1547//used function Set_Cur_Calibration.
1548void Random_Calibration(t_tree *tree)
1549{
1550  int rnd, num_comb;
1551  t_cal *calib;
1552
1553  calib = tree -> rates -> calib;
1554
1555  num_comb =  Number_Of_Comb(calib);
1556
1557  srand(time(NULL));
1558  rnd = rand()%(num_comb);
1559 
1560  Set_Current_Calibration(rnd, tree);
1561  TIMES_Set_All_Node_Priors(tree); 
1562
1563}
1564
1565
1566
1567//////////////////////////////////////////////////////////////
1568//////////////////////////////////////////////////////////////
1569//Variable curr_nd_for_cal is a vector of node numbers, the length of that vector is a number of calibrations.
1570//Function randomly updates that vector by randomly changing one node and setting times limits with respect
1571//to a new vector. 
1572int RND_Calibration_And_Node_Number(t_tree *tree)
1573{
1574  int i, j, tot_num_cal, cal_num, node_ind, node_num, *curr_nd_for_cal;
1575  phydbl *t_prior_min, *t_prior_max; //*times_partial_proba;
1576  short int *t_has_prior;
1577  t_cal *cal;
1578
1579  tot_num_cal = tree -> rates -> tot_num_cal;
1580  t_prior_min = tree -> rates -> t_prior_min;
1581  t_prior_max = tree -> rates -> t_prior_max;
1582  t_has_prior = tree -> rates -> t_has_prior;
1583  //times_partial_proba = tree -> rates -> times_partial_proba;
1584  curr_nd_for_cal = tree -> rates -> curr_nd_for_cal;
1585  cal = tree -> rates -> calib;
1586
1587  cal_num = rand()%(tot_num_cal - 1);
1588   
1589  i = 0;
1590  while (i != cal_num)
1591    {
1592      cal = cal -> next;
1593      i++;
1594    }
1595
1596  node_ind = rand()%(cal -> n_all_applies_to);
1597  node_num = cal -> all_applies_to[node_ind] -> num;
1598
1599  curr_nd_for_cal[cal_num] = node_num;
1600
1601  for(j = tree -> n_otu; j < 2 * tree -> n_otu - 1; j++) 
1602    {
1603      t_prior_min[j] = -BIG;
1604      t_prior_max[j] = BIG;
1605      t_has_prior[j] = NO; 
1606    }
1607
1608  while(cal -> prev) cal = cal -> prev;
1609 
1610  i = 0;
1611  do
1612    {
1613      t_prior_min[curr_nd_for_cal[i]] = cal -> lower;
1614      t_prior_max[curr_nd_for_cal[i]] = cal -> upper;
1615      t_has_prior[curr_nd_for_cal[i]] = YES; 
1616      i++;
1617      if(cal->next) cal = cal -> next;
1618      else break;   
1619    }
1620  while(cal);
1621
1622  while(cal -> prev) cal = cal -> prev;
1623
1624  TIMES_Set_All_Node_Priors(tree);
1625
1626  return(node_num);
1627}
1628
1629
1630
1631//////////////////////////////////////////////////////////////
1632//////////////////////////////////////////////////////////////
1633//Return the value uniformly distributed between two values.
1634phydbl Randomize_One_Node_Time(phydbl min, phydbl max)
1635{
1636  phydbl u;
1637
1638  u = Uni();
1639  u *= (max - min);
1640  u += min;
1641
1642  return(u);
1643}
1644
1645
1646
1647//////////////////////////////////////////////////////////////
1648//////////////////////////////////////////////////////////////
1649//Calculates the Hastings ratio for the analysis. Used in case of
1650//calibration conditional jump. NOT THE RIGHT ONE TO USE!
1651void Lk_Hastings_Ratio_Times(t_node *a, t_node *d, phydbl *tot_prob, t_tree *tree)
1652{
1653  phydbl t_low, t_up;
1654  phydbl *t_prior_min, *t_prior_max, *nd_t;
1655
1656  t_prior_min = tree -> rates -> t_prior_min;
1657  t_prior_max = tree -> rates -> t_prior_max;
1658  nd_t = tree -> rates -> nd_t;
1659
1660  if(d -> tax) return;
1661  else
1662    { 
1663      t_low = MAX(t_prior_min[d -> num], nd_t[d -> anc -> num]);
1664      t_up = MIN(t_prior_max[d -> num], MIN(nd_t[d -> v[1] -> num], nd_t[d -> v[2] -> num]));
1665
1666      (*tot_prob) += LOG(1) - LOG(t_up - t_low);
1667
1668      int i;
1669      For(i,3) 
1670        if((d -> v[i] != d -> anc) && (d -> b[i] != tree -> e_root))
1671          { 
1672              Lk_Hastings_Ratio_Times(d, d -> v[i], tot_prob, tree);
1673          }
1674    } 
1675}
1676
1677
1678
1679//////////////////////////////////////////////////////////////
1680//////////////////////////////////////////////////////////////
1681//Updates nodes which are below a randomized node in case if new proposed time
1682//for that node is below the current value.
1683void Update_Descendent_Cond_Jump(t_node *a, t_node *d, phydbl *L_Hast_ratio, t_tree *tree)
1684{
1685  int result = TRUE;
1686  phydbl t_low, t_up;
1687  phydbl *t_prior_min, *t_prior_max, *nd_t;
1688
1689  t_prior_min = tree -> rates -> t_prior_min;
1690  t_prior_max = tree -> rates -> t_prior_max;
1691  nd_t = tree -> rates -> nd_t;
1692
1693  Check_Node_Time(tree -> n_root, tree -> n_root -> v[1], &result, tree); 
1694  Check_Node_Time(tree -> n_root, tree -> n_root -> v[2], &result, tree);
1695
1696  if(d -> tax) return;
1697  else
1698  {
1699    if(result != TRUE)
1700      { 
1701        int i;
1702        t_low = MAX(nd_t[a -> num], t_prior_min[d -> num]);
1703        if(t_low < MIN(nd_t[d -> v[1] -> num], nd_t[d -> v[2] -> num])) t_up  = MIN(t_prior_max[d -> num], MIN(nd_t[d -> v[1] -> num], nd_t[d -> v[2] -> num])); 
1704        else t_up  = t_prior_max[d -> num]; 
1705        nd_t[d -> num] = Randomize_One_Node_Time(t_low, t_up);
1706        (*L_Hast_ratio) += LOG(1) - LOG(t_up - t_low);
1707        For(i,3) 
1708          if((d -> v[i] != d -> anc) && (d -> b[i] != tree -> e_root)) 
1709            Update_Descendent_Cond_Jump(d, d -> v[i], L_Hast_ratio, tree);
1710      }
1711    else return;
1712  }
1713}
1714
1715
1716
1717//////////////////////////////////////////////////////////////
1718//////////////////////////////////////////////////////////////
1719//Updates nodes which are above a randomized node in case if new proposed time
1720//for that node is above the current value.
1721void Update_Ancestor_Cond_Jump(t_node *d, phydbl *L_Hast_ratio, t_tree *tree)
1722{
1723  int result = TRUE;
1724  phydbl t_low, t_up;
1725  phydbl *t_prior_min, *t_prior_max, *nd_t;
1726
1727  t_prior_min = tree -> rates -> t_prior_min;
1728  t_prior_max = tree -> rates -> t_prior_max;
1729  nd_t = tree -> rates -> nd_t;
1730 
1731  Check_Node_Time(tree -> n_root, tree -> n_root -> v[1], &result, tree); 
1732  Check_Node_Time(tree -> n_root, tree -> n_root -> v[2], &result, tree);
1733
1734  if(result != TRUE) 
1735    {     
1736      if(d == tree -> n_root)
1737        {
1738                     
1739          t_low = t_prior_min[d -> num];
1740          t_up  = MIN(t_prior_max[d -> num], MIN(nd_t[d -> v[1] -> num], nd_t[d -> v[2] -> num])); 
1741          nd_t[d -> num] = Randomize_One_Node_Time(t_low, t_up);
1742          (*L_Hast_ratio) += LOG(1) - LOG(t_up - t_low);
1743          return;
1744        }
1745      else
1746        { 
1747          t_up  = MIN(t_prior_max[d -> num], MIN(nd_t[d -> v[1] -> num], nd_t[d -> v[2] -> num]));           
1748          if(nd_t[d -> anc -> num] > t_up) t_low =  t_prior_min[d -> num];
1749          else  t_low  = MAX(t_prior_min[d -> num], nd_t[d -> anc -> num]); 
1750          nd_t[d -> num] = Randomize_One_Node_Time(t_low, t_up);
1751          (*L_Hast_ratio) += LOG(1) - LOG(t_up - t_low);
1752          Update_Ancestor_Cond_Jump(d -> anc, L_Hast_ratio, tree); 
1753        }
1754     
1755    }
1756  else return;
1757}
1758
1759//////////////////////////////////////////////////////////////
1760//////////////////////////////////////////////////////////////
1761//when made a calibration conditional jump, updates node times
1762//with respect to the new calibration which was made with respect
1763//to the randomly chosen node, the root is fixed. Updates only those nodes
1764//that are not within new intervals. Traverse up and down.
1765void Update_Times_RND_Node_Ancestor_Descendant(int rnd_node, phydbl *L_Hast_ratio, t_tree *tree)
1766{
1767  int i;
1768  phydbl *t_prior_min, *t_prior_max, *nd_t;
1769  phydbl new_time_rnd_node = 0.0;
1770
1771  t_prior_min = tree -> rates -> t_prior_min;
1772  t_prior_max = tree -> rates -> t_prior_max;
1773  nd_t = tree -> rates -> nd_t;
1774
1775  new_time_rnd_node = Randomize_One_Node_Time(t_prior_min[rnd_node], t_prior_max[rnd_node]);
1776
1777  nd_t[rnd_node] = new_time_rnd_node; 
1778
1779  Update_Ancestor_Cond_Jump(tree -> a_nodes[rnd_node] -> anc, L_Hast_ratio, tree);
1780  For(i,3) 
1781    if((tree -> a_nodes[rnd_node] -> v[i] != tree -> a_nodes[rnd_node] -> anc) && (tree -> a_nodes[rnd_node] -> b[i] != tree -> e_root)) 
1782      Update_Descendent_Cond_Jump(tree -> a_nodes[rnd_node], tree -> a_nodes[rnd_node] -> v[i], L_Hast_ratio, tree);
1783 
1784}
1785
1786
1787
1788//////////////////////////////////////////////////////////////
1789//////////////////////////////////////////////////////////////
1790//when made a calibration conditional jump, updates node times
1791//with respect to the new calibration which was made with respect
1792//to the randomly chosen node, starting from the root down to the tips.
1793//Updates only those nodes that are not within new intervals.
1794void Update_Times_Down_Tree(t_node *a, t_node *d, phydbl *L_Hastings_ratio, t_tree *tree)
1795{
1796  int i; 
1797  phydbl *t_prior_min, *t_prior_max, *nd_t, t_low, t_up;
1798
1799  t_prior_min = tree -> rates -> t_prior_min;
1800  t_prior_max = tree -> rates -> t_prior_max;
1801  nd_t = tree -> rates -> nd_t;
1802
1803
1804  t_low = MAX(t_prior_min[d -> num], nd_t[a -> num]);
1805  t_up  = t_prior_max[d -> num]; 
1806
1807  //printf("\n. [1] Node number: [%d] \n", d -> num);
1808  if(d -> tax) return;
1809  else
1810    {   
1811      if(nd_t[d -> num] > t_up || nd_t[d -> num] < t_low)
1812        { 
1813          //printf("\n. [2] Node number: [%d] \n", d -> num);
1814          //(*L_Hastings_ratio) += (LOG(1) - LOG(t_up - t_low));
1815          (*L_Hastings_ratio) += (- LOG(t_up - t_low));
1816          nd_t[d -> num] = Randomize_One_Node_Time(t_low, t_up);
1817          /* t_prior_min[d -> num] = t_low;  */
1818          /* t_prior_max[d -> num] = t_up;   */
1819        }   
1820     
1821      For(i,3) 
1822        if((d -> v[i] != d -> anc) && (d -> b[i] != tree -> e_root)) 
1823          Update_Times_Down_Tree(d, d -> v[i], L_Hastings_ratio, tree);
1824    }
1825}
1826
1827//////////////////////////////////////////////////////////////
1828//////////////////////////////////////////////////////////////
1829
1830xml_node *XML_Search_Node_Attribute_Value_Clade(char *attr_name, char *value, int skip, xml_node *node)
1831{
1832  xml_node *match;
1833
1834  //printf("\n. Node name [%s] Attr name [%s] Attr value [%s] \n", node -> name, attr_name, value);
1835  if(!node)
1836    {
1837      PhyML_Printf("\n== node: %p attr: %p",node,node?node->attr:NULL);
1838      PhyML_Printf("\n== Err in file %s at line %d\n",__FILE__,__LINE__);
1839      Exit("\n");         
1840    }
1841 
1842  match = NULL;
1843  if(skip == NO && node -> attr)
1844    {
1845      xml_attr *attr;
1846      attr = node -> attr;
1847      do
1848        {
1849          if(!strcmp(attr -> name, attr_name) && !strcmp(attr -> value, value))
1850            {
1851              match = node;
1852              break;
1853            }
1854          attr = attr->next;
1855          if(!attr) break;
1856        }
1857      while(1);
1858    }
1859  if(match) return(match);
1860  if(node -> next && !match)
1861    {
1862      match = XML_Search_Node_Attribute_Value_Clade(attr_name, value, NO, node -> next);
1863      return match;
1864    }
1865  return NULL;
1866}
1867
1868//////////////////////////////////////////////////////////////
1869//////////////////////////////////////////////////////////////
1870char **XML_Reading_Clade(xml_node *n_clade, t_tree *tree)
1871{
1872  int i, n_otu;
1873  char **clade;
1874
1875  i = 0;
1876  n_otu = tree -> n_otu;
1877
1878  clade = (char **)mCalloc(n_otu, sizeof(char *));
1879  if(n_clade)
1880    {
1881      do
1882        {
1883          clade[i] =  n_clade -> attr -> value; 
1884          i++;
1885          if(n_clade -> next) n_clade = n_clade -> next;
1886          else break;
1887        }
1888      while(n_clade);
1889    }
1890  else
1891    {
1892      PhyML_Printf("==Clade is empty. \n");
1893      PhyML_Printf("\n. Err in file %s at line %d\n",__FILE__,__LINE__);
1894      Exit("\n");
1895    }
1896
1897  return(clade);                         
1898}
1899
1900//////////////////////////////////////////////////////////////
1901//////////////////////////////////////////////////////////////
1902int XML_Number_Of_Taxa_In_Clade(xml_node *n_clade)
1903{
1904  int clade_size = 0;
1905  if(n_clade)
1906    {
1907      do
1908        {
1909          clade_size++; 
1910          if(n_clade -> next) n_clade = n_clade -> next;
1911          else break;
1912        }
1913      while(n_clade);
1914    }
1915  else
1916    {
1917      PhyML_Printf("==Clade is empty. \n");
1918      PhyML_Printf("\n. Err in file %s at line %d\n",__FILE__,__LINE__);
1919      Exit("\n");
1920    }
1921  return(clade_size);
1922}
1923
1924
1925
1926//////////////////////////////////////////////////////////////
1927//////////////////////////////////////////////////////////////
1928void TIMES_Set_All_Node_Priors_S(int *result, t_tree *tree)
1929{
1930  int i;
1931  phydbl min_prior;
1932
1933
1934  /* Set all t_prior_max values */
1935  TIMES_Set_All_Node_Priors_Bottom_Up_S(tree->n_root,tree->n_root->v[2], result, tree);
1936  TIMES_Set_All_Node_Priors_Bottom_Up_S(tree->n_root,tree->n_root->v[1], result, tree);
1937
1938  tree->rates->t_prior_max[tree->n_root->num] = 
1939    MIN(tree->rates->t_prior_max[tree->n_root->num],
1940        MIN(tree->rates->t_prior_max[tree->n_root->v[2]->num],
1941            tree->rates->t_prior_max[tree->n_root->v[1]->num]));
1942
1943
1944  /* Set all t_prior_min values */
1945  if(!tree->rates->t_has_prior[tree->n_root->num])
1946    {
1947      min_prior = 1.E+10;
1948      For(i,2*tree->n_otu-2)
1949        {
1950          if(tree->rates->t_has_prior[i])
1951            {
1952              if(tree->rates->t_prior_min[i] < min_prior)
1953                min_prior = tree->rates->t_prior_min[i];
1954            }
1955        }
1956      tree->rates->t_prior_min[tree->n_root->num] = 2.0 * min_prior;
1957      /* tree->rates->t_prior_min[tree->n_root->num] = 10. * min_prior; */
1958    }
1959 
1960  if(tree->rates->t_prior_min[tree->n_root->num] > 0.0)
1961    {
1962      /* PhyML_Printf("\n== Failed to set the lower bound for the root node."); */
1963      /* PhyML_Printf("\n== Make sure at least one of the calibration interval"); */
1964      /* PhyML_Printf("\n== provides a lower bound."); */
1965      /* Exit("\n"); */
1966      *result = FALSE;
1967    }
1968
1969
1970  TIMES_Set_All_Node_Priors_Top_Down_S(tree->n_root,tree->n_root->v[2], result, tree);
1971  TIMES_Set_All_Node_Priors_Top_Down_S(tree->n_root,tree->n_root->v[1], result, tree);
1972
1973  /* Get_Node_Ranks(tree); */
1974  /* TIMES_Set_Floor(tree); */
1975}
1976
1977//////////////////////////////////////////////////////////////
1978//////////////////////////////////////////////////////////////
1979
1980
1981void TIMES_Set_All_Node_Priors_Bottom_Up_S(t_node *a, t_node *d, int *result, t_tree *tree)
1982{
1983  int i;
1984  phydbl t_sup;
1985
1986  if(d->tax) return;
1987  else 
1988    {
1989      t_node *v1, *v2; /* the two sons of d */
1990
1991      For(i,3)
1992        {
1993          if((d->v[i] != a) && (d->b[i] != tree->e_root))
1994            {
1995              TIMES_Set_All_Node_Priors_Bottom_Up_S(d,d->v[i], result, tree);         
1996            }
1997        }
1998     
1999      v1 = v2 = NULL;
2000      For(i,3) if((d->v[i] != a) && (d->b[i] != tree->e_root)) 
2001        {
2002          if(!v1) v1 = d->v[i]; 
2003          else    v2 = d->v[i];
2004        }
2005     
2006      if(tree->rates->t_has_prior[d->num] == YES)
2007        {
2008          t_sup = MIN(tree->rates->t_prior_max[d->num],
2009                      MIN(tree->rates->t_prior_max[v1->num],
2010                          tree->rates->t_prior_max[v2->num]));
2011
2012          tree->rates->t_prior_max[d->num] = t_sup;
2013
2014          if(tree->rates->t_prior_max[d->num] < tree->rates->t_prior_min[d->num])
2015            {
2016              /* PhyML_Printf("\n. prior_min=%f prior_max=%f",tree->rates->t_prior_min[d->num],tree->rates->t_prior_max[d->num]); */
2017              /* PhyML_Printf("\n. Inconsistency in the prior settings detected at t_node %d",d->num); */
2018              /* PhyML_Printf("\n. Err in file %s at line %d\n\n",__FILE__,__LINE__); */
2019              /* Warn_And_Exit("\n"); */
2020              *result = FALSE;
2021              /* return; */
2022            }
2023        }
2024      else
2025        {
2026          tree->rates->t_prior_max[d->num] = 
2027            MIN(tree->rates->t_prior_max[v1->num],
2028                tree->rates->t_prior_max[v2->num]);
2029        }
2030    }
2031}
2032
2033//////////////////////////////////////////////////////////////
2034//////////////////////////////////////////////////////////////
2035
2036
2037void TIMES_Set_All_Node_Priors_Top_Down_S(t_node *a, t_node *d, int *result, t_tree *tree)
2038{
2039  if(d->tax) return;
2040  else
2041    {
2042      int i;     
2043     
2044      if(tree->rates->t_has_prior[d->num] == YES)
2045        {
2046          tree->rates->t_prior_min[d->num] = MAX(tree->rates->t_prior_min[d->num],tree->rates->t_prior_min[a->num]);
2047         
2048          if(tree->rates->t_prior_max[d->num] < tree->rates->t_prior_min[d->num])
2049            {
2050              /* PhyML_Printf("\n. prior_min=%f prior_max=%f",tree->rates->t_prior_min[d->num],tree->rates->t_prior_max[d->num]); */
2051              /* PhyML_Printf("\n. Inconsistency in the prior settings detected at t_node %d",d->num); */
2052              /* PhyML_Printf("\n. Err in file %s at line %d\n\n",__FILE__,__LINE__); */
2053              /* Warn_And_Exit("\n"); */
2054              *result = FALSE;
2055              /* return; */
2056            }
2057        }
2058      else
2059        {
2060          tree->rates->t_prior_min[d->num] = tree->rates->t_prior_min[a->num];
2061        }
2062           
2063      For(i,3)
2064        {
2065          if((d->v[i] != a) && (d->b[i] != tree->e_root))
2066            {
2067              TIMES_Set_All_Node_Priors_Top_Down_S(d,d->v[i], result, tree);
2068            }
2069        }
2070    }
2071}
2072
2073//////////////////////////////////////////////////////////////
2074//////////////////////////////////////////////////////////////
2075
Note: See TracBrowser for help on using the repository browser.