]> git.donarmstrong.com Git - samtools.git/blob - bcftools/prob1.c
Fixed errors introduced by the 8c15f916dabce475febdf508a9cc0c708c5a7747 commit.
[samtools.git] / bcftools / prob1.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <errno.h>
6 #include <assert.h>
7 #include "prob1.h"
8
9 #include "kseq.h"
10 KSTREAM_INIT(gzFile, gzread, 16384)
11
12 #define MC_MAX_EM_ITER 16
13 #define MC_EM_EPS 1e-5
14 #define MC_DEF_INDEL 0.15
15
16 unsigned char seq_nt4_table[256] = {
17         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
18         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
19         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4 /*'-'*/, 4, 4,
20         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
21         4, 0, 4, 1,  4, 4, 4, 2,  4, 4, 4, 4,  4, 4, 4, 4, 
22         4, 4, 4, 4,  3, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
23         4, 0, 4, 1,  4, 4, 4, 2,  4, 4, 4, 4,  4, 4, 4, 4, 
24         4, 4, 4, 4,  3, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
25         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
26         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
27         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
28         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
29         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
30         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
31         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
32         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4
33 };
34
35 struct __bcf_p1aux_t {
36         int n, M, n1, is_indel;
37         uint8_t *ploidy; // haploid or diploid ONLY
38         double *q2p, *pdg; // pdg -> P(D|g)
39         double *phi, *phi_indel;
40         double *z, *zswap; // aux for afs
41         double *z1, *z2, *phi1, *phi2; // only calculated when n1 is set
42         double **hg; // hypergeometric distribution
43         double *lf; // log factorial
44         double t, t1, t2;
45         double *afs, *afs1; // afs: accumulative AFS; afs1: site posterior distribution
46         const uint8_t *PL; // point to PL
47         int PL_len;
48 };
49
50 void bcf_p1_indel_prior(bcf_p1aux_t *ma, double x)
51 {
52         int i;
53         for (i = 0; i < ma->M; ++i)
54                 ma->phi_indel[i] = ma->phi[i] * x;
55         ma->phi_indel[ma->M] = 1. - ma->phi[ma->M] * x;
56 }
57
58 static void init_prior(int type, double theta, int M, double *phi)
59 {
60         int i;
61         if (type == MC_PTYPE_COND2) {
62                 for (i = 0; i <= M; ++i)
63                         phi[i] = 2. * (i + 1) / (M + 1) / (M + 2);
64         } else if (type == MC_PTYPE_FLAT) {
65                 for (i = 0; i <= M; ++i)
66                         phi[i] = 1. / (M + 1);
67         } else {
68                 double sum;
69                 for (i = 0, sum = 0.; i < M; ++i)
70                         sum += (phi[i] = theta / (M - i));
71                 phi[M] = 1. - sum;
72         }
73 }
74
75 void bcf_p1_init_prior(bcf_p1aux_t *ma, int type, double theta)
76 {
77         init_prior(type, theta, ma->M, ma->phi);
78         bcf_p1_indel_prior(ma, MC_DEF_INDEL);
79 }
80
81 void bcf_p1_init_subprior(bcf_p1aux_t *ma, int type, double theta)
82 {
83         if (ma->n1 <= 0 || ma->n1 >= ma->M) return;
84         init_prior(type, theta, 2*ma->n1, ma->phi1);
85         init_prior(type, theta, 2*(ma->n - ma->n1), ma->phi2);
86 }
87
88 int bcf_p1_read_prior(bcf_p1aux_t *ma, const char *fn)
89 {
90         gzFile fp;
91         kstring_t s;
92         kstream_t *ks;
93         long double sum;
94         int dret, k;
95         memset(&s, 0, sizeof(kstring_t));
96         fp = strcmp(fn, "-")? gzopen(fn, "r") : gzdopen(fileno(stdin), "r");
97         ks = ks_init(fp);
98         memset(ma->phi, 0, sizeof(double) * (ma->M + 1));
99         while (ks_getuntil(ks, '\n', &s, &dret) >= 0) {
100                 if (strstr(s.s, "[afs] ") == s.s) {
101                         char *p = s.s + 6;
102                         for (k = 0; k <= ma->M; ++k) {
103                                 int x;
104                                 double y;
105                                 x = strtol(p, &p, 10);
106                                 if (x != k && (errno == EINVAL || errno == ERANGE)) return -1;
107                                 ++p;
108                                 y = strtod(p, &p);
109                                 if (y == 0. && (errno == EINVAL || errno == ERANGE)) return -1;
110                                 ma->phi[ma->M - k] += y;
111                         }
112                 }
113         }
114         ks_destroy(ks);
115         gzclose(fp);
116         free(s.s);
117         for (sum = 0., k = 0; k <= ma->M; ++k) sum += ma->phi[k];
118         fprintf(stderr, "[prior]");
119         for (k = 0; k <= ma->M; ++k) ma->phi[k] /= sum;
120         for (k = 0; k <= ma->M; ++k) fprintf(stderr, " %d:%.3lg", k, ma->phi[ma->M - k]);
121         fputc('\n', stderr);
122         for (sum = 0., k = 1; k < ma->M; ++k) sum += ma->phi[ma->M - k] * (2.* k * (ma->M - k) / ma->M / (ma->M - 1));
123         fprintf(stderr, "[%s] heterozygosity=%lf, ", __func__, (double)sum);
124         for (sum = 0., k = 1; k <= ma->M; ++k) sum += k * ma->phi[ma->M - k] / ma->M;
125         fprintf(stderr, "theta=%lf\n", (double)sum);
126         bcf_p1_indel_prior(ma, MC_DEF_INDEL);
127         return 0;
128 }
129
130 bcf_p1aux_t *bcf_p1_init(int n, uint8_t *ploidy)
131 {
132         bcf_p1aux_t *ma;
133         int i;
134         ma = calloc(1, sizeof(bcf_p1aux_t));
135         ma->n1 = -1;
136         ma->n = n; ma->M = 2 * n;
137         if (ploidy) {
138                 ma->ploidy = malloc(n);
139                 memcpy(ma->ploidy, ploidy, n);
140                 for (i = 0, ma->M = 0; i < n; ++i) ma->M += ploidy[i];
141                 if (ma->M == 2 * n) {
142                         free(ma->ploidy);
143                         ma->ploidy = 0;
144                 }
145         }
146         ma->q2p = calloc(256, sizeof(double));
147         ma->pdg = calloc(3 * ma->n, sizeof(double));
148         ma->phi = calloc(ma->M + 1, sizeof(double));
149         ma->phi_indel = calloc(ma->M + 1, sizeof(double));
150         ma->phi1 = calloc(ma->M + 1, sizeof(double));
151         ma->phi2 = calloc(ma->M + 1, sizeof(double));
152         ma->z = calloc(ma->M + 1, sizeof(double));
153         ma->zswap = calloc(ma->M + 1, sizeof(double));
154         ma->z1 = calloc(ma->M + 1, sizeof(double)); // actually we do not need this large
155         ma->z2 = calloc(ma->M + 1, sizeof(double));
156         ma->afs = calloc(ma->M + 1, sizeof(double));
157         ma->afs1 = calloc(ma->M + 1, sizeof(double));
158         ma->lf = calloc(ma->M + 1, sizeof(double));
159         for (i = 0; i < 256; ++i)
160                 ma->q2p[i] = pow(10., -i / 10.);
161         for (i = 0; i <= ma->M; ++i) ma->lf[i] = lgamma(i + 1);
162         bcf_p1_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior
163         return ma;
164 }
165
166 int bcf_p1_set_n1(bcf_p1aux_t *b, int n1)
167 {
168         if (n1 == 0 || n1 >= b->n) return -1;
169         if (b->M != b->n * 2) {
170                 fprintf(stderr, "[%s] unable to set `n1' when there are haploid samples.\n", __func__);
171                 return -1;
172         }
173         b->n1 = n1;
174         return 0;
175 }
176
177 void bcf_p1_destroy(bcf_p1aux_t *ma)
178 {
179         if (ma) {
180                 int k;
181                 free(ma->lf);
182                 if (ma->hg && ma->n1 > 0) {
183                         for (k = 0; k <= 2*ma->n1; ++k) free(ma->hg[k]);
184                         free(ma->hg);
185                 }
186                 free(ma->ploidy); free(ma->q2p); free(ma->pdg);
187                 free(ma->phi); free(ma->phi_indel); free(ma->phi1); free(ma->phi2);
188                 free(ma->z); free(ma->zswap); free(ma->z1); free(ma->z2);
189                 free(ma->afs); free(ma->afs1);
190                 free(ma);
191         }
192 }
193
194 static int cal_pdg(const bcf1_t *b, bcf_p1aux_t *ma)
195 {
196         int i, j;
197     int n = (b->n_alleles+1)*b->n_alleles/2;
198     double *lk = alloca(n * sizeof(long));
199     memset(lk, 0, sizeof(double) * n);
200         for (j = 0; j < ma->n; ++j) {
201                 const uint8_t *pi = ma->PL + j * ma->PL_len;
202                 double *pdg = ma->pdg + j * 3;
203                 pdg[0] = ma->q2p[pi[2]]; pdg[1] = ma->q2p[pi[1]]; pdg[2] = ma->q2p[pi[0]];
204         for (i=0; i<n; i++) lk[i] += pi[i];
205     }
206
207     double norm=lk[0]; 
208     for (i=1; i<n; i++) if (lk[i]<norm) norm=lk[i];
209     #if DBG
210     for (i=0,j=0; i<b->n_alleles; i++)
211     {
212         int k; for (k=0; k<=i; k++) printf("%.0f\t", lk[j++]);
213         printf("\n");
214     }
215     #endif
216     for (i=0; i<n; i++) lk[i] = pow(10,-0.1*(lk[i]-norm));
217
218     // Find out the most likely alleles. In contrast to the original version,
219     //  ALT alleles may not be printed when they are more likely than REF but
220     //  significantly less likely than the most likely ALT. The only criterion
221     //  is the LK ratio now. To obtain behaviour similar to the original one,
222     //  use the pref variant below.
223     double pmax=0; //,pref=0;
224     n = ma->is_indel ? b->n_alleles : b->n_alleles-1;
225     for (i=0; i<n; i++)
226     {
227         double pr=0;
228         int k=i*(i+1)/2;
229         for (j=0; j<=i; j++) { pr+=lk[k]; k++; }
230         for (j=i+1; j<b->n_alleles; j++) { k=j*(j+1)/2+i; pr+=lk[k]; }
231         #if DBG
232         printf("%d\t%e\n", i,pr);
233         #endif
234         if (pmax<pr) pmax=pr;
235         // if (i==0) pref=pr;
236         // if (pr<pref && pr/pmax < 1e-4) break;
237         if (pr/pmax < 1e-4) break;   // Assuming the alleles are sorted by the lk
238     }
239         return i-1;
240 }
241
242 int bcf_p1_call_gt(const bcf_p1aux_t *ma, double f0, int k)
243 {
244         double sum, g[3];
245         double max, f3[3], *pdg = ma->pdg + k * 3;
246         int q, i, max_i, ploidy;
247         ploidy = ma->ploidy? ma->ploidy[k] : 2;
248         if (ploidy == 2) {
249                 f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
250         } else {
251                 f3[0] = 1. - f0; f3[1] = 0; f3[2] = f0;
252         }
253         for (i = 0, sum = 0.; i < 3; ++i)
254                 sum += (g[i] = pdg[i] * f3[i]);
255         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
256                 g[i] /= sum;
257                 if (g[i] > max) max = g[i], max_i = i;
258         }
259         max = 1. - max;
260         if (max < 1e-308) max = 1e-308;
261         q = (int)(-4.343 * log(max) + .499);
262         if (q > 99) q = 99;
263         return q<<2|max_i;
264 }
265
266 #define TINY 1e-20
267
268 static void mc_cal_y_core(bcf_p1aux_t *ma, int beg)
269 {
270         double *z[2], *tmp, *pdg;
271         int _j, last_min, last_max;
272         assert(beg == 0 || ma->M == ma->n*2);
273         z[0] = ma->z;
274         z[1] = ma->zswap;
275         pdg = ma->pdg;
276         memset(z[0], 0, sizeof(double) * (ma->M + 1));
277         memset(z[1], 0, sizeof(double) * (ma->M + 1));
278         z[0][0] = 1.;
279         last_min = last_max = 0;
280         ma->t = 0.;
281         if (ma->M == ma->n * 2) {
282                 int M = 0;
283                 for (_j = beg; _j < ma->n; ++_j) {
284                         int k, j = _j - beg, _min = last_min, _max = last_max, M0;
285                         double p[3], sum;
286                         M0 = M; M += 2;
287                         pdg = ma->pdg + _j * 3;
288                         p[0] = pdg[0]; p[1] = 2. * pdg[1]; p[2] = pdg[2];
289                         for (; _min < _max && z[0][_min] < TINY; ++_min) z[0][_min] = z[1][_min] = 0.;
290                         for (; _max > _min && z[0][_max] < TINY; --_max) z[0][_max] = z[1][_max] = 0.;
291                         _max += 2;
292                         if (_min == 0) k = 0, z[1][k] = (M0-k+1) * (M0-k+2) * p[0] * z[0][k];
293                         if (_min <= 1) k = 1, z[1][k] = (M0-k+1) * (M0-k+2) * p[0] * z[0][k] + k*(M0-k+2) * p[1] * z[0][k-1];
294                         for (k = _min < 2? 2 : _min; k <= _max; ++k)
295                                 z[1][k] = (M0-k+1)*(M0-k+2) * p[0] * z[0][k] + k*(M0-k+2) * p[1] * z[0][k-1] + k*(k-1)* p[2] * z[0][k-2];
296                         for (k = _min, sum = 0.; k <= _max; ++k) sum += z[1][k];
297                         ma->t += log(sum / (M * (M - 1.)));
298                         for (k = _min; k <= _max; ++k) z[1][k] /= sum;
299                         if (_min >= 1) z[1][_min-1] = 0.;
300                         if (_min >= 2) z[1][_min-2] = 0.;
301                         if (j < ma->n - 1) z[1][_max+1] = z[1][_max+2] = 0.;
302                         if (_j == ma->n1 - 1) { // set pop1; ma->n1==-1 when unset
303                                 ma->t1 = ma->t;
304                                 memcpy(ma->z1, z[1], sizeof(double) * (ma->n1 * 2 + 1));
305                         }
306                         tmp = z[0]; z[0] = z[1]; z[1] = tmp;
307                         last_min = _min; last_max = _max;
308                 }
309                 //for (_j = 0; _j < last_min; ++_j) z[0][_j] = 0.; // TODO: are these necessary?
310                 //for (_j = last_max + 1; _j < ma->M; ++_j) z[0][_j] = 0.;
311         } else { // this block is very similar to the block above; these two might be merged in future
312                 int j, M = 0;
313                 for (j = 0; j < ma->n; ++j) {
314                         int k, M0, _min = last_min, _max = last_max;
315                         double p[3], sum;
316                         pdg = ma->pdg + j * 3;
317                         for (; _min < _max && z[0][_min] < TINY; ++_min) z[0][_min] = z[1][_min] = 0.;
318                         for (; _max > _min && z[0][_max] < TINY; --_max) z[0][_max] = z[1][_max] = 0.;
319                         M0 = M;
320                         M += ma->ploidy[j];
321                         if (ma->ploidy[j] == 1) {
322                                 p[0] = pdg[0]; p[1] = pdg[2];
323                                 _max++;
324                                 if (_min == 0) k = 0, z[1][k] = (M0+1-k) * p[0] * z[0][k];
325                                 for (k = _min < 1? 1 : _min; k <= _max; ++k)
326                                         z[1][k] = (M0+1-k) * p[0] * z[0][k] + k * p[1] * z[0][k-1];
327                                 for (k = _min, sum = 0.; k <= _max; ++k) sum += z[1][k];
328                                 ma->t += log(sum / M);
329                                 for (k = _min; k <= _max; ++k) z[1][k] /= sum;
330                                 if (_min >= 1) z[1][_min-1] = 0.;
331                                 if (j < ma->n - 1) z[1][_max+1] = 0.;
332                         } else if (ma->ploidy[j] == 2) {
333                                 p[0] = pdg[0]; p[1] = 2 * pdg[1]; p[2] = pdg[2];
334                                 _max += 2;
335                                 if (_min == 0) k = 0, z[1][k] = (M0-k+1) * (M0-k+2) * p[0] * z[0][k];
336                                 if (_min <= 1) k = 1, z[1][k] = (M0-k+1) * (M0-k+2) * p[0] * z[0][k] + k*(M0-k+2) * p[1] * z[0][k-1];
337                                 for (k = _min < 2? 2 : _min; k <= _max; ++k)
338                                         z[1][k] = (M0-k+1)*(M0-k+2) * p[0] * z[0][k] + k*(M0-k+2) * p[1] * z[0][k-1] + k*(k-1)* p[2] * z[0][k-2];
339                                 for (k = _min, sum = 0.; k <= _max; ++k) sum += z[1][k];
340                                 ma->t += log(sum / (M * (M - 1.)));
341                                 for (k = _min; k <= _max; ++k) z[1][k] /= sum;
342                                 if (_min >= 1) z[1][_min-1] = 0.;
343                                 if (_min >= 2) z[1][_min-2] = 0.;
344                                 if (j < ma->n - 1) z[1][_max+1] = z[1][_max+2] = 0.;
345                         }
346                         tmp = z[0]; z[0] = z[1]; z[1] = tmp;
347                         last_min = _min; last_max = _max;
348                 }
349         }
350         if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (ma->M + 1));
351 }
352
353 static void mc_cal_y(bcf_p1aux_t *ma)
354 {
355         if (ma->n1 > 0 && ma->n1 < ma->n && ma->M == ma->n * 2) { // NB: ma->n1 is ineffective when there are haploid samples
356                 int k;
357                 long double x;
358                 memset(ma->z1, 0, sizeof(double) * (2 * ma->n1 + 1));
359                 memset(ma->z2, 0, sizeof(double) * (2 * (ma->n - ma->n1) + 1));
360                 ma->t1 = ma->t2 = 0.;
361                 mc_cal_y_core(ma, ma->n1);
362                 ma->t2 = ma->t;
363                 memcpy(ma->z2, ma->z, sizeof(double) * (2 * (ma->n - ma->n1) + 1));
364                 mc_cal_y_core(ma, 0);
365                 // rescale z
366                 x = expl(ma->t - (ma->t1 + ma->t2));
367                 for (k = 0; k <= ma->M; ++k) ma->z[k] *= x;
368         } else mc_cal_y_core(ma, 0);
369 }
370
371 #define CONTRAST_TINY 1e-30
372
373 extern double kf_gammaq(double s, double z); // incomplete gamma function for chi^2 test
374
375 static inline double chi2_test(int a, int b, int c, int d)
376 {
377         double x, z;
378         x = (double)(a+b) * (c+d) * (b+d) * (a+c);
379         if (x == 0.) return 1;
380         z = a * d - b * c;
381         return kf_gammaq(.5, .5 * z * z * (a+b+c+d) / x);
382 }
383
384 // chi2=(a+b+c+d)(ad-bc)^2/[(a+b)(c+d)(a+c)(b+d)]
385 static inline double contrast2_aux(const bcf_p1aux_t *p1, double sum, int k1, int k2, double x[3])
386 {
387         double p = p1->phi[k1+k2] * p1->z1[k1] * p1->z2[k2] / sum * p1->hg[k1][k2];
388         int n1 = p1->n1, n2 = p1->n - p1->n1;
389         if (p < CONTRAST_TINY) return -1;
390         if (.5*k1/n1 < .5*k2/n2) x[1] += p;
391         else if (.5*k1/n1 > .5*k2/n2) x[2] += p;
392         else x[0] += p;
393         return p * chi2_test(k1, k2, (n1<<1) - k1, (n2<<1) - k2);
394 }
395
396 static double contrast2(bcf_p1aux_t *p1, double ret[3])
397 {
398         int k, k1, k2, k10, k20, n1, n2;
399         double sum;
400         // get n1 and n2
401         n1 = p1->n1; n2 = p1->n - p1->n1;
402         if (n1 <= 0 || n2 <= 0) return 0.;
403         if (p1->hg == 0) { // initialize the hypergeometric distribution
404                 /* NB: the hg matrix may take a lot of memory when there are many samples. There is a way
405                    to avoid precomputing this matrix, but it is slower and quite intricate. The following
406                    computation in this block can be accelerated with a similar strategy, but perhaps this
407                    is not a serious concern for now. */
408                 double tmp = lgamma(2*(n1+n2)+1) - (lgamma(2*n1+1) + lgamma(2*n2+1));
409                 p1->hg = calloc(2*n1+1, sizeof(void*));
410                 for (k1 = 0; k1 <= 2*n1; ++k1) {
411                         p1->hg[k1] = calloc(2*n2+1, sizeof(double));
412                         for (k2 = 0; k2 <= 2*n2; ++k2)
413                                 p1->hg[k1][k2] = exp(lgamma(k1+k2+1) + lgamma(p1->M-k1-k2+1) - (lgamma(k1+1) + lgamma(k2+1) + lgamma(2*n1-k1+1) + lgamma(2*n2-k2+1) + tmp));
414                 }
415         }
416         { // compute
417                 long double suml = 0;
418                 for (k = 0; k <= p1->M; ++k) suml += p1->phi[k] * p1->z[k];
419                 sum = suml;
420         }
421         { // get the max k1 and k2
422                 double max;
423                 int max_k;
424                 for (k = 0, max = 0, max_k = -1; k <= 2*n1; ++k) {
425                         double x = p1->phi1[k] * p1->z1[k];
426                         if (x > max) max = x, max_k = k;
427                 }
428                 k10 = max_k;
429                 for (k = 0, max = 0, max_k = -1; k <= 2*n2; ++k) {
430                         double x = p1->phi2[k] * p1->z2[k];
431                         if (x > max) max = x, max_k = k;
432                 }
433                 k20 = max_k;
434         }
435         { // We can do the following with one nested loop, but that is an O(N^2) thing. The following code block is much faster for large N.
436                 double x[3], y;
437                 long double z = 0., L[2];
438                 x[0] = x[1] = x[2] = 0; L[0] = L[1] = 0;
439                 for (k1 = k10; k1 >= 0; --k1) {
440                         for (k2 = k20; k2 >= 0; --k2) {
441                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
442                                 else z += y;
443                         }
444                         for (k2 = k20 + 1; k2 <= 2*n2; ++k2) {
445                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
446                                 else z += y;
447                         }
448                 }
449                 ret[0] = x[0]; ret[1] = x[1]; ret[2] = x[2];
450                 x[0] = x[1] = x[2] = 0;
451                 for (k1 = k10 + 1; k1 <= 2*n1; ++k1) {
452                         for (k2 = k20; k2 >= 0; --k2) {
453                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
454                                 else z += y;
455                         }
456                         for (k2 = k20 + 1; k2 <= 2*n2; ++k2) {
457                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
458                                 else z += y;
459                         }
460                 }
461                 ret[0] += x[0]; ret[1] += x[1]; ret[2] += x[2];
462                 if (ret[0] + ret[1] + ret[2] < 0.95) { // in case of bad things happened
463                         ret[0] = ret[1] = ret[2] = 0; L[0] = L[1] = 0;
464                         for (k1 = 0, z = 0.; k1 <= 2*n1; ++k1)
465                                 for (k2 = 0; k2 <= 2*n2; ++k2)
466                                         if ((y = contrast2_aux(p1, sum, k1, k2, ret)) >= 0) z += y;
467                         if (ret[0] + ret[1] + ret[2] < 0.95) // It seems that this may be caused by floating point errors. I do not really understand why...
468                                 z = 1.0, ret[0] = ret[1] = ret[2] = 1./3;
469                 }
470                 return (double)z;
471         }
472 }
473
474 static double mc_cal_afs(bcf_p1aux_t *ma, double *p_ref_folded, double *p_var_folded)
475 {
476         int k;
477         long double sum = 0., sum2;
478         double *phi = ma->is_indel? ma->phi_indel : ma->phi;
479         memset(ma->afs1, 0, sizeof(double) * (ma->M + 1));
480         mc_cal_y(ma);
481         // compute AFS
482         for (k = 0, sum = 0.; k <= ma->M; ++k)
483                 sum += (long double)phi[k] * ma->z[k];
484         for (k = 0; k <= ma->M; ++k) {
485                 ma->afs1[k] = phi[k] * ma->z[k] / sum;
486                 if (isnan(ma->afs1[k]) || isinf(ma->afs1[k])) return -1.;
487         }
488         // compute folded variant probability
489         for (k = 0, sum = 0.; k <= ma->M; ++k)
490                 sum += (long double)(phi[k] + phi[ma->M - k]) / 2. * ma->z[k];
491         for (k = 1, sum2 = 0.; k < ma->M; ++k)
492                 sum2 += (long double)(phi[k] + phi[ma->M - k]) / 2. * ma->z[k];
493         *p_var_folded = sum2 / sum;
494         *p_ref_folded = (phi[k] + phi[ma->M - k]) / 2. * (ma->z[ma->M] + ma->z[0]) / sum;
495         // the expected frequency
496         for (k = 0, sum = 0.; k <= ma->M; ++k) {
497                 ma->afs[k] += ma->afs1[k];
498                 sum += k * ma->afs1[k];
499         }
500         return sum / ma->M;
501 }
502
503 int bcf_p1_cal(const bcf1_t *b, int do_contrast, bcf_p1aux_t *ma, bcf_p1rst_t *rst)
504 {
505         int i, k;
506         long double sum = 0.;
507         ma->is_indel = bcf_is_indel(b);
508         rst->perm_rank = -1;
509         // set PL and PL_len
510         for (i = 0; i < b->n_gi; ++i) {
511                 if (b->gi[i].fmt == bcf_str2int("PL", 2)) {
512                         ma->PL = (uint8_t*)b->gi[i].data;
513                         ma->PL_len = b->gi[i].len;
514                         break;
515                 }
516         }
517         if (i == b->n_gi) return -1; // no PL
518         if (b->n_alleles < 2) return -1; // FIXME: find a better solution
519         // 
520         rst->rank0 = cal_pdg(b, ma);
521         rst->f_exp = mc_cal_afs(ma, &rst->p_ref_folded, &rst->p_var_folded);
522         rst->p_ref = ma->afs1[ma->M];
523         for (k = 0, sum = 0.; k < ma->M; ++k)
524                 sum += ma->afs1[k];
525         rst->p_var = (double)sum;
526         { // compute the allele count
527                 double max = -1;
528                 rst->ac = -1;
529                 for (k = 0; k <= ma->M; ++k)
530                         if (max < ma->z[k]) max = ma->z[k], rst->ac = k;
531                 rst->ac = ma->M - rst->ac;
532         }
533         // calculate f_flat and f_em
534         for (k = 0, sum = 0.; k <= ma->M; ++k)
535                 sum += (long double)ma->z[k];
536         rst->f_flat = 0.;
537         for (k = 0; k <= ma->M; ++k) {
538                 double p = ma->z[k] / sum;
539                 rst->f_flat += k * p;
540         }
541         rst->f_flat /= ma->M;
542         { // estimate equal-tail credible interval (95% level)
543                 int l, h;
544                 double p;
545                 for (i = 0, p = 0.; i <= ma->M; ++i)
546                         if (p + ma->afs1[i] > 0.025) break;
547                         else p += ma->afs1[i];
548                 l = i;
549                 for (i = ma->M, p = 0.; i >= 0; --i)
550                         if (p + ma->afs1[i] > 0.025) break;
551                         else p += ma->afs1[i];
552                 h = i;
553                 rst->cil = (double)(ma->M - h) / ma->M; rst->cih = (double)(ma->M - l) / ma->M;
554         }
555         if (ma->n1 > 0) { // compute LRT
556                 double max0, max1, max2;
557                 for (k = 0, max0 = -1; k <= ma->M; ++k)
558                         if (max0 < ma->z[k]) max0 = ma->z[k];
559                 for (k = 0, max1 = -1; k <= ma->n1 * 2; ++k)
560                         if (max1 < ma->z1[k]) max1 = ma->z1[k];
561                 for (k = 0, max2 = -1; k <= ma->M - ma->n1 * 2; ++k)
562                         if (max2 < ma->z2[k]) max2 = ma->z2[k];
563                 rst->lrt = log(max1 * max2 / max0);
564                 rst->lrt = rst->lrt < 0? 1 : kf_gammaq(.5, rst->lrt);
565         } else rst->lrt = -1.0;
566         rst->cmp[0] = rst->cmp[1] = rst->cmp[2] = rst->p_chi2 = -1.0;
567         if (do_contrast && rst->p_var > 0.5) // skip contrast2() if the locus is a strong non-variant
568                 rst->p_chi2 = contrast2(ma, rst->cmp);
569         return 0;
570 }
571
572 void bcf_p1_dump_afs(bcf_p1aux_t *ma)
573 {
574         int k;
575         fprintf(stderr, "[afs]");
576         for (k = 0; k <= ma->M; ++k)
577                 fprintf(stderr, " %d:%.3lf", k, ma->afs[ma->M - k]);
578         fprintf(stderr, "\n");
579         memset(ma->afs, 0, sizeof(double) * (ma->M + 1));
580 }