source: trunk/GDE/PHYLIP/fitch.c

Last change on this file was 19480, checked in by westram, 15 months ago
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 27.3 KB
Line 
1
2#include "phylip.h"
3#include "dist.h"
4
5/* version 3.6. (c) Copyright 1993-2002 by the University of Washington.
6   Written by Joseph Felsenstein, Akiko Fuseki, Sean Lamont, and Andrew Keeffe.
7   Permission is granted to copy and use this program provided no fee is
8   charged for it and provided that this copyright notice is not removed. */
9
10#define smoothings      4   /* number of zero-branch correction iterations */
11#define epsilonf        0.000001   /* a very small but not too small number  */
12#define delta           0.01      /* a not quite so small number */
13#define MAXNUMTREES     1000000 /* a number bigger than conceivable numtrees */
14
15
16#ifndef OLDC
17/* function prototypes */
18void   getoptions(void);
19void   allocrest(void);
20void   doinit(void);
21void   inputoptions(void);
22void   fitch_getinput(void);
23void   secondtraverse(node *, double , long *, double *);
24void   firsttraverse(node *, long *, double *);
25double evaluate(tree *);
26void   nudists(node *, node *);
27void   makedists(node *);
28
29void   makebigv(node *);
30void   correctv(node *);
31void   alter(node *, node *);
32void   nuview(node *);
33void   update(node *);
34void   smooth(node *);
35void   filltraverse(node *, node *, boolean);
36void   fillin(node *, node *, boolean);
37void   insert_(node *, node *, boolean);
38void   copynode(node *, node *);
39
40void   copy_(tree *, tree *);
41void   setuptipf(long, tree *);
42void   buildnewtip(long , tree *, long);
43void   buildsimpletree(tree *, long);
44void   addtraverse(node *, node *, boolean, long *, boolean *);
45void   re_move(node **, node **);
46void   rearrange(node *, long *, long *, boolean *);
47void   describe(node *);
48void   nodeinit(node *);
49
50void   initrav(node *);
51void   treevaluate(void);
52void   maketree(void);
53/* function prototypes */
54#endif
55
56
57
58Char infilename[FNMLNGTH], outfilename[FNMLNGTH], intreename[FNMLNGTH], outtreename[FNMLNGTH];
59long nonodes2, outgrno, nums, col, datasets, ith, njumble, jumb=0;
60long inseed;
61vector *global_x;
62intvector *reps;
63boolean minev, global, jumble, lengths, usertree, lower, upper, negallowed,
64        outgropt, replicates, trout, printdata, progress, treeprint,
65        mulsets, firstset;
66double power;
67double trweight; /* to make treeread happy */
68boolean goteof, haslengths;  /* ditto ... */
69boolean first; /* ditto ... */
70
71longer seed;
72long *enterorder;
73tree curtree, priortree, bestree, bestree2;
74Char global_ch;
75char *progname;
76
77
78
79void getoptions()
80{
81  /* interactively set options */
82  long inseed0=0, loopcount;
83  Char ch;
84  boolean done=false;
85
86  putchar('\n');
87  minev = false;
88  global = false;
89  jumble = false;
90  njumble = 1;
91  lengths = false;
92  lower = false;
93  negallowed = false;
94  outgrno = 1;
95  outgropt = false;
96  power = 2.0;
97  replicates = false;
98  trout = true;
99  upper = false;
100  usertree = false;
101  printdata = false;
102  progress = true;
103  treeprint = true;
104  loopcount = 0;
105  do {
106    cleerhome();
107    printf("\nFitch-Margoliash method version %s\n\n",VERSION);
108    printf("Settings for this run:\n");
109    printf("  D      Method (F-M, Minimum Evolution)?  %s\n",
110             (minev ? "Minimum Evolution" : "Fitch-Margoliash"));
111    printf("  U                 Search for best tree?  %s\n",
112           (usertree ? "No, use user trees in input file" : "Yes"));
113    if (usertree) {
114      printf("  N          Use lengths from user trees?  %s\n",
115             (lengths ? "Yes" : "No"));
116    }
117    printf("  P                                Power?%9.5f\n",power);
118    printf("  -      Negative branch lengths allowed?  %s\n",
119           negallowed ? "Yes" : "No");
120    printf("  O                        Outgroup root?");
121    if (outgropt)
122      printf("  Yes, at species number%3ld\n", outgrno);
123    else
124      printf("  No, use as outgroup species%3ld\n", outgrno);
125    printf("  L         Lower-triangular data matrix?");
126    if (lower)
127      printf("  Yes\n");
128    else
129      printf("  No\n");
130    printf("  R         Upper-triangular data matrix?");
131    if (upper)
132      printf("  Yes\n");
133    else
134      printf("  No\n");
135    printf("  S                        Subreplicates?");
136    if (replicates)
137      printf("  Yes\n");
138    else
139      printf("  No\n");
140    if (!usertree) {
141      printf("  G                Global rearrangements?");
142      if (global)
143        printf("  Yes\n");
144      else
145        printf("  No\n");
146      printf("  J     Randomize input order of species?");
147      if (jumble)
148        printf("  Yes (seed =%8ld,%3ld times)\n", inseed0, njumble);
149      else
150        printf("  No. Use input order\n");
151    }
152    printf("  M           Analyze multiple data sets?");
153    if (mulsets)
154      printf("  Yes, %2ld sets\n", datasets);
155    else
156      printf("  No\n");
157    printf("  0   Terminal type (IBM PC, ANSI, none)?");
158    if (ibmpc)
159      printf("  IBM PC\n");
160    if (ansi)
161      printf("  ANSI\n");
162    if (!(ibmpc || ansi))
163      printf("  (none)\n");
164    printf("  1    Print out the data at start of run");
165    if (printdata)
166      printf("  Yes\n");
167    else
168      printf("  No\n");
169    printf("  2  Print indications of progress of run");
170    if (progress)
171      printf("  Yes\n");
172    else
173      printf("  No\n");
174    printf("  3                        Print out tree");
175    if (treeprint)
176      printf("  Yes\n");
177    else
178      printf("  No\n");
179    printf("  4       Write out trees onto tree file?");
180    if (trout)
181      printf("  Yes\n");
182    else
183      printf("  No\n");
184    printf(
185   "\n  Y to accept these or type the letter for one to change\n");
186#ifdef WIN32
187    phyFillScreenColor();
188#endif
189    scanf("%c%*[^\n]", &ch);
190    getchar();
191    uppercase(&ch);
192    done = (ch == 'Y');
193   if (!done) {
194      if (strchr("DJOUNPG-LRSM01234",ch) != NULL) {
195        switch (ch) {
196
197        case 'D':
198          minev = !minev;
199          if (minev && (!negallowed))
200            negallowed = true;
201          break;
202
203        case '-':
204          negallowed = !negallowed;
205          break;
206
207        case 'G':
208          global = !global;
209          break;
210
211        case 'J':
212          jumble = !jumble;
213          if (jumble)
214            initjumble(&inseed, &inseed0, seed, &njumble);
215          else njumble = 1;
216          break;
217
218        case 'L':
219          lower = !lower;
220          break;
221
222         case 'N':
223           lengths = !lengths;
224           break;
225
226        case 'O':
227          outgropt = !outgropt;
228          if (outgropt)
229            initoutgroup(&outgrno, spp);
230          break;
231
232        case 'P':
233          initpower(&power);
234          break;
235
236        case 'R':
237          upper = !upper;
238          break;
239
240        case 'S':
241          replicates = !replicates;
242          break;
243
244        case 'U':
245          usertree = !usertree;
246          break;
247
248        case 'M':
249          mulsets = !mulsets;
250          if (mulsets)
251            initdatasets(&datasets);
252          break;
253
254        case '0':
255          initterminal(&ibmpc, &ansi);
256          break;
257
258        case '1':
259          printdata = !printdata;
260          break;
261
262        case '2':
263          progress = !progress;
264          break;
265
266        case '3':
267          treeprint = !treeprint;
268          break;
269
270        case '4':
271          trout = !trout;
272          break;
273        }
274      } else
275        printf("Not a possible option!\n");
276    }
277    countup(&loopcount, 100);
278  } while (!done);
279  if (lower && upper) {
280    printf("ERROR: Data matrix cannot be both uppeR and Lower triangular\n");
281    exxit(-1);
282  }
283}  /* getoptions */
284
285
286void allocrest()
287{
288  long i;
289
290  global_x = (vector *)Malloc(spp*sizeof(vector));
291  reps = (intvector *)Malloc(spp*sizeof(intvector));
292  for (i=0;i<spp;++i){
293    global_x[i]=(vector)Malloc(nonodes2 * sizeof(double));
294    reps[i]=(intvector)Malloc(spp * sizeof(long));
295  }
296  nayme = (naym *)Malloc(spp*sizeof(naym));
297  enterorder = (long *)Malloc(spp*sizeof(long));
298}
299
300
301void doinit()
302{
303  /* initializes variables */
304
305  inputnumbers2(&spp, &nonodes2, 2);
306  getoptions();
307  alloctree(&curtree.nodep, nonodes2);
308  allocd(nonodes2, curtree.nodep);
309  allocw(nonodes2, curtree.nodep);
310  if (!usertree) {
311    alloctree(&bestree.nodep, nonodes2);
312    allocd(nonodes2, bestree.nodep);
313    allocw(nonodes2, bestree.nodep);
314    alloctree(&priortree.nodep, nonodes2);
315    allocd(nonodes2, priortree.nodep);
316    allocw(nonodes2, priortree.nodep);
317    if (njumble > 1) {
318      alloctree(&bestree2.nodep, nonodes2);
319      allocd(nonodes2, bestree2.nodep);
320      allocw(nonodes2, bestree2.nodep);
321    }
322  }
323  allocrest();
324}  /* doinit */
325
326
327void inputoptions()
328{
329  /* print options information */
330  if (!firstset)
331    samenumsp2(ith);
332  fprintf(outfile, "\nFitch-Margoliash method version %s\n\n",VERSION);
333  if (minev)
334    fprintf(outfile, "Minimum evolution method option\n\n");
335  fprintf(outfile, "                  __ __             2\n");
336  fprintf(outfile, "                  \\  \\   (Obs - Exp)\n");
337  fprintf(outfile, "Sum of squares =  /_ /_  ------------\n");
338  fprintf(outfile, "                               ");
339  if (power == (long)power)
340    fprintf(outfile, "%2ld\n", (long)power);
341  else
342    fprintf(outfile, "%4.1f\n", power);
343  fprintf(outfile, "                   i  j      Obs\n\n");
344  fprintf(outfile, "Negative branch lengths ");
345  if (!negallowed)
346    fprintf(outfile, "not ");
347  fprintf(outfile, "allowed\n\n");
348  if (global)
349    fprintf(outfile, "global optimization\n\n");
350}  /* inputoptions */
351
352
353void fitch_getinput()
354{
355  /* reads the input data */
356  inputoptions();
357}  /* fitch_getinput */
358
359
360void secondtraverse(node *q, double y, long *nx, double *sum)
361{
362  /* from each of those places go back to all others */
363   /* nx comes from firsttraverse */
364   /* sum comes from evaluate via firsttraverse */
365  double z=0.0, TEMP=0.0;
366
367  z = y + q->v;
368  if (q->tip) {
369    TEMP = q->d[(*nx) - 1] - z;
370    *sum += q->w[(*nx) - 1] * (TEMP * TEMP);
371  } else {
372    secondtraverse(q->next->back, z, nx, sum);
373    secondtraverse(q->next->next->back, z, nx,sum);
374  }
375}  /* secondtraverse */
376
377
378void firsttraverse(node *p, long *nx, double *sum)
379{
380  /* go through tree calculating branch lengths */
381  if (minev && (p != curtree.start))
382    *sum += p->v;
383  if (p->tip) {
384    if (!minev) {
385      *nx = p->index;
386      secondtraverse(p->back, 0.0, nx, sum);
387      }
388  } else {
389    firsttraverse(p->next->back, nx,sum);
390    firsttraverse(p->next->next->back, nx,sum);
391  }
392}  /* firsttraverse */
393
394
395double evaluate(tree *t)
396{
397  double sum=0.0;
398  long nx=0;
399  /* evaluate likelihood of a tree */
400  firsttraverse(t->start->back ,&nx, &sum);
401  firsttraverse(t->start, &nx, &sum);
402  if ((!minev) && replicates && (lower || upper))
403    sum /= 2;
404  t->likelihood = -sum;
405  return (-sum);
406}  /* evaluate */
407
408
409void nudists(node *x, node *y)
410{
411  /* compute distance between an interior node and tips */
412  long nq=0, nr=0, nx=0, ny=0;
413  double dil=0, djl=0, wil=0, wjl=0, vi=0, vj=0;
414  node *qprime, *rprime;
415
416  qprime = x->next;
417  rprime = qprime->next->back;
418  qprime = qprime->back;
419  ny = y->index;
420  dil = qprime->d[ny - 1];
421  djl = rprime->d[ny - 1];
422  wil = qprime->w[ny - 1];
423  wjl = rprime->w[ny - 1];
424  vi = qprime->v;
425  vj = rprime->v;
426  x->w[ny - 1] = wil + wjl;
427  if (wil + wjl <= 0.0)
428    x->d[ny - 1] = 0.0;
429  else
430    x->d[ny - 1] = ((dil - vi) * wil + (djl - vj) * wjl) / (wil + wjl);
431  nx = x->index;
432  nq = qprime->index;
433  nr = rprime->index;
434  dil = y->d[nq - 1];
435  djl = y->d[nr - 1];
436  wil = y->w[nq - 1];
437  wjl = y->w[nr - 1];
438  y->w[nx - 1] = wil + wjl;
439  if (wil + wjl <= 0.0)
440    y->d[nx - 1] = 0.0;
441  else
442    y->d[nx - 1] = ((dil - vi) * wil + (djl - vj) * wjl) / (wil + wjl);
443}  /* nudists */
444
445
446void makedists(node *p)
447{
448  /* compute distances among three neighbors of a node */
449  long i=0, nr=0, ns=0;
450  node *q, *r, *s;
451
452  r = p->back;
453  nr = r->index;
454  for (i = 1; i <= 3; i++) {
455    q = p->next;
456    s = q->back;
457    ns = s->index;
458    if (s->w[nr - 1] + r->w[ns - 1] <= 0.0)
459      p->dist = 0.0;
460    else
461      p->dist = (s->w[nr - 1] * s->d[nr - 1] + r->w[ns - 1] * r->d[ns - 1]) /
462                (s->w[nr - 1] + r->w[ns - 1]);
463    p = q;
464    r = s;
465    nr = ns;
466  }
467}  /* makedists */
468
469
470void makebigv(node *p)
471{
472  /* make new branch length */
473  long i=0;
474  node *temp, *q, *r;
475
476  q = p->next;
477  r = q->next;
478  for (i = 1; i <= 3; i++) {
479    if (p->iter) {
480      p->v = (p->dist + r->dist - q->dist) / 2.0;
481      p->back->v = p->v;
482    }
483    temp = p;
484    p = q;
485    q = r;
486    r = temp;
487  }
488}  /* makebigv */
489
490
491void correctv(node *p)
492{
493  /* iterate branch lengths if some are to be zero */
494  node *q, *r, *temp;
495  long i=0, j=0, n=0, nq=0, nr=0, ntemp=0;
496  double wq=0.0, wr=0.0;
497
498  q = p->next;
499  r = q->next;
500  n = p->back->index;
501  nq = q->back->index;
502  nr = r->back->index;
503  for (i = 1; i <= smoothings; i++) {
504    for (j = 1; j <= 3; j++) {
505      if (p->iter) {
506        wr = r->back->w[n - 1] + p->back->w[nr - 1];
507        wq = q->back->w[n - 1] + p->back->w[nq - 1];
508        if (wr + wq <= 0.0 && !negallowed)
509          p->v = 0.0;
510        else
511          p->v = ((p->dist - q->v) * wq + (r->dist - r->v) * wr) / (wr + wq);
512        if (p->v < 0 && !negallowed)
513          p->v = 0.0;
514        p->back->v = p->v;
515      }
516      temp = p;
517      p = q;
518      q = r;
519      r = temp;
520      ntemp = n;
521      n = nq;
522      nq = nr;
523      nr = ntemp;
524    }
525  }
526}  /* correctv */
527
528
529void alter(node *x, node *y)
530{
531  /* traverse updating these views */
532  nudists(x, y);
533  if (!y->tip) {
534    alter(x, y->next->back);
535    alter(x, y->next->next->back);
536  }
537}  /* alter */
538
539
540void nuview(node *p)
541{
542  /* renew information about subtrees */
543  long i=0;
544  node *q, *r, *pprime, *temp;
545
546  q = p->next;
547  r = q->next;
548  for (i = 1; i <= 3; i++) {
549    temp = p;
550    pprime = p->back;
551    alter(p, pprime);
552    p = q;
553    q = r;
554    r = temp;
555  }
556}  /* nuview */
557
558
559void update(node *p)
560{
561  /* update branch lengths around a node */
562
563  if (p->tip)
564    return;
565  makedists(p);
566  if (p->iter || p->next->iter || p->next->next->iter) {
567    makebigv(p);
568    correctv(p);
569  }
570  nuview(p);
571}  /* update */
572
573
574void smooth(node *p)
575{
576  /* go through tree getting new branch lengths and views */
577  if (p->tip)
578    return;
579  update(p);
580  smooth(p->next->back);
581  smooth(p->next->next->back);
582}  /* smooth */
583
584
585void filltraverse(node *pb, node *qb, boolean contin)
586{
587  if (qb->tip)
588    return;
589  if (contin) {
590    filltraverse(pb, qb->next->back,contin);
591    filltraverse(pb, qb->next->next->back,contin);
592    nudists(qb, pb);
593    return;
594  }
595  if (!qb->next->back->tip)
596    nudists(qb->next->back, pb);
597  if (!qb->next->next->back->tip)
598    nudists(qb->next->next->back, pb);
599}  /* filltraverse */
600
601
602void fillin(node *pa, node *qa, boolean contin)
603{
604  if (!pa->tip) {
605    fillin(pa->next->back, qa, contin);
606    fillin(pa->next->next->back, qa, contin);
607  }
608  filltraverse(pa, qa, contin);
609}  /* fillin */
610
611
612void insert_(node *p, node *q, boolean contin_)
613{
614  /* put p and q together and iterate info. on resulting tree */
615  double x=0.0, oldlike /*, dummy*/;
616  hookup(p->next->next, q->back);
617  hookup(p->next, q);
618  x = q->v / 2.0;
619  p->v = 0.0;
620  p->back->v = 0.0;
621  p->next->v = x;
622  p->next->back->v = x;
623  p->next->next->back->v = x;
624  p->next->next->v = x;
625  fillin(p->back, p, contin_);
626  /* dummy =*/ evaluate(&curtree);
627  do {
628    oldlike = curtree.likelihood;
629    smooth(p);
630    smooth(p->back);
631    /*dummy =*/ evaluate(&curtree);
632  } while (fabs(curtree.likelihood - oldlike) > delta);
633}  /* insert_ */
634
635
636void copynode(node *c, node *d)
637{
638  /* make a copy of a node */
639
640  memcpy(d->d, c->d, nonodes2*sizeof(double));
641  memcpy(d->w, c->w, nonodes2*sizeof(double));
642  d->v = c->v;
643  d->iter = c->iter;
644  d->dist = c->dist;
645  d->xcoord = c->xcoord;
646  d->ycoord = c->ycoord;
647  d->ymin = c->ymin;
648  d->ymax = c->ymax;
649}  /* copynode */
650
651
652void copy_(tree *a, tree *b)
653{
654  /* make a copy of a tree */
655  long i, j=0;
656  node *p, *q;
657
658  for (i = 0; i < spp; i++) {
659    copynode(a->nodep[i], b->nodep[i]);
660    if (a->nodep[i]->back) {
661      if (a->nodep[i]->back == a->nodep[a->nodep[i]->back->index - 1])
662        b->nodep[i]->back = b->nodep[a->nodep[i]->back->index - 1];
663      else if (a->nodep[i]->back
664                 == a->nodep[a->nodep[i]->back->index - 1]->next)
665        b->nodep[i]->back = b->nodep[a->nodep[i]->back->index - 1]->next;
666      else
667        b->nodep[i]->back
668          = b->nodep[a->nodep[i]->back->index - 1]->next->next;
669    }
670    else b->nodep[i]->back = NULL;
671  }
672  for (i = spp; i < nonodes2; i++) {
673    p = a->nodep[i];
674    q = b->nodep[i];
675    for (j = 1; j <= 3; j++) {
676      copynode(p, q);
677      if (p->back) {
678        if (p->back == a->nodep[p->back->index - 1])
679          q->back = b->nodep[p->back->index - 1];
680        else if (p->back == a->nodep[p->back->index - 1]->next)
681          q->back = b->nodep[p->back->index - 1]->next;
682        else
683          q->back = b->nodep[p->back->index - 1]->next->next;
684      }
685      else
686        q->back = NULL;
687      p = p->next;
688      q = q->next;
689    }
690  }
691  b->likelihood = a->likelihood;
692  b->start = a->start;
693}  /* copy_ */
694
695
696void setuptipf(long m, tree *t)
697{
698  /* initialize branch lengths and views in a tip */
699  long i=0;
700  intvector n=(long *)Malloc(spp * sizeof(long)); 
701  node *WITH;
702
703  WITH = t->nodep[m - 1];
704  memcpy(WITH->d, global_x[m - 1], (nonodes2 * sizeof(double)));
705  memcpy(n, reps[m - 1], (spp * sizeof(long)));
706  for (i = 0; i < spp; i++) {
707    if (i + 1 != m && n[i] > 0) {
708      if (WITH->d[i] < epsilonf)
709        WITH->d[i] = epsilonf;
710      WITH->w[i] = n[i] / exp(power * log(WITH->d[i]));
711    } else {
712      WITH->w[m - 1] = 0.0;
713      WITH->d[m - 1] = 0.0;
714    }
715  }
716  for (i = spp; i < nonodes2; i++) {
717    WITH->w[i] = 1.0;
718    WITH->d[i] = 0.0;
719  }
720  WITH->index = m;
721  if (WITH->iter) WITH->v = 0.0;
722  free(n);
723}  /* setuptipf */
724
725
726void buildnewtip(long m, tree *t, long nextsp)
727{
728  /* initialize and hook up a new tip */
729  node *p;
730  setuptipf(m, t);
731  p = t->nodep[nextsp + spp - 3];
732  hookup(t->nodep[m - 1], p);
733}  /* buildnewtip */
734
735
736void buildsimpletree(tree *t, long nextsp)
737{
738  /* make and initialize a three-species tree */
739  curtree.start=curtree.nodep[enterorder[0] - 1]; 
740  setuptipf(enterorder[0], t);
741  setuptipf(enterorder[1], t);
742  hookup(t->nodep[enterorder[0] - 1], t->nodep[enterorder[1] - 1]);
743  buildnewtip(enterorder[2], t, nextsp);
744  insert_(t->nodep[enterorder[2] - 1]->back, t->nodep[enterorder[0] - 1],
745          false);
746}  /* buildsimpletree */
747
748
749void addtraverse(node *p, node *q, boolean contin, long *numtrees,
750                        boolean *succeeded)
751{
752 /* traverse through a tree, finding best place to add p */
753  insert_(p, q, true);
754  (*numtrees)++;
755  if (evaluate(&curtree) > bestree.likelihood){
756    copy_(&curtree, &bestree);
757    (*succeeded)=true;
758  }
759  copy_(&priortree, &curtree);
760  if (!q->tip && contin) {
761    addtraverse(p, q->next->back, contin,numtrees,succeeded);
762    addtraverse(p, q->next->next->back, contin,numtrees,succeeded);
763  }
764}  /* addtraverse */
765
766
767void re_move(node **p, node **q)
768{
769  /* re_move p and record in q where it was */
770  *q = (*p)->next->back;
771  hookup(*q, (*p)->next->next->back);
772  (*p)->next->back = NULL;
773  (*p)->next->next->back = NULL;
774  update(*q);
775  update((*q)->back);
776}  /* re_move */
777
778
779void rearrange(node *p, long *numtrees, long *nextsp, boolean *succeeded)
780{
781  node *q, *r;
782  if (!p->tip && !p->back->tip) {
783    r = p->next->next;
784    re_move(&r, &q);
785    copy_(&curtree, &priortree);
786    addtraverse(r, q->next->back, (boolean)(global && ((*nextsp) == spp)),
787                numtrees,succeeded);
788    addtraverse(r, q->next->next->back, (boolean)(global && ((*nextsp) == spp)),
789                numtrees,succeeded);
790    copy_(&bestree, &curtree);
791    if (global && ((*nextsp) == spp)) {
792      putchar('.');
793      fflush(stdout);
794    }
795    if (global && ((*nextsp) == spp) && !(*succeeded)) {
796      if (r->back->tip) {
797        r = r->next->next;
798        re_move(&r, &q);
799        q = q->back;
800        copy_(&curtree, &priortree);
801        if (!q->tip) {
802          addtraverse(r, q->next->back, true, numtrees,succeeded);
803          addtraverse(r, q->next->next->back, true, numtrees,succeeded);
804        }
805        q = q->back;
806        if (!q->tip) {
807          addtraverse(r, q->next->back, true, numtrees,succeeded);
808          addtraverse(r, q->next->next->back, true, numtrees,succeeded);
809        }
810        copy_(&bestree, &curtree);
811      }
812    }
813  }
814  if (!p->tip) {
815    rearrange(p->next->back, numtrees,nextsp,succeeded);
816    rearrange(p->next->next->back, numtrees,nextsp,succeeded);
817  }
818}  /* rearrange */
819
820
821void describe(node *p)
822{
823  /* print out information for one branch */
824  long i=0;
825  node *q;
826
827  q = p->back;
828  fprintf(outfile, "%4ld          ", q->index - spp);
829  if (p->tip) {
830    for (i = 0; i < nmlngth; i++)
831      putc(nayme[p->index - 1][i], outfile);
832  } else
833    fprintf(outfile, "%4ld      ", p->index - spp);
834  fprintf(outfile, "%15.5f\n", q->v);
835  if (!p->tip) {
836    describe(p->next->back);
837    describe(p->next->next->back);
838  }
839}  /* describe */
840
841
842void summarize()
843{
844  /* print out branch lengths etc. */
845  long i, j, totalnum;
846
847  fprintf(outfile, "\nremember:");
848  if (outgropt)
849    fprintf(outfile, " (although rooted by outgroup)");
850  fprintf(outfile, " this is an unrooted tree!\n\n");
851  if (!minev)
852    fprintf(outfile, "Sum of squares = %11.5f\n\n", -curtree.likelihood);
853  else
854    fprintf(outfile, "Sum of branch lengths = %11.5f\n\n", -curtree.likelihood);
855  if ((power == 2.0) && !minev) {
856    totalnum = 0;
857    for (i = 1; i <= nums; i++) {
858      for (j = 1; j <= nums; j++) {
859        if (i != j)
860          totalnum += reps[i - 1][j - 1];
861      }
862    }
863    fprintf(outfile, "Average percent standard deviation = ");
864    fprintf(outfile, "%11.5f\n\n",
865            100 * sqrt(-curtree.likelihood / (totalnum - 2)));
866  }
867  fprintf(outfile, "Between        And            Length\n");
868  fprintf(outfile, "-------        ---            ------\n");
869  describe(curtree.start->next->back);
870  describe(curtree.start->next->next->back);
871  describe(curtree.start->back);
872  fprintf(outfile, "\n\n");
873  if (trout) {
874    col = 0;
875    treeout(curtree.start, &col, 0.43429445222, true,
876              curtree.start);
877  }
878}  /* summarize */
879
880
881void nodeinit(node *p)
882{
883  /* initialize a node */
884  long i, j;
885
886  for (i = 1; i <= 3; i++) {
887    for (j = 0; j < nonodes2; j++) {
888      p->w[j] = 1.0;
889      p->d[j] = 0.0;
890    }
891    p = p->next;
892  }
893  if (p->iter)
894    p->v = 1.0;
895  if (p->back->iter)
896    p->back->v = 1.0;
897}  /* nodeinit */
898
899
900void initrav(node *p)
901{
902  /* traverse to initialize */
903  if (p->tip)
904    return;
905  nodeinit(p);
906  initrav(p->next->back);
907  initrav(p->next->next->back);
908}  /* initrav */
909
910
911void treevaluate()
912{
913  /* evaluate user-defined tree, iterating branch lengths */
914  long i;
915  double /*dummy,*/ oldlike;
916
917  for (i = 1; i <= spp; i++)
918    setuptipf(i, &curtree);
919  initrav(curtree.start);
920  if (curtree.start->back != NULL) {
921    initrav(curtree.start->back);
922    /*dummy =*/ evaluate(&curtree);
923    do {
924      oldlike = curtree.likelihood;
925      smooth(curtree.start);
926      /*dummy =*/ evaluate(&curtree);
927    } while (fabs(curtree.likelihood - oldlike) > delta);
928  }
929  /*dummy =*/ evaluate(&curtree);
930}  /* treevaluate */
931
932
933void maketree()
934{
935  /* contruct the tree */
936  long nextsp,numtrees;
937  boolean succeeded=false;
938  long i, j, which;
939
940  if (usertree) {
941    inputdata(replicates, printdata, lower, upper, global_x, reps);
942    setuptree(&curtree, nonodes2);
943    for (which = 1; which <= spp; which++)
944      setuptipf(which, &curtree);
945    if (eoln(infile))
946      scan_eoln(infile);
947    openfile(&intree,INTREE,"input tree file","r",progname,intreename);
948    numtrees = countsemic(&intree);
949    if (numtrees > MAXNUMTREES) {
950      printf("\nERROR: number of input trees is read incorrectly from %s\n",
951        intreename);
952      exxit(-1);
953    }
954    if (treeprint) {
955      fprintf(outfile, "User-defined tree");
956      if (numtrees > 1)
957        putc('s', outfile);
958      fprintf(outfile, ":\n\n");
959    }
960    first = true;
961    which = 1;
962    while (which <= numtrees) {
963      treeread2 (intree, &curtree.start, curtree.nodep,
964        lengths, &trweight, &goteof, &haslengths, &spp);
965      nums = spp;
966      curtree.start = curtree.nodep[outgrno - 1]->back;
967      treevaluate();
968      printree(curtree.start, treeprint, false, false);
969      summarize();
970      which++;
971    }
972    FClose(intree);
973  } else {
974    if (jumb == 1) {
975      inputdata(replicates, printdata, lower, upper, global_x, reps);
976      setuptree(&curtree, nonodes2);
977      setuptree(&priortree, nonodes2);
978      setuptree(&bestree, nonodes2);
979      if (njumble > 1) setuptree(&bestree2, nonodes2);
980    }
981    for (i = 1; i <= spp; i++)
982      enterorder[i - 1] = i;
983    if (jumble)
984      randumize(seed, enterorder);
985    nextsp = 3;
986    buildsimpletree(&curtree, nextsp);
987    curtree.start = curtree.nodep[enterorder[0] - 1]->back;
988    if (jumb == 1) numtrees = 1;
989    nextsp = 4;
990    if (progress) {
991      printf("Adding species:\n");
992      writename(0, 3, enterorder);
993#ifdef WIN32
994      phyFillScreenColor();
995#endif
996    }
997    while (nextsp <= spp) {
998      nums = nextsp;
999      buildnewtip(enterorder[nextsp - 1], &curtree, nextsp);
1000      copy_(&curtree, &priortree);
1001      bestree.likelihood = -99999.0;
1002      addtraverse(curtree.nodep[enterorder[nextsp - 1] - 1]->back,
1003                  curtree.start, true, &numtrees,&succeeded);
1004      copy_(&bestree, &curtree);
1005      if (progress) {
1006        writename(nextsp  - 1, 1, enterorder);
1007#ifdef WIN32
1008        phyFillScreenColor();
1009#endif
1010      }
1011      if (global && nextsp == spp) {
1012        if (progress) {
1013          printf("Doing global rearrangements\n");
1014          printf("  !");
1015          for (j = 1; j <= (spp - 2); j++)
1016            putchar('-');
1017          printf("!\n");
1018          printf("   ");
1019        }
1020      }
1021      succeeded = true;
1022      while (succeeded) {
1023        succeeded = false;
1024        rearrange(curtree.start,
1025                  &numtrees,&nextsp,&succeeded);
1026        if (global && ((nextsp) == spp) && progress)
1027          printf("\n   ");
1028      }
1029      if (global && nextsp == spp) {
1030        putc('\n', outfile);
1031        if (progress)
1032          putchar('\n');
1033      }
1034      if (njumble > 1) {
1035        if (jumb == 1 && nextsp == spp)
1036          copy_(&bestree, &bestree2);
1037        else if (nextsp == spp) {
1038          if (bestree2.likelihood < bestree.likelihood)
1039            copy_(&bestree, &bestree2);
1040        }
1041      }
1042      if (nextsp == spp && jumb == njumble) {
1043        if (njumble > 1) copy_(&bestree2, &curtree);
1044        curtree.start = curtree.nodep[outgrno - 1]->back;
1045        printree(curtree.start, treeprint, true, false);
1046        summarize();
1047      }
1048      nextsp++;
1049    }
1050  }
1051  if (jumb == njumble && progress) {
1052    printf("\nOutput written to file \"%s\"\n\n", outfilename);
1053    if (trout) {
1054      printf("Tree also written onto file \"%s\"\n", outtreename);
1055      putchar('\n');
1056    }
1057  }
1058}  /* maketree */
1059
1060
1061int main(int argc, Char *argv[])
1062{
1063  int i;
1064#ifdef MAC
1065  argc = 1;                /* macsetup("Fitch","");        */
1066  argv[0]="Fitch";
1067#endif
1068  init(argc,argv);
1069  progname = argv[0];
1070  openfile(&infile,INFILE,"input file","r",argv[0],infilename);
1071  openfile(&outfile,OUTFILE,"output file","w",argv[0],outfilename);
1072
1073  ibmpc = IBMCRT;
1074  ansi = ANSICRT;
1075  mulsets = false;
1076  datasets = 1;
1077  firstset = true;
1078  doinit();
1079  if (trout)
1080    openfile(&outtree,OUTTREE,"output tree file","w",argv[0],outtreename);
1081  for (i=0;i<spp;++i){
1082    enterorder[i]=0;}
1083  for (ith = 1; ith <= datasets; ith++) {
1084    if (datasets > 1) {
1085      fprintf(outfile, "Data set # %ld:\n\n",ith);
1086      if (progress)
1087        printf("\nData set # %ld:\n\n",ith);
1088    }
1089    fitch_getinput();
1090    for (jumb = 1; jumb <= njumble; jumb++)
1091        maketree();
1092    firstset = false;
1093    if (eoln(infile) && (ith < datasets))
1094      scan_eoln(infile);
1095  }
1096  if (trout)
1097    FClose(outtree);
1098  FClose(outfile);
1099  FClose(infile);
1100#ifdef MAC
1101  fixmacfile(outfilename);
1102  fixmacfile(outtreename);
1103#endif
1104  printf("Done.\n\n");
1105#ifdef WIN32
1106  phyRestoreConsoleAttributes();
1107#endif
1108  return 0;
1109}
Note: See TracBrowser for help on using the repository browser.