source: branches/profile/GDE/PHYLIP/dnaml.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: 64.9 KB
Line 
1
2#include "phylip.h"
3#include "seq.h"
4
5/* version 3.6. (c) Copyright 1993-2002 by the University of Washington.
6   Written by Joseph Felsenstein, Akiko Fuseki, Sean Lamont, Andrew Keeffe,
7   Dan Fineman, and Patrick Colacurcio.
8   Permission is granted to copy and use this program provided no fee is
9   charged for it and provided that this copyright notice is not removed. */
10
11typedef struct valrec {
12  double rat, ratxi, ratxv, orig_zz, z1, y1, z1zz, z1yy, xiz1, xiy1xv;
13  double *ww, *zz, *wwzz, *vvzz; 
14} valrec;
15
16typedef long vall[maxcategs];
17typedef double contribarr[maxcategs];
18
19#ifndef OLDC
20/* function prototypes */
21void   dnamlcopy(tree *, tree *, long, long);
22void   getoptions(void);
23void   allocrest(void);
24void   doinit(void);
25void   inputoptions(void);
26void   makeweights(void);
27void   getinput(void);
28void   inittable_for_usertree(FILE *);
29void   inittable(void);
30double evaluate(node *, boolean);
31
32void   alloc_nvd (long, nuview_data *);
33void   free_nvd (nuview_data *);
34void   nuview(node *);
35void   slopecurv(node *, double, double *, double *, double *);
36void   makenewv(node *);
37void   update(node *);
38void   smooth(node *);
39void   insert_(node *, node *, boolean);
40void   dnaml_re_move(node **, node **);
41void   buildnewtip(long, tree *);
42
43void   buildsimpletree(tree *);
44void   addtraverse(node *, node *, boolean);
45void   rearrange(node *, node *);
46void   initdnamlnode(node **, node **, node *, long, long, long *, long *,
47                 initops, pointarray, pointarray, Char *, Char *, FILE *);
48void   dnaml_coordinates(node *, double, long *, double *);
49void   dnaml_printree(void);
50void   sigma(node *, double *, double *, double *);
51void   describe(node *);
52void   reconstr(node *, long);
53void   rectrav(node *, long, long);
54void   summarize(void);
55void   dnaml_treeout(node *);
56void   inittravtree(node *);
57void   treevaluate(void);
58void   maketree(void);
59void   clean_up(void);
60void   reallocsites(void);
61/* function prototypes */
62#endif
63
64
65extern sequence y;
66
67double fracchange;
68long rcategs;
69boolean haslengths;
70
71Char infilename[FNMLNGTH], outfilename[FNMLNGTH], intreename[FNMLNGTH], outtreename[FNMLNGTH],
72     catfilename[FNMLNGTH], weightfilename[FNMLNGTH];
73double *rate, *rrate, *probcat;
74long nonodes2, sites, weightsum, categs, datasets, ith, njumble, jumb;
75long parens;
76boolean  freqsfrom, global, jumble, weights, trout, usertree, reconsider,
77  ctgry, rctgry, auto_, hypstate, ttr, progress, mulsets, justwts,
78  firstset, improve, smoothit, polishing, lngths, gama, invar;
79tree curtree, bestree, bestree2, priortree;
80node *qwhere, *grbg;
81double xi, xv, ttratio, ttratio0, freqa, freqc, freqg, freqt, freqr,
82  freqy, freqar, freqcy, freqgr, freqty, 
83  cv, alpha, lambda, invarfrac, bestyet;
84long *enterorder, inseed, inseed0;
85steptr aliasweight;
86contribarr *contribution, like, nulike, clai;
87double **term, **slopeterm, **curveterm;
88longer seed;
89Char* progname;
90char basechar[16]="acmgrsvtwyhkdbn";
91
92/* Local variables for maketree, propagated globally for c version: */
93long k, nextsp, numtrees, maxwhich, mx, mx0, mx1;
94double dummy, maxlogl;
95boolean succeeded, smoothed;
96double **l0gf;
97double *l0gl;
98valrec ***tbl;
99Char ch, ch2;
100long col;
101vall *mp=NULL;
102
103
104void dnamlcopy(tree *a, tree *b, long nonodes, long categs)
105{
106  /* used in dnaml & dnamlk */
107  long i, j;
108  node *p, *q;
109
110  for (i = 0; i < spp; i++) {
111    copynode(a->nodep[i], b->nodep[i], categs);
112    if (a->nodep[i]->back) {
113      if (a->nodep[i]->back == a->nodep[a->nodep[i]->back->index - 1])
114        b->nodep[i]->back = b->nodep[a->nodep[i]->back->index - 1];
115      else if (a->nodep[i]->back == a->nodep[a->nodep[i]->back->index - 1]->next)
116        b->nodep[i]->back = b->nodep[a->nodep[i]->back->index - 1]->next;
117      else
118        b->nodep[i]->back = b->nodep[a->nodep[i]->back->index - 1]->next->next;
119    }
120    else b->nodep[i]->back = NULL;
121  }
122  for (i = spp; i < nonodes; i++) {
123    p = a->nodep[i];
124    q = b->nodep[i];
125    for (j = 1; j <= 3; j++) {
126      copynode(p, q, categs);
127      if (p->back) {
128        if (p->back == a->nodep[p->back->index - 1])
129          q->back = b->nodep[p->back->index - 1];
130        else if (p->back == a->nodep[p->back->index - 1]->next)
131          q->back = b->nodep[p->back->index - 1]->next;
132        else
133          q->back = b->nodep[p->back->index - 1]->next->next;
134      }
135      else
136        q->back = NULL;
137      p = p->next;
138      q = q->next;
139    }
140  }
141  b->likelihood = a->likelihood;
142  b->start = a->start;               /* start used in dnaml only */
143  b->root = a->root;                 /* root used in dnamlk only */
144}  /* dnamlcopy plc*/
145
146
147void getoptions()
148{
149  /* interactively set options */
150  long i, loopcount, loopcount2;
151  Char ch;
152  boolean didchangecat, didchangercat;
153  double probsum;
154
155  fprintf(outfile, "\nNucleic acid sequence Maximum Likelihood");
156  fprintf(outfile, " method, version %s\n\n",VERSION);
157  putchar('\n');
158  ctgry = false;
159  didchangecat = false;
160  rctgry = false;
161  didchangercat = false;
162  categs = 1;
163  rcategs = 1;
164  auto_ = false;
165  freqsfrom = true;
166  gama = false;
167  global = false;
168  hypstate = false;
169  improve = false;
170  invar = false;
171  jumble = false;
172  njumble = 1;
173  lngths = false;
174  lambda = 1.0;
175  outgrno = 1;
176  outgropt = false;
177  reconsider = false;
178  trout = true;
179  ttratio = 2.0;
180  ttr = false;
181  usertree = false;
182  weights = false;
183  printdata = false;
184  progress = true;
185  treeprint = true;
186  interleaved = true;
187  loopcount = 0;
188  for (;;){
189    cleerhome();
190    printf("Nucleic acid sequence Maximum Likelihood");
191    printf(" method, version %s\n\n",VERSION);
192    printf("Settings for this run:\n");
193    printf("  U                 Search for best tree?  %s\n",
194           (usertree ? "No, use user trees in input file" : "Yes"));
195    if (usertree) {
196      printf("  L          Use lengths from user trees?  %s\n",
197               (lngths ? "Yes" : "No"));
198    }
199    printf("  T        Transition/transversion ratio:%8.4f\n",
200           (ttr ? ttratio : 2.0));
201    printf("  F       Use empirical base frequencies?  %s\n",
202           (freqsfrom ? "Yes" : "No"));
203    printf("  C                One category of sites?");
204    if (!ctgry || categs == 1)
205      printf("  Yes\n");
206    else
207      printf("  %ld categories of sites\n", categs);
208    printf("  R           Rate variation among sites?");
209    if (!rctgry)
210      printf("  constant rate\n");
211    else {
212      if (gama)
213        printf("  Gamma distributed rates\n");
214      else {
215        if (invar)
216          printf("  Gamma+Invariant sites\n");
217        else
218          printf("  user-defined HMM of rates\n");
219      }
220      printf("  A   Rates at adjacent sites correlated?");
221      if (!auto_)
222        printf("  No, they are independent\n");
223      else
224        printf("  Yes, mean block length =%6.1f\n", 1.0 / lambda);
225    }
226    printf("  W                       Sites weighted?  %s\n",
227           (weights ? "Yes" : "No"));
228    if ((!usertree) || reconsider) {
229      printf("  S        Speedier but rougher analysis?  %s\n",
230             (improve ? "No, not rough" : "Yes"));
231      printf("  G                Global rearrangements?  %s\n",
232             (global ? "Yes" : "No"));
233    }
234    if (!usertree) {
235      printf("  J   Randomize input order of sequences?");
236      if (jumble)
237        printf("  Yes (seed =%8ld,%3ld times)\n", inseed0, njumble);
238      else
239        printf("  No. Use input order\n");
240    }
241    else
242       printf("  V    Rearrange starting with user tree?  %s\n",
243             (reconsider ? "Yes" : "No"));
244    printf("  O                        Outgroup root?  %s%3ld\n",
245           (outgropt ? "Yes, at sequence number" :
246                       "No, use as outgroup species"),outgrno);
247    printf("  M           Analyze multiple data sets?");
248    if (mulsets)
249      printf("  Yes, %2ld %s\n", datasets,
250               (justwts ? "sets of weights" : "data sets"));
251    else
252      printf("  No\n");
253    printf("  I          Input sequences interleaved?  %s\n",
254           (interleaved ? "Yes" : "No, sequential"));
255    printf("  0   Terminal type (IBM PC, ANSI, none)?  %s\n",
256           (ibmpc ? "IBM PC" : ansi  ? "ANSI" : "(none)"));
257    printf("  1    Print out the data at start of run  %s\n",
258           (printdata ? "Yes" : "No"));
259    printf("  2  Print indications of progress of run  %s\n",
260           (progress ? "Yes" : "No"));
261    printf("  3                        Print out tree  %s\n",
262           (treeprint ? "Yes" : "No"));
263    printf("  4       Write out trees onto tree file?  %s\n",
264           (trout ? "Yes" : "No"));
265    printf("  5   Reconstruct hypothetical sequences?  %s\n",
266           (hypstate ? "Yes" : "No"));
267    printf("\n  Y to accept these or type the letter for one to change\n");
268#ifdef WIN32
269    phyFillScreenColor();
270#endif
271    scanf("%c%*[^\n]", &ch);
272    getchar();
273    if (ch == '\n')
274      ch = ' ';
275    uppercase(&ch);
276    if (ch == 'Y')
277      break;
278    if (strchr("ULTFCRAWSGJVOMI012345",ch) != NULL){
279      switch (ch) {
280
281      case 'F':
282        freqsfrom = !freqsfrom;
283        if (!freqsfrom) {
284          initfreqs(&freqa, &freqc, &freqg, &freqt);
285        }
286        break;
287       
288      case 'C':
289        ctgry = !ctgry;
290        if (ctgry) {
291          printf("\nSitewise user-assigned categories:\n\n");
292          initcatn(&categs);
293          if (rate){
294            free(rate);
295          }
296          rate    = (double *) Malloc(categs * sizeof(double));
297          didchangecat = true;
298          initcategs(categs, rate);
299        }
300        break;
301
302      case 'R':
303        if (!rctgry) {
304          rctgry = true;
305          gama = true;
306        } else {
307          if (gama) {
308            gama = false;
309            invar = true;
310          } else {
311            if (invar)
312              invar = false;
313            else
314              rctgry = false;
315          }
316        }
317        break;
318       
319      case 'A':
320        auto_ = !auto_;
321        if (auto_)
322          initlambda(&lambda);
323        break;
324       
325      case 'W':
326        weights = !weights;
327        break;
328
329      case 'S':
330        improve = !improve;
331        break;
332
333      case 'G':
334        global = !global;
335        break;
336       
337      case 'J':
338        jumble = !jumble;
339        if (jumble)
340          initjumble(&inseed, &inseed0, seed, &njumble);
341        else njumble = 1;
342        break;
343       
344      case 'L':
345        lngths = !lngths;
346        break;
347       
348      case 'O':
349        outgropt = !outgropt;
350        if (outgropt)
351          initoutgroup(&outgrno, spp);
352        break;
353       
354      case 'T':
355        ttr = !ttr;
356        if (ttr) {
357          initratio(&ttratio);
358        }
359        break;
360       
361      case 'U':
362        usertree = !usertree;
363        break;
364
365      case 'V':
366        reconsider = !reconsider;
367        break;
368
369      case 'M':
370        mulsets = !mulsets;
371        if (mulsets) {
372          printf("Multiple data sets or multiple weights?");
373          loopcount2 = 0;
374          do {
375            printf(" (type D or W)\n");
376            scanf("%c%*[^\n]", &ch2);
377            getchar();
378            if (ch2 == '\n')
379              ch2 = ' ';
380            uppercase(&ch2);
381            countup(&loopcount2, 10);
382          } while ((ch2 != 'W') && (ch2 != 'D'));
383          justwts = (ch2 == 'W');
384          if (justwts)
385            justweights(&datasets);
386          else
387            initdatasets(&datasets);
388          if (!jumble) {
389            jumble = true;
390            initjumble(&inseed, &inseed0, seed, &njumble);
391          }
392        }
393        break;
394       
395      case 'I':
396        interleaved = !interleaved;
397        break;
398       
399      case '0':
400        initterminal(&ibmpc, &ansi);
401        break;
402       
403      case '1':
404        printdata = !printdata;
405        break;
406       
407      case '2':
408        progress = !progress;
409        break;
410       
411      case '3':
412        treeprint = !treeprint;
413        break;
414       
415      case '4':
416        trout = !trout;
417        break;
418
419      case '5':
420        hypstate = !hypstate;
421        break;
422      }
423    } else
424      printf("Not a possible option!\n");
425    countup(&loopcount, 100);
426  }
427  if (gama || invar) {
428    loopcount = 0;
429    do {
430      printf(
431"\nCoefficient of variation of substitution rate among sites (must be positive)\n");
432      printf(
433  " In gamma distribution parameters, this is 1/(square root of alpha)\n");
434#ifdef WIN32
435      phyFillScreenColor();
436#endif
437      scanf("%lf%*[^\n]", &cv);
438      getchar();
439      countup(&loopcount, 10);
440    } while (cv <= 0.0);
441    alpha = 1.0 / (cv * cv);
442  }
443  if (!rctgry)
444    auto_ = false;
445  if (rctgry) {
446    printf("\nRates in HMM");
447    if (invar)
448      printf(" (including one for invariant sites)");
449    printf(":\n");
450    initcatn(&rcategs);
451    if (probcat){
452      free(probcat);
453      free(rrate);
454    }
455    probcat = (double *) Malloc(rcategs * sizeof(double));
456    rrate   = (double *) Malloc(rcategs * sizeof(double));
457    didchangercat = true;
458    if (gama)
459      initgammacat(rcategs, alpha, rrate, probcat); 
460    else {
461      if (invar) {
462        loopcount = 0;
463        do {
464          printf("Fraction of invariant sites?\n");
465          scanf("%lf%*[^\n]", &invarfrac);
466          getchar();
467          countup (&loopcount, 10);
468        } while ((invarfrac <= 0.0) || (invarfrac >= 1.0));
469        initgammacat(rcategs-1, alpha, rrate, probcat); 
470        for (i = 0; i < rcategs-1; i++)
471          probcat[i] = probcat[i]*(1.0-invarfrac);
472        probcat[rcategs-1] = invarfrac;
473        rrate[rcategs-1] = 0.0;
474      } else {
475        initcategs(rcategs, rrate);
476        initprobcat(rcategs, &probsum, probcat);
477      }
478    }
479  }
480  if (!didchangercat){
481    rrate      = (double *) Malloc(rcategs*sizeof(double));
482    probcat    = (double *) Malloc(rcategs*sizeof(double));
483    rrate[0]   = 1.0;
484    probcat[0] = 1.0;
485  }
486  if (!didchangecat){
487    rate       = (double *) Malloc(categs*sizeof(double));
488    rate[0]    = 1.0;
489  }
490}  /* getoptions */
491
492
493void reallocsites(void)
494{
495  long i;
496
497  for (i=0; i < spp; i++) {
498    free(y[i]);
499    y[i] = (Char *) Malloc(sites*sizeof(Char));
500  }
501  free(category);
502  free(weight);
503  free(alias);
504  free(ally);
505  free(location);
506  free(aliasweight);
507
508  category    = (long *) Malloc(sites*sizeof(long));
509  weight      = (long *) Malloc(sites*sizeof(long));
510  alias       = (long *) Malloc(sites*sizeof(long));
511  ally        = (long *) Malloc(sites*sizeof(long));
512  location    = (long *) Malloc(sites*sizeof(long));
513  aliasweight = (long *) Malloc(sites*sizeof(long));
514
515}
516
517void allocrest()
518{
519  long i;
520
521  y = (Char **) Malloc(spp*sizeof(Char *));
522  for (i = 0; i < spp; i++)
523    y[i] = (Char *) Malloc(sites*sizeof(Char));
524  nayme       = (naym *) Malloc(spp*sizeof(naym));;
525  enterorder  = (long *) Malloc(spp*sizeof(long));
526  category    = (long *) Malloc(sites*sizeof(long));
527  weight      = (long *) Malloc(sites*sizeof(long));
528  alias       = (long *) Malloc(sites*sizeof(long));
529  ally        = (long *) Malloc(sites*sizeof(long));
530  location    = (long *) Malloc(sites*sizeof(long));
531  aliasweight = (long *) Malloc(sites*sizeof(long));
532}  /* allocrest */
533
534
535void doinit()
536{ /* initializes variables */
537
538  inputnumbers(&spp, &sites, &nonodes2, 2);
539  getoptions();
540  if (printdata)
541    fprintf(outfile, "%2ld species, %3ld  sites\n", spp, sites);
542  alloctree(&curtree.nodep, nonodes2, usertree);
543  allocrest();
544  if (usertree && !reconsider)
545    return;
546  alloctree(&bestree.nodep, nonodes2, 0);
547  alloctree(&priortree.nodep, nonodes2, 0);
548  if (njumble <= 1)
549    return;
550  alloctree(&bestree2.nodep, nonodes2, 0);
551}  /* doinit */
552
553
554void inputoptions()
555{
556  long i;
557
558  if (!firstset && !justwts) {
559    samenumsp(&sites, ith);
560    reallocsites();
561  }
562  for (i = 0; i < sites; i++)
563    category[i] = 1;
564  for (i = 0; i < sites; i++)
565    weight[i] = 1;
566 
567  if (justwts || weights)
568    inputweights(sites, weight, &weights);
569  weightsum = 0;
570  for (i = 0; i < sites; i++)
571    weightsum += weight[i];
572  if (ctgry && categs > 1) {
573    inputcategs(0, sites, category, categs, "DnaML");
574    if (printdata)
575      printcategs(outfile, sites, category, "Site categories");
576  }
577  if (weights && printdata)
578    printweights(outfile, 0, sites, weight, "Sites");
579}  /* inputoptions */
580
581
582void makeweights()
583{
584  /* make up weights vector to avoid duplicate computations */
585  long i;
586
587  for (i = 1; i <= sites; i++) {
588    alias[i - 1] = i;
589    ally[i - 1] = 0;
590    aliasweight[i - 1] = weight[i - 1];
591    location[i - 1] = 0;
592  }
593  sitesort2   (sites, aliasweight);
594  sitecombine2(sites, aliasweight);
595  sitescrunch2(sites, 1, 2, aliasweight);
596  for (i = 1; i <= sites; i++) {
597    if (aliasweight[i - 1] > 0)
598      endsite = i;
599  }
600  for (i = 1; i <= endsite; i++) {
601    location[alias[i - 1] - 1] = i;
602    ally[alias[i - 1] - 1] = alias[i - 1];
603  }
604  term = (double **) Malloc( endsite * sizeof(double *));
605  for (i = 0; i < endsite; i++)
606     term[i] = (double *) Malloc( rcategs * sizeof(double));
607  slopeterm = (double **) Malloc( endsite * sizeof(double *));
608  for (i = 0; i < endsite; i++)
609     slopeterm[i] = (double *) Malloc( rcategs * sizeof(double));
610  curveterm = (double **) Malloc(endsite * sizeof(double *));
611  for (i = 0; i < endsite; i++)
612     curveterm[i] = (double *) Malloc( rcategs * sizeof(double));
613  mp = (vall *) Malloc( sites*sizeof(vall));
614  contribution = (contribarr *) Malloc( endsite*sizeof(contribarr));
615}  /* makeweights */
616
617
618void getinput()
619{
620  /* reads the input data */
621  inputoptions();
622  if (!freqsfrom)
623    getbasefreqs(freqa, freqc, freqg, freqt, &freqr, &freqy, &freqar, &freqcy,
624                 &freqgr, &freqty, &ttratio, &xi, &xv, &fracchange,
625                 freqsfrom, true);
626  if (!justwts || firstset)
627    inputdata(sites);
628  makeweights();
629  setuptree2(curtree);
630  if (!usertree || reconsider) {
631    setuptree2(bestree);
632    setuptree2(priortree);
633    if (njumble > 1)
634      setuptree2(bestree2);
635  }
636  allocx(nonodes2, rcategs, curtree.nodep, usertree);
637  if (!usertree || reconsider) {
638    allocx(nonodes2, rcategs, bestree.nodep, 0);
639    allocx(nonodes2, rcategs, priortree.nodep, 0);
640    if (njumble > 1)
641      allocx(nonodes2, rcategs, bestree2.nodep, 0);
642  }
643  makevalues2(rcategs, curtree.nodep, endsite, spp, y, alias); 
644  if (freqsfrom) {
645    empiricalfreqs(&freqa, &freqc, &freqg, &freqt, aliasweight,
646                    curtree.nodep);
647    getbasefreqs(freqa, freqc, freqg, freqt, &freqr, &freqy, &freqar, &freqcy,
648                 &freqgr, &freqty, &ttratio, &xi, &xv, &fracchange,
649                 freqsfrom, true);
650  }
651  if (!justwts || firstset)
652    fprintf(outfile, "\nTransition/transversion ratio = %10.6f\n\n", ttratio);
653}  /* getinput */
654
655
656void inittable_for_usertree(FILE *intree)
657{
658  /* If there's a user tree, then the ww/zz/wwzz/vvzz elements need
659     to be allocated appropriately. */
660  long num_comma;
661  long i, j;
662
663  /* First, figure out the largest possible furcation, i.e. the number
664     of commas plus one */
665  countcomma(&intree, &num_comma);
666  num_comma++;
667 
668  for (i = 0; i < rcategs; i++) {
669    for (j = 0; j < categs; j++) {
670      /* Free the stuff allocated assuming bifurcations */
671      free (tbl[i][j]->ww);
672      free (tbl[i][j]->zz);
673      free (tbl[i][j]->wwzz);
674      free (tbl[i][j]->vvzz);
675
676      /* Then allocate for worst-case multifurcations */
677      tbl[i][j]->ww   = (double *) Malloc( num_comma * sizeof (double));
678      tbl[i][j]->zz   = (double *) Malloc( num_comma * sizeof (double));
679      tbl[i][j]->wwzz = (double *) Malloc( num_comma * sizeof (double));
680      tbl[i][j]->vvzz = (double *) Malloc( num_comma * sizeof (double));
681    }
682  }
683}  /* inittable_for_usertree */
684
685void inittable()
686{
687  /* Define a lookup table. Precompute values and print them out in tables */
688  long i, j;
689  double sumrates;
690 
691  tbl = (valrec ***) Malloc(rcategs * sizeof(valrec **));
692  for (i = 0; i < rcategs; i++) {
693    tbl[i] = (valrec **) Malloc(categs*sizeof(valrec *));
694    for (j = 0; j < categs; j++)
695      tbl[i][j] = (valrec *) Malloc(sizeof(valrec));
696  }
697
698  for (i = 0; i < rcategs; i++) {
699    for (j = 0; j < categs; j++) {
700      tbl[i][j]->rat = rrate[i]*rate[j];
701      tbl[i][j]->ratxi = tbl[i][j]->rat * xi;
702      tbl[i][j]->ratxv = tbl[i][j]->rat * xv;
703
704      /* Allocate assuming bifurcations, will be changed later if
705         neccesarry (i.e. there's a user tree) */
706      tbl[i][j]->ww   = (double *) Malloc( 2 * sizeof (double));
707      tbl[i][j]->zz   = (double *) Malloc( 2 * sizeof (double));
708      tbl[i][j]->wwzz = (double *) Malloc( 2 * sizeof (double));
709      tbl[i][j]->vvzz = (double *) Malloc( 2 * sizeof (double));
710    }
711  }
712  if (!lngths) {            /* restandardize rates */
713    sumrates = 0.0;
714    for (i = 0; i < endsite; i++) {
715      for (j = 0; j < rcategs; j++)
716        sumrates += aliasweight[i] * probcat[j] 
717          * tbl[j][category[alias[i] - 1] - 1]->rat;
718    }
719    sumrates /= (double)sites;
720    for (i = 0; i < rcategs; i++)
721      for (j = 0; j < categs; j++) {
722        tbl[i][j]->rat /= sumrates;
723        tbl[i][j]->ratxi /= sumrates;
724        tbl[i][j]->ratxv /= sumrates;
725      }
726  }
727  if (rcategs > 1) {
728    if (gama) {
729      fprintf(outfile, "\nDiscrete approximation to gamma distributed rates\n");
730      fprintf(outfile,
731      " Coefficient of variation of rates = %f  (alpha = %f)\n",
732        cv, alpha);
733    }
734    fprintf(outfile, "\nState in HMM    Rate of change    Probability\n\n");
735    for (i = 0; i < rcategs; i++)
736      if (probcat[i] < 0.0001)
737        fprintf(outfile, "%9ld%16.3f%20.6f\n", i+1, rrate[i], probcat[i]);
738      else if (probcat[i] < 0.001)
739          fprintf(outfile, "%9ld%16.3f%19.5f\n", i+1, rrate[i], probcat[i]);
740        else if (probcat[i] < 0.01)
741            fprintf(outfile, "%9ld%16.3f%18.4f\n", i+1, rrate[i], probcat[i]);
742          else
743            fprintf(outfile, "%9ld%16.3f%17.3f\n", i+1, rrate[i], probcat[i]);
744    putc('\n', outfile);
745    if (auto_)
746      fprintf(outfile,
747     "Expected length of a patch of sites having the same rate = %8.3f\n",
748             1/lambda);
749    putc('\n', outfile);
750  }
751  if (categs > 1) {
752    fprintf(outfile, "\nSite category   Rate of change\n\n");
753    for (i = 0; i < categs; i++)
754      fprintf(outfile, "%9ld%16.3f\n", i+1, rate[i]);
755  }
756  if ((rcategs  > 1) || (categs >> 1))
757    fprintf(outfile, "\n\n");
758}  /* inittable */
759
760
761double evaluate(node *p, boolean saveit)
762{
763  contribarr tterm;
764  double sum, sum2, sumc, y, lz, y1, z1zz, z1yy, prod12, prod1, prod2, prod3,
765         sumterm, lterm;
766  long i, j, k, lai;
767  node *q;
768  sitelike x1, x2;
769
770  sum = 0.0;
771  q = p->back;
772  y = p->v;
773  lz = -y;
774  for (i = 0; i < rcategs; i++)
775    for (j = 0; j < categs; j++) {
776    tbl[i][j]->orig_zz = exp(tbl[i][j]->ratxi * lz);
777    tbl[i][j]->z1 = exp(tbl[i][j]->ratxv * lz);
778    tbl[i][j]->z1zz = tbl[i][j]->z1 * tbl[i][j]->orig_zz;
779    tbl[i][j]->z1yy = tbl[i][j]->z1 - tbl[i][j]->z1zz;
780  }
781  for (i = 0; i < endsite; i++) {
782    k = category[alias[i]-1] - 1;
783    for (j = 0; j < rcategs; j++) {
784      if (y > 0.0) {
785        y1 = 1.0 - tbl[j][k]->z1;
786        z1zz = tbl[j][k]->z1zz;
787        z1yy = tbl[j][k]->z1yy;
788      } else {
789        y1 = 0.0;
790        z1zz = 1.0;
791        z1yy = 0.0;
792      }
793      memcpy(x1, p->x[i][j], sizeof(sitelike));
794      prod1 = freqa * x1[0] + freqc * x1[(long)C - (long)A] +
795             freqg * x1[(long)G - (long)A] + freqt * x1[(long)T - (long)A];
796      memcpy(x2, q->x[i][j], sizeof(sitelike));
797      prod2 = freqa * x2[0] + freqc * x2[(long)C - (long)A] +
798             freqg * x2[(long)G - (long)A] + freqt * x2[(long)T - (long)A];
799      prod3 = (x1[0] * freqa + x1[(long)G - (long)A] * freqg) *
800              (x2[0] * freqar + x2[(long)G - (long)A] * freqgr) +
801          (x1[(long)C - (long)A] * freqc + x1[(long)T - (long)A] * freqt) *
802         (x2[(long)C - (long)A] * freqcy + x2[(long)T - (long)A] * freqty);
803      prod12 = freqa * x1[0] * x2[0] +
804               freqc * x1[(long)C - (long)A] * x2[(long)C - (long)A] +
805               freqg * x1[(long)G - (long)A] * x2[(long)G - (long)A] +
806               freqt * x1[(long)T - (long)A] * x2[(long)T - (long)A];
807      tterm[j] = z1zz * prod12 + z1yy * prod3 + y1 * prod1 * prod2;
808    }
809    sumterm = 0.0;
810    for (j = 0; j < rcategs; j++)
811      sumterm += probcat[j] * tterm[j];
812    lterm = log(sumterm);
813    for (j = 0; j < rcategs; j++)
814      clai[j] = tterm[j] / sumterm;
815    memcpy(contribution[i], clai, rcategs*sizeof(double));
816    if (saveit && !auto_ && usertree)
817      l0gf[which - 1][i] = lterm;
818    sum += aliasweight[i] * lterm;
819  }
820  for (j = 0; j < rcategs; j++)
821    like[j] = 1.0;
822  for (i = 0; i < sites; i++) {
823    sumc = 0.0;
824    for (k = 0; k < rcategs; k++)
825      sumc += probcat[k] * like[k];
826    sumc *= lambda;
827    if ((ally[i] > 0) && (location[ally[i]-1] > 0)) {
828      lai = location[ally[i] - 1];
829      memcpy(clai, contribution[lai - 1], rcategs*sizeof(double));
830      for (j = 0; j < rcategs; j++)
831        nulike[j] = ((1.0 - lambda) * like[j] + sumc) * clai[j];
832    } else {
833      for (j = 0; j < rcategs; j++)
834        nulike[j] = ((1.0 - lambda) * like[j] + sumc);
835    }
836    memcpy(like, nulike, rcategs*sizeof(double));
837  }
838  sum2 = 0.0;
839  for (i = 0; i < rcategs; i++)
840    sum2 += probcat[i] * like[i];
841  sum += log(sum2);
842  curtree.likelihood = sum;
843  if (!saveit || auto_ || !usertree)
844    return sum;
845  l0gl[which - 1] = sum;
846  if (which == 1) {
847    maxwhich = 1;
848    maxlogl = sum;
849    return sum;
850  }
851  if (sum > maxlogl) {
852    maxwhich = which;
853    maxlogl = sum;
854  }
855  return sum;
856}  /* evaluate */
857
858
859void alloc_nvd (long num_sibs, nuview_data *local_nvd)
860{
861  /* Allocate blocks of memory appropriate for the number of siblings
862     a given node has */
863  local_nvd->yy     = (double *) Malloc( num_sibs * sizeof (double));
864  local_nvd->wwzz   = (double *) Malloc( num_sibs * sizeof (double));
865  local_nvd->vvzz   = (double *) Malloc( num_sibs * sizeof (double));
866  local_nvd->vzsumr = (double *) Malloc( num_sibs * sizeof (double));
867  local_nvd->vzsumy = (double *) Malloc( num_sibs * sizeof (double));
868  local_nvd->sum    = (double *) Malloc( num_sibs * sizeof (double));
869  local_nvd->sumr   = (double *) Malloc( num_sibs * sizeof (double));
870  local_nvd->sumy   = (double *) Malloc( num_sibs * sizeof (double));
871  local_nvd->xx     = (sitelike *) Malloc( num_sibs * sizeof (sitelike));
872}  /* alloc_nvd */
873
874
875void free_nvd (nuview_data *local_nvd)
876{
877  /* The natural complement to the alloc version */
878  free (local_nvd->yy);
879  free (local_nvd->wwzz);
880  free (local_nvd->vvzz);
881  free (local_nvd->vzsumr);
882  free (local_nvd->vzsumy);
883  free (local_nvd->sum);
884  free (local_nvd->sumr);
885  free (local_nvd->sumy);
886  free (local_nvd->xx);
887}  /* free_nvd */
888
889
890void nuview(node *p)
891{
892  long i, j, k, num_sibs, sib_index;
893  nuview_data *local_nvd;
894  node *sib_ptr, *sib_back_ptr;
895  sitelike p_xx;
896  double lw;
897
898  /* Figure out how many siblings the current node has */
899  num_sibs    = count_sibs (p);
900
901  /* Recursive calls, should be called for all children */
902  sib_ptr = p;
903  for (i=0 ; i < num_sibs; i++) {
904    sib_ptr      = sib_ptr->next;
905    sib_back_ptr = sib_ptr->back;
906
907    if (!sib_back_ptr->tip &&
908        !sib_back_ptr->initialized)
909      nuview (sib_back_ptr);
910  }
911
912  /* Allocate the structure and blocks therein for variables used in
913     this function */
914  local_nvd = (nuview_data *) Malloc( sizeof (nuview_data));
915  alloc_nvd (num_sibs, local_nvd);
916
917
918  /* Loop 1: makes assignments to tbl based on some combination of
919     what's already in tbl and the children's value of v */
920  sib_ptr = p;
921  for (sib_index=0; sib_index < num_sibs; sib_index++) {
922    sib_ptr      = sib_ptr->next;
923    sib_back_ptr = sib_ptr->back;
924   
925    lw = - (sib_back_ptr->v);
926
927    for (i = 0; i < rcategs; i++)
928      for (j = 0; j < categs; j++) {
929        tbl[i][j]->ww[sib_index]   = exp(tbl[i][j]->ratxi * lw);
930        tbl[i][j]->zz[sib_index]   = exp(tbl[i][j]->ratxv * lw);
931        tbl[i][j]->wwzz[sib_index] = tbl[i][j]->ww[sib_index] * tbl[i][j]->zz[sib_index];
932        tbl[i][j]->vvzz[sib_index] = (1.0 - tbl[i][j]->ww[sib_index]) *
933          tbl[i][j]->zz[sib_index];
934      }
935  }
936 
937  /* Loop 2: */
938  for (i = 0; i < endsite; i++) {
939    k = category[alias[i]-1] - 1;
940    for (j = 0; j < rcategs; j++) {
941
942      /* Loop 2.1 */
943      sib_ptr = p;
944      for (sib_index=0; sib_index < num_sibs; sib_index++) {
945        sib_ptr         = sib_ptr->next;
946        sib_back_ptr    = sib_ptr->back;
947
948        local_nvd->wwzz[sib_index] = tbl[j][k]->wwzz[sib_index];
949        local_nvd->vvzz[sib_index] = tbl[j][k]->vvzz[sib_index];
950        local_nvd->yy[sib_index]   = 1.0 - tbl[j][k]->zz[sib_index];
951        memcpy(local_nvd->xx[sib_index],
952               sib_back_ptr->x[i][j],
953               sizeof(sitelike));
954      }
955
956      /* Loop 2.2 */
957      for (sib_index=0; sib_index < num_sibs; sib_index++) {
958        local_nvd->sum[sib_index] =
959          local_nvd->yy[sib_index] *
960          (freqa * local_nvd->xx[sib_index][(long)A] +
961           freqc * local_nvd->xx[sib_index][(long)C] +
962           freqg * local_nvd->xx[sib_index][(long)G] +
963           freqt * local_nvd->xx[sib_index][(long)T]);
964        local_nvd->sumr[sib_index] =
965          freqar * local_nvd->xx[sib_index][(long)A] +
966          freqgr * local_nvd->xx[sib_index][(long)G];
967        local_nvd->sumy[sib_index] =
968          freqcy * local_nvd->xx[sib_index][(long)C] +
969          freqty * local_nvd->xx[sib_index][(long)T];
970        local_nvd->vzsumr[sib_index] =
971          local_nvd->vvzz[sib_index] * local_nvd->sumr[sib_index];
972        local_nvd->vzsumy[sib_index] =
973          local_nvd->vvzz[sib_index] * local_nvd->sumy[sib_index];
974      }
975
976      /* Initialize to one, multiply incremental values for every
977         sibling a node has */
978      p_xx[(long)A] = 1 ;
979      p_xx[(long)C] = 1 ; 
980      p_xx[(long)G] = 1 ;
981      p_xx[(long)T] = 1 ;
982
983      for (sib_index=0; sib_index < num_sibs; sib_index++) {
984        p_xx[(long)A] *=
985          local_nvd->sum[sib_index] +
986          local_nvd->wwzz[sib_index] *
987          local_nvd->xx[sib_index][(long)A] +
988          local_nvd->vzsumr[sib_index];
989        p_xx[(long)C] *=
990          local_nvd->sum[sib_index] +
991          local_nvd->wwzz[sib_index] *
992          local_nvd->xx[sib_index][(long)C] +
993          local_nvd->vzsumy[sib_index];
994        p_xx[(long)G] *=
995          local_nvd->sum[sib_index] +
996          local_nvd->wwzz[sib_index] *
997          local_nvd->xx[sib_index][(long)G] +
998          local_nvd->vzsumr[sib_index];
999        p_xx[(long)T] *=
1000          local_nvd->sum[sib_index] +
1001          local_nvd->wwzz[sib_index] *
1002          local_nvd->xx[sib_index][(long)T] +
1003          local_nvd->vzsumy[sib_index];
1004      }
1005
1006      /* And the final point of this whole function: */
1007      memcpy(p->x[i][j], p_xx, sizeof(sitelike));
1008    }
1009  }
1010
1011  p->initialized = true;
1012
1013  free_nvd (local_nvd);
1014  free (local_nvd);
1015}  /* nuview */
1016
1017
1018void slopecurv(node *p,double y,double *like,double *slope,double *curve)
1019{
1020   /* compute log likelihood, slope and curvature at node p */
1021  long i, j, k, lai;
1022  double sum, sumc, sumterm, lterm, sumcs, sumcc, sum2, slope2, curve2,
1023        temp;
1024  double lz, zz, z1, zzs, z1s, zzc, z1c, aa, bb, cc,
1025         prod1, prod2, prod12, prod3;
1026  contribarr thelike, nulike, nuslope, nucurve,
1027    theslope, thecurve, clai, cslai, cclai;
1028  node *q;
1029  sitelike x1, x2;
1030
1031  q = p->back;
1032  sum = 0.0;
1033  lz = -y;
1034  for (i = 0; i < rcategs; i++)
1035    for (j = 0; j < categs; j++) {
1036      tbl[i][j]->orig_zz = exp(tbl[i][j]->rat * lz);
1037      tbl[i][j]->z1 = exp(tbl[i][j]->ratxv * lz);
1038    }
1039  for (i = 0; i < endsite; i++) {
1040    k = category[alias[i]-1] - 1;
1041    for (j = 0; j < rcategs; j++) {
1042      if (y > 0.0) {
1043        zz = tbl[j][k]->orig_zz;
1044        z1 = tbl[j][k]->z1;
1045      } else {
1046        zz = 1.0;
1047        z1 = 1.0;
1048      }
1049      zzs = -tbl[j][k]->rat * zz ;
1050      z1s = -tbl[j][k]->ratxv * z1 ;
1051      temp = tbl[j][k]->rat;
1052      zzc = temp * temp * zz;
1053      temp = tbl[j][k]->ratxv;
1054      z1c = temp * temp * z1;
1055      memcpy(x1, p->x[i][j], sizeof(sitelike));
1056      prod1 = freqa * x1[0] + freqc * x1[(long)C - (long)A] +
1057            freqg * x1[(long)G - (long)A] + freqt * x1[(long)T - (long)A];
1058      memcpy(x2, q->x[i][j], sizeof(sitelike));
1059      prod2 = freqa * x2[0] + freqc * x2[(long)C - (long)A] +
1060            freqg * x2[(long)G - (long)A] + freqt * x2[(long)T - (long)A];
1061      prod3 = (x1[0] * freqa + x1[(long)G - (long)A] * freqg) *
1062              (x2[0] * freqar + x2[(long)G - (long)A] * freqgr) +
1063        (x1[(long)C - (long)A] * freqc + x1[(long)T - (long)A] * freqt) *
1064        (x2[(long)C - (long)A] * freqcy + x2[(long)T - (long)A] * freqty);
1065      prod12 = freqa * x1[0] * x2[0] +
1066               freqc * x1[(long)C - (long)A] * x2[(long)C - (long)A] +
1067               freqg * x1[(long)G - (long)A] * x2[(long)G - (long)A] +
1068               freqt * x1[(long)T - (long)A] * x2[(long)T - (long)A];
1069      aa = prod12 - prod3;
1070      bb = prod3 - prod1*prod2;
1071      cc = prod1 * prod2;
1072      term[i][j] = zz * aa + z1 * bb + cc;
1073      slopeterm[i][j] = zzs * aa + z1s * bb;
1074      curveterm[i][j] = zzc * aa + z1c * bb;
1075    }
1076    sumterm = 0.0;
1077    for (j = 0; j < rcategs; j++)
1078      sumterm += probcat[j] * term[i][j];
1079    lterm = log(sumterm);
1080    for (j = 0; j < rcategs; j++) {
1081      term[i][j] = term[i][j] / sumterm;
1082      slopeterm[i][j] = slopeterm[i][j] / sumterm;
1083      curveterm[i][j] = curveterm[i][j] / sumterm; 
1084    }
1085    sum += aliasweight[i] * lterm;
1086  }
1087  for (i = 0; i < rcategs; i++) {
1088    thelike[i] = 1.0;
1089    theslope[i] = 0.0;
1090    thecurve[i] = 0.0;
1091  }
1092  for (i = 0; i < sites; i++) {
1093    sumc = 0.0;
1094    sumcs = 0.0;
1095    sumcc = 0.0;
1096    for (k = 0; k < rcategs; k++) {
1097      sumc += probcat[k] * thelike[k];
1098      sumcs += probcat[k] * theslope[k];
1099      sumcc += probcat[k] * thecurve[k];
1100    }
1101    sumc *= lambda;
1102    sumcs *= lambda;
1103    sumcc *= lambda;
1104    if ((ally[i] > 0) && (location[ally[i]-1] > 0)) {
1105      lai = location[ally[i] - 1];
1106      memcpy(clai, term[lai - 1], rcategs*sizeof(double));
1107      memcpy(cslai, slopeterm[lai - 1], rcategs*sizeof(double));
1108      memcpy(cclai, curveterm[lai - 1], rcategs*sizeof(double));
1109      if (weight[i] > 1) {
1110        for (j = 0; j < rcategs; j++) {
1111          if (clai[j] > 0.0)
1112            clai[j] = exp(weight[i]*log(clai[j]));
1113          else clai[j] = 0.0;
1114          if (cslai[j] > 0.0)
1115            cslai[j] = exp(weight[i]*log(cslai[j]));
1116          else cslai[j] = 0.0;
1117          if (cclai[j] > 0.0)
1118            cclai[j] = exp(weight[i]*log(cclai[j])); 
1119          else cclai[j] = 0.0;
1120      }
1121    }
1122      for (j = 0; j < rcategs; j++) {
1123        nulike[j] = ((1.0 - lambda) * thelike[j] + sumc) * clai[j];
1124        nuslope[j] = ((1.0 - lambda) * theslope[j] + sumcs) * clai[j]
1125                   + ((1.0 - lambda) * thelike[j] + sumc) * cslai[j];
1126        nucurve[j] = ((1.0 - lambda) * thecurve[j] + sumcc) * clai[j]
1127             + 2.0 * ((1.0 - lambda) * theslope[j] + sumcs) * cslai[j]
1128                   + ((1.0 - lambda) * thelike[j] + sumc) * cclai[j];
1129      }
1130    } else {
1131      for (j = 0; j < rcategs; j++) {
1132        nulike[j] = ((1.0 - lambda) * thelike[j] + sumc);
1133        nuslope[j] = ((1.0 - lambda) * theslope[j] + sumcs);
1134        nucurve[j] = ((1.0 - lambda) * thecurve[j] + sumcc);
1135      }
1136    }
1137    memcpy(thelike, nulike, rcategs*sizeof(double));
1138    memcpy(theslope, nuslope, rcategs*sizeof(double));
1139    memcpy(thecurve, nucurve, rcategs*sizeof(double));
1140  }
1141  sum2 = 0.0;
1142  slope2 = 0.0;
1143  curve2 = 0.0;
1144  for (i = 0; i < rcategs; i++) {
1145    sum2 += probcat[i] * thelike[i];
1146    slope2 += probcat[i] * theslope[i];
1147    curve2 += probcat[i] * thecurve[i];
1148  }
1149  sum += log(sum2);
1150  (*like) = sum;
1151  (*slope) = slope2 / sum2;
1152  (*curve) = (curve2 - slope2 * slope2 / sum2) / sum2;
1153} /* slopecurv */
1154
1155
1156void makenewv(node *p)
1157{
1158  /* Newton-Raphson algorithm improvement of a branch length */
1159  long it, ite;
1160  double y, yold=0, yorig, like, slope, curve, oldlike=0;
1161  boolean done, firsttime, better;
1162  node *q;
1163
1164  q = p->back;
1165  y = p->v;
1166  yorig = y;
1167  done = false;
1168  firsttime = true;
1169  it = 1;
1170  ite = 0;
1171  while ((it < iterations) && (ite < 20) && (!done)) {
1172    slopecurv (p, y, &like, &slope, &curve);
1173    better = false;
1174    if (firsttime) {
1175      yold = y;
1176      oldlike = like;
1177      firsttime = false;
1178      better = true;
1179    } else {
1180      if (like > oldlike) {
1181        yold = y;
1182        oldlike = like;
1183        better = true;
1184        it++;
1185      }
1186    }
1187    if (better) {
1188      y = y + slope/fabs(curve);
1189      if (yold < epsilon)
1190        yold = epsilon;
1191    } else {
1192      if (fabs(y - yold) < epsilon)
1193        ite = 20;
1194      y = (y + 19*yold) / 20.0;
1195    }
1196    ite++;
1197    done = fabs(y-yold) < epsilon;
1198  }
1199  smoothed = (fabs(yold-yorig) < epsilon) && (yorig > 1000.0*epsilon);
1200  p->v = yold;
1201  q->v = yold;
1202  curtree.likelihood = oldlike;
1203}  /* makenewv */
1204
1205
1206void update(node *p)
1207{
1208  long num_sibs, i;
1209  node *sib_ptr;
1210
1211  if (!p->tip && !p->initialized)
1212    nuview(p);
1213  if (!p->back->tip && !p->back->initialized)
1214    nuview(p->back);
1215  if ((!usertree) || (usertree && !lngths) || p->iter) {
1216    makenewv(p);
1217
1218    if (!p->tip) {
1219      num_sibs = count_sibs (p);
1220      sib_ptr  = p;
1221      for (i=0; i < num_sibs; i++) {
1222        sib_ptr              = sib_ptr->next;
1223        sib_ptr->initialized = false;
1224      }
1225    }
1226
1227    if (!p->back->tip) {
1228      num_sibs = count_sibs (p->back);
1229      sib_ptr  = p->back;
1230      for (i=0; i < num_sibs; i++) {
1231        sib_ptr              = sib_ptr->next;
1232        sib_ptr->initialized = false;
1233      }
1234    }
1235  }
1236}  /* update */
1237
1238
1239void smooth(node *p)
1240{
1241  long i, num_sibs;
1242  node *sib_ptr;
1243 
1244  smoothed = false;
1245  update (p);
1246  if (p->tip)
1247    return;
1248
1249  num_sibs = count_sibs (p);
1250  sib_ptr  = p;
1251 
1252  for (i=0; i < num_sibs; i++) {
1253    sib_ptr = sib_ptr->next;
1254   
1255    if (polishing || (smoothit && !smoothed)) {
1256      smooth(sib_ptr->back);
1257      p->initialized = false;
1258      sib_ptr->initialized = false;
1259    }
1260  }
1261}  /* smooth */
1262
1263
1264void insert_(node *p, node *q, boolean dooinit)
1265{
1266  /* Insert q near p */
1267  long i, j, num_sibs;
1268  node *r, *sib_ptr;
1269
1270  r = p->next->next;
1271  hookup(r, q->back);
1272  hookup(p->next, q);
1273  q->v = 0.5 * q->v;
1274  q->back->v = q->v;
1275  r->v = q->v;
1276  r->back->v = r->v;
1277  p->initialized = false;
1278  p->next->initialized = false;
1279  p->next->next->initialized = false;
1280  if (dooinit) {
1281    inittrav(p);
1282    inittrav(q);
1283    inittrav(q->back);
1284  }
1285  i = 1;
1286  while (i <= smoothings) {
1287    smooth (p);
1288    if (!smoothit) {
1289      if (!p->tip) {
1290        num_sibs = count_sibs (p);
1291        sib_ptr  = p;
1292        for (j=0; j < num_sibs; j++) {
1293          smooth (sib_ptr->next->back);
1294          sib_ptr = sib_ptr->next;
1295        }
1296      }
1297    }
1298    else 
1299      smooth(p->back);
1300    i++;
1301  }
1302}  /* insert_ */
1303
1304
1305void dnaml_re_move(node **p, node **q)
1306{
1307  /* remove p and record in q where it was */
1308  long i;
1309
1310  /** assumes bifurcations */
1311  *q = (*p)->next->back;
1312  hookup(*q, (*p)->next->next->back);
1313  (*p)->next->back = NULL;
1314  (*p)->next->next->back = NULL;
1315  inittrav((*q));
1316  inittrav((*q)->back);
1317  i = 1;
1318  while (i <= smoothings) {
1319    smooth(*q);
1320    if (smoothit)
1321      smooth((*q)->back);
1322    i++;
1323  }
1324}  /* dnaml_re_move */
1325
1326
1327void buildnewtip(long m, tree *tr)
1328{
1329  node *p;
1330
1331  p = tr->nodep[nextsp + spp - 3];
1332  hookup(tr->nodep[m - 1], p);
1333  p->v = initialv;
1334  p->back->v = initialv;
1335}  /* buildnewtip */
1336
1337
1338void buildsimpletree(tree *tr)
1339{
1340  hookup(tr->nodep[enterorder[0] - 1], tr->nodep[enterorder[1] - 1]);
1341  tr->nodep[enterorder[0] - 1]->v = 0.1;
1342  tr->nodep[enterorder[0] - 1]->back->v = 0.1;
1343  tr->nodep[enterorder[1] - 1]->v = 0.1;
1344  tr->nodep[enterorder[1] - 1]->back->v = 0.1;
1345  buildnewtip(enterorder[2], tr);
1346  insert_(tr->nodep[enterorder[2] - 1]->back,
1347          tr->nodep[enterorder[0] - 1], false);
1348}  /* buildsimpletree2 */
1349
1350
1351void addtraverse(node *p, node *q, boolean contin)
1352{
1353  /* try adding p at q, proceed recursively through tree */
1354  long i, num_sibs;
1355  double like, vsave = 0;
1356  node *qback = NULL, *sib_ptr;
1357
1358  if (!smoothit) {
1359    vsave = q->v;
1360    qback = q->back;
1361  }
1362  insert_(p, q, false);
1363  like = evaluate(p, false);
1364  if (like > bestyet || bestyet == UNDEFINED) {
1365    bestyet = like;
1366    if (smoothit)
1367      dnamlcopy(&curtree, &bestree, nonodes2, rcategs);
1368    else
1369      qwhere = q;
1370    succeeded = true;
1371  }
1372  if (smoothit)
1373    dnamlcopy(&priortree, &curtree, nonodes2, rcategs);
1374  else {
1375    hookup (q, qback);
1376    q->v = vsave;
1377    q->back->v = vsave;
1378    curtree.likelihood = bestyet;
1379  }
1380  if (!q->tip && contin) {
1381    num_sibs = count_sibs (q);
1382    if (q == curtree.start)
1383      num_sibs++;
1384    sib_ptr  = q;
1385    for (i=0; i < num_sibs; i++) {
1386      addtraverse(p, sib_ptr->next->back, contin);
1387      sib_ptr = sib_ptr->next;
1388    }
1389  }
1390
1391}  /* addtraverse */
1392
1393
1394void rearrange(node *p, node *pp)
1395{
1396  /* rearranges the tree, globally or locally moving pp around near p */
1397  long i, num_sibs;
1398  double v3 = 0, v4 = 0, v5 = 0;
1399  node *q, *r, *sib_ptr;
1400
1401  if (!p->tip && !p->back->tip) {
1402    curtree.likelihood = bestyet;
1403    if (p->back->next != pp)
1404      r = p->back->next;
1405    else
1406      r = p->back->next->next;
1407    /* assumes bifurcations? */
1408    if (!smoothit) {
1409      v3 = r->v;
1410      v4 = r->next->v;
1411      v5 = r->next->next->v;
1412    }
1413    else
1414      dnamlcopy(&curtree, &bestree, nonodes2, rcategs);
1415    dnaml_re_move(&r, &q);
1416    if (smoothit)
1417      dnamlcopy(&curtree, &priortree, nonodes2, rcategs);
1418    else
1419      qwhere = q;
1420    num_sibs = count_sibs (p);
1421    sib_ptr  = p;
1422    for (i=0; i < num_sibs; i++) {
1423      sib_ptr = sib_ptr->next;
1424      addtraverse(r, sib_ptr->back, (boolean)(global && (nextsp == spp)));
1425    }
1426    if (global && nextsp == spp && !succeeded) {
1427      p = p->back;
1428      if (!p->tip) {
1429        num_sibs = count_sibs (p);
1430        sib_ptr  = p;
1431        for (i=0; i < num_sibs; i++) {
1432          sib_ptr = sib_ptr->next;
1433          addtraverse(r, sib_ptr->back,(boolean)(global && (nextsp == spp)));
1434        }
1435      }
1436      p = p->back;
1437    }
1438    if (smoothit)
1439      dnamlcopy(&bestree, &curtree, nonodes2, rcategs);
1440    else {
1441      insert_(r, qwhere, true);
1442      if (qwhere == q) {
1443        r->v = v3;
1444        r->back->v = v3;
1445        r->next->v = v4;
1446        r->next->back->v = v4;
1447        r->next->next->v = v5;
1448        r->next->next->back->v = v5;
1449        curtree.likelihood = bestyet;
1450      }
1451      else {
1452        smoothit = true;
1453        for (i = 1; i<=smoothings; i++) {
1454          smooth (r);
1455          smooth (r->back);
1456        }
1457        smoothit = false;
1458        dnamlcopy(&curtree, &bestree, nonodes2, rcategs);
1459      }
1460    }
1461    if (global && nextsp == spp && progress) {
1462      putchar('.');
1463      fflush(stdout);
1464    }
1465  }
1466  if (!p->tip) {
1467    num_sibs = count_sibs (p);
1468    if (p == curtree.start)
1469      num_sibs++;
1470    sib_ptr  = p;
1471    for (i=0; i < num_sibs; i++) {
1472      sib_ptr = sib_ptr->next;
1473      rearrange(sib_ptr->back, p);
1474    }
1475  }
1476}  /* rearrange */
1477
1478
1479void initdnamlnode(node **p, node **grbg, node *q, long len, long nodei,
1480                        long *ntips, long *parens, initops whichinit,
1481                        pointarray treenode, pointarray nodep, Char *str, Char *ch,
1482                        FILE *intree)
1483{
1484  /* initializes a node */
1485  boolean minusread;
1486  double valyew, divisor;
1487
1488  switch (whichinit) {
1489  case bottom:
1490    gnu(grbg, p);
1491    (*p)->index = nodei;
1492    (*p)->tip = false;
1493    malloc_pheno((*p), endsite, rcategs);
1494    nodep[(*p)->index - 1] = (*p);
1495    break;
1496  case nonbottom:
1497    gnu(grbg, p);
1498    malloc_pheno(*p, endsite, rcategs);
1499    (*p)->index = nodei;
1500    break;
1501  case tip:
1502    match_names_to_data (str, nodep, p, spp);
1503    break;
1504  case iter:
1505    (*p)->initialized = false;
1506    (*p)->v = initialv;
1507    (*p)->iter = true;
1508    if ((*p)->back != NULL)
1509      (*p)->back->iter = true;
1510    break;
1511  case length:
1512    processlength(&valyew, &divisor, ch, &minusread, intree, parens);
1513    (*p)->v = valyew / divisor / fracchange;
1514    (*p)->iter = false;
1515    if ((*p)->back != NULL) {
1516      (*p)->back->v = (*p)->v;
1517      (*p)->back->iter = false;
1518    }
1519    break;
1520  case hsnolength:
1521    haslengths = false;
1522    break;
1523  default:        /* cases hslength, treewt, unittrwt */
1524    break;        /* should never occur                                */
1525  }
1526} /* initdnamlnode */
1527
1528
1529void dnaml_coordinates(node *p, double lengthsum, long *tipy,
1530                        double *tipmax)
1531{
1532  /* establishes coordinates of nodes */
1533  node *q, *first, *last;
1534  double xx;
1535
1536  if (p->tip) {
1537    p->xcoord = (long)(over * lengthsum + 0.5);
1538    p->ycoord = (*tipy);
1539    p->ymin = (*tipy);
1540    p->ymax = (*tipy);
1541    (*tipy) += down;
1542    if (lengthsum > (*tipmax))
1543      (*tipmax) = lengthsum;
1544    return;
1545  }
1546  q = p->next;
1547  do {
1548    xx = fracchange * q->v;
1549    if (xx > 100.0)
1550      xx = 100.0;
1551    dnaml_coordinates(q->back, lengthsum + xx, tipy,tipmax);
1552    q = q->next;
1553  } while ((p == curtree.start || p != q) &&
1554           (p != curtree.start || p->next != q));
1555  first = p->next->back;
1556  q = p;
1557  while (q->next != p)
1558    q = q->next;
1559  last = q->back;
1560  p->xcoord = (long)(over * lengthsum + 0.5);
1561  if (p == curtree.start)
1562    p->ycoord = p->next->next->back->ycoord;
1563  else
1564    p->ycoord = (first->ycoord + last->ycoord) / 2;
1565  p->ymin = first->ymin;
1566  p->ymax = last->ymax;
1567}  /* dnaml_coordinates */
1568
1569
1570void dnaml_printree()
1571{
1572  /* prints out diagram of the tree2 */
1573  long tipy;
1574  double scale, tipmax;
1575  long i;
1576
1577  if (!treeprint)
1578    return;
1579  putc('\n', outfile);
1580  tipy = 1;
1581  tipmax = 0.0;
1582  dnaml_coordinates(curtree.start, 0.0, &tipy, &tipmax);
1583  scale = 1.0 / (long)(tipmax + 1.000);
1584  for (i = 1; i <= (tipy - down); i++)
1585    drawline2(i, scale, curtree);
1586  putc('\n', outfile);
1587}  /* dnaml_printree */
1588
1589
1590void sigma(node *p, double *sumlr, double *s1, double *s2)
1591{
1592  /* compute standard deviation */
1593  double tt, aa, like, slope, curv;
1594
1595  slopecurv (p, p->v, &like, &slope, &curv);
1596  tt = p->v;
1597  p->v = epsilon;
1598  p->back->v = epsilon;
1599  aa = evaluate(p, false);
1600  p->v = tt;
1601  p->back->v = tt;
1602  (*sumlr) = evaluate(p, false) - aa;
1603  if (curv < -epsilon) {
1604    (*s1) = p->v + (-slope - sqrt(slope * slope -  3.841 * curv)) / curv;
1605    (*s2) = p->v + (-slope + sqrt(slope * slope -  3.841 * curv)) / curv;
1606  }
1607  else {
1608    (*s1) = -1.0;
1609    (*s2) = -1.0;
1610  }
1611}  /* sigma */
1612
1613
1614void describe(node *p)
1615{
1616  /* print out information for one branch */
1617  long i, num_sibs;
1618  node *q, *sib_ptr;
1619  double sumlr, sigma1, sigma2;
1620
1621  if (!p->tip && !p->initialized)
1622    nuview(p);
1623  if (!p->back->tip && !p->back->initialized)
1624    nuview(p->back);
1625  q = p->back;
1626  if (q->tip) {
1627    fprintf(outfile, " ");
1628    for (i = 0; i < nmlngth; i++)
1629      putc(nayme[q->index-1][i], outfile);
1630    fprintf(outfile, "    ");
1631  } else
1632    fprintf(outfile, "  %4ld          ", q->index - spp);
1633  if (p->tip) {
1634    for (i = 0; i < nmlngth; i++)
1635      putc(nayme[p->index-1][i], outfile);
1636  } else
1637    fprintf(outfile, "%4ld      ", p->index - spp);
1638  fprintf(outfile, "%15.5f", q->v * fracchange);
1639  if (!usertree || (usertree && !lngths) || p->iter) {
1640    sigma(q, &sumlr, &sigma1, &sigma2);
1641    if (sigma1 <= sigma2)
1642      fprintf(outfile, "     (     zero,    infinity)");
1643    else {
1644      fprintf(outfile, "     (");
1645      if (sigma2 <= 0.0)
1646        fprintf(outfile, "     zero");
1647      else
1648        fprintf(outfile, "%9.5f", sigma2 * fracchange);
1649      fprintf(outfile, ",%12.5f", sigma1 * fracchange);
1650      putc(')', outfile);
1651      }
1652    if (sumlr > 1.9205)
1653      fprintf(outfile, " *");
1654    if (sumlr > 2.995)
1655      putc('*', outfile);
1656    }
1657  putc('\n', outfile);
1658  if (!p->tip) {
1659    num_sibs = count_sibs (p);
1660    sib_ptr  = p;
1661    for (i=0; i < num_sibs; i++) {
1662      sib_ptr = sib_ptr->next;
1663      describe(sib_ptr->back);
1664    }
1665  }
1666}  /* describe */
1667
1668
1669void reconstr(node *p, long n)
1670{
1671  /* reconstruct and print out base at site n+1 at node p */
1672  long i, j, k, m, first, second, num_sibs;
1673  double f, sum, xx[4];
1674  node *q;
1675
1676  if ((ally[n] == 0) || (location[ally[n]-1] == 0))
1677    putc('.', outfile);
1678  else {
1679    j = location[ally[n]-1] - 1;
1680    for (i = 0; i < 4; i++) {
1681      f = p->x[j][mx-1][i];
1682      num_sibs = count_sibs(p);
1683      q = p;
1684      for (k = 0; k < num_sibs; k++) {
1685        q = q->next;
1686        f *= q->x[j][mx-1][i];
1687      }
1688      f = sqrt(f);
1689      xx[i] = f;
1690    }
1691    xx[0] *= freqa;
1692    xx[1] *= freqc;
1693    xx[2] *= freqg;
1694    xx[3] *= freqt;
1695    sum = xx[0]+xx[1]+xx[2]+xx[3];
1696    for (i = 0; i < 4; i++)
1697      xx[i] /= sum;
1698    first = 0;
1699    for (i = 1; i < 4; i++)
1700      if (xx [i] > xx[first])
1701        first = i;
1702    if (first == 0)
1703      second = 1;
1704    else
1705      second = 0;
1706    for (i = 0; i < 4; i++)
1707      if ((i != first) && (xx[i] > xx[second]))
1708        second = i;
1709    m = 1 << first;
1710    if (xx[first] < 0.4999995)
1711      m = m + (1 << second);
1712    if (xx[first] > 0.95)
1713      putc(toupper(basechar[m - 1]), outfile);
1714    else
1715      putc(basechar[m - 1], outfile);
1716    if (rctgry && rcategs > 1)
1717      mx = mp[n][mx - 1];   
1718    else
1719      mx = 1;
1720  }
1721} /* reconstr */
1722
1723
1724void rectrav(node *p, long m, long n)
1725{
1726  /* print out segment of reconstructed sequence for one branch */
1727  long i;
1728
1729  putc(' ', outfile);
1730  if (p->tip) {
1731    for (i = 0; i < nmlngth; i++)
1732      putc(nayme[p->index-1][i], outfile);
1733  } else
1734    fprintf(outfile, "%4ld      ", p->index - spp);
1735  fprintf(outfile, "  ");
1736  mx = mx0;
1737  for (i = m; i <= n; i++) {
1738    if ((i % 10 == 0) && (i != m))
1739      putc(' ', outfile);
1740    if (p->tip)
1741      putc(y[p->index-1][i], outfile);
1742    else
1743      reconstr(p, i);
1744  }
1745  putc('\n', outfile);
1746  if (!p->tip) {
1747    rectrav(p->next->back, m, n);
1748    rectrav(p->next->next->back, m, n);
1749  }
1750  mx1 = mx;
1751}  /* rectrav */
1752
1753
1754void summarize()
1755{
1756  /* print out branch length information and node numbers */
1757  long i, j, mm, num_sibs;
1758  double mode, sum;
1759  double like[maxcategs], nulike[maxcategs];
1760  double **marginal;
1761  node   *sib_ptr;
1762
1763  if (!treeprint)
1764    return;
1765  fprintf(outfile, "\nremember: ");
1766  if (outgropt)
1767    fprintf(outfile, "(although rooted by outgroup) ");
1768  fprintf(outfile, "this is an unrooted tree!\n\n");
1769  fprintf(outfile, "Ln Likelihood = %11.5f\n", curtree.likelihood);
1770  fprintf(outfile, "\n Between        And            Length");
1771  if (!(usertree && lngths && haslengths))
1772    fprintf(outfile, "      Approx. Confidence Limits");
1773  fprintf(outfile, "\n");
1774  fprintf(outfile, " -------        ---            ------");
1775  if (!(usertree && lngths && haslengths))
1776    fprintf(outfile, "      ------- ---------- ------");
1777  fprintf(outfile, "\n\n");
1778  for (i = spp; i < nonodes2; i++) {
1779    /* So this works with arbitrary multifurcations */
1780    if (curtree.nodep[i]) {
1781      num_sibs = count_sibs (curtree.nodep[i]);
1782      sib_ptr  = curtree.nodep[i];
1783      for (j = 0; j < num_sibs; j++) {
1784        sib_ptr->initialized = false;
1785        sib_ptr              = sib_ptr->next;
1786      }
1787    }
1788  }
1789
1790  describe(curtree.start->back);
1791
1792  /* So this works with arbitrary multifurcations */
1793  num_sibs = count_sibs (curtree.start);
1794  sib_ptr  = curtree.start;
1795  for (i=0; i < num_sibs; i++) {
1796    sib_ptr = sib_ptr->next;
1797    describe(sib_ptr->back);
1798  }
1799
1800  fprintf(outfile, "\n");
1801  if (!(usertree && lngths && haslengths)) {
1802    fprintf(outfile, "     *  = significantly positive, P < 0.05\n");
1803    fprintf(outfile, "     ** = significantly positive, P < 0.01\n\n");
1804  }
1805  dummy = evaluate(curtree.start, false);
1806  if (rctgry && rcategs > 1) {
1807    for (i = 0; i < rcategs; i++)
1808      like[i] = 1.0;
1809    for (i = sites - 1; i >= 0; i--) {
1810      sum = 0.0;
1811      for (j = 0; j < rcategs; j++) {
1812        nulike[j] = (1.0 - lambda + lambda * probcat[j]) * like[j];
1813        mp[i][j] = j + 1;
1814        for (k = 1; k <= rcategs; k++) {
1815          if (k != j + 1) {
1816            if (lambda * probcat[k - 1] * like[k - 1] > nulike[j]) {
1817              nulike[j] = lambda * probcat[k - 1] * like[k - 1];
1818              mp[i][j] = k;
1819            }
1820          }
1821        }
1822        if ((ally[i] > 0) && (location[ally[i]-1] > 0))
1823          nulike[j] *= contribution[location[ally[i] - 1] - 1][j];
1824        sum += nulike[j];
1825      }
1826      for (j = 0; j < rcategs; j++)
1827        nulike[j] /= sum;
1828      memcpy(like, nulike, rcategs * sizeof(double));
1829    }
1830    mode = 0.0;
1831    mx = 1;
1832    for (i = 1; i <= rcategs; i++) {
1833      if (probcat[i - 1] * like[i - 1] > mode) {
1834        mx = i;
1835        mode = probcat[i - 1] * like[i - 1];
1836      }
1837    }
1838    mx0 = mx;
1839    fprintf(outfile,
1840 "Combination of categories that contributes the most to the likelihood:\n\n");
1841    for (i = 1; i <= nmlngth + 3; i++)
1842      putc(' ', outfile);
1843    for (i = 1; i <= sites; i++) {
1844      fprintf(outfile, "%ld", mx);
1845      if (i % 10 == 0)
1846        putc(' ', outfile);
1847      if (i % 60 == 0 && i != sites) {
1848        putc('\n', outfile);
1849        for (j = 1; j <= nmlngth + 3; j++)
1850          putc(' ', outfile);
1851      }
1852      mx = mp[i - 1][mx - 1];
1853    }
1854    fprintf(outfile, "\n\n");
1855    marginal = (double **) Malloc( sites*sizeof(double *));
1856    for (i = 0; i < sites; i++)
1857      marginal[i] = (double *) Malloc( rcategs*sizeof(double));
1858    for (i = 0; i < rcategs; i++)
1859      like[i] = 1.0;
1860    for (i = sites - 1; i >= 0; i--) {
1861      sum = 0.0;
1862      for (j = 0; j < rcategs; j++) {
1863        nulike[j] = (1.0 - lambda + lambda * probcat[j]) * like[j];
1864        for (k = 1; k <= rcategs; k++) {
1865          if (k != j + 1)
1866              nulike[j] += lambda * probcat[k - 1] * like[k - 1];
1867        }
1868        if ((ally[i] > 0) && (location[ally[i]-1] > 0))
1869          nulike[j] *= contribution[location[ally[i] - 1] - 1][j];
1870        sum += nulike[j];
1871      }
1872      for (j = 0; j < rcategs; j++) {
1873        nulike[j] /= sum;
1874        marginal[i][j] = nulike[j];
1875      }
1876      memcpy(like, nulike, rcategs * sizeof(double));
1877    }
1878    for (i = 0; i < rcategs; i++)
1879      like[i] = 1.0;
1880    for (i = 0; i < sites; i++) {
1881      sum = 0.0;
1882      for (j = 0; j < rcategs; j++) {
1883        nulike[j] = (1.0 - lambda + lambda * probcat[j]) * like[j];
1884        for (k = 1; k <= rcategs; k++) {
1885          if (k != j + 1)
1886              nulike[j] += lambda * probcat[k - 1] * like[k - 1];
1887        }
1888        marginal[i][j] *= like[j] * probcat[j];
1889        sum += nulike[j];
1890      }
1891      for (j = 0; j < rcategs; j++)
1892        nulike[j] /= sum;
1893      memcpy(like, nulike, rcategs * sizeof(double));
1894      sum = 0.0;
1895      for (j = 0; j < rcategs; j++)
1896        sum += marginal[i][j];
1897      for (j = 0; j < rcategs; j++)
1898        marginal[i][j] /= sum;
1899    }
1900    fprintf(outfile, "Most probable category at each site if > 0.95 probability (\".\" otherwise)\n\n");
1901    for (i = 1; i <= nmlngth + 3; i++)
1902      putc(' ', outfile);
1903    for (i = 0; i < sites; i++) {
1904      sum = 0.0;
1905      for (j = 0; j < rcategs; j++)
1906        if (marginal[i][j] > sum) {
1907          sum = marginal[i][j];
1908          mm = j;
1909        }
1910      if (sum >= 0.95)
1911        fprintf(outfile, "%ld", mm+1);
1912      else
1913        putc('.', outfile);
1914      if ((i+1) % 60 == 0) {
1915        if (i != 0) {
1916          putc('\n', outfile);
1917          for (j = 1; j <= nmlngth + 3; j++)
1918            putc(' ', outfile);
1919        }
1920      }
1921      else if ((i+1) % 10 == 0)
1922        putc(' ', outfile);
1923    }
1924    putc('\n', outfile);
1925    for (i = 0; i < sites; i++)
1926      free(marginal[i]);
1927    free(marginal);
1928  }
1929  putc('\n', outfile);
1930  if (hypstate) {
1931    fprintf(outfile, "Probable sequences at interior nodes:\n\n");
1932    fprintf(outfile, "  node       ");
1933    for (i = 0; (i < 13) && (i < ((sites + (sites-1)/10 - 39) / 2)); i++)
1934      putc(' ', outfile);
1935    fprintf(outfile, "Reconstructed sequence (caps if > 0.95)\n\n");
1936    if (!rctgry || (rcategs == 1))
1937      mx0 = 1;
1938    for (i = 0; i < sites; i += 60) {
1939      k = i + 59;
1940      if (k >= sites)
1941        k = sites - 1;
1942      rectrav(curtree.start, i, k);
1943      rectrav(curtree.start->back, i, k);
1944      putc('\n', outfile);
1945      mx0 = mx1;
1946    }
1947  }
1948}  /* summarize */
1949
1950
1951void dnaml_treeout(node *p)
1952{
1953  /* write out file with representation of final tree2 */
1954  /* Only works for bifurcations! */
1955  long i, n, w;
1956  Char c;
1957  double x;
1958
1959  if (p->tip) {
1960    n = 0;
1961    for (i = 1; i <= nmlngth; i++) {
1962      if (nayme[p->index-1][i - 1] != ' ')
1963        n = i;
1964    }
1965    for (i = 0; i < n; i++) {
1966      c = nayme[p->index-1][i];
1967      if (c == ' ')
1968        c = '_';
1969      putc(c, outtree);
1970    }
1971    col += n;
1972  } else {
1973    putc('(', outtree);
1974    col++;
1975    dnaml_treeout(p->next->back);
1976    putc(',', outtree);
1977    col++;
1978    if (col > 45) {
1979      putc('\n', outtree);
1980      col = 0;
1981    }
1982    dnaml_treeout(p->next->next->back);
1983    if (p == curtree.start) {
1984      putc(',', outtree);
1985      col++;
1986      if (col > 45) {
1987        putc('\n', outtree);
1988        col = 0;
1989      }
1990      dnaml_treeout(p->back);
1991    }
1992    putc(')', outtree);
1993    col++;
1994  }
1995  x = p->v * fracchange;
1996  if (x > 0.0)
1997    w = (long)(0.43429448222 * log(x));
1998  else if (x == 0.0)
1999    w = 0;
2000  else
2001    w = (long)(0.43429448222 * log(-x)) + 1;
2002  if (w < 0)
2003    w = 0;
2004  if (p == curtree.start)
2005    fprintf(outtree, ";\n");
2006  else {
2007    fprintf(outtree, ":%*.5f", (int)(w + 7), x);
2008    col += w + 8;
2009  }
2010}  /* dnaml_treeout */
2011
2012
2013void inittravtree(node *p)
2014{
2015  /* traverse tree to set initialized and v to initial values */
2016
2017  p->initialized = false;
2018  p->back->initialized = false;
2019  if (!p->tip) {
2020    inittravtree(p->next->back);
2021    inittravtree(p->next->next->back);
2022  }
2023} /* inittravtree */
2024
2025
2026void treevaluate()
2027{
2028  /* evaluate a user tree */
2029  long i;
2030
2031  inittravtree(curtree.start);
2032  polishing = true;
2033  smoothit = true;
2034  for (i = 1; i <= smoothings * 4; i++)
2035    smooth (curtree.start);
2036  dummy = evaluate(curtree.start, true);
2037}  /* treevaluate */
2038
2039
2040void maketree()
2041{
2042  long i, j;
2043  boolean dummy_first, goteof;
2044  pointarray dummy_treenode=NULL;
2045  long nextnode;
2046  node *root, *q, *r;
2047
2048  inittable();
2049
2050  if (usertree) {
2051    openfile(&intree,INTREE,"input tree file", "r",progname,intreename);
2052    inittable_for_usertree (intree);
2053    numtrees = countsemic(&intree);
2054    if (numtrees > 2)
2055      initseed(&inseed, &inseed0, seed);
2056    l0gl = (double *) Malloc(numtrees * sizeof(double));
2057    l0gf = (double **) Malloc(numtrees * sizeof(double *));
2058    for (i=0; i < numtrees; ++i)
2059      l0gf[i] = (double *) Malloc(endsite * sizeof(double));
2060    if (treeprint) {
2061      fprintf(outfile, "User-defined tree");
2062      if (numtrees > 1)
2063        putc('s', outfile);
2064      fprintf(outfile, ":\n\n");
2065    }
2066    which = 1;
2067
2068    /* This taken out of tree read, used to be [spp-1], but referring
2069       to [0] produces output identical to what the pre-modified dnaml
2070       produced. */
2071
2072    while (which <= numtrees) {
2073     
2074      /* These initializations required each time through the loop
2075         since multiple trees require re-initialization */
2076      haslengths = true;
2077      nextnode         = 0;
2078      dummy_first      = true;
2079      goteof           = false;
2080
2081      treeread(intree, &root, dummy_treenode, &goteof, &dummy_first,
2082               curtree.nodep, &nextnode,
2083               &haslengths, &grbg, initdnamlnode);
2084      q = root;
2085      r = root;
2086      while (!(q->next == root))
2087        q = q->next;
2088      q->next = root->next;
2089      root = q;
2090      chuck(&grbg, r);
2091      curtree.nodep[spp] = q;
2092      if (goteof && (which <= numtrees)) {
2093        /* if we hit the end of the file prematurely */
2094        printf ("\n");
2095        printf ("ERROR: trees missing at end of file.\n");
2096        printf ("\tExpected number of trees:\t\t%ld\n", numtrees);
2097        printf ("\tNumber of trees actually in file:\t%ld.\n\n", which - 1);
2098        exxit(-1);
2099      }
2100
2101      curtree.start = curtree.nodep[0]->back;
2102     
2103      treevaluate();
2104      if (reconsider) {
2105        bestyet = UNDEFINED;
2106        succeeded = true;
2107        while (succeeded) {
2108          succeeded = false;
2109          rearrange(curtree.start, curtree.start->back);
2110        }
2111        treevaluate();
2112      }
2113      dnaml_printree();
2114      summarize();
2115      if (trout) {
2116        col = 0;
2117        dnaml_treeout(curtree.start);
2118      }
2119      which++;
2120    }
2121    FClose(intree);
2122    putc('\n', outfile);
2123    if (!auto_ && numtrees > 1 && weightsum > 1 )
2124      standev2(numtrees, maxwhich, 0, endsite-1, maxlogl,
2125               l0gl, l0gf, aliasweight, seed);
2126  } else {
2127    /* If there's no input user tree, */
2128    for (i = 1; i <= spp; i++)
2129      enterorder[i - 1] = i;
2130    if (jumble)
2131      randumize(seed, enterorder);
2132    if (progress) {
2133      printf("\nAdding species:\n");
2134      writename(0, 3, enterorder);
2135#ifdef WIN32
2136      phyFillScreenColor();
2137#endif
2138    }
2139    nextsp = 3;
2140    polishing = false;
2141    buildsimpletree(&curtree);
2142    curtree.start = curtree.nodep[enterorder[0] - 1]->back;
2143    smoothit = improve;
2144    nextsp = 4;
2145    while (nextsp <= spp) {
2146      buildnewtip(enterorder[nextsp - 1], &curtree);
2147      bestyet = UNDEFINED;
2148      if (smoothit)
2149        dnamlcopy(&curtree, &priortree, nonodes2, rcategs);
2150      addtraverse(curtree.nodep[enterorder[nextsp - 1] - 1]->back,
2151                  curtree.start, true);
2152      if (smoothit)
2153        dnamlcopy(&bestree, &curtree, nonodes2, rcategs);
2154      else {
2155        insert_(curtree.nodep[enterorder[nextsp - 1] - 1]->back, qwhere, true);
2156        smoothit = true;
2157        for (i = 1; i<=smoothings; i++) {
2158          smooth (curtree.start);
2159          smooth (curtree.start->back);
2160        }
2161        smoothit = false;
2162        dnamlcopy(&curtree, &bestree, nonodes2, rcategs);
2163        bestyet = curtree.likelihood;
2164      }
2165      if (progress) {
2166        writename(nextsp - 1, 1, enterorder);
2167#ifdef WIN32
2168        phyFillScreenColor();
2169#endif
2170      }
2171      if (global && nextsp == spp && progress) {
2172        printf("Doing global rearrangements\n");
2173        printf("  !");
2174        for (j = 1; j <= (spp - 3); j++)
2175          putchar('-');
2176        printf("!\n");
2177      }
2178      succeeded = true;
2179      while (succeeded) {
2180        succeeded = false;
2181        if (global && nextsp == spp && progress) {
2182          printf("   ");
2183          fflush(stdout);
2184        }
2185        rearrange(curtree.start, curtree.start->back);
2186        if (global && nextsp == spp && progress)
2187          putchar('\n');
2188      }
2189      for (i = spp; i < nextsp + spp - 2; i++) {
2190        curtree.nodep[i]->initialized = false;
2191        curtree.nodep[i]->next->initialized = false;
2192        curtree.nodep[i]->next->next->initialized = false;
2193      }
2194      if (!smoothit) {
2195        smoothit = true;
2196        for (i = 1; i<=smoothings; i++) {
2197          smooth (curtree.start);
2198          smooth (curtree.start->back);
2199        }
2200        smoothit = false;
2201        dnamlcopy(&curtree, &bestree, nonodes2, rcategs);
2202        bestyet = curtree.likelihood;
2203      }
2204      nextsp++;
2205    }
2206    if (global && progress) {
2207      putchar('\n');
2208      fflush(stdout);
2209    }
2210    if (njumble > 1) {
2211      if (jumb == 1)
2212        dnamlcopy(&bestree, &bestree2, nonodes2, rcategs);
2213      else
2214        if (bestree2.likelihood < bestree.likelihood)
2215          dnamlcopy(&bestree, &bestree2, nonodes2, rcategs);
2216    }
2217    if (jumb == njumble) {
2218      if (njumble > 1)
2219        dnamlcopy(&bestree2, &curtree, nonodes2, rcategs);
2220      curtree.start = curtree.nodep[outgrno - 1]->back;
2221      for (i = 0; i < nonodes2; i++) {
2222        if (i < spp)
2223          curtree.nodep[i]->initialized = false;
2224        else {
2225          curtree.nodep[i]->initialized = false;
2226          curtree.nodep[i]->next->initialized = false;
2227          curtree.nodep[i]->next->next->initialized = false;
2228        } 
2229      }
2230      treevaluate();
2231      dnaml_printree();
2232      summarize();
2233      if (trout) {
2234        col = 0;
2235        dnaml_treeout(curtree.start);
2236      }
2237    }
2238  }
2239  if (usertree) {
2240    free(l0gl);
2241    for (i=0; i < numtrees; i++)
2242      free(l0gf[i]);
2243    free(l0gf);
2244  }
2245  for (i = 0; i < rcategs; i++)
2246    for (j = 0; j < categs; j++)
2247      free(tbl[i][j]);
2248  for (i = 0; i < rcategs; i++)
2249    free(tbl[i]);
2250  free(tbl);
2251  if (jumb < njumble)
2252    return;
2253  free(contribution);
2254  free(mp);
2255  for (i=0; i < endsite; i++)
2256     free(term[i]);
2257  free(term);
2258  for (i=0; i < endsite; i++)
2259     free(slopeterm[i]);
2260  free(slopeterm);
2261  for (i=0; i < endsite; i++)
2262     free(curveterm[i]);
2263  free(curveterm);
2264  free_all_x_in_array (nonodes2, curtree.nodep);
2265  if (!usertree || reconsider) {
2266    free_all_x_in_array (nonodes2, bestree.nodep);
2267    free_all_x_in_array (nonodes2, priortree.nodep);
2268    if (njumble > 1)
2269      free_all_x_in_array (nonodes2, bestree2.nodep);
2270  }
2271  if (progress) {
2272    printf("\n\nOutput written to file \"%s\"\n\n", outfilename);
2273    if (trout)
2274      printf("Tree also written onto file \"%s\"\n", outtreename);
2275    putchar('\n');
2276  }
2277}  /* maketree */
2278
2279
2280void clean_up()
2281{
2282  /* Free and/or close stuff */
2283  long i;
2284
2285    free (rrate);
2286    free (probcat);
2287    free (rate);
2288    /* Seems to require freeing every time... */
2289    for (i = 0; i < spp; i++) {
2290      free (y[i]);
2291    }
2292    free (y);
2293    free (nayme);
2294    free (enterorder);
2295    free (category);
2296    free (weight);
2297    free (alias);
2298    free (ally);
2299    free (location);
2300    free (aliasweight);
2301
2302#if 0          /* ???? debug ???? */
2303  freetree2(curtree.nodep, nonodes2);
2304
2305  if (! (usertree && !reconsider)) {
2306    freetree2(bestree.nodep, nonodes2);
2307    freetree2(priortree.nodep, nonodes2);
2308  }
2309 
2310  if (! (njumble <= 1))
2311    freetree2(bestree2.nodep, nonodes2);
2312#endif
2313  FClose(infile);
2314  FClose(outfile);
2315  FClose(outtree);
2316#ifdef MAC
2317  fixmacfile(outfilename);
2318  fixmacfile(outtreename);
2319#endif
2320}   /* clean_up */
2321
2322
2323int main(int argc, Char *argv[])
2324{  /* DNA Maximum Likelihood */
2325#ifdef MAC
2326  argc = 1;                /* macsetup("DnaML","");        */
2327  argv[0] = "DnaML";
2328#endif
2329  init(argc,argv);
2330  progname = argv[0];
2331  openfile(&infile,INFILE,"input file","r",argv[0],infilename);
2332  openfile(&outfile,OUTFILE,"output file","w",argv[0],outfilename);
2333  mulsets = false;
2334  datasets = 1;
2335  firstset = true;
2336  ibmpc = IBMCRT;
2337  ansi = ANSICRT;
2338  grbg = NULL;
2339  doinit();
2340  ttratio0 = ttratio;
2341  if (ctgry)
2342    openfile(&catfile,CATFILE,"categories file","r",argv[0],catfilename);
2343  if (weights || justwts)
2344    openfile(&weightfile,WEIGHTFILE,"weights file","r",argv[0],weightfilename);
2345  if (trout)
2346    openfile(&outtree,OUTTREE,"output tree file","w",argv[0],outtreename);
2347  for (ith = 1; ith <= datasets; ith++) {
2348    if (datasets > 1) {
2349      fprintf(outfile, "Data set # %ld:\n", ith);
2350      printf("\nData set # %ld:\n", ith);
2351    }
2352    ttratio = ttratio0;
2353    getinput();
2354    if (ith == 1)
2355      firstset = false;
2356    for (jumb = 1; jumb <= njumble; jumb++)
2357      maketree();
2358  }
2359
2360  clean_up();
2361  printf("Done.\n\n");
2362#ifdef WIN32
2363  phyRestoreConsoleAttributes();
2364#endif
2365  return 0;
2366}  /* DNA Maximum Likelihood */
2367 
Note: See TracBrowser for help on using the repository browser.