1 | #include "phylip.h" |
---|
2 | |
---|
3 | /* version 3.52c. (c) Copyright 1993 by Joseph Felsenstein. |
---|
4 | Written by Joseph Felsenstein, Akiko Fuseki, Sean Lamont, and Andrew Keeffe. |
---|
5 | Permission is granted to copy and use this program provided no fee is |
---|
6 | charged for it and provided that this copyright notice is not removed. */ |
---|
7 | |
---|
8 | #define maxcutter 8 /* maximum number of bases in a site */ |
---|
9 | #define maxtrees 8 /* maximum number of user trees */ |
---|
10 | #define smoothings 10 /* number of passes in smoothing algorithm */ |
---|
11 | #define iterations 10 /* number of iterates of EM for each branch */ |
---|
12 | #define nmlngth 10 /* number of characters max. in species name */ |
---|
13 | |
---|
14 | #define epsilon 0.00001 /* used in update */ |
---|
15 | #define extrap0 100.0 /* extrapolation factor to speed iteration */ |
---|
16 | #define initialv 0.1 /* starting value of branch length */ |
---|
17 | |
---|
18 | #define ibmpc0 false |
---|
19 | #define ansi0 true |
---|
20 | #define vt520 false |
---|
21 | #define down 2 |
---|
22 | #define over 60 |
---|
23 | #define numsp2 (numsp * 2 - 2) |
---|
24 | |
---|
25 | |
---|
26 | typedef double sitelike[maxcutter + 1]; |
---|
27 | typedef sitelike *phenotype; |
---|
28 | typedef Char **sequence; |
---|
29 | typedef Char naym[nmlngth]; |
---|
30 | typedef short longer[6]; |
---|
31 | typedef double **transmatrix; |
---|
32 | typedef transmatrix *transptr; |
---|
33 | |
---|
34 | typedef struct node { |
---|
35 | struct node *next, *back; |
---|
36 | boolean tip, iter, initialized; |
---|
37 | short branchnum, number; |
---|
38 | phenotype x; |
---|
39 | naym nayme; |
---|
40 | double v; |
---|
41 | short xcoord, ycoord, ymin, ymax; |
---|
42 | } node; |
---|
43 | |
---|
44 | typedef struct tree { |
---|
45 | node **nodep; |
---|
46 | transptr trans, transprod; |
---|
47 | double likelihood; |
---|
48 | node *start; |
---|
49 | } tree; |
---|
50 | |
---|
51 | /* Local variables for maketree, propagated globally for C version: */ |
---|
52 | short nextsp, numtrees, which, maxwhich, col; |
---|
53 | double maxlogl; |
---|
54 | boolean succeeded, smoothed; |
---|
55 | transmatrix tempmatrix, tempprod; |
---|
56 | sitelike pi; |
---|
57 | double l0gl[maxtrees]; |
---|
58 | double *l0gf[maxtrees]; |
---|
59 | Char ch; |
---|
60 | |
---|
61 | /* variables declared in the other segment */ |
---|
62 | extern boolean trunc8,usertree,trout,lengths,treeprint,usertree,global, |
---|
63 | progress,outgropt,jumble; |
---|
64 | extern short endsite,enzymes,weightsum,numsp,outgrno,sitelength, |
---|
65 | jumb,njumble; |
---|
66 | extern tree curtree,priortree,bestree,bestree2; |
---|
67 | extern short *weight, *alias, *aliasweight, *enterorder; |
---|
68 | extern FILE *infile, *outfile, *treefile; |
---|
69 | extern double extrapol; |
---|
70 | extern longer seed; |
---|
71 | |
---|
72 | |
---|
73 | double randum(seed) |
---|
74 | short *seed; |
---|
75 | { |
---|
76 | /* random number generator -- slow but machine independent */ |
---|
77 | short i, j, k, sum; |
---|
78 | longer mult, newseed; |
---|
79 | double x; |
---|
80 | |
---|
81 | mult[0] = 13; |
---|
82 | mult[1] = 24; |
---|
83 | mult[2] = 22; |
---|
84 | mult[3] = 6; |
---|
85 | for (i = 0; i <= 5; i++) |
---|
86 | newseed[i] = 0; |
---|
87 | for (i = 0; i <= 5; i++) { |
---|
88 | sum = newseed[i]; |
---|
89 | k = i; |
---|
90 | if (i > 3) |
---|
91 | k = 3; |
---|
92 | for (j = 0; j <= k; j++) |
---|
93 | sum += mult[j] * seed[i - j]; |
---|
94 | newseed[i] = sum; |
---|
95 | for (j = i; j <= 4; j++) { |
---|
96 | newseed[j + 1] += newseed[j] / 64; |
---|
97 | newseed[j] &= 63; |
---|
98 | } |
---|
99 | } |
---|
100 | memcpy(seed, newseed, sizeof(longer)); |
---|
101 | seed[5] &= 3; |
---|
102 | x = 0.0; |
---|
103 | for (i = 0; i <= 5; i++) |
---|
104 | x = x / 64.0 + seed[i]; |
---|
105 | x /= 4.0; |
---|
106 | return x; |
---|
107 | } /* randum */ |
---|
108 | |
---|
109 | |
---|
110 | void copymatrix(tomat,frommat) |
---|
111 | transmatrix frommat,tomat; |
---|
112 | { |
---|
113 | int i,j; |
---|
114 | for (i=0;i<=sitelength;++i){ |
---|
115 | for (j=0;j<=sitelength;++j) |
---|
116 | tomat[i][j] = frommat[i][j]; |
---|
117 | } |
---|
118 | } |
---|
119 | |
---|
120 | void setuppi() |
---|
121 | { |
---|
122 | /* set up equilibrium probabilities of being a given |
---|
123 | number of bases away from a restriction site */ |
---|
124 | short i; |
---|
125 | double sum; |
---|
126 | |
---|
127 | pi[0] = 1.0; |
---|
128 | sum = pi[0]; |
---|
129 | for (i = 1; i <= sitelength; i++) { |
---|
130 | pi[i] = 3 * pi[i - 1] * (sitelength - i + 1) / i; |
---|
131 | sum += pi[i]; |
---|
132 | } |
---|
133 | for (i = 0; i <= sitelength; i++) |
---|
134 | pi[i] /= sum; |
---|
135 | } /* setuppi */ |
---|
136 | |
---|
137 | void maketrans(p) |
---|
138 | double p; |
---|
139 | { |
---|
140 | /* make transition matrix, product matrix with change |
---|
141 | probability p. Put the results in tempmatrix, tempprod */ |
---|
142 | short i, j, k, m1, m2; |
---|
143 | double sump, sumn, pover3, pijk, nijk, f; |
---|
144 | double binom1[maxcutter + 1], binom2[maxcutter + 1]; |
---|
145 | |
---|
146 | pover3 = p / 3; |
---|
147 | f = 2 * pover3 / (1 - pover3); |
---|
148 | for (i = 0; i <= sitelength; i++) { |
---|
149 | if (p > 1.0 - epsilon) |
---|
150 | p = 1.0 - epsilon; |
---|
151 | binom1[0] = exp((sitelength - i) * log(1 - p)); |
---|
152 | for (k = 1; k <= sitelength - i; k++) |
---|
153 | binom1[k] = binom1[k - 1] * (p / (1 - p)) * (sitelength - i - k + 1) / k; |
---|
154 | binom2[0] = exp(i * log(1 - pover3)); |
---|
155 | for (k = 1; k <= i; k++) |
---|
156 | binom2[k] = binom2[k - 1] * (pover3 / (1 - pover3)) * (i - k + 1) / k; |
---|
157 | for (j = 0; j <= sitelength; ++j) { |
---|
158 | sump = 0.0; |
---|
159 | sumn = 0.0; |
---|
160 | if (i - j > 0) |
---|
161 | m1 = i - j; |
---|
162 | else |
---|
163 | m1 = 0; |
---|
164 | if (sitelength - j < i) |
---|
165 | m2 = sitelength - j; |
---|
166 | else |
---|
167 | m2 = i; |
---|
168 | for (k = m1; k <= m2; k++) { |
---|
169 | nijk = j - i + k * 2 + f * (i - k); |
---|
170 | pijk = binom1[j - i + k] * binom2[k]; |
---|
171 | sump += pijk; |
---|
172 | sumn += pijk * nijk; |
---|
173 | } |
---|
174 | tempmatrix[i][j] = sump; |
---|
175 | tempprod[i][j] = sumn; |
---|
176 | } |
---|
177 | } |
---|
178 | } /* maketrans */ |
---|
179 | |
---|
180 | void branchtrans(i, p) |
---|
181 | short i; |
---|
182 | double p; |
---|
183 | { |
---|
184 | /* make branch transition matrix, product matrices for branch |
---|
185 | i with probability of change p*/ |
---|
186 | maketrans(p); |
---|
187 | copymatrix(curtree.trans[i - 1], tempmatrix); |
---|
188 | copymatrix(curtree.transprod[i - 1], tempprod); |
---|
189 | } /* branchtrans */ |
---|
190 | |
---|
191 | double evaluate(tr) |
---|
192 | tree *tr; |
---|
193 | { |
---|
194 | /* evaluates the likelihood, using info. at one branch */ |
---|
195 | double sum, sum2, y, liketerm, like0, lnlike0, term; |
---|
196 | short i, j, k; |
---|
197 | node *p, *q; |
---|
198 | sitelike x1, x2; |
---|
199 | |
---|
200 | sum = 0.0; |
---|
201 | p = tr->start; |
---|
202 | q = p->back; |
---|
203 | y = p->v; |
---|
204 | maketrans(y); |
---|
205 | memcpy(x1, p->x[0], sizeof(sitelike)); |
---|
206 | memcpy(x2, q->x[0], sizeof(sitelike)); |
---|
207 | if (trunc8) { |
---|
208 | like0 = 0.0; |
---|
209 | for (j = 0; j <= sitelength; j++) { |
---|
210 | liketerm = pi[j] * x1[j]; |
---|
211 | for (k = 0; k <= sitelength; k++) |
---|
212 | like0 += liketerm * tempmatrix[j][k] * x2[k]; |
---|
213 | } |
---|
214 | lnlike0 = log(enzymes * (1.0 - like0)); |
---|
215 | } |
---|
216 | for (i = 1; i <= endsite; i++) { |
---|
217 | memcpy(x1, p->x[i], sizeof(sitelike)); |
---|
218 | memcpy(x2, q->x[i], sizeof(sitelike)); |
---|
219 | sum2 = 0.0; |
---|
220 | for (j = 0; j <= sitelength; j++) { |
---|
221 | liketerm = pi[j] * x1[j]; |
---|
222 | for (k = 0; k <= sitelength; k++) |
---|
223 | sum2 += liketerm * tempmatrix[j][k] * x2[k]; |
---|
224 | } |
---|
225 | term = log(sum2); |
---|
226 | if (trunc8) |
---|
227 | term -= lnlike0; |
---|
228 | if (usertree && which <= maxtrees) |
---|
229 | l0gf[which - 1][i - 1] = term; |
---|
230 | sum += weight[i] * term; |
---|
231 | } |
---|
232 | if (usertree && which <= maxtrees) { |
---|
233 | l0gl[which - 1] = sum; |
---|
234 | if (which == 1) { |
---|
235 | maxwhich = 1; |
---|
236 | maxlogl = sum; |
---|
237 | } else if (sum > maxlogl) { |
---|
238 | maxwhich = which; |
---|
239 | maxlogl = sum; |
---|
240 | } |
---|
241 | } |
---|
242 | tr->likelihood = sum; |
---|
243 | return sum; |
---|
244 | } /* evaluate */ |
---|
245 | |
---|
246 | void nuview(p) |
---|
247 | node *p; |
---|
248 | { |
---|
249 | /* recompute fractional likelihoods for one part of tree */ |
---|
250 | short i, j, k, lowlim; |
---|
251 | double sumq, sumr; |
---|
252 | node *q, *r; |
---|
253 | sitelike xq, xr, xp; |
---|
254 | transmatrix tempq, tempr; |
---|
255 | tempq = (transmatrix)Malloc((sitelength+1) * sizeof(double *)); |
---|
256 | tempr = (transmatrix)Malloc((sitelength+1) * sizeof(double *)); |
---|
257 | for (i=0;i<=sitelength;++i){ |
---|
258 | tempq[i] = (double *)Malloc((sitelength+1) * sizeof (double)); |
---|
259 | tempr[i] = (double *)Malloc((sitelength+1) * sizeof (double)); |
---|
260 | } |
---|
261 | if (trunc8) |
---|
262 | lowlim = 0; |
---|
263 | else |
---|
264 | lowlim = 1; |
---|
265 | q = p->next->back; |
---|
266 | r = p->next->next->back; |
---|
267 | copymatrix(tempq,curtree.trans[q->branchnum - 1]); |
---|
268 | copymatrix(tempr,curtree.trans[r->branchnum - 1]); |
---|
269 | for (i = lowlim; i <= endsite; i++) { |
---|
270 | memcpy(xq, q->x[i], sizeof(sitelike)); |
---|
271 | memcpy(xr, r->x[i], sizeof(sitelike)); |
---|
272 | for (j = 0; j <= sitelength; j++) { |
---|
273 | sumq = 0.0; |
---|
274 | sumr = 0.0; |
---|
275 | for (k = 0; k <= sitelength; k++) { |
---|
276 | sumq += tempq[j][k] * xq[k]; |
---|
277 | sumr += tempr[j][k] * xr[k]; |
---|
278 | } |
---|
279 | xp[j] = sumq * sumr; |
---|
280 | } |
---|
281 | memcpy(p->x[i], xp, sizeof(sitelike)); |
---|
282 | } |
---|
283 | for (i=0;i<=sitelength;++i){ |
---|
284 | free(tempq[i]); |
---|
285 | free(tempr[i]); |
---|
286 | } |
---|
287 | free(tempq); |
---|
288 | free(tempr); |
---|
289 | } /* nuview */ |
---|
290 | |
---|
291 | |
---|
292 | void makenewv(p) |
---|
293 | node *p; |
---|
294 | { |
---|
295 | /* EM algorithm improvement of a branch length */ |
---|
296 | short i, j, k, lowlim, it, ite; |
---|
297 | double sum, sumlike, sumprod, liketerm, liket, y, yold, yorig, like0, |
---|
298 | prod0, like, oldlike, olderlike, extrap; |
---|
299 | boolean done; |
---|
300 | node *q; |
---|
301 | sitelike xx1, xx2; |
---|
302 | |
---|
303 | extrap = extrapol; |
---|
304 | q = p->back; |
---|
305 | y = p->v; |
---|
306 | yorig = y; |
---|
307 | if (trunc8) |
---|
308 | lowlim = 0; |
---|
309 | else |
---|
310 | lowlim = 1; |
---|
311 | done = false; |
---|
312 | oldlike = 0.0; |
---|
313 | olderlike = 0.0; |
---|
314 | it = 1; |
---|
315 | ite = 1; |
---|
316 | while ((it < iterations) && (!done)) { |
---|
317 | like = 0.0; |
---|
318 | maketrans(y); |
---|
319 | yold = y; |
---|
320 | sumlike = 0.0; |
---|
321 | for (i = lowlim; i <= endsite; i++) { |
---|
322 | memcpy(xx1, p->x[i], sizeof(sitelike)); |
---|
323 | memcpy(xx2, q->x[i], sizeof(sitelike)); |
---|
324 | sum = 0.0; |
---|
325 | sumprod = 0.0; |
---|
326 | for (j = 0; j <= sitelength; j++) { |
---|
327 | liket = xx1[j] * pi[j]; |
---|
328 | for (k = 0; k <= sitelength; k++) { |
---|
329 | liketerm = liket * xx2[k]; |
---|
330 | sumprod += tempprod[j][k] * liketerm; |
---|
331 | sum += tempmatrix[j][k] * liketerm; |
---|
332 | } |
---|
333 | } |
---|
334 | if (i == 0) { |
---|
335 | like0 = sum; |
---|
336 | prod0 = sumprod; |
---|
337 | } else { |
---|
338 | sumlike += weight[i] * sumprod / sum; |
---|
339 | like += weight[i] * log(sum); |
---|
340 | } |
---|
341 | } |
---|
342 | if (trunc8 && fabs(like0 - 1.0) > 1.0e-10) |
---|
343 | like -= weightsum * log(enzymes * (1.0 - like0)); |
---|
344 | y = sumlike / (sitelength * weightsum); |
---|
345 | if (trunc8) |
---|
346 | y = (1.0 - like0) * y + prod0 / sitelength; |
---|
347 | if (ite >= 3 && like > oldlike && like - oldlike < oldlike - olderlike) { |
---|
348 | extrap = (oldlike - olderlike) / (2.0 * oldlike - olderlike - like) * |
---|
349 | extrapol; |
---|
350 | ite = 1; |
---|
351 | } else |
---|
352 | extrap = extrapol; |
---|
353 | if (extrap < 0.0) |
---|
354 | extrap = extrapol; |
---|
355 | olderlike = oldlike; |
---|
356 | oldlike = like; |
---|
357 | y = extrap * (y - yold) + yold; |
---|
358 | if (y < epsilon) |
---|
359 | y = 10.0 * epsilon; |
---|
360 | if (y >= 0.75) |
---|
361 | y = 0.75; |
---|
362 | done = fabs(y-yorig) < epsilon; |
---|
363 | it++; |
---|
364 | ite++; |
---|
365 | } |
---|
366 | smoothed = (smoothed && done); |
---|
367 | p->v = y; |
---|
368 | q->v = y; |
---|
369 | branchtrans(p->branchnum, y); |
---|
370 | } /* makenewv */ |
---|
371 | |
---|
372 | void update(p) |
---|
373 | node *p; |
---|
374 | { |
---|
375 | /* improve branch length and views for one branch */ |
---|
376 | |
---|
377 | if (!p->tip) |
---|
378 | nuview(p); |
---|
379 | if (!p->back->tip) |
---|
380 | nuview(p->back); |
---|
381 | if (p->iter) |
---|
382 | makenewv(p); |
---|
383 | } /* update */ |
---|
384 | |
---|
385 | void smooth(p) |
---|
386 | node *p; |
---|
387 | { |
---|
388 | /* update nodes throughout the tree, recursively */ |
---|
389 | update(p); |
---|
390 | if (!p->tip) { |
---|
391 | smooth(p->next->back); |
---|
392 | smooth(p->next->next->back); |
---|
393 | update (p); |
---|
394 | } |
---|
395 | } /* smooth */ |
---|
396 | |
---|
397 | void hookup(p, q) |
---|
398 | node *p, *q; |
---|
399 | { |
---|
400 | /* connect two nodes */ |
---|
401 | p->back = q; |
---|
402 | q->back = p; |
---|
403 | } /* hookup */ |
---|
404 | |
---|
405 | void insert_(p, q) |
---|
406 | node *p, *q; |
---|
407 | { |
---|
408 | /* insert a subtree into a branch, improve lengths in tree */ |
---|
409 | short i, m, n; |
---|
410 | node *r; |
---|
411 | |
---|
412 | r = p->next->next; |
---|
413 | hookup(r, q->back); |
---|
414 | hookup(p->next, q); |
---|
415 | if (q->v >= 0.75) |
---|
416 | q->v = 0.75; |
---|
417 | else |
---|
418 | q->v = 0.75 * (1 - sqrt(1 - 1.333333 * q->v)); |
---|
419 | q->back->v = q->v; |
---|
420 | r->v = q->v; |
---|
421 | r->back->v = r->v; |
---|
422 | if (q->branchnum == q->number) { |
---|
423 | m = q->branchnum; |
---|
424 | n = r->number; |
---|
425 | } else { |
---|
426 | m = r->number; |
---|
427 | n = q->branchnum; |
---|
428 | } |
---|
429 | q->branchnum = m; |
---|
430 | q->back->branchnum = m; |
---|
431 | r->branchnum = n; |
---|
432 | r->back->branchnum = n; |
---|
433 | branchtrans(q->branchnum, q->v); |
---|
434 | branchtrans(r->branchnum, r->v); |
---|
435 | smoothed = false; |
---|
436 | i = 1; |
---|
437 | while (i < smoothings && !smoothed) { |
---|
438 | smoothed = true; |
---|
439 | smooth(p); |
---|
440 | smooth(p->back); |
---|
441 | i++; |
---|
442 | } |
---|
443 | } /* insert */ |
---|
444 | |
---|
445 | void re_move(p, q) |
---|
446 | node **p, **q; |
---|
447 | { |
---|
448 | /* remove p and record in q where it was */ |
---|
449 | *q = (*p)->next->back; |
---|
450 | hookup(*q, (*p)->next->next->back); |
---|
451 | (*q)->back->branchnum = (*q)->branchnum; |
---|
452 | branchtrans((*q)->branchnum, 0.75*(1 - (1 - 1.333333*(*q)->v) |
---|
453 | * (1 - 1.333333*(*p)->next->v))); |
---|
454 | (*p)->next->back = NULL; |
---|
455 | (*p)->next->next->back = NULL; |
---|
456 | update(*q); |
---|
457 | update((*q)->back); |
---|
458 | } /* re_move */ |
---|
459 | |
---|
460 | void copynode(c, d) |
---|
461 | node *c, *d; |
---|
462 | { |
---|
463 | /* copy a node */ |
---|
464 | |
---|
465 | d->branchnum = c->branchnum; |
---|
466 | memcpy(d->x, c->x, (endsite+1)*sizeof(sitelike)); |
---|
467 | memcpy(d->nayme, c->nayme, sizeof(naym)); |
---|
468 | d->v = c->v; |
---|
469 | d->iter = c->iter; |
---|
470 | d->xcoord = c->xcoord; |
---|
471 | d->ycoord = c->ycoord; |
---|
472 | d->ymin = c->ymin; |
---|
473 | d->ymax = c->ymax; |
---|
474 | } /* copynode */ |
---|
475 | |
---|
476 | void copy_(a, b) |
---|
477 | tree *a, *b; |
---|
478 | { |
---|
479 | /* copy a tree */ |
---|
480 | short i,j; |
---|
481 | node *p, *q; |
---|
482 | |
---|
483 | for (i = 0; i < numsp; i++) { |
---|
484 | copynode(a->nodep[i], b->nodep[i]); |
---|
485 | if (a->nodep[i]->back) { |
---|
486 | if (a->nodep[i]->back == a->nodep[a->nodep[i]->back->number - 1]) |
---|
487 | b->nodep[i]->back = b->nodep[a->nodep[i]->back->number - 1]; |
---|
488 | else if (a->nodep[i]->back |
---|
489 | == a->nodep[a->nodep[i]->back->number - 1]->next) |
---|
490 | b->nodep[i]->back = b->nodep[a->nodep[i]->back->number - 1]->next; |
---|
491 | else |
---|
492 | b->nodep[i]->back |
---|
493 | = b->nodep[a->nodep[i]->back->number - 1]->next->next; |
---|
494 | } |
---|
495 | else b->nodep[i]->back = NULL; |
---|
496 | } |
---|
497 | for (i = numsp; i < numsp2; i++) { |
---|
498 | p = a->nodep[i]; |
---|
499 | q = b->nodep[i]; |
---|
500 | for (j = 1; j <= 3; j++) { |
---|
501 | copynode(p, q); |
---|
502 | if (p->back) { |
---|
503 | if (p->back == a->nodep[p->back->number - 1]) |
---|
504 | q->back = b->nodep[p->back->number - 1]; |
---|
505 | else if (p->back == a->nodep[p->back->number - 1]->next) |
---|
506 | q->back = b->nodep[p->back->number - 1]->next; |
---|
507 | else |
---|
508 | q->back = b->nodep[p->back->number - 1]->next->next; |
---|
509 | } |
---|
510 | else |
---|
511 | q->back = NULL; |
---|
512 | p = p->next; |
---|
513 | q = q->next; |
---|
514 | } |
---|
515 | } |
---|
516 | b->likelihood = a->likelihood; |
---|
517 | for (i=0;i<numsp2;++i){ |
---|
518 | copymatrix(b->trans[i],a->trans[i]); |
---|
519 | copymatrix(b->transprod[i],a->transprod[i]); |
---|
520 | } |
---|
521 | b->start = a->start; |
---|
522 | } /* copy */ |
---|
523 | |
---|
524 | void buildnewtip(m, tr) |
---|
525 | short m; |
---|
526 | tree *tr; |
---|
527 | { |
---|
528 | /* set up a new tip and interior node it is connected to */ |
---|
529 | node *p; |
---|
530 | |
---|
531 | p = tr->nodep[nextsp + numsp - 3]; |
---|
532 | hookup(tr->nodep[m - 1], p); |
---|
533 | p->v = initialv; |
---|
534 | p->back->v = initialv; |
---|
535 | branchtrans(m, initialv); |
---|
536 | p->branchnum = m; |
---|
537 | p->next->branchnum = p->number; |
---|
538 | p->next->next->branchnum = p->number; |
---|
539 | p->back->branchnum = m; |
---|
540 | } /* buildnewtip */ |
---|
541 | |
---|
542 | void buildsimpletree(tr) |
---|
543 | tree *tr; |
---|
544 | { |
---|
545 | /* set up and adjust branch lengths of a three-species tree */ |
---|
546 | hookup(tr->nodep[enterorder[0] - 1], tr->nodep[enterorder[1] - 1]); |
---|
547 | tr->nodep[enterorder[0] - 1]->v = initialv; |
---|
548 | tr->nodep[enterorder[1] - 1]->v = initialv; |
---|
549 | branchtrans(enterorder[1], initialv); |
---|
550 | tr->nodep[enterorder[0] - 1]->branchnum = 2; |
---|
551 | tr->nodep[enterorder[1] - 1]->branchnum = 2; |
---|
552 | buildnewtip(enterorder[2], tr); |
---|
553 | insert_(tr->nodep[enterorder[2] - 1]->back, tr->nodep[enterorder[1] - 1]); |
---|
554 | } /* buildsimpletree */ |
---|
555 | |
---|
556 | void addtraverse(p, q, contin) |
---|
557 | node *p, *q; |
---|
558 | boolean contin; |
---|
559 | { |
---|
560 | int i; |
---|
561 | double like; |
---|
562 | /* try adding p at q, proceed recursively through tree */ |
---|
563 | insert_(p, q); |
---|
564 | numtrees++; |
---|
565 | like = evaluate(&curtree); |
---|
566 | if (like > bestree.likelihood) { |
---|
567 | copy_(&curtree, &bestree); |
---|
568 | bestree.likelihood = like; |
---|
569 | succeeded = true; |
---|
570 | } |
---|
571 | copy_(&priortree, &curtree); |
---|
572 | if (!q->tip && contin) { |
---|
573 | addtraverse(p, q->next->back, contin); |
---|
574 | addtraverse(p, q->next->next->back, contin); |
---|
575 | } |
---|
576 | } /* addtraverse */ |
---|
577 | |
---|
578 | void rearrange(p) |
---|
579 | node *p; |
---|
580 | { |
---|
581 | /* rearranges the tree, globally or locally */ |
---|
582 | node *q, *r; |
---|
583 | |
---|
584 | if (!p->tip && !p->back->tip) { |
---|
585 | r = p->next->next; |
---|
586 | re_move(&r, &q); |
---|
587 | copy_(&curtree, &priortree); |
---|
588 | addtraverse(r, q->next->back, global && nextsp == numsp); |
---|
589 | addtraverse(r, q->next->next->back, global && nextsp == numsp); |
---|
590 | copy_(&bestree, &curtree); |
---|
591 | if (global && nextsp == numsp && progress) |
---|
592 | putchar('.'); |
---|
593 | if (global && nextsp == numsp && !succeeded) { |
---|
594 | if (r->back->tip) { |
---|
595 | r = r->next->next; |
---|
596 | re_move(&r, &q); |
---|
597 | q = q->back; |
---|
598 | if (!q->tip) { |
---|
599 | copy_(&curtree, &priortree); |
---|
600 | addtraverse(r, q->next->back, true); |
---|
601 | copy_(&curtree, &priortree); |
---|
602 | addtraverse(r, q->next->next->back, true); |
---|
603 | } |
---|
604 | q = q->back; |
---|
605 | if (!q->tip) { |
---|
606 | copy_(&curtree, &priortree); |
---|
607 | addtraverse(r, q->next->back, true); |
---|
608 | copy_(&curtree, &priortree); |
---|
609 | addtraverse(r, q->next->next->back, true); |
---|
610 | } |
---|
611 | copy_(&bestree, &curtree); |
---|
612 | } |
---|
613 | } |
---|
614 | } |
---|
615 | if (!p->tip) { |
---|
616 | rearrange(p->next->back); |
---|
617 | rearrange(p->next->next->back); |
---|
618 | } |
---|
619 | } /* rearrange */ |
---|
620 | |
---|
621 | |
---|
622 | void coordinates(p, lengthsum, tipy,tipmax,x) |
---|
623 | node *p; |
---|
624 | double lengthsum; |
---|
625 | short *tipy; |
---|
626 | double *tipmax, *x; |
---|
627 | |
---|
628 | { |
---|
629 | /* establishes coordinates of nodes */ |
---|
630 | node *q, *first, *last; |
---|
631 | |
---|
632 | if (p->tip) { |
---|
633 | p->xcoord = (short)(over * lengthsum + 0.5); |
---|
634 | p->ycoord = (*tipy); |
---|
635 | p->ymin = (*tipy); |
---|
636 | p->ymax = (*tipy); |
---|
637 | (*tipy) += down; |
---|
638 | if (lengthsum > (*tipmax)) |
---|
639 | (*tipmax) = lengthsum; |
---|
640 | return; |
---|
641 | } |
---|
642 | q = p->next; |
---|
643 | do { |
---|
644 | (*x) = -0.75 * log(1.0 - 1.333333 * q->v); |
---|
645 | coordinates(q->back, lengthsum + (*x),tipy,tipmax,x); |
---|
646 | q = q->next; |
---|
647 | } while ((p == curtree.start->back || p != q) && |
---|
648 | (p != curtree.start->back || p->next != q)); |
---|
649 | first = p->next->back; |
---|
650 | q = p; |
---|
651 | while (q->next != p) |
---|
652 | q = q->next; |
---|
653 | last = q->back; |
---|
654 | p->xcoord = (short)(over * lengthsum + 0.5); |
---|
655 | if (p == curtree.start->back) |
---|
656 | p->ycoord = p->next->next->back->ycoord; |
---|
657 | else |
---|
658 | p->ycoord = (first->ycoord + last->ycoord) / 2; |
---|
659 | p->ymin = first->ymin; |
---|
660 | p->ymax = last->ymax; |
---|
661 | } /* coordinates */ |
---|
662 | |
---|
663 | void drawline(i, scale) |
---|
664 | short i; |
---|
665 | double scale; |
---|
666 | { |
---|
667 | /* draws one row of the tree diagram by moving up tree */ |
---|
668 | node *p, *q; |
---|
669 | short n, j; |
---|
670 | boolean extra; |
---|
671 | node *r, *first, *last; |
---|
672 | boolean done; |
---|
673 | |
---|
674 | p = curtree.start->back; |
---|
675 | q = curtree.start->back; |
---|
676 | extra = false; |
---|
677 | if (i == p->ycoord && p == curtree.start->back) { |
---|
678 | if (p->number - numsp >= 10) |
---|
679 | fprintf(outfile, "-%2hd", p->number - numsp); |
---|
680 | else |
---|
681 | fprintf(outfile, "--%hd", p->number - numsp); |
---|
682 | extra = true; |
---|
683 | } else |
---|
684 | fprintf(outfile, " "); |
---|
685 | do { |
---|
686 | if (!p->tip) { |
---|
687 | r = p->next; |
---|
688 | done = false; |
---|
689 | do { |
---|
690 | if (i >= r->back->ymin && i <= r->back->ymax) { |
---|
691 | q = r->back; |
---|
692 | done = true; |
---|
693 | } |
---|
694 | r = r->next; |
---|
695 | } while (!(done || p != curtree.start->back && r == p || |
---|
696 | p == curtree.start->back && r == p->next)); |
---|
697 | first = p->next->back; |
---|
698 | r = p; |
---|
699 | while (r->next != p) |
---|
700 | r = r->next; |
---|
701 | last = r->back; |
---|
702 | if (p == curtree.start->back) |
---|
703 | last = p->back; |
---|
704 | } |
---|
705 | done = (p->tip || p == q); |
---|
706 | n = (short)(scale * (q->xcoord - p->xcoord) + 0.5); |
---|
707 | if (n < 3 && !q->tip) |
---|
708 | n = 3; |
---|
709 | if (extra) { |
---|
710 | n--; |
---|
711 | extra = false; |
---|
712 | } |
---|
713 | if (q->ycoord == i && !done) { |
---|
714 | if (p->ycoord != q->ycoord) |
---|
715 | putc('+', outfile); |
---|
716 | else |
---|
717 | putc('-', outfile); |
---|
718 | if (!q->tip) { |
---|
719 | for (j = 1; j <= n - 2; j++) |
---|
720 | putc('-', outfile); |
---|
721 | if (q->number - numsp >= 10) |
---|
722 | fprintf(outfile, "%2hd", q->number - numsp); |
---|
723 | else |
---|
724 | fprintf(outfile, "-%hd", q->number - numsp); |
---|
725 | extra = true; |
---|
726 | } else { |
---|
727 | for (j = 1; j < n; j++) |
---|
728 | putc('-', outfile); |
---|
729 | } |
---|
730 | } else if (!p->tip) { |
---|
731 | if (last->ycoord > i && first->ycoord < i && i != p->ycoord) { |
---|
732 | putc('!', outfile); |
---|
733 | for (j = 1; j < n; j++) |
---|
734 | putc(' ', outfile); |
---|
735 | } else { |
---|
736 | for (j = 1; j <= n; j++) |
---|
737 | putc(' ', outfile); |
---|
738 | } |
---|
739 | } else { |
---|
740 | for (j = 1; j <= n; j++) |
---|
741 | putc(' ', outfile); |
---|
742 | } |
---|
743 | if (q != p) |
---|
744 | p = q; |
---|
745 | } while (!done); |
---|
746 | if (p->ycoord == i && p->tip) { |
---|
747 | for (j = 0; j < nmlngth; j++) |
---|
748 | putc(p->nayme[j], outfile); |
---|
749 | } |
---|
750 | putc('\n', outfile); |
---|
751 | } /* drawline */ |
---|
752 | |
---|
753 | void printree() |
---|
754 | { |
---|
755 | /* prints out diagram of the tree */ |
---|
756 | short tipy,i; |
---|
757 | double scale, tipmax, x; |
---|
758 | |
---|
759 | putc('\n', outfile); |
---|
760 | if (!treeprint) |
---|
761 | return; |
---|
762 | putc('\n', outfile); |
---|
763 | tipy = 1; |
---|
764 | tipmax = 0.0; |
---|
765 | coordinates(curtree.start->back, 0.0, &tipy,&tipmax,&x); |
---|
766 | scale = 1.0 / (tipmax + 1.000); |
---|
767 | for (i = 1; i <= tipy - down; i++) |
---|
768 | drawline(i, scale); |
---|
769 | putc('\n', outfile); |
---|
770 | } /* printree */ |
---|
771 | |
---|
772 | |
---|
773 | double sigma(q,sumlr) |
---|
774 | node *q; |
---|
775 | double *sumlr; |
---|
776 | { |
---|
777 | /* get 1.95996 * approximate standard error of branch length */ |
---|
778 | double sump, sumr, sums, sumc, p, pover3, pijk, Qjk, liketerm, f; |
---|
779 | double slopef,curvef; |
---|
780 | short i, j, k, m1, m2; |
---|
781 | double binom1[maxcutter + 1], binom2[maxcutter + 1]; |
---|
782 | transmatrix Prob, slopeP, curveP; |
---|
783 | node *r; |
---|
784 | sitelike x1, x2; |
---|
785 | double TEMP, TEMP1; |
---|
786 | Prob = (transmatrix)Malloc((sitelength+1) * sizeof(double *)); |
---|
787 | slopeP = (transmatrix)Malloc((sitelength+1) * sizeof(double *)); |
---|
788 | curveP = (transmatrix)Malloc((sitelength+1) * sizeof(double *)); |
---|
789 | for (i=0;i<=sitelength;++i){ |
---|
790 | Prob[i] = (double *)Malloc((sitelength+1) * sizeof(double)); |
---|
791 | slopeP[i] = (double *)Malloc((sitelength+1) * sizeof(double)); |
---|
792 | curveP[i] = (double *)Malloc((sitelength+1) * sizeof(double)); |
---|
793 | } |
---|
794 | p = q->v; |
---|
795 | pover3 = p / 3; |
---|
796 | for (i = 0; i <= sitelength; i++) { |
---|
797 | binom1[0] = exp((sitelength - i) * log(1 - p)); |
---|
798 | for (k = 1; k <= (sitelength - i); k++) |
---|
799 | binom1[k] = binom1[k - 1] * (p / (1 - p)) * (sitelength - i - k + 1) / k; |
---|
800 | binom2[0] = exp(i * log(1 - pover3)); |
---|
801 | for (k = 1; k <= i; k++) |
---|
802 | binom2[k] = binom2[k - 1] * (pover3 / (1 - pover3)) * (i - k + 1) / k; |
---|
803 | for (j = 0; j <= sitelength; j++) { |
---|
804 | sump = 0.0; |
---|
805 | sums = 0.0; |
---|
806 | sumc = 0.0; |
---|
807 | if (i - j > 0) |
---|
808 | m1 = i - j; |
---|
809 | else |
---|
810 | m1 = 0; |
---|
811 | if (sitelength - j < i) |
---|
812 | m2 = sitelength - j; |
---|
813 | else |
---|
814 | m2 = i; |
---|
815 | for (k = m1; k <= m2; k++) { |
---|
816 | pijk = binom1[j - i + k] * binom2[k]; |
---|
817 | sump += pijk; |
---|
818 | sums += pijk * ((j - i + k * 2) / p + (j - sitelength + k) / (1 - p) + |
---|
819 | (k - i) / (3 - p)); |
---|
820 | TEMP = 1 - p; |
---|
821 | TEMP1 = 3 - p; |
---|
822 | sumc += pijk * ((i - j - k * 2) / (p * p) + (j - sitelength + k) / |
---|
823 | (TEMP * TEMP) + (k - i) / (TEMP1 * TEMP1)); |
---|
824 | } |
---|
825 | Prob[i][j] = sump; |
---|
826 | slopeP[i][j] = sums; |
---|
827 | curveP[i][j] = sumc; |
---|
828 | } |
---|
829 | } |
---|
830 | (*sumlr) = 0.0; |
---|
831 | sumc = 0.0; |
---|
832 | sums = 0.0; |
---|
833 | r = q->back; |
---|
834 | for (i = 1; i <= endsite; i++) { |
---|
835 | f = 0.0; |
---|
836 | slopef = 0.0; |
---|
837 | curvef = 0.0; |
---|
838 | sumr = 0.0; |
---|
839 | memcpy(x1, q->x[i], sizeof(sitelike)); |
---|
840 | memcpy(x2, r->x[i], sizeof(sitelike)); |
---|
841 | for (j = 0; j <= sitelength; j++) { |
---|
842 | liketerm = pi[j] * x1[j]; |
---|
843 | sumr += liketerm * x2[j]; |
---|
844 | for (k = 0; k <= sitelength; k++) { |
---|
845 | Qjk = liketerm * x2[k]; |
---|
846 | f += Qjk * Prob[j][k]; |
---|
847 | slopef += Qjk * slopeP[j][k]; |
---|
848 | curvef += Qjk * curveP[j][k]; |
---|
849 | } |
---|
850 | } |
---|
851 | (*sumlr) += weight[i] * log(f / sumr); |
---|
852 | sums += weight[i] * slopef / f; |
---|
853 | TEMP = slopef / f; |
---|
854 | sumc += weight[i] * (curvef / f - TEMP * TEMP); |
---|
855 | } |
---|
856 | if (trunc8) { |
---|
857 | f = 0.0; |
---|
858 | slopef = 0.0; |
---|
859 | curvef = 0.0; |
---|
860 | sumr = 0.0; |
---|
861 | memcpy(x1, q->x[0], sizeof(sitelike)); |
---|
862 | memcpy(x2, r->x[0], sizeof(sitelike)); |
---|
863 | for (j = 0; j <= sitelength; j++) { |
---|
864 | liketerm = pi[j] * x1[j]; |
---|
865 | sumr += liketerm * x2[j]; |
---|
866 | for (k = 0; k <= sitelength; k++) { |
---|
867 | Qjk = liketerm * x2[k]; |
---|
868 | f += Qjk * Prob[j][k]; |
---|
869 | slopef += Qjk * slopeP[j][k]; |
---|
870 | curvef += Qjk * curveP[j][k]; |
---|
871 | } |
---|
872 | } |
---|
873 | (*sumlr) += weightsum * log((1.0 - sumr) / (1.0 - f)); |
---|
874 | sums += weightsum * slopef / (1.0 - f); |
---|
875 | TEMP = slopef / (1.0 - f); |
---|
876 | sumc += weightsum * (curvef / (1.0 - f) + TEMP * TEMP); |
---|
877 | } |
---|
878 | for (i=0;i<sitelength;++i){ |
---|
879 | free(Prob[i]); |
---|
880 | free(slopeP[i]); |
---|
881 | free(curveP[i]); |
---|
882 | } |
---|
883 | free(Prob); |
---|
884 | free(slopeP); |
---|
885 | free(curveP); |
---|
886 | if (sumc < -1.0e-6) |
---|
887 | return ((-sums - sqrt(sums * sums - 3.841 * sumc)) / sumc); |
---|
888 | else |
---|
889 | return -1.0; |
---|
890 | } /* sigma */ |
---|
891 | |
---|
892 | void describe(p) |
---|
893 | node *p; |
---|
894 | { |
---|
895 | /* print out information on one branch */ |
---|
896 | double sumlr; |
---|
897 | short i; |
---|
898 | node *q; |
---|
899 | double s; |
---|
900 | |
---|
901 | q = p->back; |
---|
902 | fprintf(outfile, "%4hd ", q->number - numsp); |
---|
903 | fprintf(outfile, " "); |
---|
904 | if (p->tip) { |
---|
905 | for (i = 0; i < nmlngth; i++) |
---|
906 | putc(p->nayme[i], outfile); |
---|
907 | } else |
---|
908 | fprintf(outfile, "%4hd ", p->number - numsp); |
---|
909 | if (q->v >= 0.75) |
---|
910 | fprintf(outfile, " infinity"); |
---|
911 | else |
---|
912 | fprintf(outfile, "%13.5f", -0.75 * log(1 - 1.333333 * q->v)); |
---|
913 | if (p->iter) { |
---|
914 | s = sigma(q, &sumlr); |
---|
915 | if (s < 0.0) |
---|
916 | fprintf(outfile, " ( zero, infinity)"); |
---|
917 | else { |
---|
918 | fprintf(outfile, " ("); |
---|
919 | if (q->v - s <= 0.0) |
---|
920 | fprintf(outfile, " zero"); |
---|
921 | else |
---|
922 | fprintf(outfile, "%9.5f", -0.75 * log(1 - 1.333333 * (q->v - s))); |
---|
923 | putc(',', outfile); |
---|
924 | if (q->v + s >= 0.75) |
---|
925 | fprintf(outfile, " infinity"); |
---|
926 | else |
---|
927 | fprintf(outfile, "%12.5f", -0.75 * log(1 - 1.333333 * (q->v + s))); |
---|
928 | putc(')', outfile); |
---|
929 | } |
---|
930 | if (sumlr > 1.9205) |
---|
931 | fprintf(outfile, " *"); |
---|
932 | if (sumlr > 2.995) |
---|
933 | putc('*', outfile); |
---|
934 | } |
---|
935 | else |
---|
936 | fprintf(outfile, " (not varied)"); |
---|
937 | putc('\n', outfile); |
---|
938 | if (!p->tip) { |
---|
939 | describe(p->next->back); |
---|
940 | describe(p->next->next->back); |
---|
941 | } |
---|
942 | } /* describe */ |
---|
943 | |
---|
944 | void summarize() |
---|
945 | { |
---|
946 | /* print out information on branches of tree */ |
---|
947 | |
---|
948 | fprintf(outfile, "\nremember: "); |
---|
949 | if (outgropt) |
---|
950 | fprintf(outfile, "(although rooted by outgroup) "); |
---|
951 | fprintf(outfile, "this is an unrooted tree!\n\n"); |
---|
952 | fprintf(outfile, "Ln Likelihood = %11.5f\n\n", curtree.likelihood); |
---|
953 | if (!usertree) |
---|
954 | fprintf(outfile, "Examined %4hd trees\n", numtrees); |
---|
955 | fprintf(outfile, " \n"); |
---|
956 | fprintf(outfile, " Between And Length"); |
---|
957 | fprintf(outfile, " Approx. Confidence Limits\n"); |
---|
958 | fprintf(outfile, " ------- --- ------"); |
---|
959 | fprintf(outfile, " ------- ---------- ------\n"); |
---|
960 | describe(curtree.start->back->next->back); |
---|
961 | describe(curtree.start->back->next->next->back); |
---|
962 | describe(curtree.start); |
---|
963 | fprintf(outfile, "\n * = significantly positive, P < 0.05\n"); |
---|
964 | fprintf(outfile, " ** = significantly positive, P < 0.01\n\n\n"); |
---|
965 | } /* summarize */ |
---|
966 | |
---|
967 | void treeout(p) |
---|
968 | node *p; |
---|
969 | { |
---|
970 | /* write out file with representation of final tree */ |
---|
971 | short i, n, w; |
---|
972 | Char c; |
---|
973 | double x; |
---|
974 | |
---|
975 | if (p->tip) { |
---|
976 | n = 0; |
---|
977 | for (i = 1; i <= nmlngth; i++) { |
---|
978 | if (p->nayme[i - 1] != ' ') |
---|
979 | n = i; |
---|
980 | } |
---|
981 | for (i = 0; i < n; i++) { |
---|
982 | c = p->nayme[i]; |
---|
983 | if (c == ' ') |
---|
984 | c = '_'; |
---|
985 | putc(c, treefile); |
---|
986 | } |
---|
987 | col += n; |
---|
988 | } else { |
---|
989 | putc('(', treefile); |
---|
990 | col++; |
---|
991 | treeout(p->next->back); |
---|
992 | putc(',', treefile); |
---|
993 | col++; |
---|
994 | if (col > 45) { |
---|
995 | putc('\n', treefile); |
---|
996 | col = 0; |
---|
997 | } |
---|
998 | treeout(p->next->next->back); |
---|
999 | if (p == curtree.start->back) { |
---|
1000 | putc(',', treefile); |
---|
1001 | col++; |
---|
1002 | if (col > 45) { |
---|
1003 | putc('\n', treefile); |
---|
1004 | col = 0; |
---|
1005 | } |
---|
1006 | treeout(p->back); |
---|
1007 | } |
---|
1008 | putc(')', treefile); |
---|
1009 | col++; |
---|
1010 | } |
---|
1011 | if (p->v >= 0.75) |
---|
1012 | x = -1.0; |
---|
1013 | else |
---|
1014 | x = -0.75 * log(1 - 1.333333 * p->v); |
---|
1015 | if (x > 0.0) |
---|
1016 | w = (short)(0.43429448222 * log(x)); |
---|
1017 | else if (x == 0.0) |
---|
1018 | w = 0; |
---|
1019 | else |
---|
1020 | w = (short)(0.43429448222 * log(-x)) + 1; |
---|
1021 | if (w < 0) |
---|
1022 | w = 0; |
---|
1023 | if (p == curtree.start->back) |
---|
1024 | fprintf(treefile, ";\n"); |
---|
1025 | else { |
---|
1026 | fprintf(treefile, ":%*.5f", (int)(w + 7), x); |
---|
1027 | col += w + 8; |
---|
1028 | } |
---|
1029 | } /* treeout */ |
---|
1030 | |
---|
1031 | |
---|
1032 | void getch(c) |
---|
1033 | Char *c; |
---|
1034 | { |
---|
1035 | /* get next nonblank character */ |
---|
1036 | do { |
---|
1037 | if (eoln(infile)) { |
---|
1038 | fscanf(infile, "%*[^\n]"); |
---|
1039 | getc(infile); |
---|
1040 | } |
---|
1041 | *c = getc(infile); |
---|
1042 | if (*c == '\n') |
---|
1043 | *c = ' '; |
---|
1044 | } while (*c == ' '); |
---|
1045 | } /* getch */ |
---|
1046 | |
---|
1047 | void findch(c,lparens,rparens) |
---|
1048 | Char c; |
---|
1049 | short *lparens,*rparens; |
---|
1050 | { |
---|
1051 | /* read forward to next occurrence of character c */ |
---|
1052 | boolean done; |
---|
1053 | |
---|
1054 | done = false; |
---|
1055 | while (!done) { |
---|
1056 | if (c == ',') { |
---|
1057 | if (ch == '(' || ch == ')' || ch == ':' || ch == ';') { |
---|
1058 | printf("\nERROR IN USER TREE: UNMATCHED PARENTHESIS"); |
---|
1059 | printf("OR MISSING COMMA\n OR NOT TRIFURCATED BASE\n"); |
---|
1060 | exit(-1); |
---|
1061 | } else if (ch == ',') |
---|
1062 | done = true; |
---|
1063 | } else if (c == ')') { |
---|
1064 | if (ch == '(' || ch == ',' || ch == ':' || ch == ';') { |
---|
1065 | printf( "\nERROR IN USER TREE:"); |
---|
1066 | printf(" UNMATCHED PARENTHESIS OR NOT BIFURCATED NODE\n"); |
---|
1067 | exit(-1); |
---|
1068 | } else if (ch == ')') { |
---|
1069 | (*rparens)++; |
---|
1070 | if ((*lparens) > 0 && (*lparens) == (*rparens)) { |
---|
1071 | if ((*lparens) == numsp - 2) { |
---|
1072 | if (eoln(infile)) { |
---|
1073 | fscanf(infile, "%*[^\n]"); |
---|
1074 | getc(infile); |
---|
1075 | } |
---|
1076 | ch = getc(infile); |
---|
1077 | if (ch == '\n') |
---|
1078 | ch = ' '; |
---|
1079 | if (ch != ';') { |
---|
1080 | printf("\nERROR IN USER TREE:"); |
---|
1081 | printf(" UNMATCHED PARENTHESIS OR MISSING SEMICOLON\n"); |
---|
1082 | exit(-1); |
---|
1083 | } |
---|
1084 | } |
---|
1085 | } |
---|
1086 | done = true; |
---|
1087 | } |
---|
1088 | } |
---|
1089 | if (done && ch != ')') |
---|
1090 | continue; |
---|
1091 | if (eoln(infile)) { |
---|
1092 | fscanf(infile, "%*[^\n]"); |
---|
1093 | getc(infile); |
---|
1094 | } |
---|
1095 | ch = getc(infile); |
---|
1096 | if (ch == '\n') |
---|
1097 | ch = ' '; |
---|
1098 | } |
---|
1099 | } /* findch */ |
---|
1100 | |
---|
1101 | void processlength(p) |
---|
1102 | node *p; |
---|
1103 | { |
---|
1104 | short digit, ordzero; |
---|
1105 | double valyew, divisor; |
---|
1106 | boolean pointread; |
---|
1107 | |
---|
1108 | ordzero = '0'; |
---|
1109 | pointread = false; |
---|
1110 | valyew = 0.0; |
---|
1111 | divisor = 1.0; |
---|
1112 | getch(&ch); |
---|
1113 | digit = ch - ordzero; |
---|
1114 | while (((unsigned short)digit <= 9) | ch == '.') { |
---|
1115 | if (ch == '.') |
---|
1116 | pointread = true; |
---|
1117 | else { |
---|
1118 | valyew = valyew * 10.0 + digit; |
---|
1119 | if (pointread) |
---|
1120 | divisor *= 10.0; |
---|
1121 | } |
---|
1122 | getch(&ch); |
---|
1123 | digit = ch - ordzero; |
---|
1124 | } |
---|
1125 | if (lengths) { |
---|
1126 | p->v = 0.75*(1.0-exp(-1.333333 * valyew / divisor)); |
---|
1127 | p->back->v = p->v; |
---|
1128 | p->iter = false; |
---|
1129 | p->back->iter = false; |
---|
1130 | } |
---|
1131 | } /* processlength */ |
---|
1132 | |
---|
1133 | |
---|
1134 | void addelement(p, nextnode,lparens,rparens,names,nolengths) |
---|
1135 | node *p; |
---|
1136 | short *nextnode,*lparens,*rparens; |
---|
1137 | boolean *names, *nolengths; |
---|
1138 | |
---|
1139 | { |
---|
1140 | /* add one node to user-defined tree */ |
---|
1141 | node *q; |
---|
1142 | short i, j, n; |
---|
1143 | boolean found; |
---|
1144 | Char str[nmlngth]; |
---|
1145 | |
---|
1146 | do { |
---|
1147 | if (eoln(infile)) { |
---|
1148 | fscanf(infile, "%*[^\n]"); |
---|
1149 | getc(infile); |
---|
1150 | } |
---|
1151 | ch = getc(infile); |
---|
1152 | if (ch == '\n') |
---|
1153 | ch = ' '; |
---|
1154 | } while (ch == ' '); |
---|
1155 | if (ch == '(') { |
---|
1156 | (*lparens)++; |
---|
1157 | if ((*lparens) > ( numsp - 2 )) { |
---|
1158 | printf("\nERROR IN USER TREE: TOO MANY LEFT PARENTHESES\n"); |
---|
1159 | exit(-1); |
---|
1160 | } else { |
---|
1161 | (*nextnode)++; |
---|
1162 | q = curtree.nodep[(*nextnode) - 1]; |
---|
1163 | n = (*nextnode); |
---|
1164 | hookup(p, q); |
---|
1165 | addelement(q->next, nextnode,lparens,rparens,names,nolengths); |
---|
1166 | findch(',',lparens,rparens); |
---|
1167 | addelement(q->next->next, nextnode,lparens,rparens,names,nolengths); |
---|
1168 | if (*nolengths && lengths) |
---|
1169 | printf("NO LENGTHS FOUND IN INPUT FILE WITH LENGTH OPTION CHOSEN"); |
---|
1170 | findch(')',lparens,rparens); |
---|
1171 | for (i = 0; i <= endsite; i++) { |
---|
1172 | for (j = 0; j <= sitelength; j++) { |
---|
1173 | q->x[i][j] = 1.0; |
---|
1174 | q->next->x[i][j] = 1.0; |
---|
1175 | q->next->next->x[i][j] = 1.0; |
---|
1176 | } |
---|
1177 | } |
---|
1178 | } |
---|
1179 | } else { |
---|
1180 | for (i = 0; i < nmlngth; i++) |
---|
1181 | str[i] = ' '; |
---|
1182 | n = 1; |
---|
1183 | do { |
---|
1184 | if (ch == '_') |
---|
1185 | ch = ' '; |
---|
1186 | str[n - 1] = ch; |
---|
1187 | if (eoln(infile)) { |
---|
1188 | fscanf(infile, "%*[^\n]"); |
---|
1189 | getc(infile); |
---|
1190 | } |
---|
1191 | ch = getc(infile); |
---|
1192 | if (ch == '\n') |
---|
1193 | ch = ' '; |
---|
1194 | n++; |
---|
1195 | } while (ch != ':' && ch != ',' && ch != ')' && n <= nmlngth); |
---|
1196 | n = 1; |
---|
1197 | do { |
---|
1198 | found = true; |
---|
1199 | for (i = 0; i < nmlngth; i++) |
---|
1200 | found = (found && str[i] == curtree.nodep[n - 1]->nayme[i]); |
---|
1201 | if (found) { |
---|
1202 | q = curtree.nodep[n - 1]; |
---|
1203 | if (names[n - 1] == false) |
---|
1204 | names[n - 1] = true; |
---|
1205 | else { |
---|
1206 | printf("\nERROR IN USER TREE: DUPLICATE NAME FOUND -- "); |
---|
1207 | for (i = 0; i < nmlngth; i++) |
---|
1208 | putchar(curtree.nodep[n - 1]->nayme[i]); |
---|
1209 | putchar('\n'); |
---|
1210 | exit(-1); |
---|
1211 | } |
---|
1212 | } else |
---|
1213 | n++; |
---|
1214 | } while (!(n > numsp || found)); |
---|
1215 | if (n > numsp) { |
---|
1216 | printf("Cannot find species: "); |
---|
1217 | for (i = 0; i < nmlngth; i++) |
---|
1218 | putchar(str[i]); |
---|
1219 | putchar('\n'); |
---|
1220 | } |
---|
1221 | hookup(p, q); |
---|
1222 | if (curtree.start->number > n) |
---|
1223 | curtree.start = curtree.nodep[n - 1]; |
---|
1224 | } |
---|
1225 | p->branchnum = q->number; |
---|
1226 | q->branchnum = q->number; |
---|
1227 | q->v = initialv; |
---|
1228 | p->v = initialv; |
---|
1229 | branchtrans(q->branchnum, q->v); |
---|
1230 | if (ch == ':') { |
---|
1231 | processlength(p); |
---|
1232 | *nolengths = false; |
---|
1233 | } |
---|
1234 | } /* addelement */ |
---|
1235 | |
---|
1236 | void treeread() |
---|
1237 | { |
---|
1238 | /* read a user-defined tree */ |
---|
1239 | /* Local variables for treeread: */ |
---|
1240 | short nextnode, lparens, rparens; |
---|
1241 | boolean *names, nolengths; |
---|
1242 | node *p; |
---|
1243 | short i, j; |
---|
1244 | |
---|
1245 | curtree.start = curtree.nodep[numsp - 1]; |
---|
1246 | do { |
---|
1247 | ch = getc(infile); |
---|
1248 | if (ch == '\n') |
---|
1249 | ch = ' '; |
---|
1250 | } while (ch == ' '); |
---|
1251 | if (ch != '(') |
---|
1252 | return; |
---|
1253 | names = (boolean *)Malloc(numsp*sizeof(boolean)); |
---|
1254 | for (i = 0; i < numsp; i++) |
---|
1255 | names[i] = false; |
---|
1256 | lparens = 1; |
---|
1257 | rparens = 0; |
---|
1258 | nolengths = true; |
---|
1259 | nextnode = numsp + 1; |
---|
1260 | p = curtree.nodep[nextnode - 1]; |
---|
1261 | for (i = 1; i <= 2; i++) { |
---|
1262 | addelement(p, &nextnode,&lparens,&rparens,names,&nolengths); |
---|
1263 | p = p->next; |
---|
1264 | findch(',',&lparens,&rparens); |
---|
1265 | } |
---|
1266 | addelement(p, &nextnode,&lparens,&rparens,names,&nolengths); |
---|
1267 | findch(')',&lparens,&rparens); |
---|
1268 | for (i = 0; i <= endsite; i++) { |
---|
1269 | for (j = 0; j <= sitelength; j++) { |
---|
1270 | p->x[i][j] = 1.0; |
---|
1271 | p->next->x[i][j] = 1.0; |
---|
1272 | p->next->next->x[i][j] = 1.0; |
---|
1273 | } |
---|
1274 | } |
---|
1275 | fscanf(infile, "%*[^\n]"); |
---|
1276 | getc(infile); |
---|
1277 | free(names); |
---|
1278 | } /* treeread */ |
---|
1279 | |
---|
1280 | |
---|
1281 | void travinit(p) |
---|
1282 | node *p; |
---|
1283 | { |
---|
1284 | /* traverse to set up initial values */ |
---|
1285 | if (p == NULL) |
---|
1286 | return; |
---|
1287 | if (p->tip) |
---|
1288 | return; |
---|
1289 | if (p->initialized) |
---|
1290 | return; |
---|
1291 | travinit(p->next->back); |
---|
1292 | travinit(p->next->next->back); |
---|
1293 | nuview(p); |
---|
1294 | p->initialized = true; |
---|
1295 | } /* travinit */ |
---|
1296 | |
---|
1297 | |
---|
1298 | void travsp(p) |
---|
1299 | node *p; |
---|
1300 | { |
---|
1301 | /* traverse to find tips */ |
---|
1302 | if (p == curtree.start) |
---|
1303 | travinit(p); |
---|
1304 | if (p->tip) |
---|
1305 | travinit(p->back); |
---|
1306 | else { |
---|
1307 | travsp(p->next->back); |
---|
1308 | travsp(p->next->next->back); |
---|
1309 | } |
---|
1310 | } /* travsp */ |
---|
1311 | |
---|
1312 | |
---|
1313 | void treevaluate() |
---|
1314 | { |
---|
1315 | /* find maximum likelihood branch lengths of user tree */ |
---|
1316 | short i; |
---|
1317 | double dummy; |
---|
1318 | |
---|
1319 | for (i = 1; i <= numsp; i++) |
---|
1320 | curtree.nodep[i-1]->initialized = false; |
---|
1321 | for (i = numsp+1; i <= numsp2; i++) { |
---|
1322 | curtree.nodep[i-1]->initialized = false; |
---|
1323 | curtree.nodep[i-1]->next->initialized = false; |
---|
1324 | curtree.nodep[i-1]->next->next->initialized = false; |
---|
1325 | } |
---|
1326 | for (i = 1; i <= smoothings * 4; i++) |
---|
1327 | smooth(curtree.start->back); |
---|
1328 | dummy = evaluate(&curtree); |
---|
1329 | } /* treevaluate */ |
---|
1330 | |
---|
1331 | |
---|
1332 | void maketree() |
---|
1333 | { |
---|
1334 | /* construct and rearrange tree */ |
---|
1335 | short i, j, k, num; |
---|
1336 | double sum, sum2, sd; |
---|
1337 | double TEMP; |
---|
1338 | node *p; |
---|
1339 | |
---|
1340 | tempmatrix = (transmatrix)Malloc((sitelength+1) * sizeof(double *)); |
---|
1341 | for (i=0;i<=sitelength;++i) |
---|
1342 | tempmatrix[i] = (double *)Malloc((sitelength+1) * sizeof (double)); |
---|
1343 | tempprod = (transmatrix)Malloc((sitelength+1) * sizeof(double *)); |
---|
1344 | for (i=0;i<=sitelength;++i) |
---|
1345 | tempprod[i] = (double *)Malloc((sitelength+1) * sizeof (double)); |
---|
1346 | |
---|
1347 | for (i = 0; i < maxtrees; i++) |
---|
1348 | l0gf[i] = (double *)Malloc(endsite*sizeof(double)); |
---|
1349 | setuppi(); |
---|
1350 | if (usertree) { |
---|
1351 | fscanf(infile, "%hd%*[^\n]", &numtrees); |
---|
1352 | getc(infile); |
---|
1353 | if (treeprint) { |
---|
1354 | fprintf(outfile, "User-defined tree"); |
---|
1355 | if (numtrees > 1) |
---|
1356 | putc('s', outfile); |
---|
1357 | fprintf(outfile, ":\n\n"); |
---|
1358 | } |
---|
1359 | which = 1; |
---|
1360 | while (which <= numtrees) { |
---|
1361 | treeread(); |
---|
1362 | treevaluate(); |
---|
1363 | printree(); |
---|
1364 | summarize(); |
---|
1365 | if (trout) { |
---|
1366 | col = 0; |
---|
1367 | treeout(curtree.start->back); |
---|
1368 | } |
---|
1369 | which++; |
---|
1370 | } |
---|
1371 | putc('\n', outfile); |
---|
1372 | if (numtrees > 1 && weightsum > 1) { |
---|
1373 | fprintf(outfile, "Tree Ln L Diff Ln L Its S.D."); |
---|
1374 | fprintf(outfile, " Significantly worse?\n\n"); |
---|
1375 | if (numtrees > maxtrees) |
---|
1376 | num = maxtrees; |
---|
1377 | else |
---|
1378 | num = numtrees; |
---|
1379 | for (i = 1; i <= num; i++) { |
---|
1380 | fprintf(outfile, "%3hd%12.5f", i, l0gl[i - 1]); |
---|
1381 | if (maxwhich == i) |
---|
1382 | fprintf(outfile, " <------ best\n"); |
---|
1383 | else { |
---|
1384 | sum = 0.0; |
---|
1385 | sum2 = 0.0; |
---|
1386 | for (j = 1; j <= endsite; j++) { |
---|
1387 | sum += weight[j] * (l0gf[maxwhich - 1][j - 1]-l0gf[i - 1][j - 1]); |
---|
1388 | TEMP = l0gf[maxwhich - 1][j - 1] - l0gf[i - 1][j - 1]; |
---|
1389 | sum2 += weight[j] * (TEMP * TEMP); |
---|
1390 | } |
---|
1391 | sd = sqrt(weightsum / (weightsum - 1.0) |
---|
1392 | * (sum2 - sum * sum / weightsum)); |
---|
1393 | fprintf(outfile, "%12.5f%12.4f", l0gl[i - 1] - maxlogl, sd); |
---|
1394 | if (sum > 1.95996 * sd) |
---|
1395 | fprintf(outfile, " Yes\n"); |
---|
1396 | else |
---|
1397 | fprintf(outfile, " No\n"); |
---|
1398 | } |
---|
1399 | } |
---|
1400 | fprintf(outfile, "\n\n"); |
---|
1401 | } |
---|
1402 | } else { |
---|
1403 | for (i = 1; i <= numsp; i++) |
---|
1404 | enterorder[i - 1] = i; |
---|
1405 | if (jumble) { |
---|
1406 | for (i = 0; i < numsp; i++) { |
---|
1407 | j = (short)(randum(seed) * numsp) + 1; |
---|
1408 | k = enterorder[j - 1]; |
---|
1409 | enterorder[j - 1] = enterorder[i]; |
---|
1410 | enterorder[i] = k; |
---|
1411 | } |
---|
1412 | } |
---|
1413 | if (progress) { |
---|
1414 | printf("\nAdding species:\n"); |
---|
1415 | printf(" "); |
---|
1416 | for (i = 0; i < nmlngth; i++) |
---|
1417 | putchar(curtree.nodep[enterorder[0] - 1]->nayme[i]); |
---|
1418 | printf("\n "); |
---|
1419 | for (i = 0; i < nmlngth; i++) |
---|
1420 | putchar(curtree.nodep[enterorder[1] - 1]->nayme[i]); |
---|
1421 | printf("\n "); |
---|
1422 | for (i = 0; i < nmlngth; i++) |
---|
1423 | putchar(curtree.nodep[enterorder[2] - 1]->nayme[i]); |
---|
1424 | putchar('\n'); |
---|
1425 | } |
---|
1426 | nextsp = 3; |
---|
1427 | buildsimpletree(&curtree); |
---|
1428 | curtree.start = curtree.nodep[enterorder[0] - 1]; |
---|
1429 | if (jumb == 1) numtrees = 1; |
---|
1430 | nextsp = 4; |
---|
1431 | while (nextsp <= numsp) { |
---|
1432 | buildnewtip(enterorder[nextsp - 1], &curtree); |
---|
1433 | bestree.likelihood = -999999.0; |
---|
1434 | copy_(&curtree, &priortree); |
---|
1435 | addtraverse(curtree.nodep[enterorder[nextsp - 1] - 1]->back, |
---|
1436 | curtree.start->back, true); |
---|
1437 | copy_(&bestree, &curtree); |
---|
1438 | if (progress) { |
---|
1439 | printf(" "); |
---|
1440 | for (j = 0; j < nmlngth; j++) |
---|
1441 | putchar(curtree.nodep[enterorder[nextsp - 1] - 1]->nayme[j]); |
---|
1442 | putchar('\n'); |
---|
1443 | } |
---|
1444 | if (global && nextsp == numsp) { |
---|
1445 | if (progress) { |
---|
1446 | printf("Doing global rearrangements\n"); |
---|
1447 | printf(" "); |
---|
1448 | } |
---|
1449 | } |
---|
1450 | succeeded = true; |
---|
1451 | while (succeeded) { |
---|
1452 | succeeded = false; |
---|
1453 | rearrange(curtree.start->back); |
---|
1454 | } |
---|
1455 | if (njumble > 1) { |
---|
1456 | if (jumb == 1 && nextsp == numsp) |
---|
1457 | copy_(&bestree, &bestree2); |
---|
1458 | else if (nextsp == numsp) { |
---|
1459 | if (bestree2.likelihood < bestree.likelihood) |
---|
1460 | copy_(&bestree, &bestree2); |
---|
1461 | } |
---|
1462 | } |
---|
1463 | if (nextsp == numsp && jumb == njumble) { |
---|
1464 | if (njumble > 1) copy_(&bestree2, &curtree); |
---|
1465 | curtree.start = curtree.nodep[outgrno - 1]; |
---|
1466 | printree(); |
---|
1467 | summarize(); |
---|
1468 | if (trout && nextsp == numsp) { |
---|
1469 | col = 0; |
---|
1470 | treeout(curtree.start->back); |
---|
1471 | } |
---|
1472 | } |
---|
1473 | nextsp++; |
---|
1474 | } |
---|
1475 | } |
---|
1476 | if (jumb == njumble) { |
---|
1477 | if (progress) { |
---|
1478 | printf("\nOutput written to output file\n\n"); |
---|
1479 | if (trout) |
---|
1480 | printf("Tree also written onto file\n"); |
---|
1481 | putchar('\n'); |
---|
1482 | } |
---|
1483 | for (i = 0; i < numsp; i++) |
---|
1484 | free(curtree.nodep[i]->x); |
---|
1485 | for (i = numsp; i < numsp2; i++) { |
---|
1486 | p = curtree.nodep[i]; |
---|
1487 | for (j = 1; j <= 3; j++) { |
---|
1488 | free(p->x); |
---|
1489 | p = p->next; |
---|
1490 | } |
---|
1491 | } |
---|
1492 | if (!usertree) { |
---|
1493 | for (i = 0; i < numsp; i++) |
---|
1494 | free(priortree.nodep[i]->x); |
---|
1495 | for (i = numsp; i < numsp2; i++) { |
---|
1496 | p = priortree.nodep[i]; |
---|
1497 | for (j = 1; j <= 3; j++) { |
---|
1498 | free(p->x); |
---|
1499 | p = p->next; |
---|
1500 | } |
---|
1501 | } |
---|
1502 | for (i = 0; i < numsp; i++) |
---|
1503 | free(bestree.nodep[i]->x); |
---|
1504 | for (i = numsp; i < numsp2; i++) { |
---|
1505 | p = bestree.nodep[i]; |
---|
1506 | for (j = 1; j <= 3; j++) { |
---|
1507 | free(p->x); |
---|
1508 | p = p->next; |
---|
1509 | } |
---|
1510 | } |
---|
1511 | if (njumble > 1) { |
---|
1512 | for (i = 0; i < numsp; i++) |
---|
1513 | free(bestree2.nodep[i]->x); |
---|
1514 | for (i = numsp; i < numsp2; i++) { |
---|
1515 | p = bestree2.nodep[i]; |
---|
1516 | for (j = 1; j <= 3; j++) { |
---|
1517 | free(p->x); |
---|
1518 | p = p->next; |
---|
1519 | } |
---|
1520 | } |
---|
1521 | } |
---|
1522 | } |
---|
1523 | for (i = 0; i < maxtrees; i++) |
---|
1524 | free(l0gf[i]); |
---|
1525 | } |
---|
1526 | } /* maketree */ |
---|
1527 | |
---|
1528 | |
---|