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