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