source: tags/arb-6.0/GDE/PHYLIP/fitch.c

Last change on this file was 2175, checked in by westram, 20 years ago

upgrade to PHYLIP 3.6

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