source: trunk/GDE/PHYLIP/discrete.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: 79.0 KB
Line 
1#include "phylip.h"
2#include "discrete.h"
3
4/* version 3.6. (c) Copyright 1993-2000 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
9long nonodes, endsite, outgrno, nextree, which;
10boolean interleaved, printdata, outgropt, treeprint, dotdiff;
11steptr weight, category, alias, location, ally;
12sequence y, convtab;
13
14
15void inputdata(long chars)
16{
17  /* input the names and sequences for each species */
18  /* used by pars */
19  long i, j, k, l;
20  long basesread=0, basesnew=0, nsymbol=0, convsymboli=0;
21  Char charstate;
22  boolean allread, done, found;
23
24  if (printdata)
25    headings(chars, "Sequences", "---------");
26  basesread = 0;
27  allread = false;
28  while (!(allread)) {
29    allread = true;
30    if (eoln(infile))
31      scan_eoln(infile);
32    i = 1;
33    while (i <= spp) {
34      if ((interleaved && basesread == 0) || !interleaved)
35        initname(i - 1);
36      j = (interleaved) ? basesread : 0;
37      done = false;
38      while (!done && !eoff(infile)) {
39        if (interleaved)
40          done = true;
41        while (j < chars && !(eoln(infile) || eoff(infile))) {
42          charstate = gettc(infile);
43          if (charstate == '\n')
44            charstate = ' ';
45          if (charstate == ' ')
46            continue;
47          if ((strchr("!\"#$%&'()*+,-./0123456789:;<=>?@\
48                       ABCDEFGHIJKLMNOPQRSTUVWXYZ[\\]^_`\
49                       abcdefghijklmnopqrstuvwxyz{|}~",charstate)) == NULL){
50            printf(
51                  "\n\nERROR: Bad base: %c at position %ld of species %ld\n\n",
52                   charstate, j+1, i);
53            exxit(-1);
54          }
55          j++;
56          y[i - 1][j - 1] = charstate;
57        }
58        if (interleaved)
59          continue;
60        if (j < chars) 
61          scan_eoln(infile);
62        else if (j == chars)
63          done = true;
64      }
65      if (interleaved && i == 1)
66        basesnew = j;
67
68      scan_eoln(infile);
69     
70      if ((interleaved && j != basesnew) ||
71          (!interleaved && j != chars)) {
72        printf("\n\nERROR: Sequences out of alignment at position %ld\n\n", j);
73        exxit(-1);
74      }
75      i++;
76    }
77    if (interleaved) {
78      basesread = basesnew;
79      allread = (basesread == chars);
80    } else
81      allread = (i > spp);
82  }
83  if (printdata) {
84    for (i = 1; i <= ((chars - 1) / 60 + 1); i++) {
85      for (j = 1; j <= spp; j++) {
86        for (k = 0; k < nmlngth; k++)
87          putc(nayme[j - 1][k], outfile);
88        fprintf(outfile, "   ");
89        l = i * 60;
90        if (l > chars)
91          l = chars;
92        for (k = (i - 1) * 60 + 1; k <= l; k++) {
93          if (dotdiff && (j > 1 && y[j - 1][k - 1] == y[0][k - 1]))
94            charstate = '.';
95          else
96            charstate = y[j - 1][k - 1];
97          putc(charstate, outfile);
98          if (k % 10 == 0 && k % 60 != 0)
99            putc(' ', outfile);
100        }
101        putc('\n', outfile);
102      }
103      putc('\n', outfile);
104    }
105    putc('\n', outfile);
106  }
107  for (i = 1; i <= chars; i++) {
108    nsymbol = 0;
109    for (j = 1; j <= spp; j++) {
110      if ((nsymbol == 0) && (y[j - 1][i - 1] != '?')) {
111        nsymbol = 1;
112        convsymboli = 1;
113        convtab[0][i-1] = y[j-1][i-1];
114      } else if (y[j - 1][i - 1] != '?'){
115        found = false;
116        for (k = 1; k <= nsymbol; k++) {
117          if (convtab[k - 1][i - 1] == y[j - 1][i - 1]) {
118            found = true;
119            convsymboli = k;
120          }
121        }
122        if (!found) {
123          nsymbol++;
124          convtab[nsymbol-1][i - 1] = y[j - 1][i - 1];
125          convsymboli = nsymbol;
126        } 
127      }
128      if (nsymbol <= 8) {
129        if (y[j - 1][i - 1] != '?')
130          y[j - 1][i - 1] = (Char)('0' + (convsymboli - 1));
131      } else {
132        printf(
133          "\n\nERROR: More than maximum of 8 symbols in column %ld\n\n", i);
134        exxit(-1);
135      }
136    }
137  }
138}  /* inputdata */
139
140
141void alloctree(pointarray *treenode, long local_nonodes, boolean usertree)
142{
143  /* allocate treenode dynamically */
144  /* used in pars */
145  long i, j;
146  node *p, *q;
147
148  *treenode = (pointarray)Malloc(local_nonodes*sizeof(node *));
149  for (i = 0; i < spp; i++) {
150    (*treenode)[i] = (node *)Malloc(sizeof(node));
151    (*treenode)[i]->tip = true;
152    (*treenode)[i]->index = i+1;
153    (*treenode)[i]->iter = true;
154    (*treenode)[i]->branchnum = i+1;
155    (*treenode)[i]->initialized = true;
156  }
157  if (!usertree)
158    for (i = spp; i < local_nonodes; i++) {
159      q = NULL;
160      for (j = 1; j <= 3; j++) {
161        p = (node *)Malloc(sizeof(node));
162        p->tip = false;
163        p->index = i+1;
164        p->iter = true;
165        p->branchnum = i+1;
166        p->initialized = false;
167        p->next = q;
168        q = p;
169      }
170      p->next->next->next = p;
171      (*treenode)[i] = p;
172    }
173} /* alloctree */
174
175
176void setuptree(pointarray treenode, long local_nonodes, boolean usertree)
177{
178  /* initialize treenodes */
179  long i;
180  node *p;
181
182  for (i = 1; i <= local_nonodes; i++) {
183    if (i <= spp || !usertree) {
184      treenode[i-1]->back = NULL;
185      treenode[i-1]->tip = (i <= spp);
186      treenode[i-1]->index = i;
187      treenode[i-1]->numdesc = 0;
188      treenode[i-1]->iter = true;
189      treenode[i-1]->initialized = true;
190    }
191  }
192  if (!usertree) {
193    for (i = spp + 1; i <= local_nonodes; i++) {
194      p = treenode[i-1]->next;
195      while (p != treenode[i-1]) {
196        p->back = NULL;
197        p->tip = false;
198        p->index = i;
199        p->numdesc = 0;
200        p->iter = true;
201        p->initialized = false;
202        p = p->next;
203      }
204    }
205  }
206} /* setuptree */
207
208
209void alloctip(node *p, long *zeros, unsigned char *zeros2)
210{ /* allocate a tip node */
211  /* used by pars */
212
213  p->numsteps = (steptr)Malloc(endsite*sizeof(long));
214  p->oldnumsteps = (steptr)Malloc(endsite*sizeof(long));
215  p->discbase = (discbaseptr)Malloc(endsite*sizeof(unsigned char));
216  p->olddiscbase = (discbaseptr)Malloc(endsite*sizeof(unsigned char));
217  memcpy(p->discbase, zeros2, endsite*sizeof(unsigned char));
218  memcpy(p->numsteps, zeros, endsite*sizeof(long));
219  memcpy(p->olddiscbase, zeros2, endsite*sizeof(unsigned char));
220  memcpy(p->oldnumsteps, zeros, endsite*sizeof(long));
221}  /* alloctip */
222
223
224void sitesort(long chars, steptr local_weight)
225{
226  /* Shell sort keeping sites, weights in same order */
227  /* used in pars */
228  long gap, i, j, jj, jg, k, itemp;
229  boolean flip, tied;
230
231  gap = chars / 2;
232  while (gap > 0) {
233    for (i = gap + 1; i <= chars; i++) {
234      j = i - gap;
235      flip = true;
236      while (j > 0 && flip) {
237        jj = alias[j - 1];
238        jg = alias[j + gap - 1];
239        tied = true;
240        k = 1;
241        while (k <= spp && tied) {
242          flip = (y[k - 1][jj - 1] > y[k - 1][jg - 1]);
243          tied = (tied && y[k - 1][jj - 1] == y[k - 1][jg - 1]);
244          k++;
245        }
246        if (!flip)
247          break;
248        itemp = alias[j - 1];
249        alias[j - 1] = alias[j + gap - 1];
250        alias[j + gap - 1] = itemp;
251        itemp = local_weight[j - 1];
252        local_weight[j - 1] = local_weight[j + gap - 1];
253        local_weight[j + gap - 1] = itemp;
254        j -= gap;
255      }
256    }
257    gap /= 2;
258  }
259}  /* sitesort */
260
261
262void sitecombine(long chars)
263{
264  /* combine sites that have identical patterns */
265  /* used in pars */
266  long i, j, k;
267  boolean tied;
268
269  i = 1;
270  while (i < chars) {
271    j = i + 1;
272    tied = true;
273    while (j <= chars && tied) {
274      k = 1;
275      while (k <= spp && tied) {
276        tied = (tied &&
277            y[k - 1][alias[i - 1] - 1] == y[k - 1][alias[j - 1] - 1]);
278        k++;
279      }
280      if (tied) {
281        weight[i - 1] += weight[j - 1];
282        weight[j - 1] = 0;
283        ally[alias[j - 1] - 1] = alias[i - 1];
284      }
285      j++;
286    }
287    i = j - 1;
288  }
289}  /* sitecombine */
290
291
292void sitescrunch(long chars)
293{
294  /* move so one representative of each pattern of
295     sites comes first */
296  /* used in pars */
297  long i, j, itemp;
298  boolean done, found;
299
300  done = false;
301  i = 1;
302  j = 2;
303  while (!done) {
304    if (ally[alias[i - 1] - 1] != alias[i - 1]) {
305      if (j <= i)
306        j = i + 1;
307      if (j <= chars) {
308        found = false;
309        do {
310          found = (ally[alias[j - 1] - 1] == alias[j - 1]);
311          j++;
312        } while (!(found || j > chars));
313        if (found) {
314          j--;
315          itemp = alias[i - 1];
316          alias[i - 1] = alias[j - 1];
317          alias[j - 1] = itemp;
318          itemp = weight[i - 1];
319          weight[i - 1] = weight[j - 1];
320          weight[j - 1] = itemp;
321        } else
322          done = true;
323      } else
324        done = true;
325    }
326    i++;
327    done = (done || i >= chars);
328  }
329}  /* sitescrunch */
330
331
332void makevalues(pointarray treenode, long *zeros, unsigned char *zeros2,
333                        boolean usertree)
334{
335  /* set up fractional likelihoods at tips */
336  /* used by pars */
337  long i, j;
338  unsigned char ns=0;
339  node *p;
340
341  setuptree(treenode, nonodes, usertree);
342  for (i = 0; i < spp; i++)
343    alloctip(treenode[i], zeros, zeros2);
344  if (!usertree) {
345    for (i = spp; i < nonodes; i++) {
346      p = treenode[i];
347      do {
348        allocdiscnontip(p, zeros, zeros2, endsite);
349        p = p->next;
350      } while (p != treenode[i]);
351    }
352  }
353  for (j = 0; j < endsite; j++) {
354    for (i = 0; i < spp; i++) {
355      switch (y[i][alias[j] - 1]) {
356
357      case '0':
358        ns = 1 << ZERO;
359        break;
360
361      case '1':
362        ns = 1 << ONE;
363        break;
364
365      case '2':
366        ns = 1 << TWO;
367        break;
368
369      case '3':
370        ns = 1 << THREE;
371        break;
372
373      case '4':
374        ns = 1 << FOUR;
375        break;
376
377      case '5':
378        ns = 1 << FIVE;
379        break;
380
381      case '6':
382        ns = 1 << SIX;
383        break;
384
385      case '7':
386        ns = 1 << SEVEN;
387        break;
388
389      case '?':
390        ns = (1 << ZERO) | (1 << ONE) | (1 << TWO) | (1 << THREE) |
391               (1 << FOUR) | (1 << FIVE) | (1 << SIX) | (1 << SEVEN);
392        break;
393      }
394      treenode[i]->discbase[j] = ns;
395      treenode[i]->numsteps[j] = 0;
396    }
397  }
398}  /* makevalues */
399
400
401void fillin(node *p, node *left, node *rt)
402{
403  /* sets up for each node in the tree the base sequence
404     at that point and counts the changes.  */
405  long i, j, k, n;
406  node *q;
407
408  if (!left) {
409    memcpy(p->discbase, rt->discbase, endsite*sizeof(unsigned char));
410    memcpy(p->numsteps, rt->numsteps, endsite*sizeof(long));
411    q = rt;
412  } else if (!rt) {
413    memcpy(p->discbase, left->discbase, endsite*sizeof(unsigned char));
414    memcpy(p->numsteps, left->numsteps, endsite*sizeof(long));
415    q = left;
416  } else {
417    for (i = 0; i < endsite; i++) {
418      p->discbase[i] = left->discbase[i] & rt->discbase[i];
419      p->numsteps[i] = left->numsteps[i] + rt->numsteps[i];
420      if (p->discbase[i] == 0) {
421        p->discbase[i] = left->discbase[i] | rt->discbase[i];
422        p->numsteps[i] += weight[i];
423      }
424    }
425    q = rt;
426  }
427  if (left && rt) n = 2;
428  else n = 1;
429  for (i = 0; i < endsite; i++)
430    for (j = (long)ZERO; j <= (long)SEVEN; j++)
431      p->discnumnuc[i][j] = 0;
432  for (k = 1; k <= n; k++) {
433    if (k == 2) q = left;
434    for (i = 0; i < endsite; i++) {
435      for (j = (long)ZERO; j <= (long)SEVEN; j++) {
436        if (q->discbase[i] & (1 << j))
437          p->discnumnuc[i][j]++;
438      }
439    }
440  }
441}  /* fillin */
442
443
444long getlargest(long *discnumnuc)
445{
446  /* find the largest in array numnuc */
447  long i, largest;
448
449  largest = 0;
450  for (i = (long)ZERO; i <= (long)SEVEN; i++)
451    if (discnumnuc[i] > largest)
452      largest = discnumnuc[i];
453  return largest;
454} /* getlargest */
455
456
457void multifillin(node *p, node *q, long dnumdesc)
458{
459  /* sets up for each node in the tree the base sequence
460     at that point and counts the changes according to the
461     changes in q's base */
462  long i, j, largest, descsteps;
463  unsigned char b;
464
465  memcpy(p->olddiscbase, p->discbase, endsite*sizeof(unsigned char));
466  memcpy(p->oldnumsteps, p->numsteps, endsite*sizeof(long));
467  for (i = 0; i < endsite; i++) {
468    descsteps = 0;
469    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
470      b = 1 << j;
471      if ((descsteps == 0) && (p->discbase[i] & b)) 
472        descsteps = p->numsteps[i] 
473                      - (p->numdesc - dnumdesc - p->discnumnuc[i][j])
474                        * weight[i];
475    }
476    if (dnumdesc == -1)
477      descsteps -= q->oldnumsteps[i];
478    else if (dnumdesc == 0)
479      descsteps += (q->numsteps[i] - q->oldnumsteps[i]);
480    else
481      descsteps += q->numsteps[i];
482    if (q->olddiscbase[i] != q->discbase[i]) {
483      for (j = (long)ZERO; j <= (long)SEVEN; j++) {
484        b = 1 << j;
485        if ((q->olddiscbase[i] & b) && !(q->discbase[i] & b))
486          p->discnumnuc[i][j]--;
487        else if (!(q->olddiscbase[i] & b) && (q->discbase[i] & b))
488          p->discnumnuc[i][j]++;
489      }
490    }
491    largest = getlargest(p->discnumnuc[i]);
492    if (q->olddiscbase[i] != q->discbase[i]) {
493      p->discbase[i] = 0;
494      for (j = (long)ZERO; j <= (long)SEVEN; j++) {
495        if (p->discnumnuc[i][j] == largest)
496            p->discbase[i] |= (1 << j);
497      }
498    }
499    p->numsteps[i] = (p->numdesc - largest) * weight[i] + descsteps;
500  }
501} /* multifillin */
502
503
504void sumnsteps(node *p, node *left, node *rt, long a, long b)
505{
506  /* sets up for each node in the tree the base sequence
507     at that point and counts the changes. */
508  long i;
509  unsigned char ns, rs, ls;
510
511  if (!left) {
512    memcpy(p->numsteps, rt->numsteps, endsite*sizeof(long));
513    memcpy(p->discbase, rt->discbase, endsite*sizeof(unsigned char));
514  } else if (!rt) {
515    memcpy(p->numsteps, left->numsteps, endsite*sizeof(long));
516    memcpy(p->discbase, left->discbase, endsite*sizeof(unsigned char));
517  } else 
518    for (i = a; i < b; i++) {
519      ls = left->discbase[i];
520      rs = rt->discbase[i];
521      ns = ls & rs;
522      p->numsteps[i] = left->numsteps[i] + rt->numsteps[i];
523      if (ns == 0) {
524        ns = ls | rs;
525        p->numsteps[i] += weight[i];
526      }
527      p->discbase[i] = ns;
528    }
529}  /* sumnsteps */
530
531
532void sumnsteps2(node *p, node *left, node *rt, long a, long b, long *threshwt)
533{
534  /* counts the changes at each node.  */
535  long i, steps;
536  unsigned char ns, rs, ls;
537  long term;
538
539  if (a == 0) p->sumsteps = 0.0;
540  if (!left)
541    memcpy(p->numsteps, rt->numsteps, endsite*sizeof(long));
542  else if (!rt)
543    memcpy(p->numsteps, left->numsteps, endsite*sizeof(long));
544  else 
545    for (i = a; i < b; i++) {
546      ls = left->discbase[i];
547      rs = rt->discbase[i];
548      ns = ls & rs;
549      p->numsteps[i] = left->numsteps[i] + rt->numsteps[i];
550      if (ns == 0)
551        p->numsteps[i] += weight[i];
552    }
553  for (i = a; i < b; i++) {
554    steps = p->numsteps[i];
555    if ((long)steps <= threshwt[i])
556      term = steps;
557    else
558      term = threshwt[i];
559    p->sumsteps += (double)term;
560  }
561}  /* sumnsteps2 */
562
563
564void multisumnsteps(node *p, node *q, long a, long b, long *threshwt)
565{
566  /* sets up for each node in the tree the base sequence
567     at that point and counts the changes according to the
568     changes in q's base */
569  long i, j, steps, largest, descsteps;
570  long term;
571
572  if (a == 0) p->sumsteps = 0.0;
573  for (i = a; i < b; i++) {
574    largest = getlargest(p->discnumnuc[i]);
575    descsteps = 0;
576    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
577      if ((descsteps == 0) && (p->discbase[i] & (1 << j))) 
578        descsteps = p->numsteps[i] -
579                        (p->numdesc - 1 - p->discnumnuc[i][j]) * weight[i];
580    }
581    descsteps += q->numsteps[i];
582    largest = 0;
583    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
584      if (q->discbase[i] & (1 << j))
585        p->discnumnuc[i][j]++;
586      if (p->discnumnuc[i][j] > largest)
587        largest = p->discnumnuc[i][j];
588    }
589    steps = ((p->numdesc - largest) * weight[i] + descsteps);
590    if ((long)steps <= threshwt[i])
591      term = steps;
592    else
593      term = threshwt[i];
594    p->sumsteps += (double)term;
595  }
596} /* multisumnsteps */
597
598
599void multisumnsteps2(node *p)
600{
601  /* counts the changes at each multi-way node. Sums up
602     steps of all descendants */
603  long i, j, largest;
604  node *q;
605  discbaseptr b;
606
607  for (i = 0; i < endsite; i++) {
608    p->numsteps[i] = 0;
609    q = p->next;
610    while (q != p) {
611      if (q->back) {
612        p->numsteps[i] += q->back->numsteps[i];
613        b = q->back->discbase;
614        for (j = (long)ZERO; j <= (long)SEVEN; j++)
615          if (b[i] & (1 << j))
616            p->discnumnuc[i][j]++;
617      }
618      q = q->next;
619    }
620    largest = getlargest(p->discnumnuc[i]);
621    p->numsteps[i] += ((p->numdesc - largest) * weight[i]);
622    p->discbase[i] = 0;
623    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
624      if (p->discnumnuc[i][j] == largest)
625       p->discbase[i] |= (1 << j);
626    }
627  }
628}  /* multisumnsteps2 */
629
630
631boolean alltips(node *forknode, node *p)
632{
633  /* returns true if all descendants of forknode except p are tips;
634     false otherwise.  */
635  node *q, *r;
636  boolean tips;
637
638  tips = true;
639  r = forknode;
640  q = forknode->next;
641  do {
642    if (q->back && q->back != p && !q->back->tip)
643      tips = false;
644    q = q->next;
645  } while (tips && q != r);
646  return tips;
647} /* alltips */
648
649
650void gdispose(node *p, node **grbg, pointarray treenode)
651{
652  /* go through tree throwing away nodes */
653  node *q, *r;
654
655  p->back = NULL;
656  if (p->tip)
657    return;
658  treenode[p->index - 1] = NULL;
659  q = p->next;
660  while (q != p) {
661    gdispose(q->back, grbg, treenode);
662    q->back = NULL;
663    r = q;
664    q = q->next;
665    chucktreenode(grbg, r);
666  }
667  chucktreenode(grbg, q);
668}  /* gdispose */
669
670
671void preorder(node *p, node *r, node *root, node *removing, node *adding,
672                        node *changing, long dnumdesc)
673{
674  /* recompute number of steps in preorder taking both ancestoral and
675     descendent steps into account. removing points to a node being
676     removed, if any */
677  node *q, *p1, *p2;
678
679  if (p && !p->tip && p != adding) {
680    q = p;
681    do {
682      if (p->back != r) {
683        if (p->numdesc > 2) {
684          if (changing)
685            multifillin (p, r, dnumdesc);
686          else
687            multifillin (p, r, 0);
688        } else {
689          p1 = p->next;
690          if (!removing)
691            while (!p1->back)
692              p1 = p1->next;
693          else
694            while (!p1->back || p1->back == removing)
695              p1 = p1->next;
696          p2 = p1->next;
697          if (!removing)
698            while (!p2->back)
699              p2 = p2->next;
700          else
701            while (!p2->back || p2->back == removing)
702              p2 = p2->next;
703          p1 = p1->back;
704          p2 = p2->back;
705          if (p->back == p1) p1 = NULL;
706          else if (p->back == p2) p2 = NULL;
707          memcpy(p->olddiscbase, p->discbase, endsite*sizeof(unsigned char));
708          memcpy(p->oldnumsteps, p->numsteps, endsite*sizeof(long));
709          fillin(p, p1, p2);
710        }
711      }
712      p = p->next;
713    } while (p != q);
714    q = p;
715    do {
716      preorder(p->next->back, p->next, root, removing, adding, NULL, 0);
717      p = p->next;
718    } while (p->next != q);
719  }
720} /* preorder */
721
722
723void updatenumdesc(node *p, node *root, long n)
724{
725  /* set p's numdesc to n.  If p is the root, numdesc of p's
726  descendants are set to n-1. */
727  node *q;
728
729  q = p;
730  if (p == root && n > 0) {
731    p->numdesc = n;
732    n--;
733    q = q->next;
734  }
735  do {
736    q->numdesc = n;
737    q = q->next;
738  } while (q != p);
739}
740
741
742void add(node *below, node *newtip, node *newfork, node **root, boolean recompute,
743                        pointarray treenode, node **grbg, long *zeros, unsigned char *zeros2)
744{
745  /* inserts the nodes newfork and its left descendant, newtip,
746     to the tree.  below becomes newfork's right descendant.
747     if newfork is NULL, newtip is added as below's sibling */
748  /* used in pars */
749  node *p;
750
751  if (below != treenode[below->index - 1])
752    below = treenode[below->index - 1];
753  if (newfork) {
754    if (below->back != NULL)
755      below->back->back = newfork;
756    newfork->back = below->back;
757    below->back = newfork->next->next;
758    newfork->next->next->back = below;
759    newfork->next->back = newtip;
760    newtip->back = newfork->next;
761    if (*root == below)
762      *root = newfork;
763    updatenumdesc(newfork, *root, 2);
764  } else {
765    gnudisctreenode(grbg, &p, below->index, endsite, zeros, zeros2);
766    p->back = newtip;
767    newtip->back = p;
768    p->next = below->next;
769    below->next = p;
770    updatenumdesc(below, *root, below->numdesc + 1);
771  }
772  if (!newtip->tip)
773    updatenumdesc(newtip, *root, newtip->numdesc);
774  (*root)->back = NULL;
775  if (!recompute)
776    return;
777  if (!newfork) {
778    memcpy(newtip->back->discbase, below->discbase, endsite*sizeof(unsigned char));
779    memcpy(newtip->back->numsteps, below->numsteps, endsite*sizeof(long));
780    memcpy(newtip->back->discnumnuc, below->discnumnuc, endsite*sizeof(discnucarray));
781    if (below != *root) {
782      memcpy(below->back->olddiscbase, zeros2, endsite*sizeof(unsigned char));
783      memcpy(below->back->oldnumsteps, zeros, endsite*sizeof(long));
784      multifillin(newtip->back, below->back, 1);
785    }
786    if (!newtip->tip) {
787      memcpy(newtip->back->olddiscbase, zeros2, endsite*sizeof(unsigned char));
788      memcpy(newtip->back->oldnumsteps, zeros, endsite*sizeof(long));
789      preorder(newtip, newtip->back, *root, NULL, NULL, below, 1);
790    }
791    memcpy(newtip->olddiscbase, zeros2, endsite*sizeof(unsigned char));
792    memcpy(newtip->oldnumsteps, zeros, endsite*sizeof(long));
793    preorder(below, newtip, *root, NULL, newtip, below, 1);
794    if (below != *root)
795      preorder(below->back, below, *root, NULL, NULL, NULL, 0);
796  } else {
797    fillin(newtip->back, newtip->back->next->back,
798             newtip->back->next->next->back);
799    if (!newtip->tip) {
800      memcpy(newtip->back->olddiscbase, zeros2, endsite*sizeof(unsigned char));
801      memcpy(newtip->back->oldnumsteps, zeros, endsite*sizeof(long));
802      preorder(newtip, newtip->back, *root, NULL, NULL, newfork, 1);
803    }
804    if (newfork != *root) {
805      memcpy(below->back->discbase, newfork->back->discbase, endsite*sizeof(unsigned char));
806      memcpy(below->back->numsteps, newfork->back->numsteps, endsite*sizeof(long));
807      preorder(newfork, newtip, *root, NULL, newtip, NULL, 0);
808    } else {
809      fillin(below->back, newtip, NULL);
810      fillin(newfork, newtip, below);
811      memcpy(below->back->olddiscbase, zeros2, endsite*sizeof(unsigned char));
812      memcpy(below->back->oldnumsteps, zeros, endsite*sizeof(long));
813      preorder(below, below->back, *root, NULL, NULL, newfork, 1);
814    }
815    if (newfork != *root) {
816      memcpy(newfork->olddiscbase, below->discbase, endsite*sizeof(unsigned char));
817      memcpy(newfork->oldnumsteps, below->numsteps, endsite*sizeof(long));
818      preorder(newfork->back, newfork, *root, NULL, NULL, NULL, 0);
819    }
820  }
821}  /* add */
822
823
824void findbelow(node **below, node *item, node *fork)
825{
826  /* decide which of fork's binary children is below */
827
828  if (fork->next->back == item)
829    *below = fork->next->next->back;
830  else
831    *below = fork->next->back;
832} /* findbelow */
833
834
835void re_move(node *item, node **fork, node **root, boolean recompute,
836                        pointarray treenode, node **grbg, long *zeros, unsigned char *zeros2)
837{
838  /* removes nodes item and its ancestor, fork, from the tree.
839     the new descendant of fork's ancestor is made to be
840     fork's second descendant (other than item).  Also
841     returns pointers to the deleted nodes, item and fork.
842     If item belongs to a node with more than 2 descendants,
843     fork will not be deleted */
844  /* used in pars */
845  node *p, *q, *other, *otherback = NULL;
846
847  if (item->back == NULL) {
848    *fork = NULL;
849    return;
850  }
851  *fork = treenode[item->back->index - 1];
852  if ((*fork)->numdesc == 2) {
853    updatenumdesc(*fork, *root, 0);
854    findbelow(&other, item, *fork);
855    otherback = other->back;
856    if (*root == *fork) {
857      if (other->tip)
858        *root = NULL;
859      else {
860        *root = other;
861        updatenumdesc(other, *root, other->numdesc);
862      }
863    }
864    p = item->back->next->back;
865    q = item->back->next->next->back;
866    if (p != NULL)
867      p->back = q;
868    if (q != NULL)
869      q->back = p;
870    (*fork)->back = NULL;
871    p = (*fork)->next;
872    while (p != *fork) {
873      p->back = NULL;
874      p = p->next;
875    }
876  } else {
877    updatenumdesc(*fork, *root, (*fork)->numdesc - 1);
878    p = *fork;
879    while (p->next != item->back)
880      p = p->next;
881    p->next = item->back->next;
882  }
883  if (!item->tip) {
884    updatenumdesc(item, item, item->numdesc);
885    if (recompute) {
886      memcpy(item->back->olddiscbase, item->back->discbase,
887               endsite*sizeof(unsigned char));
888      memcpy(item->back->oldnumsteps, item->back->numsteps, endsite*sizeof(long));
889      memcpy(item->back->discbase, zeros2, endsite*sizeof(unsigned char));
890      memcpy(item->back->numsteps, zeros, endsite*sizeof(long));
891      preorder(item, item->back, *root, item->back, NULL, item, -1);
892    }
893  }
894  if ((*fork)->numdesc >= 2)
895    chucktreenode(grbg, item->back);
896  item->back = NULL;
897  if (!recompute)
898    return;
899  if ((*fork)->numdesc == 0) {
900    memcpy(otherback->olddiscbase, otherback->discbase,
901             endsite*sizeof(unsigned char)); 
902    memcpy(otherback->oldnumsteps, otherback->numsteps, endsite*sizeof(long));
903    if (other == *root) {
904      memcpy(otherback->discbase, zeros2, endsite*sizeof(unsigned char));
905      memcpy(otherback->numsteps, zeros, endsite*sizeof(long));
906    } else {
907      memcpy(otherback->discbase, other->back->discbase,
908               endsite*sizeof(unsigned char));
909      memcpy(otherback->numsteps, other->back->numsteps, endsite*sizeof(long));
910    }
911    p = other->back;
912    other->back = otherback;
913    if (other == *root)
914      preorder(other, otherback, *root, otherback, NULL, other, -1);
915    else
916      preorder(other, otherback, *root, NULL, NULL, NULL, 0);
917    other->back = p;
918    if (other != *root) {
919      memcpy(other->olddiscbase,(*fork)->discbase, endsite*sizeof(unsigned char));
920      memcpy(other->oldnumsteps,(*fork)->numsteps, endsite*sizeof(long));
921      preorder(other->back, other, *root, NULL, NULL, NULL, 0);
922    }
923  } else {
924    memcpy(item->olddiscbase, item->discbase, endsite*sizeof(unsigned char));
925    memcpy(item->oldnumsteps, item->numsteps, endsite*sizeof(long));
926    memcpy(item->discbase, zeros2, endsite*sizeof(unsigned char));
927    memcpy(item->numsteps, zeros, endsite*sizeof(long));
928    preorder(*fork, item, *root, NULL, NULL, *fork, -1);
929    if (*fork != *root)
930      preorder((*fork)->back, *fork, *root, NULL, NULL, NULL, 0);
931    memcpy(item->discbase, item->olddiscbase, endsite*sizeof(unsigned char));
932    memcpy(item->numsteps, item->oldnumsteps, endsite*sizeof(long));
933  }
934}  /* re_move */
935
936
937void postorder(node *p)
938{
939  /* traverses an n-ary tree, suming up steps at a node's descendants */
940  /* used in pars */
941  node *q;
942
943  if (p->tip)
944    return;
945  q = p->next;
946  while (q != p) {
947    postorder(q->back);
948    q = q->next;
949  }
950  zerodiscnumnuc(p, endsite);
951  if (p->numdesc > 2)
952    multisumnsteps2(p);
953  else
954    fillin(p, p->next->back, p->next->next->back);
955}  /* postorder */
956
957
958void getnufork(node **nufork, node **grbg, pointarray treenode, long *zeros,
959                        unsigned char *zeros2)
960{
961  /* find a fork not used currently */
962  long i;
963
964  i = spp;
965  while (treenode[i] && treenode[i]->numdesc > 0) i++;
966  if (!treenode[i])
967    gnudisctreenode(grbg, &treenode[i], i, endsite, zeros, zeros2);
968  *nufork = treenode[i];
969} /* getnufork */
970
971
972void reroot(node *outgroup, node *root)
973{
974  /* reorients tree, putting outgroup in desired position. used if
975     the root is binary. */
976  /* used in pars */
977  node *p, *q;
978
979  if (outgroup->back->index == root->index)
980    return;
981  p = root->next;
982  q = root->next->next;
983  p->back->back = q->back;
984  q->back->back = p->back;
985  p->back = outgroup;
986  q->back = outgroup->back;
987  outgroup->back->back = q;
988  outgroup->back = p;
989}  /* reroot */
990
991
992void reroot2(node *outgroup, node *root)
993{
994  /* reorients tree, putting outgroup in desired position. */
995  /* used in pars */
996  node *p;
997
998  p = outgroup->back->next;
999  while (p->next != outgroup->back)
1000    p = p->next;
1001  root->next = outgroup->back;
1002  p->next = root;
1003}  /* reroot2 */
1004
1005
1006void reroot3(node *outgroup,node *root,node *root2,node *lastdesc,node **grbg)
1007{
1008  /* reorients tree, putting back outgroup in original position. */
1009  /* used in pars */
1010  node *p;
1011
1012  p = root->next;
1013  while (p->next != root)
1014    p = p->next;
1015  chucktreenode(grbg, root);
1016  p->next = outgroup->back;
1017  root2->next = lastdesc->next;
1018  lastdesc->next = root2;
1019}  /* reroot3 */
1020
1021
1022void savetraverse(node *p)
1023{
1024  /* sets BOOLEANs that indicate which way is down */
1025  node *q;
1026
1027  p->bottom = true;
1028  if (p->tip)
1029    return;
1030  q = p->next;
1031  while (q != p) {
1032    q->bottom = false;
1033    savetraverse(q->back);
1034    q = q->next;
1035  }
1036}  /* savetraverse */
1037
1038
1039void newindex(long i, node *p)
1040{
1041  /* assigns index i to node p */
1042
1043  while (p->index != i) {
1044    p->index = i;
1045    p = p->next;
1046  }
1047} /* newindex */
1048
1049
1050void flipindexes(long nextnode, pointarray treenode)
1051{
1052  /* flips index of nodes between nextnode and last node.  */
1053  long last;
1054  node *temp;
1055
1056  last = nonodes;
1057  while (treenode[last - 1]->numdesc == 0)
1058    last--;
1059  if (last > nextnode) {
1060    temp = treenode[nextnode - 1];
1061    treenode[nextnode - 1] = treenode[last - 1];
1062    treenode[last - 1] = temp;
1063    newindex(nextnode, treenode[nextnode - 1]);
1064    newindex(last, treenode[last - 1]);
1065  }
1066} /* flipindexes */ 
1067
1068
1069boolean parentinmulti(node *anode)
1070{
1071  /* sees if anode's parent has more than 2 children */
1072  node *p;
1073
1074  while (!anode->bottom) anode = anode->next;
1075  p = anode->back;
1076  while (!p->bottom)
1077    p = p->next;
1078  return (p->numdesc > 2);
1079} /* parentinmulti */
1080
1081
1082long sibsvisited(node *anode, long *place)
1083{
1084  /* computes the number of nodes which are visited earlier than anode among
1085  its siblings */
1086  node *p;
1087  long nvisited;
1088
1089  while (!anode->bottom) anode = anode->next;
1090  p = anode->back->next;
1091  nvisited = 0;
1092  do {
1093    if (!p->bottom && place[p->back->index - 1] != 0)
1094      nvisited++;
1095    p = p->next;
1096  } while (p != anode->back);
1097  return nvisited;
1098}  /* sibsvisited */
1099
1100
1101long smallest(node *anode, long *place)
1102{
1103  /* finds the smallest index of sibling of anode */
1104  node *p;
1105  long min;
1106
1107  while (!anode->bottom) anode = anode->next;
1108  p = anode->back->next;
1109  if (p->bottom) p = p->next;
1110  min = nonodes;
1111  do {
1112    if (p->back && place[p->back->index - 1] != 0) {
1113      if (p->back->index <= spp) {
1114        if (p->back->index < min)
1115          min = p->back->index;
1116      } else {
1117        if (place[p->back->index - 1] < min)
1118          min = place[p->back->index - 1];
1119      }
1120    }
1121    p = p->next;
1122    if (p->bottom) p = p->next;
1123  } while (p != anode->back);
1124  return min;
1125}  /* smallest */
1126
1127
1128void bintomulti(node **root, node **binroot, node **grbg, long *zeros,
1129                        unsigned char *zeros2)
1130{  /* attaches root's left child to its right child and makes
1131      the right child new root */
1132  node *left, *right, *newnode, *temp;
1133
1134  right = (*root)->next->next->back;
1135  left = (*root)->next->back;
1136  if (right->tip) {
1137    (*root)->next = right->back;
1138    (*root)->next->next = left->back;
1139    temp = left;
1140    left = right;
1141    right = temp;
1142    right->back->next = *root;
1143  }
1144  gnudisctreenode(grbg, &newnode, right->index, endsite, zeros, zeros2);
1145  newnode->next = right->next;
1146  newnode->back = left;
1147  left->back = newnode;
1148  right->next = newnode;
1149  (*root)->next->back = (*root)->next->next->back = NULL;
1150  *binroot = *root;
1151  (*binroot)->numdesc = 0;
1152  *root = right;
1153  (*root)->numdesc++;
1154  (*root)->back = NULL;
1155} /* bintomulti */
1156
1157
1158void backtobinary(node **root, node *binroot, node **grbg)
1159{ /* restores binary root */
1160  node *p;
1161
1162  binroot->next->back = (*root)->next->back;
1163  (*root)->next->back->back = binroot->next;
1164  p = (*root)->next;
1165  (*root)->next = p->next;
1166  binroot->next->next->back = *root;
1167  (*root)->back = binroot->next->next;
1168  chucktreenode(grbg, p);
1169  (*root)->numdesc--;
1170  *root = binroot;
1171  (*root)->numdesc = 2;
1172} /* backtobinary */
1173
1174
1175boolean outgrin(node *root, node *outgrnode)
1176{ /* checks if outgroup node is a child of root */
1177  node *p;
1178
1179  p = root->next;
1180  while (p != root) {
1181    if (p->back == outgrnode)
1182      return true;
1183    p = p->next;
1184  }
1185  return false;
1186} /* outgrin */
1187
1188
1189void flipnodes(node *nodea, node *nodeb)
1190{ /* flip nodes */
1191  node *backa, *backb;
1192
1193  backa = nodea->back;
1194  backb = nodeb->back;
1195  backa->back = nodeb;
1196  backb->back = nodea;
1197  nodea->back = backb;
1198  nodeb->back = backa;
1199} /* flipnodes */
1200
1201
1202void moveleft(node *root, node *outgrnode, node **flipback)
1203{ /* makes outgroup node to leftmost child of root */
1204  node *p;
1205  boolean done;
1206
1207  p = root->next;
1208  done = false;
1209  while (p != root && !done) {
1210    if (p->back == outgrnode) {
1211      *flipback = p;
1212      flipnodes(root->next->back, p->back);
1213      done = true;
1214    }
1215    p = p->next;
1216  }
1217} /* moveleft */
1218
1219
1220void savetree(node *root, long *place, pointarray treenode,
1221                        node **grbg, long *zeros, unsigned char *zeros2)
1222{  /* record in place where each species has to be
1223     added to reconstruct this tree */
1224  /* used by pars */
1225  long i, j, nextnode, nvisited;
1226  node *p, *q, *r = NULL, *root2, *lastdesc,
1227                *outgrnode, *binroot, *flipback;
1228  boolean done, newfork;
1229
1230  binroot = NULL;
1231  lastdesc = NULL;
1232  root2 = NULL;
1233  flipback = NULL;
1234  outgrnode = treenode[outgrno - 1];
1235  if (root->numdesc == 2)
1236    bintomulti(&root, &binroot, grbg, zeros, zeros2);
1237  if (outgrin(root, outgrnode)) {
1238    if (outgrnode != root->next->back)
1239      moveleft(root, outgrnode, &flipback);
1240  } else {
1241    root2 = root;
1242    lastdesc = root->next;
1243    while (lastdesc->next != root)
1244      lastdesc = lastdesc->next;
1245    lastdesc->next = root->next;
1246    gnudisctreenode(grbg, &root, outgrnode->back->index, endsite, zeros, zeros2);
1247    root->numdesc = root2->numdesc;
1248    reroot2(outgrnode, root);
1249  }
1250  savetraverse(root);
1251  nextnode = spp + 1;
1252  for (i = nextnode; i <= nonodes; i++)
1253    if (treenode[i - 1]->numdesc == 0)
1254      flipindexes(i, treenode);
1255  for (i = 0; i < nonodes; i++)
1256    place[i] = 0;
1257  place[root->index - 1] = 1;
1258  for (i = 1; i <= spp; i++) {
1259    p = treenode[i - 1];
1260    while (place[p->index - 1] == 0) {
1261      place[p->index - 1] = i;
1262      while (!p->bottom)
1263        p = p->next;
1264      r = p;
1265      p = p->back;
1266    }
1267    if (i > 1) {
1268      q = treenode[i - 1]; 
1269      newfork = true;
1270      nvisited = sibsvisited(q, place);
1271      if (nvisited == 0) {
1272        if (parentinmulti(r)) {
1273          nvisited = sibsvisited(r, place);
1274          if (nvisited == 0)
1275            place[i - 1] = place[p->index - 1];
1276          else if (nvisited == 1)
1277            place[i - 1] = smallest(r, place);
1278          else {
1279            place[i - 1] = -smallest(r, place);
1280            newfork = false;
1281          }
1282        } else
1283          place[i - 1] = place[p->index - 1];
1284      } else if (nvisited == 1) {
1285        place[i - 1] = place[p->index - 1];
1286      } else {
1287        place[i - 1] = -smallest(q, place);
1288        newfork = false;
1289      }
1290      if (newfork) {
1291        j = place[p->index - 1];
1292        done = false;
1293        while (!done) {
1294          place[p->index - 1] = nextnode;
1295          while (!p->bottom)
1296            p = p->next;
1297          p = p->back;
1298          done = (p == NULL);
1299          if (!done)
1300            done = (place[p->index - 1] != j);
1301          if (done) {
1302            nextnode++;
1303          }
1304        }
1305      }
1306    }
1307  }
1308  if (flipback)
1309    flipnodes(outgrnode, flipback->back);
1310  else {
1311    if (root2) {
1312      reroot3(outgrnode, root, root2, lastdesc, grbg);
1313      root = root2;
1314    }
1315  }
1316  if (binroot)
1317    backtobinary(&root, binroot, grbg);
1318}  /* savetree */ 
1319
1320
1321void addnsave(node *p, node *item, node *nufork, node **root, node **grbg,
1322                        boolean multf, pointarray treenode, long *place, long *zeros,
1323                        unsigned char *zeros2)
1324{  /* adds item to tree and save it.  Then removes item. */
1325  node *dummy;
1326
1327  if (!multf)
1328    add(p, item, nufork, root, false, treenode, grbg, zeros, zeros2);
1329  else
1330    add(p, item, NULL, root, false, treenode, grbg, zeros, zeros2);
1331  savetree(*root, place, treenode, grbg, zeros, zeros2);
1332  if (!multf)
1333    re_move(item, &nufork, root, false, treenode, grbg, zeros, zeros2);
1334  else
1335    re_move(item, &dummy, root, false, treenode, grbg, zeros, zeros2);
1336} /* addnsave */
1337
1338
1339void addbestever(long *pos, long *local_nextree, long maxtrees, boolean collapse,
1340                        long *place, bestelm *bestrees)
1341{ /* adds first best tree */
1342
1343  *pos = 1;
1344  *local_nextree = 1;
1345  initbestrees(bestrees, maxtrees, true);
1346  initbestrees(bestrees, maxtrees, false);
1347  addtree(*pos, local_nextree, collapse, place, bestrees);
1348} /* addbestever */
1349
1350
1351void addtiedtree(long pos, long *local_nextree, long maxtrees, boolean collapse,
1352                        long *place, bestelm *bestrees)
1353{ /* add tied tree */
1354
1355  if (*local_nextree <= maxtrees)
1356    addtree(pos, local_nextree, collapse, place, bestrees);
1357} /* addtiedtree */
1358
1359
1360void clearcollapse(pointarray treenode)
1361{
1362  /* clears collapse status at a node */
1363  long i;
1364  node *p;
1365
1366  for (i = 0; i < nonodes; i++) {
1367    treenode[i]->collapse = undefined;
1368    if (!treenode[i]->tip) {
1369      p = treenode[i]->next;
1370      while (p != treenode[i]) {
1371        p->collapse = undefined;
1372        p = p->next;
1373      }
1374    }
1375  }
1376} /* clearcollapse */
1377
1378
1379void clearbottom(pointarray treenode)
1380{
1381  /* clears boolean bottom at a node */
1382  long i;
1383  node *p;
1384
1385  for (i = 0; i < nonodes; i++) {
1386    treenode[i]->bottom = false;
1387    if (!treenode[i]->tip) {
1388      p = treenode[i]->next;
1389      while (p != treenode[i]) {
1390        p->bottom = false;
1391        p = p->next;
1392      }
1393    }
1394  }
1395} /* clearbottom */
1396
1397
1398void collabranch(node *collapfrom, node *tempfrom, node *tempto)
1399{ /* collapse branch from collapfrom */
1400  long i, j, largest, descsteps;
1401  boolean done;
1402  unsigned char b;
1403
1404  for (i = 0; i < endsite; i++) {
1405    largest = getlargest(collapfrom->discnumnuc[i]);
1406    descsteps = 0;
1407    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
1408      b = 1 << j;
1409      if ((descsteps == 0) && (collapfrom->discbase[i] & b)) 
1410        descsteps = tempfrom->oldnumsteps[i] 
1411                     - (collapfrom->numdesc - collapfrom->discnumnuc[i][j])
1412                       * weight[i];
1413    }
1414    largest = getlargest(tempto->discnumnuc[i]);
1415    done = false;
1416    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
1417      b = 1 << j;
1418      if (!done && (tempto->discbase[i] & b)) {
1419        descsteps += (tempto->numsteps[i] 
1420                      - (tempto->numdesc - collapfrom->numdesc
1421                        - tempto->discnumnuc[i][j]) * weight[i]);
1422        done = true;
1423      }
1424    }
1425    for (j = (long)ZERO; j <= (long)SEVEN; j++)
1426      tempto->discnumnuc[i][j] += collapfrom->discnumnuc[i][j];
1427    largest = getlargest(tempto->discnumnuc[i]);
1428    tempto->discbase[i] = 0;
1429    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
1430      if (tempto->discnumnuc[i][j] == largest)
1431        tempto->discbase[i] |= (1 << j);
1432    }
1433    tempto->numsteps[i] = (tempto->numdesc - largest) * weight[i] + descsteps;
1434  }
1435} /* collabranch */
1436
1437
1438boolean allcommonbases(node *a, node *b, boolean *allsame)
1439{  /* see if bases are common at all sites for nodes a and b */   
1440  long i;
1441  boolean allcommon;
1442
1443  allcommon = true;
1444  *allsame = true;
1445  for (i = 0; i < endsite; i++) {
1446    if ((a->discbase[i] & b->discbase[i]) == 0)
1447      allcommon = false;
1448    else if (a->discbase[i] != b->discbase[i])
1449      *allsame = false;
1450  }
1451  return allcommon;
1452} /* allcommonbases */
1453
1454
1455void findbottom(node *p, node **local_bottom)
1456{ /* find a node with field bottom set at node p */
1457  node *q;
1458
1459  if (p->bottom)
1460    *local_bottom = p;
1461  else {
1462    q = p->next;
1463    while(!q->bottom && q != p)
1464      q = q->next;
1465    *local_bottom = q;
1466  }
1467} /* findbottom */
1468
1469
1470boolean moresteps(node *a, node *b)
1471{  /* see if numsteps of node a exceeds those of node b */   
1472  long i;
1473
1474  for (i = 0; i < endsite; i++)
1475    if (a->numsteps[i] > b->numsteps[i])
1476      return true;
1477  return false;
1478} /* moresteps */
1479
1480
1481boolean passdown(node *desc, node *parent, node *start, node *below, node *item,
1482                        node *added, node *total, node *tempdsc, node *tempprt,
1483                        boolean multf)
1484{ /* track down to node start to see if an ancestor branch can be collapsed */
1485  node *temp;
1486  boolean done, allsame;
1487
1488  done = (parent == start);
1489  while (!done) {
1490    desc = parent;
1491    findbottom(parent->back, &parent);
1492    if (multf && start == below && parent == below)
1493      parent = added;
1494    memcpy(tempdsc->discbase, tempprt->discbase, endsite*sizeof(unsigned char));
1495    memcpy(tempdsc->numsteps, tempprt->numsteps, endsite*sizeof(long));
1496    memcpy(tempdsc->olddiscbase, desc->discbase, endsite*sizeof(unsigned char));
1497    memcpy(tempdsc->oldnumsteps, desc->numsteps, endsite*sizeof(long));
1498    memcpy(tempprt->discbase, parent->discbase, endsite*sizeof(unsigned char));
1499    memcpy(tempprt->numsteps, parent->numsteps, endsite*sizeof(long));
1500    memcpy(tempprt->discnumnuc, parent->discnumnuc, endsite*sizeof(discnucarray));
1501    tempprt->numdesc = parent->numdesc;
1502    multifillin(tempprt, tempdsc, 0);
1503    if (!allcommonbases(tempprt, parent, &allsame))
1504      return false;
1505    else if (moresteps(tempprt, parent))
1506      return false;
1507    else if (allsame)
1508      return true;
1509    if (parent == added)
1510      parent = below;
1511    done = (parent == start);
1512    if (done && ((start == item) || (!multf && start == below))) {
1513      memcpy(tempdsc->discbase, tempprt->discbase, endsite*sizeof(unsigned char));
1514      memcpy(tempdsc->numsteps, tempprt->numsteps, endsite*sizeof(long));
1515      memcpy(tempdsc->olddiscbase, start->discbase, endsite*sizeof(unsigned char));
1516      memcpy(tempdsc->oldnumsteps, start->numsteps, endsite*sizeof(long));
1517      multifillin(added, tempdsc, 0);
1518      tempprt = added;
1519    }
1520  } 
1521  temp = tempdsc;
1522  if (start == below || start == item)
1523    fillin(temp, tempprt, below->back);
1524  else
1525    fillin(temp, tempprt, added);
1526  return !moresteps(temp, total);
1527} /* passdown */
1528
1529
1530boolean trycollapdesc(node *desc, node *parent, node *start, node *below,
1531                        node *item, node *added, node *total, node *tempdsc, node *tempprt,
1532                        boolean multf,long *zeros, unsigned char *zeros2)
1533{ /* see if branch between nodes desc and parent can be collapsed */
1534  boolean allsame;
1535
1536  if (desc->numdesc == 1)
1537    return true;
1538  if (multf && start == below && parent == below)
1539    parent = added;
1540  memcpy(tempdsc->discbase, zeros2, endsite*sizeof(unsigned char));
1541  memcpy(tempdsc->numsteps, zeros, endsite*sizeof(long));
1542  memcpy(tempdsc->olddiscbase, desc->discbase, endsite*sizeof(unsigned char));
1543  memcpy(tempdsc->oldnumsteps, desc->numsteps, endsite*sizeof(long));
1544  memcpy(tempprt->discbase, parent->discbase, endsite*sizeof(unsigned char));
1545  memcpy(tempprt->numsteps, parent->numsteps, endsite*sizeof(long));
1546  memcpy(tempprt->discnumnuc, parent->discnumnuc, endsite*sizeof(discnucarray));
1547  tempprt->numdesc = parent->numdesc - 1;
1548  multifillin(tempprt, tempdsc, -1);
1549  tempprt->numdesc += desc->numdesc;
1550  collabranch(desc, tempdsc, tempprt);
1551  if (!allcommonbases(tempprt, parent, &allsame) || 
1552        moresteps(tempprt, parent)) {
1553    if (parent != added) {
1554      desc->collapse = nocollap;
1555      parent->collapse = nocollap;
1556    }
1557    return false;
1558  } else if (allsame) {
1559    if (parent != added) {
1560      desc->collapse = tocollap;
1561      parent->collapse = tocollap;
1562    }
1563    return true;
1564  }
1565  if (parent == added)
1566    parent = below;
1567  if ((start == item && parent == item) ||
1568        (!multf && start == below && parent == below)) {
1569    memcpy(tempdsc->discbase, tempprt->discbase, endsite*sizeof(unsigned char));
1570    memcpy(tempdsc->numsteps, tempprt->numsteps, endsite*sizeof(long));
1571    memcpy(tempdsc->olddiscbase, start->discbase, endsite*sizeof(unsigned char));
1572    memcpy(tempdsc->oldnumsteps, start->numsteps, endsite*sizeof(long));
1573    memcpy(tempprt->discbase, added->discbase, endsite*sizeof(unsigned char));
1574    memcpy(tempprt->numsteps, added->numsteps, endsite*sizeof(long));
1575    memcpy(tempprt->discnumnuc, added->discnumnuc, endsite*sizeof(discnucarray));
1576    tempprt->numdesc = added->numdesc;
1577    multifillin(tempprt, tempdsc, 0);
1578    if (!allcommonbases(tempprt, added, &allsame))
1579      return false;
1580    else if (moresteps(tempprt, added))
1581      return false;
1582    else if (allsame)
1583      return true;
1584  }
1585  return passdown(desc, parent, start, below, item, added, total, tempdsc,
1586                    tempprt, multf);
1587} /* trycollapdesc */
1588
1589
1590void setbottom(node *p)
1591{ /* set field bottom at node p */
1592  node *q;
1593
1594  p->bottom = true;
1595  q = p->next;
1596  do {
1597    q->bottom = false;
1598    q = q->next;
1599  } while (q != p);
1600} /* setbottom */
1601
1602
1603boolean zeroinsubtree(node *subtree, node *start, node *below, node *item,
1604                        node *added, node *total, node *tempdsc, node *tempprt,
1605                        boolean multf, node* root, long *zeros, 
1606                        unsigned char *zeros2)
1607{ /* sees if subtree contains a zero length branch */
1608  node *p;
1609
1610  if (!subtree->tip) {
1611    setbottom(subtree);
1612    p = subtree->next;
1613    do {
1614      if (p->back && !p->back->tip && 
1615         !((p->back->collapse == nocollap) && (subtree->collapse == nocollap))
1616           && (subtree->numdesc != 1)) {
1617        if ((p->back->collapse == tocollap) && (subtree->collapse == tocollap)
1618             && multf && (subtree != below))
1619          return true;
1620        /* when root->numdesc == 2
1621         * there is no mandatory step at the root,
1622         * instead of checking at the root we check around it
1623         * we only need to check p because the first if
1624         * statement already gets rid of it for the subtree */
1625        else if ((p->back->index != root->index || root->numdesc > 2) &&
1626            trycollapdesc(p->back, subtree, start, below, item, added, total, 
1627                tempdsc, tempprt, multf, zeros, zeros2))
1628          return true;
1629        else if ((p->back->index == root->index && root->numdesc == 2) &&
1630            !(root->next->back->tip) && !(root->next->next->back->tip) &&
1631            trycollapdesc(root->next->back, root->next->next->back, start,
1632                below, item, added, total, tempdsc, tempprt, multf, zeros,
1633                zeros2))
1634          return true;
1635      }
1636      p = p->next;
1637    } while (p != subtree);
1638    p = subtree->next;
1639    do {
1640      if (p->back && !p->back->tip) {
1641        if (zeroinsubtree(p->back, start, below, item, added, total, 
1642                            tempdsc, tempprt, multf, root, zeros, zeros2))
1643          return true;
1644      }
1645      p = p->next;
1646    } while (p != subtree);
1647  }
1648  return false;
1649} /* zeroinsubtree */
1650
1651
1652boolean collapsible(node *item, node *below, node *temp, node *temp1, 
1653                node *tempdsc, node *tempprt, node *added, node *total, 
1654                boolean multf, node *root, long *zeros, unsigned char *zeros2, 
1655                pointarray treenode)
1656{
1657  /* sees if any branch can be collapsed */
1658  node *belowbk;
1659  boolean allsame;
1660
1661  if (multf) {
1662    memcpy(tempdsc->discbase, item->discbase, endsite*sizeof(unsigned char));
1663    memcpy(tempdsc->numsteps, item->numsteps, endsite*sizeof(long));
1664    memcpy(tempdsc->olddiscbase, zeros2, endsite*sizeof(unsigned char));
1665    memcpy(tempdsc->oldnumsteps, zeros, endsite*sizeof(long));
1666    memcpy(added->discbase, below->discbase, endsite*sizeof(unsigned char));
1667    memcpy(added->numsteps, below->numsteps, endsite*sizeof(long));
1668    memcpy(added->discnumnuc, below->discnumnuc, endsite*sizeof(discnucarray));
1669    added->numdesc = below->numdesc + 1;
1670    multifillin(added, tempdsc, 1);
1671  } else {
1672    fillin(added, item, below);
1673    added->numdesc = 2;
1674  }
1675  fillin(total, added, below->back);
1676  clearbottom(treenode);
1677  if (below->back) {
1678    if (zeroinsubtree(below->back, below->back, below, item, added, total,
1679                        tempdsc, tempprt, multf, root, zeros, zeros2))
1680      return true;
1681  }
1682  if (multf) {
1683    if (zeroinsubtree(below, below, below, item, added, total,
1684                       tempdsc, tempprt, multf, root, zeros, zeros2))
1685      return true;
1686  } else if (!below->tip) {
1687    if (zeroinsubtree(below, below, below, item, added, total,
1688                        tempdsc, tempprt, multf, root, zeros, zeros2))
1689      return true;
1690  }
1691  if (!item->tip) {
1692    if (zeroinsubtree(item, item, below, item, added, total,
1693                        tempdsc, tempprt, multf, root, zeros, zeros2))
1694      return true;
1695  }
1696  if (multf && below->back && !below->back->tip) {
1697    memcpy(tempdsc->discbase, zeros2, endsite*sizeof(unsigned char));
1698    memcpy(tempdsc->numsteps, zeros, endsite*sizeof(long));
1699    memcpy(tempdsc->olddiscbase, added->discbase, endsite*sizeof(unsigned char));
1700    memcpy(tempdsc->oldnumsteps, added->numsteps, endsite*sizeof(long));
1701    if (below->back == treenode[below->back->index - 1])
1702      belowbk = below->back->next;
1703    else
1704      belowbk = treenode[below->back->index - 1];
1705    memcpy(tempprt->discbase, belowbk->discbase, endsite*sizeof(unsigned char));
1706    memcpy(tempprt->numsteps, belowbk->numsteps, endsite*sizeof(long));
1707    memcpy(tempprt->discnumnuc, belowbk->discnumnuc, endsite*sizeof(discnucarray));
1708    tempprt->numdesc = belowbk->numdesc - 1;
1709    multifillin(tempprt, tempdsc, -1);
1710    tempprt->numdesc += added->numdesc;
1711    collabranch(added, tempdsc, tempprt);
1712    if (!allcommonbases(tempprt, belowbk, &allsame))
1713      return false;
1714    else if (allsame && !moresteps(tempprt, belowbk))
1715      return true;
1716    else if (belowbk->back) {
1717      fillin(temp, tempprt, belowbk->back);
1718      fillin(temp1, belowbk, belowbk->back);
1719      return !moresteps(temp, temp1);
1720    }
1721  }
1722  return false;
1723} /* collapsible */
1724
1725
1726void replaceback(node **oldback, node *item, node *forknode, node **grbg,
1727                        long *zeros, unsigned char *zeros2)
1728{ /* replaces back node of item with another */
1729  node *p;
1730
1731  p = forknode;
1732  while (p->next->back != item)
1733    p = p->next;
1734  *oldback = p->next;
1735  gnudisctreenode(grbg, &p->next, forknode->index, endsite, zeros, zeros2);
1736  p->next->next = (*oldback)->next;
1737  p->next->back = (*oldback)->back;
1738  p->next->back->back = p->next;
1739  (*oldback)->next = (*oldback)->back = NULL;
1740} /* replaceback */
1741
1742
1743void putback(node *oldback, node *item, node *forknode, node **grbg)
1744{ /* restores node to back of item */
1745  node *p, *q;
1746
1747  p = forknode;
1748  while (p->next != item->back)
1749    p = p->next;
1750  q = p->next;
1751  oldback->next = p->next->next;
1752  p->next = oldback;
1753  oldback->back = item;
1754  item->back = oldback;
1755  oldback->index = forknode->index;
1756  chucktreenode(grbg, q);
1757} /* putback */
1758
1759
1760void savelocrearr(node *item, node *forknode, node *below, node *tmp, node *tmp1,
1761                        node *tmp2, node *tmp3, node *tmprm, node *tmpadd, node **root,
1762                        long maxtrees, long *local_nextree, boolean multf, boolean bestever,
1763                        boolean *saved, long *place, bestelm *bestrees, pointarray treenode,
1764                        node **grbg, long *zeros, unsigned char *zeros2)
1765{ /* saves tied or better trees during local rearrangements by removing
1766     item from forknode and adding to below */
1767  node *other, *otherback=NULL, *oldfork, *nufork, *oldback;
1768  long pos;
1769  boolean found, collapse;
1770
1771  if (forknode->numdesc == 2) {
1772    findbelow(&other, item, forknode);
1773    otherback = other->back;
1774    oldback = NULL;
1775  } else {
1776    other = NULL;
1777    replaceback(&oldback, item, forknode, grbg, zeros, zeros2);
1778  }
1779  re_move(item, &oldfork, root, false, treenode, grbg, zeros, zeros2);
1780  if (!multf)
1781    getnufork(&nufork, grbg, treenode, zeros, zeros2);
1782  else
1783    nufork = NULL;
1784  addnsave(below, item, nufork, root, grbg, multf, treenode, place,
1785             zeros, zeros2);
1786  pos = 0;
1787  findtree(&found, &pos, *local_nextree, place, bestrees);
1788  if (other) {
1789    add(other, item, oldfork, root, false, treenode, grbg, zeros, zeros2);
1790    if (otherback->back != other)
1791      flipnodes(item, other);
1792  } else
1793    add(forknode, item, NULL, root, false, treenode, grbg, zeros, zeros2);
1794  *saved = false;
1795  if (found) {
1796    if (oldback)
1797      putback(oldback, item, forknode, grbg);
1798  } else {
1799    if (oldback)
1800      chucktreenode(grbg, oldback);
1801    re_move(item, &oldfork, root, true, treenode, grbg, zeros, zeros2);
1802    collapse = collapsible(item, below, tmp, tmp1, tmp2, tmp3, tmprm,
1803                     tmpadd, multf, *root, zeros, zeros2, treenode);
1804    if (!collapse) {
1805      if (bestever)
1806        addbestever(&pos, local_nextree, maxtrees, collapse, place, bestrees);
1807      else
1808        addtiedtree(pos, local_nextree, maxtrees, collapse, place, bestrees);
1809    }
1810    if (other)
1811      add(other, item, oldfork, root, true, treenode, grbg, zeros, zeros2);
1812    else
1813      add(forknode, item, NULL, root, true, treenode, grbg, zeros, zeros2);
1814    *saved = !collapse;
1815  }
1816} /* savelocrearr */
1817
1818
1819void clearvisited(pointarray treenode)
1820{
1821  /* clears boolean visited at a node */
1822  long i;
1823  node *p;
1824
1825  for (i = 0; i < nonodes; i++) {
1826    treenode[i]->visited = false;
1827    if (!treenode[i]->tip) {
1828      p = treenode[i]->next;
1829      while (p != treenode[i]) {
1830        p->visited = false;
1831        p = p->next;
1832      }
1833    }
1834  }
1835} /* clearvisited */
1836
1837
1838void hyprint(long b1,long b2,struct LOC_hyptrav *htrav,pointarray treenode)
1839{
1840  /* print out states in sites b1 through b2 at node */
1841  long i, j, k;
1842  boolean dot, found;
1843
1844  if (htrav->bottom) {
1845    if (!outgropt)
1846      fprintf(outfile, "       ");
1847    else
1848      fprintf(outfile, "root   ");
1849  } else
1850    fprintf(outfile, "%4ld   ", htrav->r->back->index - spp);
1851  if (htrav->r->tip) {
1852    for (i = 0; i < nmlngth; i++)
1853      putc(nayme[htrav->r->index - 1][i], outfile);
1854  } else
1855    fprintf(outfile, "%4ld      ", htrav->r->index - spp);
1856  if (htrav->bottom)
1857    fprintf(outfile, "          ");
1858  else if (htrav->nonzero)
1859    fprintf(outfile, "   yes    ");
1860  else if (htrav->maybe)
1861    fprintf(outfile, "  maybe   ");
1862  else
1863    fprintf(outfile, "   no     ");
1864  for (i = b1; i <= b2; i++) {
1865    j = location[ally[i - 1] - 1];
1866    htrav->tempset = htrav->r->discbase[j - 1];
1867    htrav->anc = htrav->hypset[j - 1];
1868    if (!htrav->bottom)
1869      htrav->anc = treenode[htrav->r->back->index - 1]->discbase[j - 1];
1870    dot = dotdiff && (htrav->tempset == htrav->anc && !htrav->bottom);
1871    if (dot)
1872      putc('.', outfile); 
1873    else {
1874      found = false;
1875      k = (long)ZERO;
1876      do {
1877        if (htrav->tempset == (1 << k)) {
1878          putc(convtab[k][i - 1], outfile);
1879          found = true;
1880        }
1881        k++;
1882      } while (!found && k <= (long)SEVEN);
1883      if (!found) 
1884        putc('?', outfile);
1885    }
1886    if (i % 10 == 0)
1887      putc(' ', outfile);
1888  }
1889  putc('\n', outfile);
1890}  /* hyprint */
1891
1892
1893void gnubase(gbases **p, gbases **garbage, long local_endsite)
1894{
1895  /* this and the following are do-it-yourself garbage collectors.
1896     Make a new node or pull one off the garbage list */
1897  if (*garbage != NULL) {
1898    *p = *garbage;
1899    *garbage = (*garbage)->next;
1900  } else {
1901    *p = (gbases *)Malloc(sizeof(gbases));
1902    (*p)->discbase = (discbaseptr)Malloc(local_endsite*sizeof(unsigned char));
1903  }
1904  (*p)->next = NULL;
1905}  /* gnubase */
1906
1907
1908void chuckbase(gbases *p, gbases **garbage)
1909{
1910  /* collect garbage on p -- put it on front of garbage list */
1911  p->next = *garbage;
1912  *garbage = p;
1913}  /* chuckbase */
1914
1915
1916void hyptrav(node *r_, discbaseptr hypset_, long b1, long b2, boolean bottom_,
1917                        pointarray treenode, gbases **garbage)
1918{
1919  /*  compute, print out states at one interior node */
1920  struct LOC_hyptrav Vars;
1921  long i, j, k;
1922  long largest;
1923  gbases *ancset;
1924  discnucarray *tempnuc;
1925  node *p, *q;
1926
1927  Vars.bottom = bottom_;
1928  Vars.r = r_;
1929  Vars.hypset = hypset_;
1930  gnubase(&ancset, garbage, endsite);
1931  tempnuc = (discnucarray *)Malloc(endsite*sizeof(discnucarray));
1932  Vars.maybe = false;
1933  Vars.nonzero = false;
1934  if (!Vars.r->tip)
1935    zerodiscnumnuc(Vars.r, endsite);
1936  for (i = b1 - 1; i < b2; i++) {
1937    j = location[ally[i] - 1];
1938    Vars.anc = Vars.hypset[j - 1];
1939    if (!Vars.r->tip) {
1940      p = Vars.r->next;
1941      for (k = (long)ZERO; k <= (long)SEVEN; k++)
1942        if (Vars.anc & (1 << k))
1943          Vars.r->discnumnuc[j - 1][k]++;
1944      do {
1945        for (k = (long)ZERO; k <= (long)SEVEN; k++)
1946          if (p->back->discbase[j - 1] & (1 << k))
1947            Vars.r->discnumnuc[j - 1][k]++;
1948        p = p->next;
1949      } while (p != Vars.r);
1950      largest = getlargest(Vars.r->discnumnuc[j - 1]);
1951      Vars.tempset = 0;
1952      for (k = (long)ZERO; k <= (long)SEVEN; k++) {
1953        if (Vars.r->discnumnuc[j - 1][k] == largest)
1954          Vars.tempset |= (1 << k);
1955      }
1956      Vars.r->discbase[j - 1] = Vars.tempset;
1957    }
1958    if (!Vars.bottom)
1959      Vars.anc = treenode[Vars.r->back->index - 1]->discbase[j - 1];
1960    Vars.nonzero = (Vars.nonzero || (Vars.r->discbase[j - 1] & Vars.anc) == 0);
1961    Vars.maybe = (Vars.maybe || Vars.r->discbase[j - 1] != Vars.anc);
1962  }
1963  hyprint(b1, b2, &Vars, treenode);
1964  Vars.bottom = false;
1965  if (!Vars.r->tip) {
1966    memcpy(tempnuc, Vars.r->discnumnuc, endsite*sizeof(discnucarray));
1967    q = Vars.r->next;
1968    do {
1969      memcpy(Vars.r->discnumnuc, tempnuc, endsite*sizeof(discnucarray));
1970      for (i = b1 - 1; i < b2; i++) {
1971        j = location[ally[i] - 1];
1972        for (k = (long)ZERO; k <= (long)SEVEN; k++)
1973          if (q->back->discbase[j - 1] & (1 << k))
1974            Vars.r->discnumnuc[j - 1][k]--;
1975        largest = getlargest(Vars.r->discnumnuc[j - 1]);
1976        ancset->discbase[j - 1] = 0;
1977        for (k = (long)ZERO; k <= (long)SEVEN; k++)
1978          if (Vars.r->discnumnuc[j - 1][k] == largest)
1979            ancset->discbase[j - 1] |= (1 << k);
1980        if (!Vars.bottom)
1981          Vars.anc = ancset->discbase[j - 1];
1982      }
1983      hyptrav(q->back, ancset->discbase, b1, b2, Vars.bottom,
1984                treenode, garbage);
1985      q = q->next;
1986    } while (q != Vars.r);
1987  }
1988  chuckbase(ancset, garbage);
1989}  /* hyptrav */
1990
1991
1992void hypstates(long chars, node *root, pointarray treenode, gbases **garbage)
1993{
1994  /* fill in and describe states at interior nodes */
1995  /* used in pars */
1996  long i, n;
1997  discbaseptr nothing;
1998
1999  fprintf(outfile, "\nFrom    To     Any Steps?    State at upper node\n");
2000  fprintf(outfile, "                            ");
2001  if (dotdiff)
2002    fprintf(outfile, " ( . means same as in the node below it on tree)\n");
2003  nothing = (discbaseptr)Malloc(endsite*sizeof(unsigned char));
2004  for (i = 0; i < endsite; i++)
2005    nothing[i] = 0;
2006  for (i = 1; i <= ((chars - 1) / 40 + 1); i++) {
2007    putc('\n', outfile);
2008    n = i * 40;
2009    if (n > chars)
2010      n = chars;
2011    hyptrav(root, nothing, i * 40 - 39, n, true, treenode, garbage);
2012  }
2013  free(nothing);
2014}  /* hypstates */
2015
2016
2017void initbranchlen(node *p)
2018{
2019  node *q;
2020
2021  p->v = 0.0;
2022  if (p->back)
2023    p->back->v = 0.0;
2024  if (p->tip)
2025    return;
2026  q = p->next;
2027  while (q != p) {
2028    initbranchlen(q->back);
2029    q = q->next;
2030  }
2031  q = p->next;
2032  while (q != p) {
2033    q->v = 0.0;
2034    q = q->next;
2035  }
2036} /* initbranchlen */
2037
2038
2039void initmin(node *p, long sitei, boolean internal)
2040{
2041  long i;
2042
2043  if (internal) {
2044    for (i = (long)ZERO; i <= (long)SEVEN; i++) {
2045      p->disccumlengths[i] = 0;
2046      p->discnumreconst[i] = 1;
2047    }
2048  } else {
2049    for (i = (long)ZERO; i <= (long)SEVEN; i++) {
2050      if (p->discbase[sitei - 1] & (1 << i)) {
2051        p->disccumlengths[i] = 0;
2052        p->discnumreconst[i] = 1;
2053      } else {
2054        p->disccumlengths[i] = -1;
2055        p->discnumreconst[i] = 0;
2056      }
2057    }
2058  }
2059} /* initmin */
2060
2061
2062void initbase(node *p, long sitei)
2063{
2064  /* traverse tree to initialize base at internal nodes */
2065  node *q;
2066  long i, largest;
2067
2068  if (p->tip)
2069    return;
2070  q = p->next;
2071  while (q != p) {
2072    if (q->back) {
2073      memcpy(q->discnumnuc, p->discnumnuc, endsite*sizeof(discnucarray));
2074      for (i = (long)ZERO; i <= (long)SEVEN; i++) {
2075        if (q->back->discbase[sitei - 1] & (1 << i))
2076          q->discnumnuc[sitei - 1][i]--;
2077      }
2078      if (p->back) {
2079        for (i = (long)ZERO; i <= (long)SEVEN; i++) {
2080          if (p->back->discbase[sitei - 1] & (1 << i))
2081            q->discnumnuc[sitei - 1][i]++;
2082        }
2083      }
2084      largest = getlargest(q->discnumnuc[sitei - 1]);
2085      q->discbase[sitei - 1] = 0;
2086      for (i = (long)ZERO; i <= (long)SEVEN; i++) {
2087        if (q->discnumnuc[sitei - 1][i] == largest)
2088          q->discbase[sitei - 1] |= (1 << i);
2089      }
2090    }
2091    q = q->next;
2092  }
2093  q = p->next;
2094  while (q != p) {
2095    initbase(q->back, sitei);
2096    q = q->next;
2097  }
2098} /* initbase */
2099
2100
2101void inittreetrav(node *p, long sitei)
2102{
2103  /* traverse tree to clear boolean initialized and set up base */
2104  node *q;
2105
2106  if (p->tip) {
2107    initmin(p, sitei, false);
2108    p->initialized = true;
2109    return;
2110  }
2111  q = p->next;
2112  while (q != p) {
2113    inittreetrav(q->back, sitei);
2114    q = q->next;
2115  }
2116  initmin(p, sitei, true);
2117  p->initialized = false;
2118  q = p->next;
2119  while (q != p) {
2120    initmin(q, sitei, true);
2121    q->initialized = false;
2122    q = q->next;
2123  }
2124} /* inittreetrav */
2125
2126
2127void compmin(node *p, node *desc)
2128{
2129  /* computes minimum lengths up to p */
2130  long i, j, minn, cost, desclen, descrecon=0, maxx;
2131
2132  maxx = 10 * spp;
2133  for (i = (long)ZERO; i <= (long)SEVEN; i++) {
2134    minn = maxx;
2135    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
2136      if (i == j) 
2137        cost = 0;
2138      else
2139        cost = 1;
2140      if (desc->disccumlengths[j] == -1) {
2141        desclen = maxx;
2142      } else {
2143        desclen = desc->disccumlengths[j];
2144      }
2145      if (minn > cost + desclen) {
2146        minn = cost + desclen;
2147        descrecon = 0;
2148      }
2149      if (minn == cost + desclen) {
2150        descrecon += desc->discnumreconst[j];
2151      }
2152    }
2153    p->disccumlengths[i] += minn;
2154    p->discnumreconst[i] *= descrecon;
2155  }
2156  p->initialized = true;
2157} /* compmin */
2158
2159
2160void minpostorder(node *p, pointarray treenode)
2161{
2162  /* traverses an n-ary tree, computing minimum steps at each node */
2163  node *q;
2164
2165  if (p->tip) {
2166    return;
2167  }
2168  q = p->next;
2169  while (q != p) {
2170    if (q->back)
2171      minpostorder(q->back, treenode);
2172    q = q->next;
2173  }
2174  if (!p->initialized) {
2175    q = p->next;
2176    while (q != p) {
2177      if (q->back)
2178        compmin(p, q->back);
2179      q = q->next;
2180    }
2181  }
2182}  /* minpostorder */
2183
2184
2185void branchlength(node *subtr1, node *subtr2, double *brlen, 
2186                               pointarray treenode)
2187{
2188  /* computes a branch length between two subtrees for a given site */
2189  long i, j, minn, cost, nom, denom;
2190  node *temp;
2191
2192  if (subtr1->tip) {
2193    temp = subtr1;
2194    subtr1 = subtr2;
2195    subtr2 = temp;
2196  }
2197  if (subtr1->index == outgrno) {
2198    temp = subtr1;
2199    subtr1 = subtr2;
2200    subtr2 = temp;
2201  }
2202  minpostorder(subtr1, treenode);
2203  minpostorder(subtr2, treenode);
2204  minn = 10 * spp;
2205  nom = 0;
2206  denom = 0;
2207  for (i = (long)ZERO; i <= (long)SEVEN; i++) {
2208    for (j = (long)ZERO; j <= (long)SEVEN; j++) {
2209      if (i == j)
2210        cost = 0;
2211      else
2212        cost = 1;
2213      if (subtr1->disccumlengths[i] != -1 && (subtr2->disccumlengths[j] != -1)) {
2214        if (subtr1->disccumlengths[i] + cost + subtr2->disccumlengths[j] < minn) {
2215          minn = subtr1->disccumlengths[i] + cost + subtr2->disccumlengths[j];
2216          nom = 0;
2217          denom = 0;
2218        }
2219        if (subtr1->disccumlengths[i] + cost + subtr2->disccumlengths[j] == minn) {
2220          nom += subtr1->discnumreconst[i] * subtr2->discnumreconst[j] * cost;
2221          denom += subtr1->discnumreconst[i] * subtr2->discnumreconst[j];
2222        }
2223      }
2224    }
2225  }
2226  *brlen = (double)nom/(double)denom;
2227} /* branchlength */ 
2228
2229
2230void printbranchlengths(node *p)
2231{
2232  node *q;
2233  long i;
2234
2235  if (p->tip)
2236    return;
2237  q = p->next;
2238  do {
2239    fprintf(outfile, "%6ld      ",q->index - spp);
2240    if (q->back->tip) {
2241      for (i = 0; i < nmlngth; i++)
2242        putc(nayme[q->back->index - 1][i], outfile);
2243    } else
2244      fprintf(outfile, "%6ld    ", q->back->index - spp);
2245    fprintf(outfile, "     %.2f\n",q->v);
2246    if (q->back)
2247      printbranchlengths(q->back);
2248    q = q->next;
2249  } while (q != p);
2250} /* printbranchlengths */
2251
2252
2253void branchlentrav(node *p, node *root, long sitei, long chars,
2254                        double *brlen, pointarray treenode)
2255  {
2256  /*  traverses the tree computing tree length at each branch */
2257  node *q;
2258
2259  if (p->tip)
2260    return;
2261  if (p->index == outgrno)
2262    p = p->back;
2263  q = p->next;
2264  do {
2265    if (q->back) {
2266      branchlength(q, q->back, brlen, treenode);
2267      q->v += ((weight[sitei - 1] / 10.0) * (*brlen));
2268      q->back->v += ((weight[sitei - 1] / 10.0) * (*brlen));
2269      if (!q->back->tip)
2270        branchlentrav(q->back, root, sitei, chars, brlen, treenode);
2271    }
2272    q = q->next;
2273  } while (q != p);
2274}  /* branchlentrav */
2275
2276
2277void treelength(node *root, long chars, pointarray treenode)
2278  {
2279  /*  calls branchlentrav at each site */
2280  long sitei;
2281  double trlen;
2282
2283  initbranchlen(root);
2284  for (sitei = 1; sitei <= endsite; sitei++) {
2285    trlen = 0.0;
2286    initbase(root, sitei);
2287    inittreetrav(root, sitei);
2288    branchlentrav(root, root, sitei, chars, &trlen, treenode);
2289  }
2290} /* treelength */
2291
2292
2293void coordinates(node *p, long *tipy, double f, long *fartemp)
2294{
2295  /* establishes coordinates of nodes for display without lengths */
2296  node *q, *first, *last, *mid1 = NULL, *mid2 = NULL;
2297  long numbranches, numb2;
2298
2299  if (p->tip) {
2300    p->xcoord = 0;
2301    p->ycoord = *tipy;
2302    p->ymin = *tipy;
2303    p->ymax = *tipy;
2304    (*tipy) += down;
2305    return;
2306  }
2307  numbranches = 0;
2308  q = p->next;
2309  do {
2310    coordinates(q->back, tipy, f, fartemp);
2311    numbranches += 1;
2312    q = q->next;
2313  } while (p != q);
2314  first = p->next->back;
2315  q = p->next;
2316  while (q->next != p) 
2317    q = q->next;
2318  last = q->back;
2319  numb2 = 1;
2320  q = p->next;
2321  while (q != p) {
2322    if (numb2 == (numbranches + 1)/2)
2323      mid1 = q->back;
2324    if (numb2 == (numbranches/2 + 1))
2325      mid2 = q->back;
2326    numb2 += 1;
2327    q = q->next;
2328  }
2329  p->xcoord = (long)((double)(last->ymax - first->ymin) * f);
2330  p->ycoord = (long)((mid1->ycoord + mid2->ycoord) / 2);
2331  p->ymin = first->ymin;
2332  p->ymax = last->ymax;
2333  if (p->xcoord > *fartemp)
2334    *fartemp = p->xcoord;
2335}  /* coordinates */
2336
2337
2338void drawline(long i, double scale, node *root)
2339{
2340  /* draws one row of the tree diagram by moving up tree */
2341  node *p, *q, *r, *first = NULL, *last = NULL;
2342  long n, j;
2343  boolean extra, done, noplus;
2344
2345  p = root;
2346  q = root;
2347  extra = false;
2348  noplus = false;
2349  if (i == (long)p->ycoord && p == root) {
2350    if (p->index - spp >= 10)
2351      fprintf(outfile, " %2ld", p->index - spp);
2352    else
2353      fprintf(outfile, "  %ld", p->index - spp);
2354    extra = true;
2355    noplus = true;
2356  } else
2357    fprintf(outfile, "  ");
2358  do {
2359    if (!p->tip) {
2360      r = p->next;
2361      done = false;
2362      do {
2363        if (i >= r->back->ymin && i <= r->back->ymax) {
2364          q = r->back;
2365          done = true;
2366        }
2367        r = r->next;
2368      } while (!(done || r == p));
2369      first = p->next->back;
2370      r = p->next;
2371      while (r->next != p)
2372        r = r->next;
2373      last = r->back;
2374    }
2375    done = (p == q);
2376    n = (long)(scale * (p->xcoord - q->xcoord) + 0.5);
2377    if (n < 3 && !q->tip)
2378      n = 3;
2379    if (extra) {
2380      n--;
2381      extra = false;
2382    }
2383    if ((long)q->ycoord == i && !done) {
2384      if (noplus) {
2385        putc('-', outfile);
2386        noplus = false;
2387      }
2388      else
2389        putc('+', outfile);
2390      if (!q->tip) {
2391        for (j = 1; j <= n - 2; j++)
2392          putc('-', outfile);
2393        if (q->index - spp >= 10)
2394          fprintf(outfile, "%2ld", q->index - spp);
2395        else
2396          fprintf(outfile, "-%ld", q->index - spp);
2397        extra = true;
2398        noplus = true;
2399      } else {
2400        for (j = 1; j < n; j++)
2401          putc('-', outfile);
2402      }
2403    } else if (!p->tip) {
2404      if ((long)last->ycoord > i && (long)first->ycoord < i
2405            && i != (long)p->ycoord) {
2406        putc('!', outfile);
2407        for (j = 1; j < n; j++)
2408          putc(' ', outfile);
2409      } else {
2410        for (j = 1; j <= n; j++)
2411          putc(' ', outfile);
2412      }
2413      noplus = false;
2414    } else {
2415      for (j = 1; j <= n; j++)
2416        putc(' ', outfile);
2417      noplus = false;
2418    }
2419    if (p != q)
2420      p = q;
2421  } while (!done);
2422  if ((long)p->ycoord == i && p->tip) {
2423    for (j = 0; j < nmlngth; j++)
2424      putc(nayme[p->index - 1][j], outfile);
2425  }
2426  putc('\n', outfile);
2427}  /* drawline */
2428
2429
2430void printree(node *root, double f)
2431{
2432  /* prints out diagram of the tree */
2433  /* used in pars */
2434  long i, tipy, dummy;
2435  double scale;
2436
2437  putc('\n', outfile);
2438  if (!treeprint)
2439    return;
2440  putc('\n', outfile);
2441  tipy = 1;
2442  dummy = 0;
2443  coordinates(root, &tipy, f, &dummy);
2444  scale = 1.5;
2445  putc('\n', outfile);
2446  for (i = 1; i <= (tipy - down); i++)
2447    drawline(i, scale, root);
2448  fprintf(outfile, "\n  remember:");
2449  if (outgropt)
2450    fprintf(outfile, " (although rooted by outgroup)");
2451  fprintf(outfile, " this is an unrooted tree!\n\n");
2452}  /* printree */
2453
2454
2455void writesteps(long chars, boolean weights, steptr oldweight, node *root)
2456{
2457  /* used in pars */
2458  long i, j, k, l;
2459
2460  putc('\n', outfile);
2461  if (weights)
2462    fprintf(outfile, "weighted ");
2463  fprintf(outfile, "steps in each site:\n");
2464  fprintf(outfile, "      ");
2465  for (i = 0; i <= 9; i++)
2466    fprintf(outfile, "%4ld", i);
2467  fprintf(outfile, "\n     *------------------------------------");
2468  fprintf(outfile, "-----\n");
2469  for (i = 0; i <= (chars / 10); i++) {
2470    fprintf(outfile, "%5ld", i * 10);
2471    putc('|', outfile);
2472    for (j = 0; j <= 9; j++) {
2473      k = i * 10 + j;
2474      if (k == 0 || k > chars)
2475        fprintf(outfile, "    ");
2476      else {
2477        l = location[ally[k - 1] - 1];
2478        if (oldweight[k - 1] > 0)
2479          fprintf(outfile, "%4ld",
2480                  oldweight[k - 1] *
2481                  (root->numsteps[l - 1] / weight[l - 1]));
2482        else
2483          fprintf(outfile, "   0");
2484      }
2485    }
2486    putc('\n', outfile);
2487  }
2488} /* writesteps */
2489
2490
2491void treeout(node *p, long local_nextree, long *col, node *root)
2492{
2493  /* write out file with representation of final tree */
2494  /* used in pars */
2495  node *q;
2496  long i, n;
2497  Char c;
2498
2499  if (p->tip) {
2500    n = 0;
2501    for (i = 1; i <= nmlngth; i++) {
2502      if (nayme[p->index - 1][i - 1] != ' ')
2503        n = i;
2504    }
2505    for (i = 0; i < n; i++) {
2506      c = nayme[p->index - 1][i];
2507      if (c == ' ')
2508        c = '_';
2509      putc(c, outtree);
2510    }
2511    *col += n;
2512  } else {
2513    putc('(', outtree);
2514    (*col)++;
2515    q = p->next;
2516    while (q != p) {
2517      treeout(q->back, local_nextree, col, root);
2518      q = q->next;
2519      if (q == p)
2520        break;
2521      putc(',', outtree);
2522      (*col)++;
2523      if (*col > 60) {
2524        putc('\n', outtree);
2525        *col = 0;
2526      }
2527    }
2528    putc(')', outtree);
2529    (*col)++;
2530  }
2531  if (p != root)
2532    return;
2533  if (local_nextree > 2)
2534    fprintf(outtree, "[%6.4f];\n", 1.0 / (local_nextree - 1));
2535  else
2536    fprintf(outtree, ";\n");
2537}  /* treeout */
2538
2539
2540void treeout3(node *p, long local_nextree, long *col, node *root)
2541{
2542  /* write out file with representation of final tree */
2543  /* used in dnapars -- writes branch lengths */
2544  node *q;
2545  long i, n, w;
2546  double x;
2547  Char c;
2548
2549  if (p->tip) {
2550    n = 0;
2551    for (i = 1; i <= nmlngth; i++) {
2552      if (nayme[p->index - 1][i - 1] != ' ')
2553        n = i;
2554    }
2555    for (i = 0; i < n; i++) {
2556      c = nayme[p->index - 1][i];
2557      if (c == ' ')
2558        c = '_';
2559      putc(c, outtree);
2560    }
2561    *col += n;
2562  } else {
2563    putc('(', outtree);
2564    (*col)++;
2565    q = p->next;
2566    while (q != p) {
2567      treeout3(q->back, local_nextree, col, root);
2568      q = q->next;
2569      if (q == p)
2570        break;
2571      putc(',', outtree);
2572      (*col)++;
2573      if (*col > 60) {
2574        putc('\n', outtree);
2575        *col = 0;
2576      }
2577    }
2578    putc(')', outtree);
2579    (*col)++;
2580  }
2581  x = p->v;
2582  if (x > 0.0)
2583    w = (long)(0.43429448222 * log(x));
2584  else if (x == 0.0)
2585    w = 0;
2586  else
2587    w = (long)(0.43429448222 * log(-x)) + 1;
2588  if (w < 0)
2589    w = 0;
2590  if (p != root) {
2591    fprintf(outtree, ":%*.2f", (int)(w + 4), x);
2592    col += w + 8;
2593  }
2594  if (p != root)
2595    return;
2596  if (local_nextree > 2)
2597    fprintf(outtree, "[%6.4f];\n", 1.0 / (local_nextree - 1));
2598  else
2599    fprintf(outtree, ";\n");
2600}  /* treeout3 */
2601
2602
2603void drawline3(long i, double scale, node *start)
2604{
2605  /* draws one row of the tree diagram by moving up tree */
2606  /* used in pars */
2607  node *p, *q;
2608  long n, j;
2609  boolean extra;
2610  node *r, *first = NULL, *last = NULL;
2611  boolean done;
2612
2613  p = start;
2614  q = start;
2615  extra = false;
2616  if (i == (long)p->ycoord) {
2617    if (p->index - spp >= 10)
2618      fprintf(outfile, " %2ld", p->index - spp);
2619    else
2620      fprintf(outfile, "  %ld", p->index - spp);
2621    extra = true;
2622  } else
2623    fprintf(outfile, "  ");
2624  do {
2625    if (!p->tip) {
2626      r = p->next;
2627      done = false;
2628      do {
2629        if (i >= r->back->ymin && i <= r->back->ymax) {
2630          q = r->back;
2631          done = true;
2632        }
2633        r = r->next;
2634      } while (!(done || (r == p))); 
2635      first = p->next->back;
2636      r = p;
2637      while (r->next != p)
2638        r = r->next;
2639      last = r->back;
2640    }
2641    done = (p->tip || p == q);
2642    n = (long)(scale * (q->xcoord - p->xcoord) + 0.5);
2643    if (n < 3 && !q->tip)
2644      n = 3;
2645    if (extra) {
2646      n--;
2647      extra = false;
2648    }
2649    if ((long)q->ycoord == i && !done) {
2650      if ((long)p->ycoord != (long)q->ycoord)
2651        putc('+', outfile);
2652      else
2653        putc('-', outfile);
2654      if (!q->tip) {
2655        for (j = 1; j <= n - 2; j++)
2656          putc('-', outfile);
2657        if (q->index - spp >= 10)
2658          fprintf(outfile, "%2ld", q->index - spp);
2659        else
2660          fprintf(outfile, "-%ld", q->index - spp);
2661        extra = true;
2662      } else {
2663        for (j = 1; j < n; j++)
2664          putc('-', outfile);
2665      }
2666    } else if (!p->tip) {
2667      if ((long)last->ycoord > i && (long)first->ycoord < i &&
2668          (i != (long)p->ycoord || p == start)) {
2669        putc('|', outfile);
2670        for (j = 1; j < n; j++)
2671          putc(' ', outfile);
2672      } else {
2673        for (j = 1; j <= n; j++)
2674          putc(' ', outfile);
2675      }
2676    } else {
2677      for (j = 1; j <= n; j++)
2678        putc(' ', outfile);
2679    }
2680    if (q != p)
2681      p = q;
2682  } while (!done);
2683  if ((long)p->ycoord == i && p->tip) {
2684    for (j = 0; j < nmlngth; j++)
2685      putc(nayme[p->index-1][j], outfile);
2686  }
2687  putc('\n', outfile);
2688}  /* drawline3 */
2689
2690
2691void standev(long chars, long numtrees, long minwhich, double minsteps,
2692                        double *nsteps, long **fsteps, longer seed)
2693{  /* compute and write standard deviation of user trees */
2694   /* used in pars */
2695  long i, j, k;
2696  double wt, sumw, sum, sum2, sd;
2697  double temp;
2698  double **covar, *P, *f;
2699
2700#define SAMPLES 1000
2701#define MAXSHIMOTREES 1000
2702/* ????? if numtrees too big for Shimo, truncate */
2703  if (numtrees == 2) {
2704    fprintf(outfile, "Kishino-Hasegawa-Templeton test\n\n");
2705    fprintf(outfile, "Tree    Steps   Diff Steps   Its S.D.");
2706    fprintf(outfile, "   Significantly worse?\n\n");
2707    which = 1;
2708    while (which <= numtrees) {
2709      fprintf(outfile, "%3ld%10.1f", which, nsteps[which - 1] / 10);
2710      if (minwhich == which)
2711        fprintf(outfile, "  <------ best\n");
2712      else {
2713        sumw = 0.0;
2714        sum = 0.0;
2715        sum2 = 0.0;
2716        for (i = 0; i < chars; i++) {
2717          if (weight[i] > 0) {
2718            wt = weight[i] / 10.0;
2719            sumw += wt;
2720            temp = (fsteps[which - 1][i] - fsteps[minwhich - 1][i]) / 10.0;
2721            sum += temp;
2722            sum2 += temp * temp / wt;
2723          }
2724        }
2725        temp = sum / sumw;
2726        sd = sqrt(sumw / (sumw - 1.0) * (sum2 - temp * temp));
2727        fprintf(outfile, "%10.1f%12.4f",
2728                (nsteps[which - 1] - minsteps) / 10, sd);
2729        if (sum > 1.95996 * sd)
2730          fprintf(outfile, "           Yes\n");
2731        else
2732          fprintf(outfile, "           No\n");
2733      }
2734      which++;
2735    }
2736    fprintf(outfile, "\n\n");
2737  } else {           /* Shimodaira-Hasegawa test using normal approximation */
2738    fprintf(outfile, "Shimodaira-Hasegawa test\n\n");
2739    covar = (double **)Malloc(numtrees*sizeof(double *)); 
2740    for (i = 0; i < numtrees; i++)
2741      covar[i] = (double *)Malloc(numtrees*sizeof(double)); 
2742    sumw = 0.0;
2743    for (i = 0; i < chars; i++)
2744      sumw += weight[i];
2745    for (i = 0; i < numtrees; i++) {        /* compute covariances of trees */
2746      sum = nsteps[i]/(10.0*sumw);
2747      for (j = 0; j <=i; j++) {
2748        sum2 = nsteps[j]/(10.0*sumw);
2749        temp = 0.0;
2750        for (k = 0; k < chars; k++) {
2751          wt = weight[k]/10.0;
2752          if (weight[k] > 0) {
2753            temp = temp + wt*(fsteps[i][k]/(10.0*wt)-sum)
2754                            *(fsteps[j][k]/(10.0*wt)-sum2);
2755          }
2756        }
2757        covar[i][j] = temp;
2758        if (i != j)
2759          covar[j][i] = temp;
2760      }
2761    }
2762    for (i = 0; i < numtrees; i++) { /* in-place Cholesky decomposition
2763                                        of trees x trees covariance matrix */
2764      sum = 0.0;
2765      for (j = 0; j <= i-1; j++)
2766        sum = sum + covar[i][j] * covar[i][j];
2767      if (covar[i][i]-sum <= 0.0)
2768        temp = 0.0;
2769      else
2770        temp = sqrt(covar[i][i] - sum);
2771      covar[i][i] = temp;
2772      for (j = i+1; j < numtrees; j++) {
2773        sum = 0.0;
2774        for (k = 0; k < i; k++)
2775          sum = sum + covar[i][k] * covar[j][k];
2776        if (fabs(temp) < 1.0E-23)
2777          covar[j][i] = 0.0;
2778        else
2779          covar[j][i] = (covar[j][i] - sum)/temp;
2780      }
2781    }
2782    f = (double *)Malloc(numtrees*sizeof(double)); /* vector of P's of trees */
2783    P = (double *)Malloc(numtrees*sizeof(double)); /* vector of P's of trees */
2784    for (i = 0; i < numtrees; i++)
2785      P[i] = 0.0;
2786    sum2 = nsteps[0]/10.0;               /* sum2 will be smallest # of steps */
2787    for (i = 1; i < numtrees; i++)
2788      if (sum2 > nsteps[i]/10.0)
2789        sum2 = nsteps[i]/10.0;
2790    for (i = 1; i < SAMPLES; i++) {           /* loop over resampled trees */
2791      for (j = 0; j < numtrees; j++) {        /* draw vectors */
2792        sum = 0.0;
2793        for (k = 0; k <= j; k++)
2794          sum += normrand(seed)*covar[j][k];
2795        f[j] = sum;
2796      }
2797      sum = f[1];
2798      for (j = 1; j < numtrees; j++)          /* get min of vector */
2799        if (f[j] < sum)
2800          sum = f[j];
2801      for (j = 0; j < numtrees; j++)          /* accumulate P's */
2802        if (nsteps[j]/10.0-sum2 < f[j] - sum)
2803          P[j] += 1.0/SAMPLES;
2804    }
2805    fprintf(outfile, "Tree    Steps   Diff Steps   P value");
2806    fprintf(outfile, "   Significantly worse?\n\n");
2807    for (i = 0; i < numtrees; i++) {
2808      fprintf(outfile, "%3ld%10.1f", i+1, nsteps[i]/10);
2809      if ((minwhich-1) == i)
2810        fprintf(outfile, "  <------ best\n");
2811      else {
2812        fprintf(outfile, "  %9.1f %10.3f", nsteps[i]/10.0-sum2, P[i]);
2813        if (P[i] < 0.05)
2814          fprintf(outfile, "           Yes\n");
2815        else
2816          fprintf(outfile, "           No\n");
2817      }
2818    }
2819  fprintf(outfile, "\n");
2820  free(P);             /* free the variables we Malloc'ed */
2821  free(f);
2822  for (i = 0; i < numtrees; i++)
2823    free(covar[i]);
2824  free(covar);
2825  }
2826}  /* standev */
2827
2828
2829void freetip(node *anode)
2830{
2831  /* used in pars */
2832
2833  free(anode->numsteps);
2834  free(anode->oldnumsteps);
2835  free(anode->discbase);
2836  free(anode->olddiscbase);
2837}  /* freetip */
2838
2839
2840void freenontip(node *anode)
2841{
2842  /* used in pars */
2843  free(anode->numsteps);
2844  free(anode->oldnumsteps);
2845  free(anode->discbase);
2846  free(anode->olddiscbase);
2847  free(anode->discnumnuc);
2848}  /* freenontip */
2849
2850
2851void freenodes(long local_nonodes, pointarray treenode)
2852{
2853  /* used in pars */
2854  long i;
2855  node *p;
2856
2857  for (i = 0; i < spp; i++)
2858    freetip(treenode[i]);
2859  for (i = spp; i < local_nonodes; i++) {
2860    if (treenode[i] != NULL) {
2861      p = treenode[i]->next;
2862      do {
2863        freenontip(p);
2864        p = p->next;
2865      } while (p != treenode[i]);
2866      freenontip(p);
2867    }
2868  }
2869}  /* freenodes */
2870
2871
2872void freenode(node **anode)
2873{
2874  /* used in pars */
2875
2876  freenontip(*anode);
2877  free(*anode);
2878}  /* freenode */
2879
2880
2881void freetree(long local_nonodes, pointarray treenode)
2882{
2883  /* used in pars */
2884  long i;
2885  node *p, *q;
2886
2887  for (i = 0; i < spp; i++)
2888    free(treenode[i]);
2889  for (i = spp; i < local_nonodes; i++) {
2890    if (treenode[i] != NULL) {
2891      p = treenode[i]->next;
2892      do {
2893        q = p->next;
2894        free(p);
2895        p = q;
2896      } while (p != treenode[i]);
2897      free(p);
2898    }
2899  }
2900  free(treenode);
2901}  /* freetree */
2902
2903
2904void freegarbage(gbases **garbage)
2905{
2906  /* used in pars */
2907  gbases *p;
2908
2909  while (*garbage) {
2910    p = *garbage;
2911    *garbage = (*garbage)->next;
2912    free(p->discbase);
2913    free(p);
2914  }
2915}  /* freegarbage */
2916
2917
2918void freegrbg(node **grbg)
2919{
2920  /* used in pars */
2921  node *p;
2922
2923  while (*grbg) {
2924    p = *grbg;
2925    *grbg = (*grbg)->next;
2926    freenontip(p);
2927    free(p);
2928  }
2929} /*freegrbg */
Note: See TracBrowser for help on using the repository browser.