1 | #include "phylip.h" |
---|
2 | #include "seq.h" |
---|
3 | |
---|
4 | /* version 3.6. (c) Copyright 1993-2002 by the University of Washington. |
---|
5 | Written by Joseph Felsenstein, Akiko Fuseki, Sean Lamont, and Andrew Keeffe. |
---|
6 | Permission is granted to copy and use this program provided no fee is |
---|
7 | charged for it and provided that this copyright notice is not removed. */ |
---|
8 | |
---|
9 | #define maxtrees 100 /* maximum number of trees to be printed out */ |
---|
10 | #define often 100 /* how often to notify how many trees examined */ |
---|
11 | #define many 1000 /* how many multiples of howoften before stop */ |
---|
12 | |
---|
13 | typedef node **pointptr; |
---|
14 | typedef long *treenumbers; |
---|
15 | typedef double *valptr; |
---|
16 | typedef long *placeptr; |
---|
17 | |
---|
18 | #ifndef OLDC |
---|
19 | /* function prototypes */ |
---|
20 | void getoptions(void); |
---|
21 | void allocrest(void); |
---|
22 | void doinit(void); |
---|
23 | void makeweights(void); |
---|
24 | void doinput(void); |
---|
25 | void supplement(node *); |
---|
26 | void evaluate(node *); |
---|
27 | void addtraverse(node *, node *, node *, long *, long *, valptr, placeptr); |
---|
28 | void addit(long ); |
---|
29 | void dnapenny_reroot(node *); |
---|
30 | void describe(void); |
---|
31 | void maketree(void); |
---|
32 | void reallocchars(void); |
---|
33 | /* function prototypes */ |
---|
34 | #endif |
---|
35 | |
---|
36 | |
---|
37 | extern sequence y; |
---|
38 | |
---|
39 | Char infilename[FNMLNGTH],outfilename[FNMLNGTH],outtreename[FNMLNGTH], weightfilename[FNMLNGTH]; |
---|
40 | node *root; |
---|
41 | long *zeros=NULL; |
---|
42 | long chars, howmany, howoften, col, msets, ith; |
---|
43 | boolean weights, thresh, simple, trout, progress, stepbox, ancseq, |
---|
44 | mulsets, firstset, justwts; |
---|
45 | double threshold; |
---|
46 | steptr oldweight; |
---|
47 | pointptr treenode; /* pointers to all nodes in tree */ |
---|
48 | double fracdone, fracinc; |
---|
49 | boolean *added; |
---|
50 | gbases *garbage; |
---|
51 | node **grbg; |
---|
52 | Char basechar[]="ACMGRSVTWYHKDBNO???????????????"; |
---|
53 | |
---|
54 | /* Variables for maketree, propagated globally for C version: */ |
---|
55 | long examined, mults; |
---|
56 | boolean firsttime, done, recompute; |
---|
57 | double like, bestyet; |
---|
58 | treenumbers *bestorders, *bestrees; |
---|
59 | treenumbers current, order; |
---|
60 | long *threshwt; |
---|
61 | baseptr nothing; |
---|
62 | node *temp, *temp1; |
---|
63 | long suppno[] = |
---|
64 | { 0,0,0,0,0,1,1,1,1,1,1,1,1,1,1,2,2,2,2,2,2,2,2,2,2,3,3,3,3,3,4}; |
---|
65 | |
---|
66 | long suppset[] = /* this was previously a function. */ |
---|
67 | { /* in C, it doesn't need to be. */ |
---|
68 | 1 << ((long)A), |
---|
69 | 1 << ((long)C), |
---|
70 | 1 << ((long)G), |
---|
71 | 1 << ((long)T), |
---|
72 | 1 << ((long)O), |
---|
73 | (1 << ((long)A)) | (1 << ((long)C)), |
---|
74 | (1 << ((long)A)) | (1 << ((long)G)), |
---|
75 | (1 << ((long)A)) | (1 << ((long)T)), |
---|
76 | (1 << ((long)A)) | (1 << ((long)O)), |
---|
77 | (1 << ((long)C)) | (1 << ((long)G)), |
---|
78 | (1 << ((long)C)) | (1 << ((long)T)), |
---|
79 | (1 << ((long)C)) | (1 << ((long)O)), |
---|
80 | (1 << ((long)G)) | (1 << ((long)T)), |
---|
81 | (1 << ((long)G)) | (1 << ((long)O)), |
---|
82 | (1 << ((long)T)) | (1 << ((long)O)), |
---|
83 | (1 << ((long)A)) | (1 << ((long)C)) | (1 << ((long)G)), |
---|
84 | (1 << ((long)A)) | (1 << ((long)C)) | (1 << ((long)T)), |
---|
85 | (1 << ((long)A)) | (1 << ((long)C)) | (1 << ((long)O)), |
---|
86 | (1 << ((long)A)) | (1 << ((long)G)) | (1 << ((long)T)), |
---|
87 | (1 << ((long)A)) | (1 << ((long)G)) | (1 << ((long)O)), |
---|
88 | (1 << ((long)A)) | (1 << ((long)T)) | (1 << ((long)O)), |
---|
89 | (1 << ((long)C)) | (1 << ((long)G)) | (1 << ((long)T)), |
---|
90 | (1 << ((long)C)) | (1 << ((long)G)) | (1 << ((long)O)), |
---|
91 | (1 << ((long)C)) | (1 << ((long)T)) | (1 << ((long)O)), |
---|
92 | (1 << ((long)G)) | (1 << ((long)T)) | (1 << ((long)O)), |
---|
93 | (1 << ((long)A))|(1 << ((long)C))|(1 << ((long)G))|(1 << ((long)T)), |
---|
94 | (1 << ((long)A))|(1 << ((long)C))|(1 << ((long)G))|(1 << ((long)O)), |
---|
95 | (1 << ((long)A))|(1 << ((long)C))|(1 << ((long)T))|(1 << ((long)O)), |
---|
96 | (1 << ((long)A))|(1 << ((long)G))|(1 << ((long)T))|(1 << ((long)O)), |
---|
97 | (1 << ((long)C))|(1 << ((long)G))|(1 << ((long)T)) | (1 << ((long)O)), |
---|
98 | (1 << ((long)A))|(1 << ((long)C))|(1 << ((long)G)) | (1 << ((long)T)) | (1 << ((long)O))}; |
---|
99 | |
---|
100 | |
---|
101 | void getoptions() |
---|
102 | { |
---|
103 | /* interactively set options */ |
---|
104 | long loopcount, loopcount2; |
---|
105 | Char ch, ch2; |
---|
106 | |
---|
107 | fprintf(outfile, "\nPenny algorithm for DNA, version %s\n",VERSION); |
---|
108 | fprintf(outfile, " branch-and-bound to find all"); |
---|
109 | fprintf(outfile, " most parsimonious trees\n\n"); |
---|
110 | howoften = often; |
---|
111 | howmany = many; |
---|
112 | outgrno = 1; |
---|
113 | outgropt = false; |
---|
114 | simple = true; |
---|
115 | thresh = false; |
---|
116 | threshold = spp; |
---|
117 | trout = true; |
---|
118 | weights = false; |
---|
119 | justwts = false; |
---|
120 | printdata = false; |
---|
121 | progress = true; |
---|
122 | treeprint = true; |
---|
123 | stepbox = false; |
---|
124 | ancseq = false; |
---|
125 | interleaved = true; |
---|
126 | loopcount = 0; |
---|
127 | for (;;) { |
---|
128 | cleerhome(); |
---|
129 | printf("\nPenny algorithm for DNA, version %s\n",VERSION); |
---|
130 | printf(" branch-and-bound to find all most parsimonious trees\n\n"); |
---|
131 | printf("Settings for this run:\n"); |
---|
132 | printf(" H How many groups of %4ld trees:%6ld\n", howoften, howmany); |
---|
133 | printf(" F How often to report, in trees: %4ld\n", howoften); |
---|
134 | printf(" S Branch and bound is simple? %s\n", |
---|
135 | (simple ? "Yes" : "No. reconsiders order of species")); |
---|
136 | printf(" O Outgroup root? %s%3ld\n", |
---|
137 | (outgropt ? "Yes, at sequence number" : |
---|
138 | "No, use as outgroup species"),outgrno); |
---|
139 | printf(" T Use Threshold parsimony?"); |
---|
140 | if (thresh) |
---|
141 | printf(" Yes, count steps up to%4.1f per site\n", threshold); |
---|
142 | else |
---|
143 | printf(" No, use ordinary parsimony\n"); |
---|
144 | printf(" W Sites weighted? %s\n", |
---|
145 | (weights ? "Yes" : "No")); |
---|
146 | printf(" M Analyze multiple data sets?"); |
---|
147 | if (mulsets) |
---|
148 | printf(" Yes, %2ld %s\n", msets, |
---|
149 | (justwts ? "sets of weights" : "data sets")); |
---|
150 | else |
---|
151 | printf(" No\n"); |
---|
152 | printf(" I Input sequences interleaved? %s\n", |
---|
153 | (interleaved ? "Yes" : "No, sequential")); |
---|
154 | printf(" 0 Terminal type (IBM PC, ANSI, none)? %s\n", |
---|
155 | (ibmpc ? "IBM PC" : ansi ? "ANSI" : "(none)")); |
---|
156 | printf(" 1 Print out the data at start of run %s\n", |
---|
157 | (printdata ? "Yes" : "No")); |
---|
158 | printf(" 2 Print indications of progress of run %s\n", |
---|
159 | (progress ? "Yes" : "No")); |
---|
160 | printf(" 3 Print out tree %s\n", |
---|
161 | (treeprint ? "Yes" : "No")); |
---|
162 | printf(" 4 Print out steps in each site %s\n", |
---|
163 | (stepbox ? "Yes" : "No" )); |
---|
164 | printf(" 5 Print sequences at all nodes of tree %s\n", |
---|
165 | (ancseq ? "Yes" : "No")); |
---|
166 | printf(" 6 Write out trees onto tree file? %s\n", |
---|
167 | (trout ? "Yes" : "No")); |
---|
168 | if(weights && justwts){ |
---|
169 | printf("WARNING: W option and Multiple Weights options are both on. "); |
---|
170 | printf( |
---|
171 | "The W menu option is unnecessary and has no additional effect. \n"); |
---|
172 | } |
---|
173 | printf("\nAre these settings correct? (type Y or the letter for one to change)\n"); |
---|
174 | #ifdef WIN32 |
---|
175 | phyFillScreenColor(); |
---|
176 | #endif |
---|
177 | scanf("%c%*[^\n]", &ch); |
---|
178 | getchar(); |
---|
179 | if (ch == '\n') |
---|
180 | ch = ' '; |
---|
181 | uppercase(&ch); |
---|
182 | if (ch == 'Y') |
---|
183 | break; |
---|
184 | uppercase(&ch); |
---|
185 | if ((strchr("WHMSOFTI1234560",ch)) != NULL){ |
---|
186 | switch (ch) { |
---|
187 | |
---|
188 | case 'H': |
---|
189 | inithowmany(&howmany, howoften); |
---|
190 | break; |
---|
191 | |
---|
192 | case 'W': |
---|
193 | weights = !weights; |
---|
194 | break; |
---|
195 | |
---|
196 | case 'M': |
---|
197 | mulsets = !mulsets; |
---|
198 | if (mulsets){ |
---|
199 | printf("Multiple data sets or multiple weights?"); |
---|
200 | loopcount2 = 0; |
---|
201 | do { |
---|
202 | printf(" (type D or W)\n"); |
---|
203 | #ifdef WIN32 |
---|
204 | phyFillScreenColor(); |
---|
205 | #endif |
---|
206 | scanf("%c%*[^\n]", &ch2); |
---|
207 | getchar(); |
---|
208 | if (ch2 == '\n') |
---|
209 | ch2 = ' '; |
---|
210 | uppercase(&ch2); |
---|
211 | countup(&loopcount2, 10); |
---|
212 | } while ((ch2 != 'W') && (ch2 != 'D')); |
---|
213 | justwts = (ch2 == 'W'); |
---|
214 | if (justwts) |
---|
215 | justweights(&msets); |
---|
216 | else |
---|
217 | initdatasets(&msets); |
---|
218 | } |
---|
219 | break; |
---|
220 | |
---|
221 | case 'F': |
---|
222 | inithowoften(&howoften); |
---|
223 | break; |
---|
224 | |
---|
225 | case 'S': |
---|
226 | simple = !simple; |
---|
227 | break; |
---|
228 | |
---|
229 | case 'O': |
---|
230 | outgropt = !outgropt; |
---|
231 | if (outgropt) |
---|
232 | initoutgroup(&outgrno, spp); |
---|
233 | else |
---|
234 | outgrno = 1; |
---|
235 | break; |
---|
236 | |
---|
237 | case 'T': |
---|
238 | thresh = !thresh; |
---|
239 | if (thresh) |
---|
240 | initthreshold(&threshold); |
---|
241 | break; |
---|
242 | |
---|
243 | case 'I': |
---|
244 | interleaved = !interleaved; |
---|
245 | break; |
---|
246 | |
---|
247 | case '0': |
---|
248 | initterminal(&ibmpc, &ansi); |
---|
249 | break; |
---|
250 | |
---|
251 | case '1': |
---|
252 | printdata = !printdata; |
---|
253 | break; |
---|
254 | |
---|
255 | case '2': |
---|
256 | progress = !progress; |
---|
257 | break; |
---|
258 | |
---|
259 | case '3': |
---|
260 | treeprint = !treeprint; |
---|
261 | break; |
---|
262 | |
---|
263 | case '4': |
---|
264 | stepbox = !stepbox; |
---|
265 | break; |
---|
266 | |
---|
267 | case '5': |
---|
268 | ancseq = !ancseq; |
---|
269 | break; |
---|
270 | |
---|
271 | case '6': |
---|
272 | trout = !trout; |
---|
273 | break; |
---|
274 | } |
---|
275 | } else |
---|
276 | printf("Not a possible option!\n"); |
---|
277 | countup(&loopcount, 100); |
---|
278 | } |
---|
279 | } /* getoptions */ |
---|
280 | |
---|
281 | void allocrest() |
---|
282 | { |
---|
283 | long i; |
---|
284 | |
---|
285 | y = (Char **)Malloc(spp*sizeof(Char *)); |
---|
286 | for (i = 0; i < spp; i++) |
---|
287 | y[i] = (Char *)Malloc(chars*sizeof(Char)); |
---|
288 | weight = (long *)Malloc(chars*sizeof(long)); |
---|
289 | oldweight = (long *)Malloc(chars*sizeof(long)); |
---|
290 | alias = (steptr)Malloc(chars*sizeof(long)); |
---|
291 | ally = (steptr)Malloc(chars*sizeof(long)); |
---|
292 | location = (steptr)Malloc(chars*sizeof(long)); |
---|
293 | nayme = (naym *)Malloc(spp*sizeof(naym)); |
---|
294 | bestorders = (treenumbers *)Malloc(maxtrees*sizeof(treenumbers)); |
---|
295 | for (i = 1; i <= maxtrees; i++) |
---|
296 | bestorders[i - 1] = (treenumbers)Malloc(spp*sizeof(long)); |
---|
297 | bestrees = (treenumbers *)Malloc(maxtrees*sizeof(treenumbers)); |
---|
298 | for (i = 1; i <= maxtrees; i++) |
---|
299 | bestrees[i - 1] = (treenumbers)Malloc(spp*sizeof(long)); |
---|
300 | current = (treenumbers)Malloc(spp*sizeof(long)); |
---|
301 | order = (treenumbers)Malloc(spp*sizeof(long)); |
---|
302 | added = (boolean *)Malloc(nonodes*sizeof(boolean)); |
---|
303 | } /* allocrest */ |
---|
304 | |
---|
305 | void reallocchars(void) |
---|
306 | {/* The amount of chars can change between runs |
---|
307 | this function reallocates all the variables |
---|
308 | whose size depends on the amount of chars */ |
---|
309 | long i; |
---|
310 | |
---|
311 | for (i = 0; i < spp; i++) { |
---|
312 | free(y[i]); |
---|
313 | y[i] = (Char *)Malloc(chars*sizeof(Char)); |
---|
314 | } |
---|
315 | |
---|
316 | free(weight); |
---|
317 | free(oldweight); |
---|
318 | free(alias); |
---|
319 | free(ally); |
---|
320 | free(location); |
---|
321 | |
---|
322 | weight = (long *)Malloc(chars*sizeof(long)); |
---|
323 | oldweight = (long *)Malloc(chars*sizeof(long)); |
---|
324 | alias = (steptr)Malloc(chars*sizeof(long)); |
---|
325 | ally = (steptr)Malloc(chars*sizeof(long)); |
---|
326 | location = (steptr)Malloc(chars*sizeof(long)); |
---|
327 | } /* reallocchars */ |
---|
328 | |
---|
329 | void doinit() |
---|
330 | { |
---|
331 | /* initializes variables */ |
---|
332 | inputnumbers(&spp, &chars, &nonodes, 1); |
---|
333 | getoptions(); |
---|
334 | if (printdata) |
---|
335 | fprintf(outfile, "%2ld species, %3ld sites\n", spp, chars); |
---|
336 | alloctree(&treenode, nonodes, false); |
---|
337 | allocrest(); |
---|
338 | } /* doinit */ |
---|
339 | |
---|
340 | void makeweights() |
---|
341 | { |
---|
342 | /* make up weights vector to avoid duplicate computations */ |
---|
343 | long i; |
---|
344 | |
---|
345 | for (i = 1; i <= chars; i++) { |
---|
346 | alias[i - 1] = i; |
---|
347 | oldweight[i - 1] = weight[i - 1]; |
---|
348 | ally[i - 1] = i; |
---|
349 | } |
---|
350 | sitesort(chars, weight); |
---|
351 | sitecombine(chars); |
---|
352 | sitescrunch(chars); |
---|
353 | endsite = 0; |
---|
354 | for (i = 1; i <= chars; i++) { |
---|
355 | if (ally[i - 1] == i) |
---|
356 | endsite++; |
---|
357 | } |
---|
358 | for (i = 1; i <= endsite; i++) |
---|
359 | location[alias[i - 1] - 1] = i; |
---|
360 | if (!thresh) |
---|
361 | threshold = spp; |
---|
362 | threshwt = (long *)Malloc(endsite*sizeof(long)); |
---|
363 | for (i = 0; i < endsite; i++) { |
---|
364 | weight[i] *= 10; |
---|
365 | threshwt[i] = (long)(threshold * weight[i] + 0.5); |
---|
366 | } |
---|
367 | if ( zeros != NULL ) |
---|
368 | free(zeros); |
---|
369 | zeros = (long *)Malloc(endsite*sizeof(long)); /*in makeweights()*/ |
---|
370 | for (i = 0; i < endsite; i++) |
---|
371 | zeros[i] = 0; |
---|
372 | } /* makeweights */ |
---|
373 | |
---|
374 | |
---|
375 | void doinput() |
---|
376 | { |
---|
377 | /* reads the input data */ |
---|
378 | long i; |
---|
379 | |
---|
380 | if (justwts) { |
---|
381 | if (firstset) |
---|
382 | inputdata(chars); |
---|
383 | for (i = 0; i < chars; i++) |
---|
384 | weight[i] = 1; |
---|
385 | inputweights(chars, weight, &weights); |
---|
386 | if (justwts) { |
---|
387 | fprintf(outfile, "\n\nWeights set # %ld:\n\n", ith); |
---|
388 | if (progress) |
---|
389 | printf("\nWeights set # %ld:\n\n", ith); |
---|
390 | } |
---|
391 | if (printdata) |
---|
392 | printweights(outfile, 0, chars, weight, "Sites"); |
---|
393 | } else { |
---|
394 | if (!firstset){ |
---|
395 | samenumsp(&chars, ith); |
---|
396 | reallocchars(); |
---|
397 | } |
---|
398 | inputdata(chars); |
---|
399 | for (i = 0; i < chars; i++) |
---|
400 | weight[i] = 1; |
---|
401 | if (weights) { |
---|
402 | inputweights(chars, weight, &weights); |
---|
403 | if (printdata) |
---|
404 | printweights(outfile, 0, chars, weight, "Sites"); |
---|
405 | } |
---|
406 | } |
---|
407 | |
---|
408 | makeweights(); |
---|
409 | makevalues(treenode, zeros, false); |
---|
410 | alloctemp(&temp, zeros, endsite); |
---|
411 | alloctemp(&temp1, zeros, endsite); |
---|
412 | } /* doinput */ |
---|
413 | |
---|
414 | |
---|
415 | void supplement(node *r) |
---|
416 | { |
---|
417 | /* determine minimum number of steps more which will |
---|
418 | be added when rest of species are put in tree */ |
---|
419 | long i, j, k, sum, sumall=0, sumadded=0; |
---|
420 | boolean doneadded, allhave, addedhave, has; |
---|
421 | long supps; |
---|
422 | |
---|
423 | for (i = 0; i < endsite; i++) { |
---|
424 | sum = 3; |
---|
425 | j = 1; |
---|
426 | doneadded = false; |
---|
427 | do { |
---|
428 | allhave = true; |
---|
429 | addedhave = true; |
---|
430 | supps = suppset[j-1]; |
---|
431 | for (k = 0; k < spp; k++) { |
---|
432 | has = ((treenode[k]->base[i] & supps) != 0); |
---|
433 | if (added[k] && !doneadded) |
---|
434 | addedhave = (addedhave && has); |
---|
435 | allhave = (allhave && has); |
---|
436 | } |
---|
437 | if (allhave) |
---|
438 | sumall = suppno[j - 1]; |
---|
439 | if (addedhave) |
---|
440 | sumadded = suppno[j - 1]; |
---|
441 | doneadded = (doneadded || addedhave); |
---|
442 | j++; |
---|
443 | } while (!(j > 31 || (allhave && doneadded))); |
---|
444 | if (addedhave && allhave) |
---|
445 | sum = sumall - sumadded; |
---|
446 | r->numsteps[i] += sum * weight[i]; |
---|
447 | } |
---|
448 | } /* supplement */ |
---|
449 | |
---|
450 | |
---|
451 | void evaluate(node *r) |
---|
452 | { |
---|
453 | /* determines the number of steps needed for a tree. this is |
---|
454 | the minimum number of steps needed to evolve sequences on |
---|
455 | this tree */ |
---|
456 | long i, steps; |
---|
457 | double sum; |
---|
458 | |
---|
459 | sum = 0.0; |
---|
460 | supplement(r); |
---|
461 | for (i = 0; i < endsite; i++) { |
---|
462 | steps = r->numsteps[i]; |
---|
463 | if ((long)steps <= threshwt[i]) |
---|
464 | sum += steps; |
---|
465 | else |
---|
466 | sum += threshwt[i]; |
---|
467 | } |
---|
468 | if (examined == 0 && mults == 0) |
---|
469 | bestyet = -1.0; |
---|
470 | like = sum; |
---|
471 | } /* evaluate */ |
---|
472 | |
---|
473 | |
---|
474 | void addtraverse(node *p, node *item, node *fork, long *m, long *n, |
---|
475 | valptr valyew, placeptr place) |
---|
476 | { |
---|
477 | /* traverse all places to add item */ |
---|
478 | if (done) |
---|
479 | return; |
---|
480 | if (*m <= 2 || (p != root && p != root->next->back)) { |
---|
481 | if (p == root) |
---|
482 | fillin(temp, item, p); |
---|
483 | else { |
---|
484 | fillin(temp1, item, p); |
---|
485 | fillin(temp, temp1, p->back); |
---|
486 | } |
---|
487 | (*n)++; |
---|
488 | evaluate(temp); |
---|
489 | examined++; |
---|
490 | if (examined == howoften) { |
---|
491 | examined = 0; |
---|
492 | mults++; |
---|
493 | if (mults == howmany) |
---|
494 | done = true; |
---|
495 | if (progress) { |
---|
496 | printf("%7ld", mults); |
---|
497 | if (bestyet >= 0) |
---|
498 | printf("%16.1f", bestyet / 10.0); |
---|
499 | else |
---|
500 | printf(" - "); |
---|
501 | printf("%17ld%20.2f\n", nextree - 1, fracdone * 100); |
---|
502 | #ifdef WIN32 |
---|
503 | phyFillScreenColor(); |
---|
504 | #endif |
---|
505 | } |
---|
506 | } |
---|
507 | valyew[(*n) - 1] = like; |
---|
508 | place[(*n) - 1] = p->index; |
---|
509 | } |
---|
510 | if (!p->tip) { |
---|
511 | addtraverse(p->next->back, item, fork, m,n,valyew,place); |
---|
512 | addtraverse(p->next->next->back, item, fork,m,n,valyew,place); |
---|
513 | } |
---|
514 | } /* addtraverse */ |
---|
515 | |
---|
516 | |
---|
517 | void addit(long m) |
---|
518 | { |
---|
519 | /* adds the species one by one, recursively */ |
---|
520 | long n; |
---|
521 | valptr valyew; |
---|
522 | placeptr place; |
---|
523 | long i, j, n1, besttoadd=0; |
---|
524 | valptr bestval; |
---|
525 | placeptr bestplace; |
---|
526 | double oldfrac, oldfdone, sum, bestsum; |
---|
527 | |
---|
528 | valyew = (valptr)Malloc(nonodes*sizeof(double)); |
---|
529 | bestval = (valptr)Malloc(nonodes*sizeof(double)); |
---|
530 | place = (placeptr)Malloc(nonodes*sizeof(long)); |
---|
531 | bestplace = (placeptr)Malloc(nonodes*sizeof(long)); |
---|
532 | if (simple && !firsttime) { |
---|
533 | n = 0; |
---|
534 | added[order[m - 1] - 1] = true; |
---|
535 | addtraverse(root, treenode[order[m - 1] - 1], |
---|
536 | treenode[spp + m - 2], &m,&n,valyew,place); |
---|
537 | besttoadd = order[m - 1]; |
---|
538 | memcpy(bestplace, place, nonodes*sizeof(long)); |
---|
539 | memcpy(bestval, valyew, nonodes*sizeof(double)); |
---|
540 | } else { |
---|
541 | bestsum = -1.0; |
---|
542 | for (i = 1; i <= spp; i++) { |
---|
543 | if (!added[i - 1]) { |
---|
544 | n = 0; |
---|
545 | added[i - 1] = true; |
---|
546 | addtraverse(root, treenode[i - 1], treenode[spp + m - 2], |
---|
547 | &m,&n,valyew,place); |
---|
548 | added[i - 1] = false; |
---|
549 | sum = 0.0; |
---|
550 | for (j = 0; j < n; j++) |
---|
551 | sum += valyew[j]; |
---|
552 | if (sum > bestsum) { |
---|
553 | bestsum = sum; |
---|
554 | besttoadd = i; |
---|
555 | memcpy(bestplace, place, nonodes*sizeof(long)); |
---|
556 | memcpy(bestval, valyew, nonodes*sizeof(double)); |
---|
557 | } |
---|
558 | } |
---|
559 | } |
---|
560 | } |
---|
561 | order[m - 1] = besttoadd; |
---|
562 | memcpy(place, bestplace, nonodes*sizeof(long)); |
---|
563 | memcpy(valyew, bestval, nonodes*sizeof(double)); |
---|
564 | shellsort(valyew, place, n); |
---|
565 | oldfrac = fracinc; |
---|
566 | oldfdone = fracdone; |
---|
567 | n1 = 0; |
---|
568 | for (i = 0; i < n; i++) { |
---|
569 | if (valyew[i] <= bestyet || bestyet < 0.0) |
---|
570 | n1++; |
---|
571 | } |
---|
572 | if (n1 > 0) |
---|
573 | fracinc /= n1; |
---|
574 | for (i = 0; i < n; i++) { |
---|
575 | if (valyew[i] <= bestyet || bestyet < 0.0) { |
---|
576 | current[m - 1] = place[i]; |
---|
577 | recompute = (m < spp); |
---|
578 | add(treenode[place[i] - 1], treenode[besttoadd - 1], |
---|
579 | treenode[spp + m - 2], &root, recompute, treenode, grbg, zeros); |
---|
580 | added[besttoadd - 1] = true; |
---|
581 | if (m < spp) |
---|
582 | addit(m + 1); |
---|
583 | else { |
---|
584 | if (valyew[i] < bestyet || bestyet < 0.0) { |
---|
585 | nextree = 1; |
---|
586 | bestyet = valyew[i]; |
---|
587 | } |
---|
588 | if (nextree <= maxtrees) { |
---|
589 | memcpy(bestorders[nextree - 1], order, |
---|
590 | spp*sizeof(long)); |
---|
591 | memcpy(bestrees[nextree - 1], current, |
---|
592 | spp*sizeof(long)); |
---|
593 | } |
---|
594 | nextree++; |
---|
595 | firsttime = false; |
---|
596 | } |
---|
597 | recompute = (m < spp); |
---|
598 | re_move(treenode[besttoadd - 1], &treenode[spp + m - 2], &root, |
---|
599 | recompute, treenode, grbg, zeros); |
---|
600 | added[besttoadd - 1] = false; |
---|
601 | } |
---|
602 | fracdone += fracinc; |
---|
603 | } |
---|
604 | fracinc = oldfrac; |
---|
605 | fracdone = oldfdone; |
---|
606 | free(valyew); |
---|
607 | free(bestval); |
---|
608 | free(place); |
---|
609 | free(bestplace); |
---|
610 | } /* addit */ |
---|
611 | |
---|
612 | |
---|
613 | void dnapenny_reroot(node *outgroup) |
---|
614 | { |
---|
615 | /* reorients tree, putting outgroup in desired position. */ |
---|
616 | node *p, *q, *newbottom, *oldbottom; |
---|
617 | |
---|
618 | if (outgroup->back->index == root->index) |
---|
619 | return; |
---|
620 | newbottom = outgroup->back; |
---|
621 | p = treenode[newbottom->index - 1]->back; |
---|
622 | while (p->index != root->index) { |
---|
623 | oldbottom = treenode[p->index - 1]; |
---|
624 | treenode[p->index - 1] = p; |
---|
625 | p = oldbottom->back; |
---|
626 | } |
---|
627 | p = root->next; |
---|
628 | q = root->next->next; |
---|
629 | p->back->back = q->back; |
---|
630 | q->back->back = p->back; |
---|
631 | p->back = outgroup; |
---|
632 | q->back = outgroup->back; |
---|
633 | outgroup->back->back = root->next->next; |
---|
634 | outgroup->back = root->next; |
---|
635 | treenode[newbottom->index - 1] = newbottom; |
---|
636 | } /* dnapenny_reroot */ |
---|
637 | |
---|
638 | |
---|
639 | void describe() |
---|
640 | { |
---|
641 | /* prints ancestors, steps and table of numbers of steps in |
---|
642 | each site */ |
---|
643 | |
---|
644 | if (stepbox) |
---|
645 | writesteps(chars, weights, oldweight, root); |
---|
646 | if (ancseq) { |
---|
647 | hypstates(chars, root, treenode, &garbage, basechar); |
---|
648 | putc('\n', outfile); |
---|
649 | } |
---|
650 | putc('\n', outfile); |
---|
651 | if (trout) { |
---|
652 | col = 0; |
---|
653 | treeout(root, nextree, &col, root); |
---|
654 | } |
---|
655 | } /* describe */ |
---|
656 | |
---|
657 | |
---|
658 | void maketree() |
---|
659 | { |
---|
660 | /* tree construction recursively by branch and bound */ |
---|
661 | long i, j, k; |
---|
662 | node *dummy; |
---|
663 | |
---|
664 | if (progress) { |
---|
665 | printf("\nHow many\n"); |
---|
666 | printf("trees looked Approximate\n"); |
---|
667 | printf("at so far Length of How many percentage\n"); |
---|
668 | printf("(multiples longest tree trees this long searched\n"); |
---|
669 | printf("of %4ld): found so far found so far so far\n", |
---|
670 | howoften); |
---|
671 | printf("---------- ------------ ------------ ------------\n"); |
---|
672 | } |
---|
673 | #ifdef WIN32 |
---|
674 | phyFillScreenColor(); |
---|
675 | #endif |
---|
676 | done = false; |
---|
677 | mults = 0; |
---|
678 | examined = 0; |
---|
679 | nextree = 1; |
---|
680 | root = treenode[0]; |
---|
681 | firsttime = true; |
---|
682 | for (i = 0; i < spp; i++) |
---|
683 | added[i] = false; |
---|
684 | added[0] = true; |
---|
685 | order[0] = 1; |
---|
686 | k = 2; |
---|
687 | fracdone = 0.0; |
---|
688 | fracinc = 1.0; |
---|
689 | bestyet = -1.0; |
---|
690 | recompute = true; |
---|
691 | addit(k); |
---|
692 | if (done) { |
---|
693 | if (progress) { |
---|
694 | printf("Search broken off! Not guaranteed to\n"); |
---|
695 | printf(" have found the most parsimonious trees.\n"); |
---|
696 | } |
---|
697 | if (treeprint) { |
---|
698 | fprintf(outfile, "Search broken off! Not guaranteed to\n"); |
---|
699 | fprintf(outfile, " have found the most parsimonious\n"); |
---|
700 | fprintf(outfile, " trees, but here is what we found:\n"); |
---|
701 | } |
---|
702 | } |
---|
703 | if (treeprint) { |
---|
704 | fprintf(outfile, "\nrequires a total of %18.3f\n\n", bestyet / 10.0); |
---|
705 | if (nextree == 2) |
---|
706 | fprintf(outfile, "One most parsimonious tree found:\n"); |
---|
707 | else |
---|
708 | fprintf(outfile, "%6ld trees in all found\n", nextree - 1); |
---|
709 | } |
---|
710 | if (nextree > maxtrees + 1) { |
---|
711 | if (treeprint) |
---|
712 | fprintf(outfile, "here are the first%4ld of them\n", (long)maxtrees); |
---|
713 | nextree = maxtrees + 1; |
---|
714 | } |
---|
715 | if (treeprint) |
---|
716 | putc('\n', outfile); |
---|
717 | for (i = 0; i < spp; i++) |
---|
718 | added[i] = true; |
---|
719 | for (i = 0; i <= nextree - 2; i++) { |
---|
720 | root = treenode[0]; |
---|
721 | for (j = k; j <= spp; j++) |
---|
722 | add(treenode[bestrees[i][j - 1] - 1], |
---|
723 | treenode[bestorders[i][j - 1] - 1], treenode[spp + j - 2], |
---|
724 | &root, recompute, treenode, grbg, zeros); |
---|
725 | dnapenny_reroot(treenode[outgrno - 1]); |
---|
726 | postorder(root); |
---|
727 | evaluate(root); |
---|
728 | printree(root, 1.0); |
---|
729 | describe(); |
---|
730 | for (j = k - 1; j < spp; j++) |
---|
731 | re_move(treenode[bestorders[i][j] - 1], &dummy, &root, |
---|
732 | recompute, treenode, grbg, zeros); |
---|
733 | } |
---|
734 | if (progress) { |
---|
735 | printf("\nOutput written to file \"%s\"\n\n", outfilename); |
---|
736 | if (trout) |
---|
737 | printf("Trees also written onto file \"%s\"\n\n", outtreename); |
---|
738 | } |
---|
739 | freetemp(&temp); |
---|
740 | freetemp(&temp1); |
---|
741 | if (ancseq) |
---|
742 | freegarbage(&garbage); |
---|
743 | } /* maketree */ |
---|
744 | |
---|
745 | |
---|
746 | int main(int argc, Char *argv[]) |
---|
747 | { /* Penny's branch-and-bound method for DNA sequences */ |
---|
748 | #ifdef MAC |
---|
749 | argc = 1; /* macsetup("Dnapenny",""); */ |
---|
750 | argv[0] = "Dnapenny"; |
---|
751 | #endif |
---|
752 | init(argc, argv); |
---|
753 | |
---|
754 | /* Reads in the number of species, number of characters, |
---|
755 | options and data. Then finds all most parsimonious trees */ |
---|
756 | openfile(&infile,INFILE,"input file", "r",argv[0],infilename); |
---|
757 | openfile(&outfile,OUTFILE,"output file", "w",argv[0],outfilename); |
---|
758 | |
---|
759 | ibmpc = IBMCRT; |
---|
760 | ansi = ANSICRT; |
---|
761 | mulsets = false; |
---|
762 | garbage = NULL; |
---|
763 | msets = 1; |
---|
764 | firstset = true; |
---|
765 | doinit(); |
---|
766 | if (weights || justwts) |
---|
767 | openfile(&weightfile,WEIGHTFILE,"weights file","r",argv[0],weightfilename); |
---|
768 | if (trout) |
---|
769 | openfile(&outtree,OUTTREE,"output tree file", "w",argv[0],outtreename); |
---|
770 | for (ith = 1; ith <= msets; ith++) { |
---|
771 | doinput(); |
---|
772 | if (ith == 1) |
---|
773 | firstset = false; |
---|
774 | if (msets > 1 && !justwts) { |
---|
775 | fprintf(outfile, "\nData set # %ld:\n",ith); |
---|
776 | if (progress) |
---|
777 | printf("\nData set # %ld:\n",ith); |
---|
778 | } |
---|
779 | maketree(); |
---|
780 | free(threshwt); |
---|
781 | freenodes(nonodes,treenode); |
---|
782 | } |
---|
783 | FClose(infile); |
---|
784 | FClose(outfile); |
---|
785 | FClose(outtree); |
---|
786 | #ifdef MAC |
---|
787 | fixmacfile(outfilename); |
---|
788 | fixmacfile(outtreename); |
---|
789 | #endif |
---|
790 | #ifdef WIN32 |
---|
791 | phyRestoreConsoleAttributes(); |
---|
792 | #endif |
---|
793 | return 0; |
---|
794 | } /* Penny's branch-and-bound method for DNA sequences */ |
---|