-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathtools.c
6940 lines (6153 loc) · 262 KB
/
tools.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/* tools.c
*/
#include "paml.h"
#ifdef USE_GSL
#include <gsl/gsl_blas.h>
#include <gsl/gsl_linalg.h>
#endif
/************************
sequences
*************************/
char BASEs[] = "TCAGUYRMKSWHBVD-N?";
char *EquateBASE[] = { "T","C","A","G", "T", "TC","AG","CA","TG","CG","TA",
"TCA","TCG","CAG","TAG", "TCAG","TCAG","TCAG" };
char CODONs[256][4];
char AAs[] = "ARNDCQEGHILKMFPSTWYV-*?X";
char nChara[256], CharaMap[256][64];
char AA3Str[] = { "AlaArgAsnAspCysGlnGluGlyHisIleLeuLysMetPheProSerThrTrpTyrVal***" };
char BINs[] = "TC";
int GeneticCode[][64] =
{ {13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,-1,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9, 9,12,16,16,16,16, 2, 2,11,11,15,15, 1, 1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 0:universal */
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,17,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9,12,12,16,16,16,16, 2, 2,11,11,15,15,-1,-1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 1:vertebrate mt.*/
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,17,17,
16,16,16,16,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9,12,12,16,16,16,16, 2, 2,11,11,15,15, 1, 1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 2:yeast mt. */
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,17,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9, 9,12,16,16,16,16, 2, 2,11,11,15,15, 1, 1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 3:mold mt. */
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,17,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9,12,12,16,16,16,16, 2, 2,11,11,15,15,15,15,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 4:invertebrate mt. */
{13,13,10,10,15,15,15,15,18,18, 5, 5, 4, 4,-1,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9, 9,12,16,16,16,16, 2, 2,11,11,15,15, 1, 1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 5:ciliate nuclear*/
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,17,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9, 9,12,16,16,16,16, 2, 2, 2,11,15,15,15,15,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 6:echinoderm mt.*/
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4, 4,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9, 9,12,16,16,16,16, 2, 2,11,11,15,15, 1, 1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 7:euplotid mt. */
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,-1,17,
10,10,10,15,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9, 9,12,16,16,16,16, 2, 2,11,11,15,15, 1, 1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7},
/* 8:alternative yeast nu.*/
{13,13,10,10,15,15,15,15,18,18,-1,-1, 4, 4,17,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9,12,12,16,16,16,16, 2, 2,11,11,15,15, 7, 7,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 9:ascidian mt. */
{13,13,10,10,15,15,15,15,18,18,-1, 5, 4, 4,-1,17,
10,10,10,10,14,14,14,14, 8, 8, 5, 5, 1, 1, 1, 1,
9, 9, 9,12,16,16,16,16, 2, 2,11,11,15,15, 1, 1,
19,19,19,19, 0, 0, 0, 0, 3, 3, 6, 6, 7, 7, 7, 7}, /* 10:blepharisma nu.*/
{ 1, 1, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3, 4, 4, 4, 4,
5, 5, 5, 5, 6, 6, 6, 6, 7, 7, 7, 7, 8, 8, 8, 8,
9, 9, 9, 9,10,10,10,10,11,11,11,11,12,12,12,12,
13,13,13,13,14,14,14,14,15,15,15,15,16,16,16,16} /* 11:Ziheng's regular code */
}; /* GeneticCode[icode][#codon] */
int noisy = 0, Iround = 0, NFunCall = 0, NEigenQ, NPMatUVRoot;
double SIZEp = 0;
int blankline(char *str)
{
char *p = str;
while (*p) if (isalnum(*p++)) return(0);
return(1);
}
int PopEmptyLines(FILE* fseq, int lline, char line[])
{
/* pop out empty lines in the sequence data file.
returns -1 if EOF.
*/
char *eqdel = ".-?", *p;
int i;
for (i = 0; ; i++) {
p = fgets(line, lline, fseq);
if (p == NULL) return(-1);
while (*p)
if (*p == eqdel[0] || *p == eqdel[1] || *p == eqdel[2] || isalpha(*p))
/*
if (*p==eqdel[0] || *p==eqdel[1] || *p==eqdel[2] || isalnum(*p))
*/
return(0);
else p++;
}
}
int picksite(char *z, int l, int begin, int gap, char *result)
{
/* pick every gap-th site, e.g., the third codon position for example.
*/
int il = begin;
for (il = 0, z += begin; il < l; il += gap, z += gap) *result++ = *z;
return(0);
}
int CodeChara(char b, int seqtype)
{
/* This codes nucleotides or amino acids into 0, 1, 2, ...
*/
int i, n = (seqtype <= 1 ? 4 : (seqtype == 2 ? 20 : 2));
char *pch = (seqtype <= 1 ? BASEs : (seqtype == 2 ? AAs : BINs));
if (seqtype <= 1)
switch (b) {
case 'T': case 'U': return(0);
case 'C': return(1);
case 'A': return(2);
case 'G': return(3);
}
else
for (i = 0; i < n; i++)
if (b == pch[i]) return (i);
if (noisy >= 9) printf("\nwarning: strange character '%c' ", b);
return (-1);
}
int dnamaker(char z[], int ls, double pi[])
{
/* sequences z[] are coded 0,1,2,3
*/
int i, j;
double p[4], r, smallv = 1e-5;
xtoy(pi, p, 4);
for (i = 1; i < 4; i++) p[i] += p[i - 1];
if (fabs(p[3] - 1) > smallv)
zerror("sum pi != 1..");
for (i = 0; i < ls; i++) {
for (j = 0, r = rndu(); j < 4; j++)
if (r < p[j]) break;
z[i] = (char)j;
}
return (0);
}
int transform(char *z, int ls, int direction, int seqtype)
{
/* direction==1 from TCAG to 0123, ==0 from 0123 to TCGA.
*/
int il, status = 0;
char *p;
char *pch = (seqtype <= 1 ? BASEs : (seqtype == 2 ? AAs : BINs));
if (direction)
for (il = 0, p = z; il < ls; il++, p++) {
if ((*p = (char)CodeChara(*p, seqtype)) == (char)(-1)) status = -1;
}
else
for (il = 0, p = z; il < ls; il++, p++) *p = pch[(int)(*p)];
return (status);
}
int f_mono_di(FILE *fout, char *z, int ls, int iring,
double fb1[], double fb2[], double CondP[])
{
/* get mono- di- nucleitide frequencies.
*/
int i, j, il;
char *s;
double t1, t2;
t1 = 1. / (double)ls;
t2 = 1. / (double)(ls - 1 + iring);
for (i = 0; i < 4; fb1[i++] = 0.0) for (j = 0; j < 4; fb2[i * 4 + j++] = 0.0);
for (il = 0, s = z; il < ls - 1; il++, s++) {
fb1[*s - 1] += t1;
fb2[(*s - 1) * 4 + *(s + 1) - 1] += t2;
}
fb1[*s - 1] += t1;
if (iring) fb2[(*s - 1) * 4 + z[0] - 1] += t2;
for (i = 0; i < 4; i++) for (j = 0; j < 4; j++) CondP[i * 4 + j] = fb2[i * 4 + j] / fb1[i];
fprintf(fout, "\nmono-\n");
for (i = 0; i < 4; i++) fprintf(fout, "%12.4f", fb1[i]);
fprintf(fout, "\n\ndi- & conditional P\n");
for (i = 0; i < 4; i++) {
for (j = 0; j < 4; j++) fprintf(fout, "%9.4f%7.4f ", fb2[i * 4 + j], CondP[i * 4 + j]);
fprintf(fout, "\n");
}
fprintf(fout, "\n");
return (0);
}
int PickExtreme(FILE *fout, char *z, int ls, int iring, int lfrag, int *ffrag)
{
/* picking up (lfrag)-tuples with extreme frequencies.
*/
char *pz = z;
int i, j, isf, n = (1 << 2 * lfrag), lvirt = ls - (lfrag - 1)*(1 - iring);
double fb1[4], fb2[4 * 4], p_2[4 * 4];
double prob1, prob2, ne1, ne2, u1, u2, ualpha = 2.0;
int ib[10];
f_mono_di(fout, z, ls, iring, fb1, fb2, p_2);
if (iring) {
zerror("change PickExtreme()");
for (i = 0; i < lfrag - 1; i++) z[ls + i] = z[i];
z[ls + i] = (char)0;
}
printf("\ncounting %d tuple frequencies", lfrag);
for (i = 0; i < n; i++) ffrag[i] = 0;
for (i = 0; i < lvirt; i++, pz++) {
for (j = 0, isf = 0; j < lfrag; j++) isf = isf * 4 + (int)pz[j] - 1;
ffrag[isf] ++;
}
/* analyze */
for (i = 0; i < n; i++) {
for (j = 0, isf = i; j < lfrag; ib[lfrag - 1 - j] = isf % 4, isf = isf / 4, j++);
for (j = 0, prob1 = 1.0; j < lfrag; prob1 *= fb1[ib[j++]]);
for (j = 0, prob2 = fb1[ib[0]]; j < lfrag - 1; j++)
prob2 *= p_2[ib[j] * 4 + ib[j + 1]];
ne1 = (double)lvirt * prob1;
ne2 = (double)lvirt * prob2;
if (ne1 <= 0.0) ne1 = 0.5;
if (ne2 <= 0.0) ne2 = 0.5;
u1 = ((double)ffrag[i] - ne1) / sqrt(ne1);
u2 = ((double)ffrag[i] - ne2) / sqrt(ne2);
if (fabs(u1) > ualpha /* && fabs(u2)>ualpha */) {
fprintf(fout, "\n");
for (i = 0; i < lfrag; i++) fprintf(fout, "%1c", BASEs[ib[j]]);
fprintf(fout, "%6d %8.1f%7.2f %8.1f%7.2f ", ffrag[i], ne1, u1, ne2, u2);
if (u1 < -ualpha && u2 < -ualpha) fprintf(fout, " %c", '-');
else if (u1 > ualpha && u2 > ualpha) fprintf(fout, " %c", '+');
else if (u1*u2<0 && fabs(u1) > ualpha && fabs(u2) > ualpha)
fprintf(fout, " %c", '?');
else
fprintf(fout, " %c", ' ');
}
}
return (0);
}
int zztox(int n31, int l, char *z1, char *z2, double *x)
{
/* x[n31][4][4] */
double t = 1. / (double)(l / n31);
int i, ib[2];
int il;
zero(x, n31 * 16);
for (i = 0; i < n31; i++) {
for (il = 0; il < l; il += n31) {
ib[0] = z1[il + i] - 1;
ib[1] = z2[il + i] - 1;
x[i * 16 + ib[0] * 4 + ib[1]] += t;
}
/*
fprintf (f1, "\nThe difference matrix X %6d\tin %6d\n", i+1,n31);
for (j=0; j<4; j++) {
for (k=0; k<4; k++) fprintf(f1, "%10.2f", x[i][j][k]);
fputc ('\n', f1);
}
*/
}
return (0);
}
int testXMat(double x[])
{
/* test whether X matrix is acceptable (0) or not (-1) */
int it = 0, i, j;
double t;
for (i = 0, t = 0; i < 4; i++) for (j = 0; j < 4; j++) {
if (x[i * 4 + j] < 0 || x[i * 4 + j]>1) it = -1;
t += x[i * 4 + j];
}
if (fabs(t - 1) > 1e-4) it = -1;
return(it);
}
int difcodonNG(char codon1[], char codon2[], double *SynSite, double *AsynSite,
double *SynDif, double *AsynDif, int transfed, int icode)
{
/* # of synonymous and non-synonymous sites and differences.
Nei, M. and T. Gojobori (1986)
returns the number of differences between two codons.
The two codons (codon1 & codon2) do not contain ambiguity characters.
dmark[i] (=0,1,2) is the i_th different codon position, with i=0,1,ndiff
step[j] (=0,1,2) is the codon position to be changed at step j (j=0,1,ndiff)
b[i][j] (=0,1,2,3) is the nucleotide at position j (0,1,2) in codon i (0,1)
I made some arbitrary decisions when the two codons have ambiguity characters
20 September 2002.
*/
int i, j, k, i1, i2, iy[2] = { 0 }, iaa[2], ic[2];
int ndiff, npath, nstop, sdpath, ndpath, dmark[3], step[3], b[2][3], bt1[3], bt2[3];
int by[3] = { 16, 4, 1 };
char str[4] = "";
for (i = 0, *SynSite = 0, nstop = 0; i < 2; i++) {
for (j = 0; j < 3; j++) {
if (transfed) b[i][j] = (i ? codon1[j] : codon2[j]);
else b[i][j] = (int)CodeChara((char)(i ? codon1[j] : codon2[j]), 0);
iy[i] += by[j] * b[i][j];
if (b[i][j] < 0 || b[i][j]>3) {
if (noisy >= 9)
printf("\nwarning ambiguity in difcodonNG: %s %s", codon1, codon2);
*SynSite = 0.5; *AsynSite = 2.5;
*SynDif = (codon1[2] != codon2[2]) / 2;
*AsynDif = *SynDif + (codon1[0] != codon2[0]) + (codon1[1] != codon2[1]);
return((int)(*SynDif + *AsynDif));
}
}
iaa[i] = GeneticCode[icode][iy[i]];
if (iaa[i] == -1) {
printf("\nNG86: stop codon %s.\n", getcodon(str, iy[i]));
exit(-1);
}
for (j = 0; j < 3; j++)
for (k = 0; k < 4; k++) {
if (k == b[i][j]) continue;
i1 = GeneticCode[icode][iy[i] + (k - b[i][j])*by[j]];
if (i1 == -1)
nstop++;
else if (i1 == iaa[i])
(*SynSite)++;
}
}
*SynSite *= 3 / 18.; /* 2 codons, 2*9 possibilities. */
*AsynSite = 3 * (1 - nstop / 18.) - *SynSite;
#if 0 /* MEGA 1.1 */
* AsynSite = 3 - *SynSite;
#endif
ndiff = 0; *SynDif = *AsynDif = 0;
for (k = 0; k < 3; k++) dmark[k] = -1;
for (k = 0; k < 3; k++)
if (b[0][k] - b[1][k]) dmark[ndiff++] = k;
if (ndiff == 0) return(0);
npath = 1;
nstop = 0;
if (ndiff > 1)
npath = (ndiff == 2 ? 2 : 6);
if (ndiff == 1) {
if (iaa[0] == iaa[1]) (*SynDif)++;
else (*AsynDif)++;
}
else { /* ndiff=2 or 3 */
for (k = 0; k < npath; k++) {
for (i1 = 0; i1 < 3; i1++)
step[i1] = -1;
if (ndiff == 2) {
step[0] = dmark[k];
step[1] = dmark[1 - k];
}
else {
step[0] = k / 2; step[1] = k % 2;
if (step[0] <= step[1]) step[1]++;
step[2] = 3 - step[0] - step[1];
}
for (i1 = 0; i1 < 3; i1++)
bt1[i1] = bt2[i1] = b[0][i1];
sdpath = ndpath = 0; /* mutations for each path */
for (i1 = 0; i1 < ndiff; i1++) { /* mutation steps for each path */
bt2[step[i1]] = b[1][step[i1]];
for (i2 = 0, ic[0] = ic[1] = 0; i2 < 3; i2++) {
ic[0] += bt1[i2] * by[i2];
ic[1] += bt2[i2] * by[i2];
}
for (i2 = 0; i2 < 2; i2++) iaa[i2] = GeneticCode[icode][ic[i2]];
if (iaa[1] == -1) {
nstop++; sdpath = ndpath = 0; break;
}
if (iaa[0] == iaa[1]) sdpath++;
else ndpath++;
for (i2 = 0; i2 < 3; i2++)
bt1[i2] = bt2[i2];
}
*SynDif += (double)sdpath;
*AsynDif += (double)ndpath;
}
}
if (npath == nstop) {
puts("NG86: All paths are through stop codons..");
if (ndiff == 2) { *SynDif = 0; *AsynDif = 2; }
else { *SynDif = 1; *AsynDif = 2; }
}
else {
*SynDif /= (double)(npath - nstop); *AsynDif /= (double)(npath - nstop);
}
return (ndiff);
}
int difcodonLWL85(char codon1[], char codon2[], double sites[3], double sdiff[3],
double vdiff[3], int transfed, int icode)
{
/* This partitions codon sites according to degeneracy, that is sites[3] has
L0, L2, L4, averaged between the two codons. It also compares the two codons
and add the differences to the transition and transversion differences for use
in LWL85 and similar methods.
The two codons (codon1 & codon2) should not contain ambiguity characters.
c[0] & c[1] are the two codons, coded 0, 1, ..., 63.
b[][] has the nucleotides, coded 0123 for TCAG.
*/
int b[2][3], by[3] = { 16, 4, 1 }, i, j, ifold[2], c[2], ct, aa[2], ibase, nsame;
char str[4] = "";
for (i = 0; i < 3; i++) sites[i] = sdiff[i] = vdiff[i] = 0;
/* check the two codons and code them */
for (i = 0; i < 2; i++) {
for (j = 0, c[i] = 0; j < 3; j++) {
if (transfed) b[i][j] = (i ? codon1[j] : codon2[j]);
else b[i][j] = (int)CodeChara((char)(i ? codon1[j] : codon2[j]), 0);
c[i] += b[i][j] * by[j];
if (b[i][j] < 0 || b[i][j]>3) {
if (noisy >= 9)
printf("\nwarning ambiguity in difcodonLWL85: %s %s", codon1, codon2);
return(0);
}
}
aa[i] = GeneticCode[icode][c[i]];
if (aa[i] == -1) {
printf("\nLWL85: stop codon %s.\n", getcodon(str, c[i]));
exit(-1);
}
}
for (j = 0; j < 3; j++) { /* three codon positions */
for (i = 0; i < 2; i++) { /* two codons */
for (ibase = 0, nsame = 0; ibase < 4; ibase++) { /* change codon i pos j into ibase */
ct = c[i] + (ibase - b[i][j])*by[j];
if (ibase != b[i][j] && aa[i] == GeneticCode[icode][ct]) nsame++;
}
if (nsame == 0) ifold[i] = 0; /* codon i pos j is 0-fold */
else if (nsame == 1 || nsame == 2) ifold[i] = 1; /* codon i pos j is 2-fold */
else ifold[i] = 2; /* codon i pos j is 4-fold */
sites[ifold[i]] += .5;
}
if (b[0][j] == b[1][j]) continue;
if (b[0][j] + b[1][j] == 1 || b[0][j] + b[1][j] == 5) { /* pos j has a transition */
sdiff[ifold[0]] += .5; sdiff[ifold[1]] += .5;
}
else { /* pos j has a transversion */
vdiff[ifold[0]] += .5; vdiff[ifold[1]] += .5;
}
}
return (0);
}
int testTransP(double P[], int n)
{
int i, j, status = 0;
double sum, smallv = 1e-10;
for (i = 0; i < n; i++) {
for (j = 0, sum = 0; j < n; sum += P[i*n + j++])
if (P[i*n + j] < -smallv) status = -1;
if (fabs(sum - 1) > smallv && status == 0) {
printf("\nrow sum (#%2d) = 1 = %10.6f", i + 1, sum);
status = -1;
}
}
return (status);
}
double testDetailedBalance(double P[], double pi[], int n)
{
/* this calculates maxdiff for the detailed balance check. maxdiff should be close
to 0 if the detailed balance condition holds.
*/
int i, j;
double maxdiff = 0, d;
for (i = 0; i < n; i++) {
for (j = 0; j < n; j++) {
d = fabs(pi[i] * P[i * n + j] - pi[j] * P[j * n + i]);
if (d > maxdiff) maxdiff = d;
}
}
return (maxdiff);
}
int PMatUVRoot(double P[], double t, int n, double U[], double V[], double Root[])
{
/* P(t) = U * exp{Root*t} * V
*/
int i, j, k;
double exptm1, uexpt, * pP, smallp = 0;
NPMatUVRoot++;
if (t < -0.01) printf("\nt = %.5f in PMatUVRoot", t);
if (t < 1e-100) {
identity(P, n);
return(0);
}
memset(P, 0, n*n * sizeof(double));
for (k = 0; k < n; k++) {
for (i = 0, pP = P, exptm1 = expm1(t*Root[k]); i < n; i++)
for (j = 0, uexpt = U[i*n + k] * exptm1; j < n; j++)
*pP++ += uexpt*V[k*n + j];
}
for (i = 0; i < n; i++) P[i*n+i] ++;
for (i = 0; i < n * n; i++)
if (P[i] < smallp) P[i] = 0;
#if (DEBUG>=5)
if (testTransP(P, n)) {
printf("\nP(%.6f) err in PMatUVRoot.\n", t);
exit(-1);
}
#endif
return (0);
}
int PMatQRev(double Q[], double pi[], double t, int n, double space[])
{
/* This calculates P(t) = exp(Q*t), where Q is the rate matrix for a
time-reversible Markov process.
Q[] or P[] has the rate matrix as input, and P(t) in return.
space[n*n*2+n*2]
*/
double *U = space, *V = U + n*n, *Root = V + n*n, *spacesqrtpi = Root + n;
eigenQREV(Q, pi, n, Root, U, V, spacesqrtpi);
PMatUVRoot(Q, t, n, U, V, Root);
return(0);
}
void pijJC69(double pij[2], double t)
{
double exptm1;
if (t < -1e-6)
printf("\nt = %.5f in pijJC69", t);
exptm1 = expm1(-4 * t / 3.);
pij[0] = 1. + 3/4.0 * exptm1;
pij[1] = -exptm1 / 4;
}
int PMatK80(double P[], double t, double kappa)
{
/* PMat for JC69 and K80
*/
int i, j;
double e1, e2;
if (t < -1e-6)
printf("\nt = %.5f in PMatK80", t);
e1 = expm1(-4 * t / (kappa + 2));
if (fabs(kappa - 1) < 1e-20) {
for (i = 0; i < 4; i++)
for (j = 0; j < 4; j++)
if (i == j) P[i * 4 + j] = 1. + 3 / 4.0 * e1;
else P[i * 4 + j] = -e1 / 4;
}
else {
e2 = expm1(-2 * t*(kappa + 1) / (kappa + 2));
for (i = 0; i < 4; i++)
P[i * 4 + i] = 1 + (e1 + 2 * e2) / 4;
P[0 * 4 + 1] = P[1 * 4 + 0] = P[2 * 4 + 3] = P[3 * 4 + 2] = (e1 - 2 * e2) / 4;
P[0 * 4 + 2] = P[0 * 4 + 3] = P[2 * 4 + 0] = P[3 * 4 + 0] =
P[1 * 4 + 2] = P[1 * 4 + 3] = P[2 * 4 + 1] = P[3 * 4 + 1] = -e1 / 4;
}
return (0);
}
int PMatT92(double P[], double t, double kappa, double pGC)
{
/* PMat for Tamura'92
t is branch lnegth, number of changes per site.
*/
double e1, e2;
t /= (pGC*(1 - pGC)*kappa + .5);
if (t < -1e-4) printf("\nt = %.5f in PMatT92", t);
e1 = expm1(-t);
e2 = expm1(-(kappa + 1)*t / 2);
P[0 * 4 + 0] = P[2 * 4 + 2] = 1 + 0.5*(1 - pGC) * e1 + pGC*e2;
P[1 * 4 + 1] = P[3 * 4 + 3] = 1 + pGC / 2 * e1 + (1 - pGC)*e2;
P[1 * 4 + 0] = P[3 * 4 + 2] = (1 - pGC) / 2 * e1 - (1 - pGC)*e2;
P[0 * 4 + 1] = P[2 * 4 + 3] = pGC / 2 * e1 - pGC*e2;
P[0 * 4 + 2] = P[2 * 4 + 0] = P[3 * 4 + 0] = P[1 * 4 + 2] = -(1 - pGC) / 2 * e1;
P[1 * 4 + 3] = P[3 * 4 + 1] = P[0 * 4 + 3] = P[2 * 4 + 1] = -pGC / 2 * e1;
return (0);
}
int PMatTN93(double P[], double a1t, double a2t, double bt, double pi[])
{
double T = pi[0], C = pi[1], A = pi[2], G = pi[3], Y = T + C, R = A + G;
double e1, e2, e3, smallv = -1e-6;
if (noisy && (a1t < smallv || a2t < smallv || bt < smallv))
printf("\nat=%12.6f %12.6f bt=%12.6f", a1t, a2t, bt);
e1 = expm1(-bt);
e2 = expm1(-(R*a2t + Y*bt));
e3 = expm1(-(Y*a1t + R*bt));
P[0 * 4 + 0] = 1 + (R*T*e1 + C*e3) / Y;
P[0 * 4 + 1] = (R*e1 - e3)*C / Y;
P[0 * 4 + 2] = -A*e1;
P[0 * 4 + 3] = -G*e1;
P[1 * 4 + 0] = (R*e1 - e3)*T / Y;
P[1 * 4 + 1] = 1 + (R*C*e1 + T*e3) / Y;
P[1 * 4 + 2] = -A*e1;
P[1 * 4 + 3] = -G*e1;
P[2 * 4 + 0] = -T*e1;
P[2 * 4 + 1] = -C*e1;
P[2 * 4 + 2] = 1 + Y*A / R*e1 + G / R*e2;
P[2 * 4 + 3] = Y*G / R*e1 - G / R*e2;
P[3 * 4 + 0] = -T*e1;
P[3 * 4 + 1] = -C*e1;
P[3 * 4 + 2] = Y*A / R*e1 - A / R*e2;
P[3 * 4 + 3] = 1 + Y*G / R*e1 + A / R*e2;
return(0);
}
int EvolveHKY85(char source[], char target[], int ls, double t,
double rates[], double pi[4], double kappa, int isHKY85)
{
/* isHKY85=1 if HKY85, =0 if F84
Use NULL for rates if rates are identical among sites.
*/
int i, j, h, n = 4;
double TransP[16], a1t, a2t, bt, r, Y = pi[0] + pi[1], R = pi[2] + pi[3];
if (isHKY85) a1t = a2t = kappa;
else { a1t = 1 + kappa / Y; a2t = 1 + kappa / R; }
bt = t / (2 * (pi[0] * pi[1] * a1t + pi[2] * pi[3] * a2t) + 2 * Y*R);
a1t *= bt; a2t *= bt;
for (h = 0; h < ls; h++) {
if (h == 0 || (rates && rates[h] != rates[h - 1])) {
r = (rates ? rates[h] : 1);
PMatTN93(TransP, a1t*r, a2t*r, bt*r, pi);
for (i = 0; i < n; i++) {
for (j = 1; j < n; j++) TransP[i*n + j] += TransP[i*n + j - 1];
if (fabs(TransP[i*n + n - 1] - 1) > 1e-5) zerror("TransP err");
}
}
for (j = 0, i = source[h], r = rndu(); j < n - 1; j++) if (r < TransP[i*n + j]) break;
target[h] = (char)j;
}
return (0);
}
int Rates4Sites(double rates[], double alpha, int ncatG, int ls, int cdf, double space[])
{
/* Rates for sites from the gamma (ncatG=0) or discrete-gamma (ncatG>1).
Rates are converted into the c.d.f. if cdf=1, which is useful for
simulation under JC69-like models.
space[ncatG*5]
*/
int h, ir, j, K = ncatG, *Lalias = (int*)(space + 3 * K), *counts = (int*)(space + 4 * K);
double *rK = space, *freqK = space + K, *Falias = space + 2 * K;
if (alpha == 0) {
if (rates)
for (h = 0; h < ls; h++) rates[h] = 1;
}
else {
if (K > 1) {
DiscreteGamma(freqK, rK, alpha, alpha, K, DGammaUseMedian);
MultiNomialAliasSetTable(K, freqK, Falias, Lalias, space + 5 * K);
MultiNomialAlias(ls, K, Falias, Lalias, counts);
for (ir = 0, h = 0; ir < K; ir++)
for (j = 0; j < counts[ir]; j++) rates[h++] = rK[ir];
}
else
for (h = 0; h < ls; h++) rates[h] = rndgamma(alpha) / alpha;
if (cdf) {
for (h = 1; h < ls; h++) rates[h] += rates[h - 1];
abyx(1 / rates[ls - 1], rates, ls); /* this rescaling may not be needed. */
}
}
return (0);
}
char *getcodon(char codon[], int icodon)
{
/* id : (0,63) */
if (icodon < 0 || icodon>63) {
printf("\ncodon %d\n", icodon);
zerror("getcodon.");
}
codon[0] = BASEs[icodon / 16];
codon[1] = BASEs[(icodon % 16) / 4];
codon[2] = BASEs[icodon % 4];
codon[3] = 0;
return (codon);
}
char *getAAstr(char *AAstr, int iaa)
{
/* iaa (0,20) with 20 meaning termination */
if (iaa < 0 || iaa>20) zerror("getAAstr: iaa err. \n");
strncpy(AAstr, AA3Str + iaa * 3, 3);
return (AAstr);
}
int NucListall(char b, int *nb, int ib[4])
{
/* Resolve an ambiguity nucleotide b into all possibilities.
nb is number of bases and ib (0,1,2,3) list all of them.
Data are complete if (nb==1).
*/
int j, k;
k = (int)(strchr(BASEs, (int)b) - BASEs);
if (k < 0)
{
printf("NucListall: strange character %c\n", b); return(-1);
}
if (k < 4) {
*nb = 1; ib[0] = k;
}
else {
*nb = (int)strlen(EquateBASE[k]);
for (j = 0; j < *nb; j++)
ib[j] = (int)(strchr(BASEs, EquateBASE[k][j]) - BASEs);
}
return(0);
}
int Codon2AA(char codon[3], char aa[3], int icode, int *iaa)
{
/* translate a triplet codon[] into amino acid (aa[] and iaa), using
genetic code icode. This deals with ambiguity nucleotides.
*iaa=(0,...,19), 20 for stop or missing data.
Distinquish between stop codon and missing data?
naa=0: only stop codons; 1: one AA; 2: more than 1 AA.
Returns 0: if one amino acid
1: if multiple amino acids (ambiguity data)
-1: if stop codon
*/
int nb[3], ib[3][4], ic, i, i0, i1, i2, iaa0 = -1, naa = 0;
for (i = 0; i < 3; i++)
NucListall(codon[i], &nb[i], ib[i]);
for (i0 = 0; i0 < nb[0]; i0++)
for (i1 = 0; i1 < nb[1]; i1++)
for (i2 = 0; i2 < nb[2]; i2++) {
ic = ib[0][i0] * 16 + ib[1][i1] * 4 + ib[2][i2];
*iaa = GeneticCode[icode][ic];
if (*iaa == -1) continue;
if (naa == 0) { iaa0 = *iaa; naa++; }
else if (*iaa != iaa0) naa = 2;
}
if (naa == 0) {
printf("stop codon %c%c%c\n", codon[0], codon[1], codon[2]);
*iaa = 20;
}
else if (naa == 2) *iaa = 20;
else *iaa = iaa0;
memcpy(aa, AA3Str + *iaa * 3, 3*sizeof(char));
return(naa == 1 ? 0 : (naa == 0 ? -1 : 1));
}
int DNA2protein(char dna[], char protein[], int lc, int icode)
{
/* translate a DNA into a protein, using genetic code icode, with lc codons.
dna[] and protein[] can be the same string.
*/
int h, iaa, k;
char aa3[4];
for (h = 0; h < lc; h++) {
k = Codon2AA(dna + h * 3, aa3, icode, &iaa);
if (k == -1) printf(" stop codon at %d out of %d\n", h + 1, lc);
protein[h] = AAs[iaa];
}
return(0);
}
int printcu(FILE *fout, double fcodon[], int icode)
{
/* output codon usage table and other related statistics
space[20+1+3*5]
Outputs the genetic code table if fcodon==NULL
*/
int wc = 8, wd = 0; /* wc: for codon, wd: decimal */
int it, i, j, k, iaa;
double faa[21], fb3x4[3 * 5]; /* chi34, Ic, lc, */
char *word = "|-", aa3[4] = " ", codon[4] = " ", ss3[4][4], *noodle;
/*
static double aawt[] = { 89.1, 174.2, 132.1, 133.1, 121.2, 146.2,
147.1, 75.1, 155.2, 131.2, 131.2, 146.2, 149.2, 165.2, 115.1,
105.1, 119.1, 204.2, 181.2, 117.1 };
*/
if (fcodon) { zero(faa, 21); zero(fb3x4, 12); }
else wc = 0;
for (i = 0; i < 4; i++) strcpy(ss3[i], "\0\0\0");
noodle = strc(4 * (10 + 2 + wc) - 2, word[1]);
fprintf(fout, "\n%s\n", noodle);
for (i = 0; i < 4; i++) {
for (j = 0; j < 4; j++) {
for (k = 0; k < 4; k++) {
it = i * 16 + k * 4 + j;
iaa = GeneticCode[icode][it];
if (iaa == -1) iaa = 20;
getcodon(codon, it); getAAstr(aa3, iaa);
if (!strcmp(ss3[k], aa3) && j > 0)
fprintf(fout, " ");
else {
fprintf(fout, "%s %c", aa3, (iaa < 20 ? AAs[iaa] : '*'));
strcpy(ss3[k], aa3);
}
fprintf(fout, " %s", codon);
if (fcodon) fprintf(fout, "%*.*f", wc, wd, fcodon[it]);
if (k < 3) fprintf(fout, " %c ", word[0]);
}
fprintf(fout, "\n");
}
fprintf(fout, "%s\n", noodle);
}
return(0);
}
int printcums(FILE *fout, int ns, double fcodons[], int icode)
{
int neach0 = 6, neach = neach0, wc = 4, wd = 0; /* wc: for codon, wd: decimal */
int iaa, it, i, j, k, i1, ngroup, igroup;
char *word = "|-", aa3[4] = " ", codon[4] = " ", ss3[4][4], *noodle;
ngroup = (ns - 1) / neach + 1;
for (igroup = 0; igroup < ngroup; igroup++) {
if (igroup == ngroup - 1)
neach = ns - neach0*igroup;
noodle = strc(4 * (10 + wc*neach) - 2, word[1]);
strcat(noodle, "\n");
fputs(noodle, fout);
for (i = 0; i < 4; i++) strcpy(ss3[i], " ");
for (i = 0; i < 4; i++) {
for (j = 0; j < 4; j++) {
for (k = 0; k < 4; k++) {
it = i * 16 + k * 4 + j;
iaa = GeneticCode[icode][it];
if (iaa == -1) iaa = 20;
getcodon(codon, it);
getAAstr(aa3, iaa);
if (!strcmp(ss3[k], aa3) && j > 0) fprintf(fout, " ");
else { fprintf(fout, "%s", aa3); strcpy(ss3[k], aa3); }
fprintf(fout, " %s", codon);
for (i1 = 0; i1 < neach; i1++)
fprintf(fout, " %*.*f", wc - 1, wd, fcodons[(igroup*neach0 + i1) * 64 + it]);
if (k < 3) fprintf(fout, " %c ", word[0]);
}
fprintf(fout, "\n");
}
fputs(noodle, fout);
}
fprintf(fout, "\n");
}
return(0);
}
int QtoPi(double Q[], double pi[], int n, double space[])
{
/* from rate matrix Q[] to pi, the stationary frequencies:
Q' * pi = 0 pi * 1 = 1
space[] is of size n*(n+1).
*/
int i, j;
double *T = space; /* T[n*(n+1)] */
for (i = 0; i < n + 1; i++) T[i] = 1;
for (i = 1; i < n; i++) {
for (j = 0; j < n; j++)
T[i*(n + 1) + j] = Q[j*n + i]; /* transpose */
T[i*(n + 1) + n] = 0.;
}
matinv(T, n, n + 1, pi);
for (i = 0; i < n; i++)
pi[i] = T[i*(n + 1) + n];
return (0);
}
int PtoPi(double P[], double pi[], int n, double space[])
{
/* from transition probability P[ij] to pi, the stationary frequencies
(P'-I) * pi = 0 pi * 1 = 1
space[] is of size n*(n+1).
*/
int i, j;
double *T = space; /* T[n*(n+1)] */
for (i = 0; i < n + 1; i++) T[i] = 1;
for (i = 1; i < n; i++) {
for (j = 0; j < n; j++)
T[i*(n + 1) + j] = P[j*n + i] - (double)(i == j); /* transpose */
T[i*(n + 1) + n] = 0;
}
matinv(T, n, n + 1, pi);
for (i = 0; i < n; i++) pi[i] = T[i*(n + 1) + n];
return (0);
}
int PtoX(double P1[], double P2[], double pi[], double X[])
{
/* from P1 & P2 to X. X = P1' diag{pi} P2
*/
int i, j, k;
for (i = 0; i < 4; i++)
for (j = 0; j < 4; j++)
for (k = 0, X[i * 4 + j] = 0.0; k < 4; k++) {
X[i * 4 + j] += pi[k] * P1[k * 4 + i] * P2[k * 4 + j];
}
return (0);
}
int ScanFastaFile(FILE *fin, int *ns, int *ls, int *aligned)
{
/* This scans a fasta alignment file to get com.ns & com.ls.
Returns -1 if the sequences are not aligned and have different lengths.
*/
int len = 0, ch, starter = '>', stop = '/'; /* both EOF and / mark the end of the file. */
char name[200], *p;
if (noisy) printf("\nprocessing fasta file");
for (*aligned = 1, *ns = -1, *ls = 0; ; ) {
ch = fgetc(fin);
if (ch == starter || ch == EOF || ch == stop) {
if (*ns >= 0) { /* process end of the sequence */
if (noisy) printf(" %7d sites", len);
if (*ns > 1 && len != *ls) {
*aligned = 0;
printf("previous sequence %s has len %d, current seq has %d\n", name, *ls, len);
}
if (len > *ls) *ls = len;
}
(*ns)++; /* next sequence */
if (ch == EOF || ch == stop) break;
/* fscanf(fin, "%s", name); */
p = name;
while ((ch = getc(fin)) != '\n' && ch != EOF) *p++ = ch;
*p = '\0';
if (noisy) printf("\nreading seq#%2d %-50s", *ns + 1, name);
len = 0;
}