source: tags/ms_r18q1/GDE/PHYLIP/dnapars.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: 47.1 KB
Line 
1#include "phylip.h"
2#include "seq.h"
3
4/* version 3.6 (c) Copyright 1993-2002 by the University of Washington.
5   Written by Joseph Felsenstein, Akiko Fuseki, Sean Lamont, and Andrew Keeffe.
6   Permission is granted to copy and use this program provided no fee is
7   charged for it and provided that this copyright notice is not removed. */
8
9#define MAXNUMTREES     1000000  /* bigger than number of user trees can be */
10
11extern sequence y;
12
13
14#ifndef OLDC
15/* function prototypes */
16void   getoptions(void);
17void   allocrest(void);
18void   reallocchars(void);
19void   doinit(void);
20void   makeweights(void);
21void   doinput(void);
22void   initdnaparsnode(node **, node **, node *, long, long, long *,
23        long *, initops, pointarray, pointarray, Char *, Char *, FILE *);
24void   evaluate(node *);
25void   tryadd(node *, node *, node *);
26void   addpreorder(node *, node *, node *);
27void   trydescendants(node *, node *, node *, node *, boolean);
28
29void   trylocal(node *, node *);
30void   trylocal2(node *, node *, node *);
31void   tryrearr(node *, boolean *);
32void   repreorder(node *, boolean *);
33void   rearrange(node **);
34void   describe(void);
35void   dnapars_coordinates(node *, double, long *, double *);
36void   dnapars_printree(void);
37void   globrearrange(void);
38void   grandrearr(void);
39
40void   maketree(void);
41void   freerest(void);
42void   load_tree(long treei);
43
44/* function prototypes */
45#endif
46
47
48Char infilename[FNMLNGTH], outfilename[FNMLNGTH], intreename[FNMLNGTH], outtreename[FNMLNGTH],
49     weightfilename[FNMLNGTH];
50char basechar[32]="ACMGRSVTWYHKDBNO???????????????";
51node *root;
52long chars, col, msets, ith, njumble, jumb, maxtrees;
53/*   chars = number of sites in actual sequences */
54long inseed, inseed0;
55double threshold;
56boolean jumble, usertree, thresh, weights, thorough, rearrfirst,
57          trout, progress, stepbox, ancseq, mulsets, justwts, firstset, mulf,
58          multf;
59steptr oldweight;
60longer seed;
61pointarray treenode;            /* pointers to all nodes in tree */
62long *enterorder;
63long *zeros;
64
65/* local variables for Pascal maketree, propagated globally for C version: */
66
67long minwhich;
68double like, minsteps, bestyet, bestlike, bstlike2;
69boolean lastrearr, recompute;
70double nsteps[maxuser];
71long **fsteps;
72node *there, *oldnufork;
73long *place;
74bestelm *bestrees;
75long *threshwt;
76baseptr nothing;
77gbases *garbage;
78node *temp, *temp1, *temp2, *tempsum, *temprm, *tempadd, *tempf, *tmp, *tmp1,
79       *tmp2, *tmp3, *tmprm, *tmpadd;
80boolean *names;
81node *grbg;
82char *progname;
83
84
85void getoptions()
86{
87  /* interactively set options */
88  long loopcount, loopcount2;
89  Char ch, ch2;
90
91  fprintf(outfile, "\nDNA parsimony algorithm, version %s\n\n",VERSION);
92  jumble = false;
93  njumble = 1;
94  outgrno = 1;
95  outgropt = false;
96  thresh = false;
97  thorough = true;
98  transvp = false;
99  rearrfirst = false;
100  maxtrees = 10000;
101  trout = true;
102  usertree = false;
103  weights = false;
104  mulsets = false;
105  printdata = false;
106  progress = true;
107  treeprint = true;
108  stepbox = false;
109  ancseq = false;
110  dotdiff = true;
111  interleaved = true;
112  loopcount = 0;
113  for (;;) {
114    cleerhome();
115    printf("\nDNA parsimony algorithm, version %s\n\n",VERSION);
116    printf("Setting for this run:\n");
117    printf("  U                 Search for best tree?  %s\n",
118           (usertree ? "No, use user trees in input file" : "Yes"));
119    if (!usertree) {
120      printf("  S                        Search option?  ");
121      if (thorough)
122        printf("More thorough search\n");
123      else if (rearrfirst)
124        printf("Rearrange on one best tree\n");
125      else
126        printf("Less thorough\n");
127      printf("  V              Number of trees to save?  %ld\n", maxtrees);
128      printf("  J   Randomize input order of sequences?");
129      if (jumble)
130        printf("  Yes (seed =%8ld,%3ld times)\n", inseed0, njumble);
131      else
132        printf("  No. Use input order\n");
133    }
134    printf("  O                        Outgroup root?");
135    if (outgropt)
136      printf("  Yes, at sequence number%3ld\n", outgrno);
137    else
138      printf("  No, use as outgroup species%3ld\n", outgrno);
139    printf("  T              Use Threshold parsimony?");
140    if (thresh)
141      printf("  Yes, count steps up to%4.1f per site\n", threshold);
142    else
143      printf("  No, use ordinary parsimony\n");
144    printf("  N           Use Transversion parsimony?");
145    if (transvp)
146      printf("  Yes, count only transversions\n");
147    else
148      printf("  No, count all steps\n");
149    printf("  W                       Sites weighted?  %s\n",
150           (weights ? "Yes" : "No"));
151    printf("  M           Analyze multiple data sets?");
152    if (mulsets)
153      printf("  Yes, %2ld %s\n", msets,
154               (justwts ? "sets of weights" : "data sets"));
155    else
156      printf("  No\n");
157    printf("  I          Input sequences interleaved?  %s\n",
158           (interleaved ? "Yes" : "No, sequential"));
159    printf("  0   Terminal type (IBM PC, ANSI, none)?  %s\n",
160           ibmpc ? "IBM PC" : ansi  ? "ANSI"  : "(none)");
161    printf("  1    Print out the data at start of run  %s\n",
162           (printdata ? "Yes" : "No"));
163    printf("  2  Print indications of progress of run  %s\n",
164           progress ? "Yes" : "No");
165    printf("  3                        Print out tree  %s\n",
166           treeprint ? "Yes" : "No");
167    printf("  4          Print out steps in each site  %s\n",
168           stepbox ? "Yes" : "No");
169    printf("  5  Print sequences at all nodes of tree  %s\n",
170           ancseq ? "Yes" : "No");
171    if (ancseq || printdata)
172      printf("  .  Use dot-differencing to display them  %s\n",
173           dotdiff ? "Yes" : "No");
174    printf("  6       Write out trees onto tree file?  %s\n",
175           trout ? "Yes" : "No");
176    printf("\n  Y to accept these or type the letter for one to change\n");
177#ifdef WIN32
178    phyFillScreenColor();
179#endif
180    scanf("%c%*[^\n]", &ch);
181    getchar();
182    if (ch == '\n')
183      ch = ' ';
184    uppercase(&ch);
185    if (ch == 'Y')
186      break;
187    if (strchr("WSVJOTNUMI12345.60",ch) != NULL){
188      switch (ch) {
189       
190      case 'J':
191        jumble = !jumble;
192        if (jumble)
193          initjumble(&inseed, &inseed0, seed, &njumble);
194        else njumble = 1;
195        break;
196       
197      case 'O':
198        outgropt = !outgropt;
199        if (outgropt)
200          initoutgroup(&outgrno, spp);
201        break;
202       
203      case 'T':
204        thresh = !thresh;
205        if (thresh)
206          initthreshold(&threshold);
207        break;
208       
209      case 'N':
210        transvp = !transvp;
211        break;
212       
213      case 'W':
214        weights = !weights;
215        break;
216
217      case 'M':
218        mulsets = !mulsets;
219        if (mulsets) {
220          printf("Multiple data sets or multiple weights?");
221          loopcount2 = 0;
222          do {
223            printf(" (type D or W)\n");
224#ifdef WIN32
225            phyFillScreenColor();
226#endif
227            scanf("%c%*[^\n]", &ch2);
228            getchar();
229            if (ch2 == '\n')
230              ch2 = ' ';
231            uppercase(&ch2);
232            countup(&loopcount2, 10);
233          } while ((ch2 != 'W') && (ch2 != 'D'));
234          justwts = (ch2 == 'W');
235          if (justwts)
236            justweights(&msets);
237          else
238            initdatasets(&msets);
239          if (!jumble) {
240            jumble = true;
241            initjumble(&inseed, &inseed0, seed, &njumble);
242          }
243        }
244        break;
245       
246      case 'U':
247        usertree = !usertree;
248        break;
249       
250      case 'S':
251        thorough = !thorough;
252        if (!thorough)
253          printf("Rearrange on just one best tree?");
254          loopcount2 = 0;
255          do {
256            printf(" (type Y or N)\n");
257#ifdef WIN32
258            phyFillScreenColor();
259#endif
260            scanf("%c%*[^\n]", &ch2);
261            getchar();
262            if (ch2 == '\n')
263              ch2 = ' ';
264            uppercase(&ch2);
265            countup(&loopcount2, 10);
266          } while ((ch2 != 'Y') && (ch2 != 'N'));
267          rearrfirst = (ch2 == 'Y');
268        break;
269
270      case 'V':
271        loopcount2 = 0;
272        do {
273          printf("type the number of trees to save\n");
274#ifdef WIN32
275          phyFillScreenColor();
276#endif
277          scanf("%ld%*[^\n]", &maxtrees);
278          if (maxtrees > MAXNUMTREES)
279            maxtrees = MAXNUMTREES;       
280          getchar();
281          countup(&loopcount2, 10);
282        } while (maxtrees < 1);
283        break;
284
285      case 'I':
286        interleaved = !interleaved;
287        break;
288       
289      case '0':
290        initterminal(&ibmpc, &ansi);
291        break;
292       
293      case '1':
294        printdata = !printdata;
295        break;
296       
297      case '2':
298        progress = !progress;
299        break;
300       
301      case '3':
302        treeprint = !treeprint;
303        break;
304       
305      case '4':
306        stepbox = !stepbox;
307        break;
308       
309      case '5':
310        ancseq = !ancseq;
311        break;
312       
313      case '.':
314        dotdiff = !dotdiff;
315        break;
316       
317      case '6':
318        trout = !trout;
319        break;
320      }
321    } else
322      printf("Not a possible option!\n");
323    countup(&loopcount, 100);
324  }
325  if (transvp)
326    fprintf(outfile, "Transversion parsimony\n\n");
327}  /* getoptions */
328
329
330void allocrest()
331{
332  long i;
333
334  y = (Char **)Malloc(spp*sizeof(Char *));
335  for (i = 0; i < spp; i++)
336    y[i] = (Char *)Malloc(chars*sizeof(Char));
337  bestrees = (bestelm *)Malloc(maxtrees*sizeof(bestelm));
338  for (i = 1; i <= maxtrees; i++)
339    bestrees[i - 1].btree = (long *)Malloc(nonodes*sizeof(long));
340  nayme = (naym *)Malloc(spp*sizeof(naym));
341  enterorder = (long *)Malloc(spp*sizeof(long));
342  place = (long *)Malloc(nonodes*sizeof(long));
343  weight = (long *)Malloc(chars*sizeof(long));
344  oldweight = (long *)Malloc(chars*sizeof(long));
345  alias = (long *)Malloc(chars*sizeof(long));
346  ally = (long *)Malloc(chars*sizeof(long));
347  location = (long *)Malloc(chars*sizeof(long));
348}  /* allocrest */
349
350
351void doinit()
352{
353  /* initializes variables */
354
355  inputnumbers(&spp, &chars, &nonodes, 1);
356  getoptions();
357  if (printdata)
358    fprintf(outfile, "%2ld species, %3ld  sites\n\n", spp, chars);
359  alloctree(&treenode, nonodes, usertree);
360}  /* doinit */
361
362
363void makeweights()
364{
365  /* make up weights vector to avoid duplicate computations */
366  long i;
367
368  for (i = 1; i <= chars; i++) {
369    alias[i - 1] = i;
370    oldweight[i - 1] = weight[i - 1];
371    ally[i - 1] = i;
372  }
373  sitesort(chars, weight);
374  sitecombine(chars);
375  sitescrunch(chars);
376  endsite = 0;
377  for (i = 1; i <= chars; i++) {
378    if (ally[i - 1] == i)
379      endsite++;
380  }
381  for (i = 1; i <= endsite; i++)
382    location[alias[i - 1] - 1] = i;
383  if (!thresh)
384    threshold = spp;
385  threshwt = (long *)Malloc(endsite*sizeof(long));
386  for (i = 0; i < endsite; i++) {
387    weight[i] *= 10;
388    threshwt[i] = (long)(threshold * weight[i] + 0.5);
389  }
390  zeros = (long *)Malloc(endsite*sizeof(long));
391  for (i = 0; i < endsite; i++)
392    zeros[i] = 0;
393}  /* makeweights */
394
395
396void doinput()
397{
398  /* reads the input data */
399  long i;
400
401  if (justwts) {
402    if (firstset)
403      inputdata(chars);
404    for (i = 0; i < chars; i++)
405      weight[i] = 1;
406    inputweights(chars, weight, &weights);
407    if (justwts) {
408      fprintf(outfile, "\n\nWeights set # %ld:\n\n", ith);
409      if (progress)
410        printf("\nWeights set # %ld:\n\n", ith);
411    }
412    if (printdata)
413      printweights(outfile, 0, chars, weight, "Sites");
414  } else {
415    if (!firstset){
416      samenumsp(&chars, ith);
417      reallocchars();
418    }
419    inputdata(chars);
420    for (i = 0; i < chars; i++)
421      weight[i] = 1;
422    if (weights) {
423      inputweights(chars, weight, &weights);
424      if (printdata)
425        printweights(outfile, 0, chars, weight, "Sites");
426    }
427  }
428  makeweights();
429  makevalues(treenode, zeros, usertree);
430  if (!usertree) {
431    allocnode(&temp, zeros, endsite);
432    allocnode(&temp1, zeros, endsite);
433    allocnode(&temp2, zeros, endsite);
434    allocnode(&tempsum, zeros, endsite);
435    allocnode(&temprm, zeros, endsite);
436    allocnode(&tempadd, zeros, endsite);
437    allocnode(&tempf, zeros, endsite);
438    allocnode(&tmp, zeros, endsite);
439    allocnode(&tmp1, zeros, endsite);
440    allocnode(&tmp2, zeros, endsite);
441    allocnode(&tmp3, zeros, endsite);
442    allocnode(&tmprm, zeros, endsite);
443    allocnode(&tmpadd, zeros, endsite);
444  }
445}  /* doinput */
446
447
448void initdnaparsnode(node **p, node **grbg, node *q, long len, long nodei,
449                        long *ntips, long *parens, initops whichinit,
450                        pointarray treenode, pointarray nodep, Char *str,
451                        Char *ch, FILE *intree)
452{
453  /* initializes a node */
454  boolean minusread;
455  double valyew, divisor;
456
457  switch (whichinit) {
458  case bottom:
459    gnutreenode(grbg, p, nodei, endsite, zeros);
460    treenode[nodei - 1] = *p;
461    break;
462  case nonbottom:
463    gnutreenode(grbg, p, nodei, endsite, zeros);
464    break;
465  case tip:
466    match_names_to_data (str, treenode, p, spp);
467    break;
468  case length:         /* if there is a length, read it and discard value */
469    processlength(&valyew, &divisor, ch, &minusread, intree, parens);
470    break;
471  default:            /*cases hslength,hsnolength,treewt,unittrwt,iter,*/
472    break;
473  }
474} /* initdnaparsnode */
475
476
477void evaluate(node *r)
478{
479  /* determines the number of steps needed for a tree. this is
480     the minimum number of steps needed to evolve sequences on
481     this tree */
482  long i, steps;
483  long term;
484  double sum;
485
486  sum = 0.0;
487  for (i = 0; i < endsite; i++) {
488    steps = r->numsteps[i];
489    if ((long)steps <= threshwt[i])
490      term = steps;
491    else
492      term = threshwt[i];
493    sum += (double)term;
494    if (usertree && which <= maxuser)
495      fsteps[which - 1][i] = term;
496  }
497  if (usertree && which <= maxuser) {
498    nsteps[which - 1] = sum;
499    if (which == 1) {
500      minwhich = 1;
501      minsteps = sum;
502    } else if (sum < minsteps) {
503      minwhich = which;
504      minsteps = sum;
505    }
506  }
507  like = -sum;
508}  /* evaluate */
509
510
511void tryadd(node *p, node *item, node *nufork)
512{
513  /* temporarily adds one fork and one tip to the tree.
514     if the location where they are added yields greater
515     "likelihood" than other locations tested up to that
516     time, then keeps that location as there */
517  long pos;
518  double belowsum, parentsum;
519  boolean found, collapse, changethere, trysave;
520
521  if (!p->tip) {
522    memcpy(temp->base, p->base, endsite*sizeof(long));
523    memcpy(temp->numsteps, p->numsteps, endsite*sizeof(long));
524    memcpy(temp->numnuc, p->numnuc, endsite*sizeof(nucarray));
525    temp->numdesc = p->numdesc + 1;
526    if (p->back) {
527      multifillin(temp, tempadd, 1);
528      sumnsteps2(tempsum, temp, p->back, 0, endsite, threshwt);
529    } else {
530      multisumnsteps(temp, tempadd, 0, endsite, threshwt);
531      tempsum->sumsteps = temp->sumsteps;
532    }
533    if (tempsum->sumsteps <= -bestyet) {
534      if (p->back)
535        sumnsteps2(tempsum, temp, p->back, endsite+1, endsite, threshwt);
536      else {
537        multisumnsteps(temp, temp1, endsite+1, endsite, threshwt);
538        tempsum->sumsteps = temp->sumsteps;
539      }
540    }
541    p->sumsteps = tempsum->sumsteps;
542  }
543  if (p == root)
544    sumnsteps2(temp, item, p, 0, endsite, threshwt);
545  else {
546    sumnsteps(temp1, item, p, 0, endsite);
547    sumnsteps2(temp, temp1, p->back, 0, endsite, threshwt);
548  }
549  if (temp->sumsteps <= -bestyet) {
550    if (p == root)
551      sumnsteps2(temp, item, p, endsite+1, endsite, threshwt);
552    else {
553      sumnsteps(temp1, item, p, endsite+1, endsite);
554      sumnsteps2(temp, temp1, p->back, endsite+1, endsite, threshwt);
555    }
556  }
557  belowsum = temp->sumsteps;
558  multf = false;
559  like = -belowsum;
560  if (!p->tip && belowsum >= p->sumsteps) {
561    multf = true;
562    like = -p->sumsteps;
563  }
564  trysave = true;
565  if (!multf && p != root) {
566    parentsum = treenode[p->back->index - 1]->sumsteps;
567    if (belowsum >= parentsum)
568      trysave = false;
569  }
570  if (lastrearr) {
571    changethere = true;
572    if (like >= bstlike2 && trysave) {
573      if (like > bstlike2)
574        found = false;
575      else {
576        addnsave(p, item, nufork, &root, &grbg, multf, treenode, place, zeros);
577        pos = 0;
578        findtree(&found, &pos, nextree, place, bestrees);
579      }
580      if (!found) {
581        collapse = collapsible(item, p, temp, temp1, temp2, tempsum, temprm,
582                     tmpadd, multf, root, zeros, treenode);
583        if (!thorough)
584          changethere = !collapse;
585        if (thorough || !collapse || like > bstlike2) {
586          if (like > bstlike2) {
587            addnsave(p, item, nufork, &root, &grbg, multf, treenode, 
588                       place, zeros);
589            bestlike = bstlike2 = like;
590            addbestever(&pos, &nextree, maxtrees, collapse, place, bestrees);
591          } else
592            addtiedtree(pos, &nextree, maxtrees, collapse, place, bestrees);
593        }
594      }
595    }
596    if (like >= bestyet) {
597      if (like > bstlike2)
598        bstlike2 = like;
599      if (changethere && trysave) {
600        bestyet = like;
601        there = p;
602        mulf = multf;
603      }
604    }
605  } else if ((like > bestyet) || (like >= bestyet && trysave)) {
606    bestyet = like;
607    there = p;
608    mulf = multf;
609  }
610}  /* tryadd */
611
612
613void addpreorder(node *p, node *item, node *nufork)
614{
615  /* traverses a n-ary tree, calling function tryadd
616     at a node before calling tryadd at its descendants */
617  node *q;
618
619  if (p == NULL)
620    return;
621  tryadd(p, item, nufork);
622  if (!p->tip) {
623    q = p->next;
624    while (q != p) {
625      addpreorder(q->back, item, nufork);
626      q = q->next;
627    }
628  }
629}  /* addpreorder */
630
631
632void trydescendants(node *item, node *forknode, node *parent,
633                        node *parentback, boolean trybelow)
634{
635  /* tries rearrangements at parent and below parent's descendants */
636  node *q, *tempblw;
637  boolean bestever=0, belowbetter, multf=0, saved, trysave;
638  double parentsum=0, belowsum;
639
640  memcpy(temp->base, parent->base, endsite*sizeof(long));
641  memcpy(temp->numsteps, parent->numsteps, endsite*sizeof(long));
642  memcpy(temp->numnuc, parent->numnuc, endsite*sizeof(nucarray));
643  temp->numdesc = parent->numdesc + 1;
644  multifillin(temp, tempadd, 1);
645  sumnsteps2(tempsum, parentback, temp, 0, endsite, threshwt);
646  belowbetter = true;
647  if (lastrearr) {
648    parentsum = tempsum->sumsteps;
649    if (-tempsum->sumsteps >= bstlike2) {
650      belowbetter = false;
651      bestever = false;
652      multf = true;
653      if (-tempsum->sumsteps > bstlike2)
654        bestever = true;
655      savelocrearr(item, forknode, parent, tmp, tmp1, tmp2, tmp3, tmprm,
656                    tmpadd, &root, maxtrees, &nextree, multf, bestever,
657                    &saved, place, bestrees, treenode, &grbg, zeros);
658      if (saved) {
659        like = bstlike2 = -tempsum->sumsteps;
660        there = parent;
661        mulf = true;
662      }
663    }
664  } else if (-tempsum->sumsteps >= like) {
665    there = parent;
666    mulf = true;
667    like = -tempsum->sumsteps;
668  }
669  if (trybelow) {
670    sumnsteps(temp, parent, tempadd, 0, endsite);
671    sumnsteps2(tempsum, temp, parentback, 0, endsite, threshwt);
672    if (lastrearr) {
673      belowsum = tempsum->sumsteps;
674      if (-tempsum->sumsteps >= bstlike2 && belowbetter && 
675            (forknode->numdesc > 2 ||
676               (forknode->numdesc == 2 && 
677                 parent->back->index != forknode->index))) {
678        trysave = false;
679        memcpy(temp->base, parentback->base, endsite*sizeof(long));
680        memcpy(temp->numsteps, parentback->numsteps, endsite*sizeof(long));
681        memcpy(temp->numnuc, parentback->numnuc, endsite*sizeof(nucarray));
682        temp->numdesc = parentback->numdesc + 1;
683        multifillin(temp, tempadd, 1);
684        sumnsteps2(tempsum, parent, temp, 0, endsite, threshwt);
685        if (-tempsum->sumsteps < bstlike2) {
686          multf = false;
687          bestever = false;
688          trysave = true;
689        }
690        if (-belowsum > bstlike2) {
691          multf = false;
692          bestever = true;
693          trysave = true;
694        }
695        if (trysave) {
696          if (treenode[parent->index - 1] != parent)
697            tempblw = parent->back;
698          else
699            tempblw = parent;
700          savelocrearr(item, forknode, tempblw, tmp, tmp1, tmp2, tmp3, tmprm,
701                         tmpadd, &root, maxtrees, &nextree, multf, bestever,
702                         &saved, place, bestrees, treenode, &grbg, zeros);
703          if (saved) {
704            like = bstlike2 = -belowsum;
705            there = tempblw;
706            mulf = false;
707          }
708        }
709      }
710    } else if (-tempsum->sumsteps > like) {
711      like = -tempsum->sumsteps;
712      if (-tempsum->sumsteps > bestyet) {
713        if (treenode[parent->index - 1] != parent)
714          tempblw = parent->back;
715        else
716          tempblw = parent;
717        there = tempblw;
718        mulf = false;
719      }
720    }
721  }
722  q = parent->next;
723  while (q != parent) {
724    if (q->back && q->back != item) {
725      memcpy(temp1->base, q->base, endsite*sizeof(long));
726      memcpy(temp1->numsteps, q->numsteps, endsite*sizeof(long));
727      memcpy(temp1->numnuc, q->numnuc, endsite*sizeof(nucarray));
728      temp1->numdesc = q->numdesc;
729      multifillin(temp1, parentback, 0);
730      if (lastrearr)
731        belowbetter = (-parentsum < bstlike2);
732      if (!q->back->tip) {
733        memcpy(temp->base, q->back->base, endsite*sizeof(long));
734        memcpy(temp->numsteps, q->back->numsteps, endsite*sizeof(long));
735        memcpy(temp->numnuc, q->back->numnuc, endsite*sizeof(nucarray));
736        temp->numdesc = q->back->numdesc + 1;
737        multifillin(temp, tempadd, 1);
738        sumnsteps2(tempsum, temp1, temp, 0, endsite, threshwt);
739        if (lastrearr) {
740          if (-tempsum->sumsteps >= bstlike2) {
741            belowbetter = false;
742            bestever = false;
743            multf = true;
744            if (-tempsum->sumsteps > bstlike2)
745              bestever = true;
746            savelocrearr(item, forknode, q->back, tmp, tmp1, tmp2, tmp3, tmprm,
747                          tmpadd, &root, maxtrees, &nextree, multf, bestever,
748                          &saved, place, bestrees, treenode, &grbg, zeros);
749            if (saved) {
750              like = bstlike2 = -tempsum->sumsteps;
751              there = q->back;
752              mulf = true;
753            }
754          }
755        } else if (-tempsum->sumsteps >= like) {
756          like = -tempsum->sumsteps;
757          there = q->back;
758          mulf = true;
759        }
760      }
761      sumnsteps(temp, q->back, tempadd, 0, endsite);
762      sumnsteps2(tempsum, temp, temp1, 0, endsite, threshwt);
763      if (lastrearr) {
764        if (-tempsum->sumsteps >= bstlike2) {
765          trysave = false;
766          multf = false;
767          if (belowbetter) {
768            bestever = false;
769            trysave = true;
770          }
771          if (-tempsum->sumsteps > bstlike2) {
772            bestever = true;
773            trysave = true;
774          }
775          if (trysave) {
776            if (treenode[q->back->index - 1] != q->back)
777              tempblw = q;
778            else
779              tempblw = q->back;
780            savelocrearr(item, forknode, tempblw, tmp, tmp1, tmp2, tmp3, tmprm,
781                        tmpadd, &root, maxtrees, &nextree, multf, bestever,
782                        &saved, place, bestrees, treenode, &grbg, zeros);
783            if (saved) {
784              like = bstlike2 = -tempsum->sumsteps;
785              there = tempblw;
786              mulf = false;
787            }
788          }
789        }
790      } else if (-tempsum->sumsteps > like) {
791        like = -tempsum->sumsteps;
792        if (-tempsum->sumsteps > bestyet) {
793          if (treenode[q->back->index - 1] != q->back)
794            tempblw = q;
795          else
796            tempblw = q->back;
797          there = tempblw;
798          mulf = false;
799        }
800      }
801    }
802    q = q->next;
803  }
804} /* trydescendants */
805
806
807void trylocal(node *item, node *forknode)
808{
809  /* rearranges below forknode, below descendants of forknode when
810     there are more than 2 descendants, then unroots the back of
811     forknode and rearranges on its descendants */
812  node *q;
813  boolean bestever, multf, saved;
814
815  memcpy(temprm->base, zeros, endsite*sizeof(long));
816  memcpy(temprm->numsteps, zeros, endsite*sizeof(long));
817  memcpy(temprm->oldbase, item->base, endsite*sizeof(long));
818  memcpy(temprm->oldnumsteps, item->numsteps, endsite*sizeof(long));
819  memcpy(tempf->base, forknode->base, endsite*sizeof(long));
820  memcpy(tempf->numsteps, forknode->numsteps, endsite*sizeof(long));
821  memcpy(tempf->numnuc, forknode->numnuc, endsite*sizeof(nucarray));
822  tempf->numdesc = forknode->numdesc - 1;
823  multifillin(tempf, temprm, -1);
824  if (!forknode->back) {
825    sumnsteps2(tempsum, tempf, tempadd, 0, endsite, threshwt);
826    if (lastrearr) {
827      if (-tempsum->sumsteps > bstlike2) {
828        bestever = true;
829        multf = false;
830        savelocrearr(item, forknode, forknode, tmp, tmp1, tmp2, tmp3, tmprm,
831                       tmpadd, &root, maxtrees, &nextree, multf, bestever,
832                       &saved, place, bestrees, treenode, &grbg, zeros);
833        if (saved) {
834          like = bstlike2 = -tempsum->sumsteps;
835          there = forknode;
836          mulf = false;
837        }
838      }
839    } else if (-tempsum->sumsteps > like) {
840      like = -tempsum->sumsteps;
841      if (-tempsum->sumsteps > bestyet) {
842        there = forknode;
843        mulf = false;
844      }
845    }
846  } else {
847    sumnsteps(temp, tempf, tempadd, 0, endsite);
848    sumnsteps2(tempsum, temp, forknode->back, 0, endsite, threshwt);
849    if (lastrearr) {
850      if (-tempsum->sumsteps > bstlike2) {
851        bestever = true;
852        multf = false;
853        savelocrearr(item, forknode, forknode, tmp, tmp1, tmp2, tmp3, tmprm,
854                       tmpadd, &root, maxtrees, &nextree, multf, bestever,
855                       &saved, place, bestrees, treenode, &grbg, zeros);
856        if (saved) {
857          like = bstlike2 = -tempsum->sumsteps;
858          there = forknode;
859          mulf = false;
860        }
861      }
862    } else if (-tempsum->sumsteps > like) {
863      like = -tempsum->sumsteps;
864      if (-tempsum->sumsteps > bestyet) {
865        there = forknode;
866        mulf = false;
867      }
868    }
869    trydescendants(item, forknode, forknode->back, tempf, false);
870  }
871  q = forknode->next;
872  while (q != forknode) {
873    if (q->back != item) {
874      memcpy(temp2->base, q->base, endsite*sizeof(long));
875      memcpy(temp2->numsteps, q->numsteps, endsite*sizeof(long));
876      memcpy(temp2->numnuc, q->numnuc, endsite*sizeof(nucarray));
877      temp2->numdesc = q->numdesc - 1;
878      multifillin(temp2, temprm, -1);
879      if (!q->back->tip) {
880        trydescendants(item, forknode, q->back, temp2, true);
881      } else {
882        sumnsteps(temp1, q->back, tempadd, 0, endsite);
883        sumnsteps2(tempsum, temp1, temp2, 0, endsite, threshwt);
884        if (lastrearr) {
885          if (-tempsum->sumsteps > bstlike2) {
886            multf = false;
887            bestever = true;
888            savelocrearr(item, forknode, q->back, tmp, tmp1, tmp2, tmp3,
889                         tmprm, tmpadd, &root, maxtrees, &nextree, multf,
890                         bestever, &saved, place, bestrees, treenode,
891                         &grbg, zeros);
892            if (saved) {
893              like = bstlike2 = -tempsum->sumsteps;
894              there = q->back;
895              mulf = false;
896            }
897          }
898        } else if ((-tempsum->sumsteps) > like) {
899          like = -tempsum->sumsteps;
900          if (-tempsum->sumsteps > bestyet) {
901            there = q->back;
902            mulf = false;
903          }
904        }
905      }
906    }
907    q = q->next;
908  }
909} /* trylocal */
910
911
912void trylocal2(node *item, node *forknode, node *other)
913{
914  /* rearranges below forknode, below descendants of forknode when
915     there are more than 2 descendants, then unroots the back of
916     forknode and rearranges on its descendants.  Used if forknode
917     has binary descendants */
918  node *q;
919  boolean bestever=0, multf, saved, belowbetter, trysave;
920
921  memcpy(tempf->base, other->base, endsite*sizeof(long));
922  memcpy(tempf->numsteps, other->numsteps, endsite*sizeof(long));
923  memcpy(tempf->oldbase, forknode->base, endsite*sizeof(long));
924  memcpy(tempf->oldnumsteps, forknode->numsteps, endsite*sizeof(long));
925  tempf->numdesc = other->numdesc;
926  if (forknode->back)
927    trydescendants(item, forknode, forknode->back, tempf, false);
928  if (!other->tip) {
929    memcpy(temp->base, other->base, endsite*sizeof(long));
930    memcpy(temp->numsteps, other->numsteps, endsite*sizeof(long));
931    memcpy(temp->numnuc, other->numnuc, endsite*sizeof(nucarray));
932    temp->numdesc = other->numdesc + 1;
933    multifillin(temp, tempadd, 1);
934    if (forknode->back)
935      sumnsteps2(tempsum, forknode->back, temp, 0, endsite, threshwt);
936    else
937      sumnsteps2(tempsum, NULL, temp, 0, endsite, threshwt);
938    belowbetter = true;
939    if (lastrearr) {
940      if (-tempsum->sumsteps >= bstlike2) {
941        belowbetter = false;
942        bestever = false;
943        multf = true;
944        if (-tempsum->sumsteps > bstlike2)
945          bestever = true;
946        savelocrearr(item, forknode, other, tmp, tmp1, tmp2, tmp3, tmprm,
947                       tmpadd, &root, maxtrees, &nextree, multf, bestever,
948                       &saved, place, bestrees, treenode, &grbg, zeros);
949        if (saved) {
950          like = bstlike2 = -tempsum->sumsteps;
951          there = other;
952          mulf = true;
953        }
954      }
955    } else if (-tempsum->sumsteps >= like) {
956      there = other;
957      mulf = true;
958      like = -tempsum->sumsteps;
959    }
960    if (forknode->back) {
961      memcpy(temprm->base, forknode->back->base, endsite*sizeof(long));
962      memcpy(temprm->numsteps, forknode->back->numsteps, endsite*sizeof(long));
963    } else {
964      memcpy(temprm->base, zeros, endsite*sizeof(long));
965      memcpy(temprm->numsteps, zeros, endsite*sizeof(long));
966    }
967    memcpy(temprm->oldbase, other->back->base, endsite*sizeof(long));
968    memcpy(temprm->oldnumsteps, other->back->numsteps, endsite*sizeof(long));
969    q = other->next;
970    while (q != other) {
971      memcpy(temp2->base, q->base, endsite*sizeof(long));
972      memcpy(temp2->numsteps, q->numsteps, endsite*sizeof(long));
973      memcpy(temp2->numnuc, q->numnuc, endsite*sizeof(nucarray));
974      if (forknode->back) {
975        temp2->numdesc = q->numdesc;
976        multifillin(temp2, temprm, 0);
977      } else {
978        temp2->numdesc = q->numdesc - 1;
979        multifillin(temp2, temprm, -1);
980      }
981      if (!q->back->tip)
982        trydescendants(item, forknode, q->back, temp2, true);
983      else {
984        sumnsteps(temp1, q->back, tempadd, 0, endsite);
985        sumnsteps2(tempsum, temp1, temp2, 0, endsite, threshwt);
986        if (lastrearr) {
987          if (-tempsum->sumsteps >= bstlike2) {
988            trysave = false;
989            multf = false;
990            if (belowbetter) {
991              bestever = false;
992              trysave = true;
993            }
994            if (-tempsum->sumsteps > bstlike2) {
995              bestever = true;
996              trysave = true;
997            }
998            if (trysave) {
999              savelocrearr(item, forknode, q->back, tmp, tmp1, tmp2, tmp3,
1000                           tmprm, tmpadd, &root, maxtrees, &nextree, multf,
1001                           bestever, &saved, place, bestrees, treenode,
1002                           &grbg, zeros);
1003              if (saved) {
1004                like = bstlike2 = -tempsum->sumsteps;
1005                there = q->back;
1006                mulf = false;
1007              }
1008            }
1009          }
1010        } else if (-tempsum->sumsteps > like) {
1011          like = -tempsum->sumsteps;
1012          if (-tempsum->sumsteps > bestyet) {
1013            there = q->back;
1014            mulf = false;
1015          }
1016        }
1017      }
1018      q = q->next;
1019    }
1020  }
1021} /* trylocal2 */
1022
1023
1024void tryrearr(node *p, boolean *success)
1025{
1026  /* evaluates one rearrangement of the tree.
1027     if the new tree has greater "likelihood" than the old
1028     one sets success = TRUE and keeps the new tree.
1029     otherwise, restores the old tree */
1030  node *forknode, *newfork, *other, *oldthere;
1031  double oldlike;
1032  boolean oldmulf;
1033
1034  if (p->back == NULL)
1035    return;
1036  forknode = treenode[p->back->index - 1]; 
1037  if (!forknode->back && forknode->numdesc <= 2 && alltips(forknode, p))
1038    return;
1039  oldlike = bestyet;
1040  like = -10.0 * spp * chars;
1041  memcpy(tempadd->base, p->base, endsite*sizeof(long));
1042  memcpy(tempadd->numsteps, p->numsteps, endsite*sizeof(long));
1043  memcpy(tempadd->oldbase, zeros, endsite*sizeof(long));
1044  memcpy(tempadd->oldnumsteps, zeros, endsite*sizeof(long));
1045  if (forknode->numdesc > 2) {
1046    oldthere = there = forknode;
1047    oldmulf = mulf = true;
1048    trylocal(p, forknode);
1049  } else {
1050    findbelow(&other, p, forknode);
1051    oldthere = there = other;
1052    oldmulf = mulf = false;
1053    trylocal2(p, forknode, other);
1054  }
1055  if ((like <= oldlike) || (there == oldthere && mulf == oldmulf))
1056    return;
1057  recompute = true;
1058  re_move(p, &forknode, &root, recompute, treenode, &grbg, zeros);
1059  if (mulf)
1060    add(there, p, NULL, &root, recompute, treenode, &grbg, zeros);
1061  else {
1062    if (forknode->numdesc > 0)
1063      getnufork(&newfork, &grbg, treenode, zeros);
1064    else
1065      newfork = forknode;
1066    add(there, p, newfork, &root, recompute, treenode, &grbg, zeros);
1067  } 
1068  if (like > oldlike) {
1069    *success = true;
1070    bestyet = like;
1071  }
1072}  /* tryrearr */
1073
1074
1075void repreorder(node *p, boolean *success)
1076{
1077  /* traverses a binary tree, calling PROCEDURE tryrearr
1078     at a node before calling tryrearr at its descendants */
1079  node *q, *this;
1080
1081  if (p == NULL)
1082    return;
1083  if (!p->visited) {
1084    tryrearr(p, success);
1085    p->visited = true;
1086  }
1087  if (!p->tip) {
1088    q = p;
1089    while (q->next != p) {
1090      this = q->next->back;
1091      repreorder(q->next->back,success);
1092      if (q->next->back == this)
1093        q = q->next;
1094    }
1095  }
1096}  /* repreorder */
1097
1098
1099void rearrange(node **r)
1100{
1101  /* traverses the tree (preorder), finding any local
1102     rearrangement which decreases the number of steps.
1103     if traversal succeeds in increasing the tree's
1104     "likelihood", PROCEDURE rearrange runs traversal again */
1105  boolean success=true;
1106
1107  while (success) {
1108    success = false;
1109    clearvisited(treenode);
1110    repreorder(*r, &success);
1111  }
1112}  /* rearrange */
1113
1114
1115void describe()
1116{
1117  /* prints ancestors, steps and table of numbers of steps in
1118     each site */
1119
1120  if (treeprint) {
1121    fprintf(outfile, "\nrequires a total of %10.3f\n", like / -10.0);
1122    fprintf(outfile, "\n  between      and       length\n");
1123    fprintf(outfile, "  -------      ---       ------\n");
1124    printbranchlengths(root);
1125  }
1126  if (stepbox)
1127    writesteps(chars, weights, oldweight, root);
1128  if (ancseq) {
1129    hypstates(chars, root, treenode, &garbage, basechar);
1130    putc('\n', outfile);
1131  }
1132  putc('\n', outfile);
1133  if (trout) {
1134    col = 0;
1135    treeout3(root, nextree, &col, root);
1136  }
1137}  /* describe */
1138
1139
1140void dnapars_coordinates(node *p, double lengthsum, long *tipy,
1141                        double *tipmax)
1142{
1143  /* establishes coordinates of nodes */
1144  node *q, *first, *last;
1145  double xx;
1146
1147  if (p == NULL)
1148    return;
1149  if (p->tip) {
1150    p->xcoord = (long)(over * lengthsum + 0.5);
1151    p->ycoord = (*tipy);
1152    p->ymin = (*tipy);
1153    p->ymax = (*tipy);
1154    (*tipy) += down;
1155    if (lengthsum > (*tipmax))
1156      (*tipmax) = lengthsum;
1157    return;
1158  }
1159  q = p->next;
1160  do {
1161    xx = q->v;
1162    if (xx > 100.0)
1163      xx = 100.0;
1164    dnapars_coordinates(q->back, lengthsum + xx, tipy,tipmax);
1165    q = q->next;
1166  } while (p != q);
1167  first = p->next->back;
1168  q = p;
1169  while (q->next != p)
1170    q = q->next;
1171  last = q->back;
1172  p->xcoord = (long)(over * lengthsum + 0.5);
1173  if ((p == root) || count_sibs(p) > 2)
1174    p->ycoord = p->next->next->back->ycoord;
1175  else
1176    p->ycoord = (first->ycoord + last->ycoord) / 2;
1177  p->ymin = first->ymin;
1178  p->ymax = last->ymax;
1179}  /* dnapars_coordinates */
1180
1181
1182void dnapars_printree()
1183{
1184  /* prints out diagram of the tree2 */
1185  long tipy;
1186  double scale, tipmax;
1187  long i;
1188 
1189  if (!treeprint)
1190    return;
1191  putc('\n', outfile);
1192  tipy = 1;
1193  tipmax = 0.0;
1194  dnapars_coordinates(root, 0.0, &tipy, &tipmax);
1195  scale = 1.0 / (long)(tipmax + 1.000);
1196  for (i = 1; i <= (tipy - down); i++)
1197    drawline3(i, scale, root);
1198  putc('\n', outfile);
1199}  /* dnapars_printree */
1200
1201
1202void globrearrange()
1203{
1204  /* does global rearrangements */
1205  long j;
1206  double gotlike;
1207  boolean frommulti;
1208  node *item, *nufork;
1209
1210  recompute = true;
1211  do {
1212    printf("   ");
1213    gotlike = bstlike2 = bestlike;
1214    for (j = 0; j < nonodes; j++) {
1215      bestyet = -10.0 * spp * chars;
1216      if (j < spp)
1217        item = treenode[enterorder[j] -1];
1218      else 
1219        item = treenode[j];
1220
1221      if ((item != root) && 
1222           ((j < spp) || ((j >= spp) && (item->numdesc > 0))) &&
1223           !((item->back->index == root->index) && (root->numdesc == 2)
1224              && alltips(root, item))) {
1225        re_move(item, &nufork, &root, recompute, treenode, &grbg, zeros);
1226        frommulti = (nufork->numdesc > 0);
1227        clearcollapse(treenode);
1228        there = root;
1229        memcpy(tempadd->base, item->base, endsite*sizeof(long));
1230        memcpy(tempadd->numsteps, item->numsteps, endsite*sizeof(long));
1231        memcpy(tempadd->oldbase, zeros, endsite*sizeof(long));
1232        memcpy(tempadd->oldnumsteps, zeros, endsite*sizeof(long));
1233        if (frommulti){
1234          oldnufork = nufork;
1235          getnufork(&nufork, &grbg, treenode, zeros);
1236        }
1237        addpreorder(root, item, nufork);
1238        if (frommulti)
1239          oldnufork = NULL;
1240        if (!mulf)
1241          add(there, item, nufork, &root, recompute, treenode, &grbg, zeros);
1242        else
1243          add(there, item, NULL, &root, recompute, treenode, &grbg, zeros);
1244      }
1245      if (progress) {
1246        if (j % ((nonodes / 72) + 1) == 0)
1247          putchar('.');
1248        fflush(stdout);
1249      }
1250    }
1251    if (progress) {
1252      putchar('\n');
1253#ifdef WIN32
1254      phyFillScreenColor();
1255#endif
1256    }
1257  } while (bestlike > gotlike);
1258} /* globrearrange */
1259
1260
1261void load_tree(long treei)
1262{ /* restores a tree from bestrees */
1263  long j, nextnode;
1264  boolean recompute = false;
1265  node *dummy;
1266
1267  for (j = spp - 1; j >= 1; j--)
1268    re_move(treenode[j], &dummy, &root, recompute, treenode, &grbg,
1269              zeros);
1270  root = treenode[0];
1271  recompute = true;
1272  add(treenode[0], treenode[1], treenode[spp], &root, recompute,
1273    treenode, &grbg, zeros);
1274  nextnode = spp + 2;
1275  for (j = 3; j <= spp; j++) {
1276    if (bestrees[treei].btree[j - 1] > 0)
1277      add(treenode[bestrees[treei].btree[j - 1] - 1], treenode[j - 1],
1278            treenode[nextnode++ - 1], &root, recompute, treenode, &grbg,
1279            zeros);
1280    else
1281      add(treenode[treenode[-bestrees[treei].btree[j-1]-1]->back->index-1],
1282            treenode[j - 1], NULL, &root, recompute, treenode, &grbg,
1283            zeros);
1284  }
1285}
1286
1287
1288void grandrearr()
1289{
1290  /* calls global rearrangement on best trees */
1291  long treei; 
1292  boolean done;
1293
1294  done = false;
1295  do {
1296    treei = findunrearranged(bestrees, nextree, true);
1297    if (treei < 0)
1298      done = true;
1299    else 
1300      bestrees[treei].gloreange = true;
1301   
1302    if (!done) {
1303      load_tree(treei);
1304      globrearrange();
1305      done = rearrfirst;
1306    }
1307  } while (!done);
1308} /* grandrearr */
1309
1310
1311void maketree()
1312{
1313  /* constructs a binary tree from the pointers in treenode.
1314     adds each node at location which yields highest "likelihood"
1315  then rearranges the tree for greatest "likelihood" */
1316  long i, j, numtrees, nextnode;
1317  boolean done, firsttree, goteof, haslengths;
1318  node *item, *nufork, *dummy;
1319  pointarray nodep;
1320
1321  if (!usertree) {
1322    for (i = 1; i <= spp; i++)
1323      enterorder[i - 1] = i;
1324    if (jumble)
1325      randumize(seed, enterorder);
1326    recompute = true;
1327    root = treenode[enterorder[0] - 1];
1328    add(treenode[enterorder[0] - 1], treenode[enterorder[1] - 1],
1329        treenode[spp], &root, recompute, treenode, &grbg, zeros);
1330    if (progress) {
1331      printf("Adding species:\n");
1332      writename(0, 2, enterorder);
1333#ifdef WIN32
1334      phyFillScreenColor();
1335#endif
1336    }
1337    lastrearr = false;
1338    oldnufork = NULL;
1339    for (i = 3; i <= spp; i++) {
1340      bestyet = -10.0 * spp * chars;
1341      item = treenode[enterorder[i - 1] - 1];
1342      getnufork(&nufork, &grbg, treenode, zeros);
1343      there = root;
1344      memcpy(tempadd->base, item->base, endsite*sizeof(long));
1345      memcpy(tempadd->numsteps, item->numsteps, endsite*sizeof(long));
1346      memcpy(tempadd->oldbase, zeros, endsite*sizeof(long));
1347      memcpy(tempadd->oldnumsteps, zeros, endsite*sizeof(long));
1348      addpreorder(root, item, nufork);
1349      if (!mulf)
1350        add(there, item, nufork, &root, recompute, treenode, &grbg, zeros);
1351      else
1352        add(there, item, NULL, &root, recompute, treenode, &grbg, zeros);
1353      like = bestyet;
1354      rearrange(&root);
1355      if (progress) {
1356        writename(i - 1, 1, enterorder);
1357#ifdef WIN32
1358        phyFillScreenColor();
1359#endif
1360      }
1361      lastrearr = (i == spp);
1362      if (lastrearr) {
1363        bestlike = bestyet;
1364        if (jumb == 1) {
1365          bstlike2 = bestlike;
1366          nextree = 1;
1367          initbestrees(bestrees, maxtrees, true);
1368          initbestrees(bestrees, maxtrees, false);
1369        }
1370        if (progress) {
1371          printf("\nDoing global rearrangements");
1372          if (rearrfirst)
1373            printf(" on the first of the trees tied for best\n");
1374          else
1375            printf(" on all trees tied for best\n");
1376          printf("  !");
1377          for (j = 0; j < nonodes; j++)
1378            if (j % ((nonodes / 72) + 1) == 0)
1379              putchar('-');
1380          printf("!\n");
1381#ifdef WIN32
1382          phyFillScreenColor();
1383#endif
1384        }
1385        globrearrange();
1386      }
1387    }
1388    done = false;
1389    while (!done && findunrearranged(bestrees, nextree, true) >= 0) {
1390      grandrearr();
1391      done = rearrfirst;
1392    }
1393    if (progress)
1394      putchar('\n'); 
1395    recompute = false;
1396    for (i = spp - 1; i >= 1; i--)
1397      re_move(treenode[i], &dummy, &root, recompute, treenode, &grbg, zeros);
1398    if (jumb == njumble) {
1399      if (thorough && (nextree > 2))
1400        reducebestrees(bestrees, &nextree);
1401      if (treeprint) {
1402        putc('\n', outfile);
1403        if (nextree == 2)
1404          fprintf(outfile, "One most parsimonious tree found:\n");
1405        else
1406          fprintf(outfile, "%6ld trees in all found\n", nextree - 1);
1407      }
1408      if (nextree > maxtrees + 1) {
1409        if (treeprint)
1410          fprintf(outfile, "here are the first %4ld of them\n", (long)maxtrees);
1411        nextree = maxtrees + 1;
1412      }
1413      if (treeprint)
1414        putc('\n', outfile);
1415      for (i = 0; i <= (nextree - 2); i++) {
1416        root = treenode[0];
1417        add(treenode[0], treenode[1], treenode[spp], &root, recompute,
1418          treenode, &grbg, zeros);
1419        nextnode = spp + 2;
1420        for (j = 3; j <= spp; j++) {
1421          if (bestrees[i].btree[j - 1] > 0)
1422            add(treenode[bestrees[i].btree[j - 1] - 1], treenode[j - 1],
1423              treenode[nextnode++ - 1], &root, recompute, treenode, &grbg,
1424              zeros);
1425          else
1426            add(treenode[treenode[-bestrees[i].btree[j - 1]-1]->back->index-1],
1427              treenode[j - 1], NULL, &root, recompute, treenode, &grbg, zeros);
1428        }
1429        reroot(treenode[outgrno - 1], root);
1430        postorder(root);
1431        evaluate(root);
1432        treelength(root, chars, treenode);
1433        dnapars_printree();
1434        describe();
1435        for (j = 1; j < spp; j++)
1436          re_move(treenode[j], &dummy, &root, recompute, treenode,
1437                    &grbg, zeros);
1438      }
1439    }
1440  } else {
1441    openfile(&intree,INTREE,"input tree", "r",progname,intreename);
1442    numtrees = countsemic(&intree);
1443    if (numtrees > 2)
1444      initseed(&inseed, &inseed0, seed);
1445    if (numtrees > MAXNUMTREES) {
1446      printf("\nERROR: number of input trees is read incorrectly from %s\n",
1447        intreename);
1448      exxit(-1);
1449    }
1450    if (treeprint) {
1451      fprintf(outfile, "User-defined tree");
1452      if (numtrees > 1)
1453        putc('s', outfile);
1454      fprintf(outfile, ":\n");
1455    }
1456    fsteps = (long **)Malloc(maxuser*sizeof(long *));
1457    for (j = 1; j <= maxuser; j++)
1458      fsteps[j - 1] = (long *)Malloc(endsite*sizeof(long));
1459    if (trout)
1460      fprintf(outtree, "%ld\n", numtrees);
1461    nodep = NULL;
1462    which = 1;
1463    while (which <= numtrees) {
1464      firsttree = true;
1465      nextnode = 0;
1466      haslengths = true;
1467      treeread(intree, &root, treenode, &goteof, &firsttree,
1468                 nodep, &nextnode, &haslengths,
1469                 &grbg, initdnaparsnode);
1470      if (treeprint)
1471        fprintf(outfile, "\n\n");
1472      if (outgropt)
1473        reroot(treenode[outgrno - 1], root);
1474      postorder(root);
1475      evaluate(root);
1476      treelength(root, chars, treenode);
1477      dnapars_printree();
1478      describe();
1479      if (which < numtrees)
1480        gdispose(root, &grbg, treenode);
1481      which++;
1482    }
1483    FClose(intree);
1484    putc('\n', outfile);
1485    if (numtrees > 1 && chars > 1 )
1486      standev(chars, numtrees, minwhich, minsteps, nsteps, fsteps, seed);
1487    for (j = 1; j <= maxuser; j++)
1488      free(fsteps[j - 1]);
1489    free(fsteps);
1490  }
1491  if (jumb == njumble) {
1492    if (progress) {
1493      printf("\nOutput written to file \"%s\"\n\n", outfilename);
1494      if (trout) {
1495        printf("Tree");
1496        if (numtrees > 1)
1497          printf("s");
1498        printf(" also written onto file \"%s\"\n\n", outtreename);
1499      }
1500    }
1501  }
1502}  /* maketree */
1503
1504
1505void reallocchars() 
1506{ /* The amount of chars can change between runs
1507     this function reallocates all the variables
1508     whose size depends on the amount of chars */
1509  long i;
1510
1511  for (i=0; i < spp; i++){
1512    free(y[i]);
1513    y[i] = (Char *)Malloc(chars*sizeof(Char));
1514  }
1515  free(weight);
1516  free(oldweight); 
1517  free(alias);
1518  free(ally);
1519  free(location);
1520
1521  weight = (long *)Malloc(chars*sizeof(long));
1522  oldweight = (long *)Malloc(chars*sizeof(long));
1523  alias = (long *)Malloc(chars*sizeof(long));
1524  ally = (long *)Malloc(chars*sizeof(long));
1525  location = (long *)Malloc(chars*sizeof(long));
1526}
1527
1528
1529void freerest()
1530{ /* free variables that are allocated each data set */
1531  long i;
1532
1533  if (!usertree) {
1534    freenode(&temp);
1535    freenode(&temp1);
1536    freenode(&temp2);
1537    freenode(&tempsum);
1538    freenode(&temprm);
1539    freenode(&tempadd);
1540    freenode(&tempf);
1541    freenode(&tmp);
1542    freenode(&tmp1);
1543    freenode(&tmp2);
1544    freenode(&tmp3);
1545    freenode(&tmprm);
1546    freenode(&tmpadd);
1547  }
1548  for (i = 0; i < spp; i++)
1549    free(y[i]);
1550  free(y);
1551  for (i = 1; i <= maxtrees; i++)
1552    free(bestrees[i - 1].btree);
1553  free(bestrees);
1554  free(nayme);
1555  free(enterorder);
1556  free(place);
1557  free(weight);
1558  free(oldweight);
1559  free(alias);
1560  free(ally);
1561  free(location);
1562  freegrbg(&grbg);
1563  if (ancseq)
1564    freegarbage(&garbage);
1565  free(threshwt);
1566  free(zeros);
1567  freenodes(nonodes, treenode);
1568}  /* freerest */
1569
1570
1571int main(int argc, Char *argv[])
1572{  /* DNA parsimony by uphill search */
1573
1574  /* reads in spp, chars, and the data. Then calls maketree to
1575     construct the tree */
1576#ifdef MAC
1577   argc = 1;        /* macsetup("Dnapars","");        */
1578   argv[0] = "Dnapars";
1579#endif
1580  init(argc, argv);
1581  progname = argv[0];
1582  openfile(&infile,INFILE,"input file", "r",argv[0],infilename);
1583  openfile(&outfile,OUTFILE,"output file", "w",argv[0],outfilename);
1584
1585  ibmpc = IBMCRT;
1586  ansi = ANSICRT;
1587  msets = 1;
1588  firstset = true;
1589  garbage = NULL;
1590  grbg = NULL;
1591  doinit();
1592  if (weights || justwts)
1593    openfile(&weightfile,WEIGHTFILE,"weights file","r",argv[0],weightfilename);
1594  if (trout)
1595    openfile(&outtree,OUTTREE,"output tree file", "w",argv[0],outtreename);
1596  for (ith = 1; ith <= msets; ith++) {
1597    if (!(justwts && !firstset))
1598      allocrest();
1599    if (msets > 1 && !justwts) {
1600    fprintf(outfile, "\nData set # %ld:\n\n", ith);
1601      if (progress)
1602        printf("\nData set # %ld:\n\n", ith);
1603    }
1604    doinput();
1605    if (ith == 1)
1606      firstset = false;
1607    for (jumb = 1; jumb <= njumble; jumb++)
1608      maketree();
1609    if (!justwts)
1610      freerest();
1611  }
1612  freetree(nonodes, treenode);
1613  FClose(infile);
1614  FClose(outfile);
1615  if (weights || justwts)
1616    FClose(weightfile);
1617  if (trout)
1618    FClose(outtree);
1619  if (usertree)
1620    FClose(intree);
1621#ifdef MAC
1622  fixmacfile(outfilename);
1623  fixmacfile(outtreename);
1624#endif
1625  if (progress)
1626    printf("Done.\n\n");
1627#ifdef WIN32
1628  phyRestoreConsoleAttributes();
1629#endif
1630  return 0;
1631}  /* DNA parsimony by uphill search */
1632
Note: See TracBrowser for help on using the repository browser.