]> git.donarmstrong.com Git - samtools.git/blob - bcftools/prob1.c
Bug fixes
[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 <limits.h>
8 #include "prob1.h"
9 #include "kstring.h"
10
11 #include "kseq.h"
12 KSTREAM_INIT(gzFile, gzread, 16384)
13
14 #define MC_MAX_EM_ITER 16
15 #define MC_EM_EPS 1e-5
16 #define MC_DEF_INDEL 0.15
17
18 unsigned char seq_nt4_table[256] = {
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, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4 /*'-'*/, 4, 4,
22         4, 4, 4, 4,  4, 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, 0, 4, 1,  4, 4, 4, 2,  4, 4, 4, 4,  4, 4, 4, 4, 
26         4, 4, 4, 4,  3, 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         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
34         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4
35 };
36
37 struct __bcf_p1aux_t {
38         int n, M, n1, is_indel;
39         uint8_t *ploidy; // haploid or diploid ONLY
40         double *q2p, *pdg; // pdg -> P(D|g)
41         double *phi, *phi_indel;
42         double *z, *zswap; // aux for afs
43         double *z1, *z2, *phi1, *phi2; // only calculated when n1 is set
44         double **hg; // hypergeometric distribution
45         double *lf; // log factorial
46         double t, t1, t2;
47         double *afs, *afs1; // afs: accumulative AFS; afs1: site posterior distribution
48         const uint8_t *PL; // point to PL
49         int PL_len;
50 };
51
52 void bcf_p1_indel_prior(bcf_p1aux_t *ma, double x)
53 {
54         int i;
55         for (i = 0; i < ma->M; ++i)
56                 ma->phi_indel[i] = ma->phi[i] * x;
57         ma->phi_indel[ma->M] = 1. - ma->phi[ma->M] * x;
58 }
59
60 static void init_prior(int type, double theta, int M, double *phi)
61 {
62         int i;
63         if (type == MC_PTYPE_COND2) {
64                 for (i = 0; i <= M; ++i)
65                         phi[i] = 2. * (i + 1) / (M + 1) / (M + 2);
66         } else if (type == MC_PTYPE_FLAT) {
67                 for (i = 0; i <= M; ++i)
68                         phi[i] = 1. / (M + 1);
69         } else {
70                 double sum;
71                 for (i = 0, sum = 0.; i < M; ++i)
72                         sum += (phi[i] = theta / (M - i));
73                 phi[M] = 1. - sum;
74         }
75 }
76
77 void bcf_p1_init_prior(bcf_p1aux_t *ma, int type, double theta)
78 {
79         init_prior(type, theta, ma->M, ma->phi);
80         bcf_p1_indel_prior(ma, MC_DEF_INDEL);
81 }
82
83 void bcf_p1_init_subprior(bcf_p1aux_t *ma, int type, double theta)
84 {
85         if (ma->n1 <= 0 || ma->n1 >= ma->M) return;
86         init_prior(type, theta, 2*ma->n1, ma->phi1);
87         init_prior(type, theta, 2*(ma->n - ma->n1), ma->phi2);
88 }
89
90 int bcf_p1_read_prior(bcf_p1aux_t *ma, const char *fn)
91 {
92         gzFile fp;
93         kstring_t s;
94         kstream_t *ks;
95         long double sum;
96         int dret, k;
97         memset(&s, 0, sizeof(kstring_t));
98         fp = strcmp(fn, "-")? gzopen(fn, "r") : gzdopen(fileno(stdin), "r");
99         ks = ks_init(fp);
100         memset(ma->phi, 0, sizeof(double) * (ma->M + 1));
101         while (ks_getuntil(ks, '\n', &s, &dret) >= 0) {
102                 if (strstr(s.s, "[afs] ") == s.s) {
103                         char *p = s.s + 6;
104                         for (k = 0; k <= ma->M; ++k) {
105                                 int x;
106                                 double y;
107                                 x = strtol(p, &p, 10);
108                                 if (x != k && (errno == EINVAL || errno == ERANGE)) return -1;
109                                 ++p;
110                                 y = strtod(p, &p);
111                                 if (y == 0. && (errno == EINVAL || errno == ERANGE)) return -1;
112                                 ma->phi[ma->M - k] += y;
113                         }
114                 }
115         }
116         ks_destroy(ks);
117         gzclose(fp);
118         free(s.s);
119         for (sum = 0., k = 0; k <= ma->M; ++k) sum += ma->phi[k];
120         fprintf(stderr, "[prior]");
121         for (k = 0; k <= ma->M; ++k) ma->phi[k] /= sum;
122         for (k = 0; k <= ma->M; ++k) fprintf(stderr, " %d:%.3lg", k, ma->phi[ma->M - k]);
123         fputc('\n', stderr);
124         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));
125         fprintf(stderr, "[%s] heterozygosity=%lf, ", __func__, (double)sum);
126         for (sum = 0., k = 1; k <= ma->M; ++k) sum += k * ma->phi[ma->M - k] / ma->M;
127         fprintf(stderr, "theta=%lf\n", (double)sum);
128         bcf_p1_indel_prior(ma, MC_DEF_INDEL);
129         return 0;
130 }
131
132 bcf_p1aux_t *bcf_p1_init(int n, uint8_t *ploidy)
133 {
134         bcf_p1aux_t *ma;
135         int i;
136         ma = calloc(1, sizeof(bcf_p1aux_t));
137         ma->n1 = -1;
138         ma->n = n; ma->M = 2 * n;
139         if (ploidy) {
140                 ma->ploidy = malloc(n);
141                 memcpy(ma->ploidy, ploidy, n);
142                 for (i = 0, ma->M = 0; i < n; ++i) ma->M += ploidy[i];
143                 if (ma->M == 2 * n) {
144                         free(ma->ploidy);
145                         ma->ploidy = 0;
146                 }
147         }
148         ma->q2p = calloc(256, sizeof(double));
149         ma->pdg = calloc(3 * ma->n, sizeof(double));
150         ma->phi = calloc(ma->M + 1, sizeof(double));
151         ma->phi_indel = calloc(ma->M + 1, sizeof(double));
152         ma->phi1 = calloc(ma->M + 1, sizeof(double));
153         ma->phi2 = calloc(ma->M + 1, sizeof(double));
154         ma->z = calloc(ma->M + 1, sizeof(double));
155         ma->zswap = calloc(ma->M + 1, sizeof(double));
156         ma->z1 = calloc(ma->M + 1, sizeof(double)); // actually we do not need this large
157         ma->z2 = calloc(ma->M + 1, sizeof(double));
158         ma->afs = calloc(ma->M + 1, sizeof(double));
159         ma->afs1 = calloc(ma->M + 1, sizeof(double));
160         ma->lf = calloc(ma->M + 1, sizeof(double));
161         for (i = 0; i < 256; ++i)
162                 ma->q2p[i] = pow(10., -i / 10.);
163         for (i = 0; i <= ma->M; ++i) ma->lf[i] = lgamma(i + 1);
164         bcf_p1_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior
165         return ma;
166 }
167
168 int bcf_p1_set_n1(bcf_p1aux_t *b, int n1)
169 {
170         if (n1 == 0 || n1 >= b->n) return -1;
171         if (b->M != b->n * 2) {
172                 fprintf(stderr, "[%s] unable to set `n1' when there are haploid samples.\n", __func__);
173                 return -1;
174         }
175         b->n1 = n1;
176         return 0;
177 }
178
179 void bcf_p1_set_ploidy(bcf1_t *b, bcf_p1aux_t *ma)
180 {
181     // bcf_p1aux_t fields are not visible outside of prob1.c, hence this wrapper.
182     // Ideally, this should set ploidy per site to allow pseudo-autosomal regions
183     b->ploidy = ma->ploidy;
184 }
185
186 void bcf_p1_destroy(bcf_p1aux_t *ma)
187 {
188         if (ma) {
189                 int k;
190                 free(ma->lf);
191                 if (ma->hg && ma->n1 > 0) {
192                         for (k = 0; k <= 2*ma->n1; ++k) free(ma->hg[k]);
193                         free(ma->hg);
194                 }
195                 free(ma->ploidy); free(ma->q2p); free(ma->pdg);
196                 free(ma->phi); free(ma->phi_indel); free(ma->phi1); free(ma->phi2);
197                 free(ma->z); free(ma->zswap); free(ma->z1); free(ma->z2);
198                 free(ma->afs); free(ma->afs1);
199                 free(ma);
200         }
201 }
202
203 extern double kf_gammap(double s, double z);
204
205 int call_multiallelic_gt(bcf1_t *b, bcf_p1aux_t *ma, double threshold)
206 {
207     int nals = 1;
208     char *p;
209     for (p=b->alt; *p; p++)
210     {
211         if ( *p=='X' || p[0]=='.' ) break;
212         if ( p[0]==',' ) nals++;
213     }
214     if ( b->alt[0] && !*p ) nals++;
215
216     if ( nals==1 ) return 1;
217
218     if ( nals>4 )
219     {
220         if ( *b->ref=='N' ) return 0;
221         fprintf(stderr,"Not ready for this, more than 4 alleles at %d: %s, %s\n", b->pos+1, b->ref,b->alt); 
222         exit(1);
223     }
224
225     // set PL and PL_len
226     uint8_t *pl = NULL;
227     int npl = 0, idp=-1;
228     int i;
229     for (i = 0; i < b->n_gi; ++i) 
230     {
231         if (b->gi[i].fmt == bcf_str2int("PL", 2)) 
232         {
233             pl  = (uint8_t*)b->gi[i].data;
234             npl = b->gi[i].len;
235         }
236         if (b->gi[i].fmt == bcf_str2int("DP", 2))  idp=i;
237     }
238     if ( !pl ) return -1;
239
240     int npdg = nals*(nals+1)/2;
241     float *pdg,*_pdg;
242     _pdg = pdg = malloc(sizeof(float)*ma->n*npdg);
243     for (i=0; i<ma->n; i++)
244     {
245         int j; 
246         float sum = 0;
247         for (j=0; j<npdg; j++)
248         {
249             _pdg[j] = pow(10,-0.1*pl[j]);
250             sum += _pdg[j];
251         }
252         if ( sum )
253             for (j=0; j<npdg; j++) _pdg[j] /= sum;
254         _pdg += npdg;
255         pl += npl;
256     }
257
258     if ((p = strstr(b->info, "QS=")) == 0) { fprintf(stderr,"INFO/QS is required with -m, exiting\n"); exit(1); }
259     float qsum[4];
260     if ( sscanf(p+3,"%f,%f,%f,%f",&qsum[0],&qsum[1],&qsum[2],&qsum[3])!=4 ) { fprintf(stderr,"Could not parse %s\n",p); exit(1); }
261
262     
263     int ia,ib,ic, max_als=0, max_als2=0;
264     float max_lk = INT_MIN, max_lk2 = INT_MIN, lk_sum = INT_MIN;
265     for (ia=0; ia<nals; ia++)
266     {
267         //if ( ia && qsum[ia]==0 ) continue;
268         float lk_tot = 0;
269         int iaa = (ia+1)*(ia+2)/2-1;
270         int isample;
271         for (isample=0; isample<ma->n; isample++)
272         {
273             float *p = pdg + isample*npdg;
274             assert( log(p[iaa]) <= 0 );
275             lk_tot += log(p[iaa]);
276         }
277         if ( max_lk<lk_tot ) { max_lk2 = max_lk; max_als2 = max_als; max_lk = lk_tot; max_als = 1<<ia; }
278         else if ( max_lk2<lk_tot ) { max_lk2 = lk_tot; max_als2 = 1<<ia; }
279         lk_sum = lk_tot>lk_sum ? lk_tot + log(1+exp(lk_sum-lk_tot)) : lk_sum + log(1+exp(lk_tot-lk_sum));
280     }
281     for (ia=0; ia<nals; ia++)
282     {
283         if ( qsum[ia]==0 ) continue;
284         //if ( ia && qsum[ia]==0 ) continue;
285         for (ib=0; ib<ia; ib++)
286         {
287             if ( qsum[ib]==0 ) continue;
288             //if ( ib && qsum[ib]==0 ) continue;
289             //if ( qsum[ia]+qsum[ib]==0 ) continue;
290             float lk_tot = 0;
291             float fa  = qsum[ia]/(qsum[ia]+qsum[ib]);
292             float fb  = qsum[ib]/(qsum[ia]+qsum[ib]);
293             float fab = 2*fa*fb; fa *= fa; fb *= fb;
294             int isample, iaa = (ia+1)*(ia+2)/2-1, ibb = (ib+1)*(ib+2)/2-1, iab = iaa - ia + ib;
295             for (isample=0; isample<ma->n; isample++)
296             {
297                 if ( b->ploidy && b->ploidy[isample]==1 ) continue;
298                 float *p = pdg + isample*npdg;
299                 assert( log(fa*p[iaa] + fb*p[ibb] + fab*p[iab]) <= 0 );
300                 lk_tot += log(fa*p[iaa] + fb*p[ibb] + fab*p[iab]);
301             }
302             if ( max_lk<lk_tot ) { max_lk2 = max_lk; max_als2 = max_als; max_lk = lk_tot; max_als = 1<<ia|1<<ib; }
303             else if ( max_lk2<lk_tot ) { max_lk2 = lk_tot; max_als2 = 1<<ia|1<<ib; }
304             lk_sum = lk_tot>lk_sum ? lk_tot + log(1+exp(lk_sum-lk_tot)) : lk_sum + log(1+exp(lk_tot-lk_sum));
305         }
306     }
307     for (ia=0; ia<nals; ia++)
308     {
309         if ( qsum[ia]==0 ) continue;
310         //if ( ia && qsum[ia]==0 ) continue;
311         for (ib=0; ib<ia; ib++)
312         {
313             if ( qsum[ib]==0 ) continue;
314             //if ( ib && qsum[ib]==0 ) continue;
315             for (ic=0; ic<ib; ic++)
316             {
317                 if ( qsum[ic]==0 ) continue;
318                 //if ( ic && qsum[ic]==0 ) continue;
319                 //if ( qsum[ia]+qsum[ib]+qsum[ic]==0 ) continue;
320                 float lk_tot = 0;
321                 float fa  = qsum[ia]/(qsum[ia]+qsum[ib]+qsum[ic]);
322                 float fb  = qsum[ib]/(qsum[ia]+qsum[ib]+qsum[ic]);
323                 float fc  = qsum[ic]/(qsum[ia]+qsum[ib]+qsum[ic]);
324                 float fab = 2*fa*fb, fac = 2*fa*fc, fbc = 2*fb*fc; fa *= fa; fb *= fb; fc *= fc;
325                 int isample, iaa = (ia+1)*(ia+2)/2-1, ibb = (ib+1)*(ib+2)/2-1, icc = (ic+1)*(ic+2)/2-1;
326                 int iab = iaa - ia + ib, iac = iaa - ia + ic, ibc = ibb - ib + ic;
327                 for (isample=0; isample<ma->n; isample++)
328                 {
329                     if ( b->ploidy && b->ploidy[isample]==1 ) continue;
330                     float *p = pdg + isample*npdg;
331                     assert( log(fa*p[iaa] + fb*p[ibb] + fc*p[icc] + fab*p[iab] + fac*p[iac] + fbc*p[ibc]) <= 0 );
332                     lk_tot += log(fa*p[iaa] + fb*p[ibb] + fc*p[icc] + fab*p[iab] + fac*p[iac] + fbc*p[ibc]);
333                 }
334                 if ( max_lk<lk_tot ) { max_lk2 = max_lk; max_als2 = max_als; max_lk = lk_tot; max_als = 1<<ia|1<<ib|1<<ic; }
335                 else if ( max_lk2<lk_tot ) { max_lk2 = lk_tot; max_als2 = 1<<ia|1<<ib|1<<ic; }
336                 lk_sum = lk_tot>lk_sum ? lk_tot + log(1+exp(lk_sum-lk_tot)) : lk_sum + log(1+exp(lk_tot-lk_sum));
337         //printf("    %f\n", lk_sum);
338             }
339         }
340     }
341
342     int n1=0, n2=0;
343     for (i=0; i<nals; i++) if ( max_als&1<<i) n1++;
344     for (i=0; i<nals; i++) if ( max_als2&1<<i) n2++;
345     if ( n2<n1 && kf_gammap(1,2.0*(max_lk-max_lk2))<threshold )
346     {
347         max_lk = max_lk2;
348         max_als = max_als2;
349     }
350
351     // Get the BCF record ready for GT and GQ
352     kstring_t s;
353     int old_n_gi = b->n_gi;
354     s.m = b->m_str; s.l = b->l_str - 1; s.s = b->str;
355     kputs(":GT:GQ", &s); kputc('\0', &s);
356     b->m_str = s.m; b->l_str = s.l; b->str = s.s;
357     bcf_sync(b);
358
359     // Call GTs
360     int isample, gts=0, ac[4] = {0,0,0,0};
361     for (isample = 0; isample < b->n_smpl; isample++) 
362     {
363         int ploidy = b->ploidy ? b->ploidy[isample] : 2;
364         float *p = pdg + isample*npdg;
365         int ia, als = 0;
366         float lk = INT_MIN, lk_sum=0;
367         for (ia=0; ia<nals; ia++)
368         {
369             if ( !(max_als&1<<ia) ) continue;
370             int iaa = (ia+1)*(ia+2)/2-1;
371             float _lk = p[iaa]*qsum[ia]*qsum[ia];
372             if ( _lk > lk ) { lk = _lk; als = ia<<3 | ia; }
373             lk_sum += _lk;
374         }
375         if ( ploidy==2 ) 
376         {
377             for (ia=0; ia<nals; ia++)
378             {
379                 if ( !(max_als&1<<ia) ) continue;
380                 int iaa = (ia+1)*(ia+2)/2-1;
381                 for (ib=0; ib<ia; ib++)
382                 {
383                     if ( !(max_als&1<<ib) ) continue;
384                     int iab = iaa - ia + ib;
385                     float _lk = 2*qsum[ia]*qsum[ib]*p[iab];
386                     if ( _lk > lk ) { lk = _lk; als = ib<<3 | ia; }
387                     lk_sum += _lk;
388                 }
389             }
390         }
391         lk = -log(1-lk/lk_sum)/0.2302585;
392         if ( idp>=0 && ((uint16_t*)b->gi[idp].data)[isample]==0 ) 
393         {
394             ((uint8_t*)b->gi[old_n_gi].data)[isample]   = als | 1<<7;
395             ((uint8_t*)b->gi[old_n_gi+1].data)[isample] = 0;
396             continue;
397         }
398         ((uint8_t*)b->gi[old_n_gi].data)[isample]   = als;
399         ((uint8_t*)b->gi[old_n_gi+1].data)[isample] = lk<100 ? (int)lk : 99;
400
401         gts |= 1<<(als>>3&7) | 1<<(als&7);
402         ac[ als>>3&7 ]++;
403         ac[ als&7 ]++;
404     }
405     bcf_fit_alt(b,max_als);
406
407     // prepare ref, alt, filter, info, format
408     memset(&s, 0, sizeof(kstring_t)); kputc('\0', &s); 
409     kputs(b->ref, &s); kputc('\0', &s);
410     kputs(b->alt, &s); kputc('\0', &s); kputc('\0', &s);
411     {
412         kstring_t info; memset(&info, 0, sizeof(kstring_t));
413         int an=0, nalts=0;
414         for (i=0; i<nals; i++)
415         {
416             an += ac[i];
417             if ( i>0 && ac[i] ) nalts++;
418         }
419         ksprintf(&info, "AN=%d;", an);
420         if ( nalts )
421         {
422             kputs("AC=", &info);
423             for (i=1; i<nals; i++)
424             {
425                 if ( !(gts&1<<i) ) continue;
426                 nalts--;
427                 ksprintf(&info,"%d", ac[i]);
428                 if ( nalts>0 ) kputc(',', &info);
429             }
430             kputc(';', &info);
431         }
432         kputs(b->info, &info); 
433         info.l -= remove_tag(info.s, "I16=", ';');
434         info.l -= remove_tag(info.s, "QS=", ';');
435         kputc('\0', &info);
436         kputs(info.s, &s); kputc('\0', &s);
437         free(info.s);
438     }
439     kputs(b->fmt, &s); kputc('\0', &s);
440     free(b->str);
441     b->m_str = s.m; b->l_str = s.l; b->str = s.s;
442     b->qual = -4.343*(log(1-exp(max_lk-lk_sum)));
443     if ( b->qual>999 ) b->qual = 999;
444     bcf_sync(b);
445
446
447     free(pdg);
448     return gts;
449 }
450
451 static int cal_pdg(const bcf1_t *b, bcf_p1aux_t *ma)
452 {
453     int i, j;
454     long *p, tmp;
455     p = alloca(b->n_alleles * sizeof(long));
456     memset(p, 0, sizeof(long) * b->n_alleles);
457     for (j = 0; j < ma->n; ++j) {
458         const uint8_t *pi = ma->PL + j * ma->PL_len;
459         double *pdg = ma->pdg + j * 3;
460         pdg[0] = ma->q2p[pi[2]]; pdg[1] = ma->q2p[pi[1]]; pdg[2] = ma->q2p[pi[0]];
461         for (i = 0; i < b->n_alleles; ++i)
462             p[i] += (int)pi[(i+1)*(i+2)/2-1];
463     }
464     for (i = 0; i < b->n_alleles; ++i) p[i] = p[i]<<4 | i;
465     for (i = 1; i < b->n_alleles; ++i) // insertion sort
466         for (j = i; j > 0 && p[j] < p[j-1]; --j)
467             tmp = p[j], p[j] = p[j-1], p[j-1] = tmp;
468     for (i = b->n_alleles - 1; i >= 0; --i)
469         if ((p[i]&0xf) == 0) break;
470     return i;
471 }
472
473
474 int bcf_p1_call_gt(const bcf_p1aux_t *ma, double f0, int k)
475 {
476         double sum, g[3];
477         double max, f3[3], *pdg = ma->pdg + k * 3;
478         int q, i, max_i, ploidy;
479         ploidy = ma->ploidy? ma->ploidy[k] : 2;
480         if (ploidy == 2) {
481                 f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
482         } else {
483                 f3[0] = 1. - f0; f3[1] = 0; f3[2] = f0;
484         }
485         for (i = 0, sum = 0.; i < 3; ++i)
486                 sum += (g[i] = pdg[i] * f3[i]);
487         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
488                 g[i] /= sum;
489                 if (g[i] > max) max = g[i], max_i = i;
490         }
491         max = 1. - max;
492         if (max < 1e-308) max = 1e-308;
493         q = (int)(-4.343 * log(max) + .499);
494         if (q > 99) q = 99;
495         return q<<2|max_i;
496 }
497
498 #define TINY 1e-20
499
500 static void mc_cal_y_core(bcf_p1aux_t *ma, int beg)
501 {
502         double *z[2], *tmp, *pdg;
503         int _j, last_min, last_max;
504         assert(beg == 0 || ma->M == ma->n*2);
505         z[0] = ma->z;
506         z[1] = ma->zswap;
507         pdg = ma->pdg;
508         memset(z[0], 0, sizeof(double) * (ma->M + 1));
509         memset(z[1], 0, sizeof(double) * (ma->M + 1));
510         z[0][0] = 1.;
511         last_min = last_max = 0;
512         ma->t = 0.;
513         if (ma->M == ma->n * 2) {
514                 int M = 0;
515                 for (_j = beg; _j < ma->n; ++_j) {
516                         int k, j = _j - beg, _min = last_min, _max = last_max, M0;
517                         double p[3], sum;
518                         M0 = M; M += 2;
519                         pdg = ma->pdg + _j * 3;
520                         p[0] = pdg[0]; p[1] = 2. * pdg[1]; p[2] = pdg[2];
521                         for (; _min < _max && z[0][_min] < TINY; ++_min) z[0][_min] = z[1][_min] = 0.;
522                         for (; _max > _min && z[0][_max] < TINY; --_max) z[0][_max] = z[1][_max] = 0.;
523                         _max += 2;
524                         if (_min == 0) k = 0, z[1][k] = (M0-k+1) * (M0-k+2) * p[0] * z[0][k];
525                         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];
526                         for (k = _min < 2? 2 : _min; k <= _max; ++k)
527                                 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];
528                         for (k = _min, sum = 0.; k <= _max; ++k) sum += z[1][k];
529                         ma->t += log(sum / (M * (M - 1.)));
530                         for (k = _min; k <= _max; ++k) z[1][k] /= sum;
531                         if (_min >= 1) z[1][_min-1] = 0.;
532                         if (_min >= 2) z[1][_min-2] = 0.;
533                         if (j < ma->n - 1) z[1][_max+1] = z[1][_max+2] = 0.;
534                         if (_j == ma->n1 - 1) { // set pop1; ma->n1==-1 when unset
535                                 ma->t1 = ma->t;
536                                 memcpy(ma->z1, z[1], sizeof(double) * (ma->n1 * 2 + 1));
537                         }
538                         tmp = z[0]; z[0] = z[1]; z[1] = tmp;
539                         last_min = _min; last_max = _max;
540                 }
541                 //for (_j = 0; _j < last_min; ++_j) z[0][_j] = 0.; // TODO: are these necessary?
542                 //for (_j = last_max + 1; _j < ma->M; ++_j) z[0][_j] = 0.;
543         } else { // this block is very similar to the block above; these two might be merged in future
544                 int j, M = 0;
545                 for (j = 0; j < ma->n; ++j) {
546                         int k, M0, _min = last_min, _max = last_max;
547                         double p[3], sum;
548                         pdg = ma->pdg + j * 3;
549                         for (; _min < _max && z[0][_min] < TINY; ++_min) z[0][_min] = z[1][_min] = 0.;
550                         for (; _max > _min && z[0][_max] < TINY; --_max) z[0][_max] = z[1][_max] = 0.;
551                         M0 = M;
552                         M += ma->ploidy[j];
553                         if (ma->ploidy[j] == 1) {
554                                 p[0] = pdg[0]; p[1] = pdg[2];
555                                 _max++;
556                                 if (_min == 0) k = 0, z[1][k] = (M0+1-k) * p[0] * z[0][k];
557                                 for (k = _min < 1? 1 : _min; k <= _max; ++k)
558                                         z[1][k] = (M0+1-k) * p[0] * z[0][k] + k * p[1] * z[0][k-1];
559                                 for (k = _min, sum = 0.; k <= _max; ++k) sum += z[1][k];
560                                 ma->t += log(sum / M);
561                                 for (k = _min; k <= _max; ++k) z[1][k] /= sum;
562                                 if (_min >= 1) z[1][_min-1] = 0.;
563                                 if (j < ma->n - 1) z[1][_max+1] = 0.;
564                         } else if (ma->ploidy[j] == 2) {
565                                 p[0] = pdg[0]; p[1] = 2 * pdg[1]; p[2] = pdg[2];
566                                 _max += 2;
567                                 if (_min == 0) k = 0, z[1][k] = (M0-k+1) * (M0-k+2) * p[0] * z[0][k];
568                                 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];
569                                 for (k = _min < 2? 2 : _min; k <= _max; ++k)
570                                         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];
571                                 for (k = _min, sum = 0.; k <= _max; ++k) sum += z[1][k];
572                                 ma->t += log(sum / (M * (M - 1.)));
573                                 for (k = _min; k <= _max; ++k) z[1][k] /= sum;
574                                 if (_min >= 1) z[1][_min-1] = 0.;
575                                 if (_min >= 2) z[1][_min-2] = 0.;
576                                 if (j < ma->n - 1) z[1][_max+1] = z[1][_max+2] = 0.;
577                         }
578                         tmp = z[0]; z[0] = z[1]; z[1] = tmp;
579                         last_min = _min; last_max = _max;
580                 }
581         }
582         if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (ma->M + 1));
583 }
584
585 static void mc_cal_y(bcf_p1aux_t *ma)
586 {
587         if (ma->n1 > 0 && ma->n1 < ma->n && ma->M == ma->n * 2) { // NB: ma->n1 is ineffective when there are haploid samples
588                 int k;
589                 long double x;
590                 memset(ma->z1, 0, sizeof(double) * (2 * ma->n1 + 1));
591                 memset(ma->z2, 0, sizeof(double) * (2 * (ma->n - ma->n1) + 1));
592                 ma->t1 = ma->t2 = 0.;
593                 mc_cal_y_core(ma, ma->n1);
594                 ma->t2 = ma->t;
595                 memcpy(ma->z2, ma->z, sizeof(double) * (2 * (ma->n - ma->n1) + 1));
596                 mc_cal_y_core(ma, 0);
597                 // rescale z
598                 x = expl(ma->t - (ma->t1 + ma->t2));
599                 for (k = 0; k <= ma->M; ++k) ma->z[k] *= x;
600         } else mc_cal_y_core(ma, 0);
601 }
602
603 #define CONTRAST_TINY 1e-30
604
605 extern double kf_gammaq(double s, double z); // incomplete gamma function for chi^2 test
606
607 static inline double chi2_test(int a, int b, int c, int d)
608 {
609         double x, z;
610         x = (double)(a+b) * (c+d) * (b+d) * (a+c);
611         if (x == 0.) return 1;
612         z = a * d - b * c;
613         return kf_gammaq(.5, .5 * z * z * (a+b+c+d) / x);
614 }
615
616 // chi2=(a+b+c+d)(ad-bc)^2/[(a+b)(c+d)(a+c)(b+d)]
617 static inline double contrast2_aux(const bcf_p1aux_t *p1, double sum, int k1, int k2, double x[3])
618 {
619         double p = p1->phi[k1+k2] * p1->z1[k1] * p1->z2[k2] / sum * p1->hg[k1][k2];
620         int n1 = p1->n1, n2 = p1->n - p1->n1;
621         if (p < CONTRAST_TINY) return -1;
622         if (.5*k1/n1 < .5*k2/n2) x[1] += p;
623         else if (.5*k1/n1 > .5*k2/n2) x[2] += p;
624         else x[0] += p;
625         return p * chi2_test(k1, k2, (n1<<1) - k1, (n2<<1) - k2);
626 }
627
628 static double contrast2(bcf_p1aux_t *p1, double ret[3])
629 {
630         int k, k1, k2, k10, k20, n1, n2;
631         double sum;
632         // get n1 and n2
633         n1 = p1->n1; n2 = p1->n - p1->n1;
634         if (n1 <= 0 || n2 <= 0) return 0.;
635         if (p1->hg == 0) { // initialize the hypergeometric distribution
636                 /* NB: the hg matrix may take a lot of memory when there are many samples. There is a way
637                    to avoid precomputing this matrix, but it is slower and quite intricate. The following
638                    computation in this block can be accelerated with a similar strategy, but perhaps this
639                    is not a serious concern for now. */
640                 double tmp = lgamma(2*(n1+n2)+1) - (lgamma(2*n1+1) + lgamma(2*n2+1));
641                 p1->hg = calloc(2*n1+1, sizeof(void*));
642                 for (k1 = 0; k1 <= 2*n1; ++k1) {
643                         p1->hg[k1] = calloc(2*n2+1, sizeof(double));
644                         for (k2 = 0; k2 <= 2*n2; ++k2)
645                                 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));
646                 }
647         }
648         { // compute
649                 long double suml = 0;
650                 for (k = 0; k <= p1->M; ++k) suml += p1->phi[k] * p1->z[k];
651                 sum = suml;
652         }
653         { // get the max k1 and k2
654                 double max;
655                 int max_k;
656                 for (k = 0, max = 0, max_k = -1; k <= 2*n1; ++k) {
657                         double x = p1->phi1[k] * p1->z1[k];
658                         if (x > max) max = x, max_k = k;
659                 }
660                 k10 = max_k;
661                 for (k = 0, max = 0, max_k = -1; k <= 2*n2; ++k) {
662                         double x = p1->phi2[k] * p1->z2[k];
663                         if (x > max) max = x, max_k = k;
664                 }
665                 k20 = max_k;
666         }
667         { // 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.
668                 double x[3], y;
669                 long double z = 0., L[2];
670                 x[0] = x[1] = x[2] = 0; L[0] = L[1] = 0;
671                 for (k1 = k10; k1 >= 0; --k1) {
672                         for (k2 = k20; k2 >= 0; --k2) {
673                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
674                                 else z += y;
675                         }
676                         for (k2 = k20 + 1; k2 <= 2*n2; ++k2) {
677                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
678                                 else z += y;
679                         }
680                 }
681                 ret[0] = x[0]; ret[1] = x[1]; ret[2] = x[2];
682                 x[0] = x[1] = x[2] = 0;
683                 for (k1 = k10 + 1; k1 <= 2*n1; ++k1) {
684                         for (k2 = k20; k2 >= 0; --k2) {
685                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
686                                 else z += y;
687                         }
688                         for (k2 = k20 + 1; k2 <= 2*n2; ++k2) {
689                                 if ((y = contrast2_aux(p1, sum, k1, k2, x)) < 0) break;
690                                 else z += y;
691                         }
692                 }
693                 ret[0] += x[0]; ret[1] += x[1]; ret[2] += x[2];
694                 if (ret[0] + ret[1] + ret[2] < 0.95) { // in case of bad things happened
695                         ret[0] = ret[1] = ret[2] = 0; L[0] = L[1] = 0;
696                         for (k1 = 0, z = 0.; k1 <= 2*n1; ++k1)
697                                 for (k2 = 0; k2 <= 2*n2; ++k2)
698                                         if ((y = contrast2_aux(p1, sum, k1, k2, ret)) >= 0) z += y;
699                         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...
700                                 z = 1.0, ret[0] = ret[1] = ret[2] = 1./3;
701                 }
702                 return (double)z;
703         }
704 }
705
706 static double mc_cal_afs(bcf_p1aux_t *ma, double *p_ref_folded, double *p_var_folded)
707 {
708         int k;
709         long double sum = 0., sum2;
710         double *phi = ma->is_indel? ma->phi_indel : ma->phi;
711         memset(ma->afs1, 0, sizeof(double) * (ma->M + 1));
712         mc_cal_y(ma);
713         // compute AFS
714         for (k = 0, sum = 0.; k <= ma->M; ++k)
715                 sum += (long double)phi[k] * ma->z[k];
716         for (k = 0; k <= ma->M; ++k) {
717                 ma->afs1[k] = phi[k] * ma->z[k] / sum;
718                 if (isnan(ma->afs1[k]) || isinf(ma->afs1[k])) return -1.;
719         }
720         // compute folded variant probability
721         for (k = 0, sum = 0.; k <= ma->M; ++k)
722                 sum += (long double)(phi[k] + phi[ma->M - k]) / 2. * ma->z[k];
723         for (k = 1, sum2 = 0.; k < ma->M; ++k)
724                 sum2 += (long double)(phi[k] + phi[ma->M - k]) / 2. * ma->z[k];
725         *p_var_folded = sum2 / sum;
726         *p_ref_folded = (phi[k] + phi[ma->M - k]) / 2. * (ma->z[ma->M] + ma->z[0]) / sum;
727         // the expected frequency
728         for (k = 0, sum = 0.; k <= ma->M; ++k) {
729                 ma->afs[k] += ma->afs1[k];
730                 sum += k * ma->afs1[k];
731         }
732         return sum / ma->M;
733 }
734
735 int bcf_p1_cal(const bcf1_t *b, int do_contrast, bcf_p1aux_t *ma, bcf_p1rst_t *rst)
736 {
737         int i, k;
738         long double sum = 0.;
739         ma->is_indel = bcf_is_indel(b);
740         rst->perm_rank = -1;
741         // set PL and PL_len
742         for (i = 0; i < b->n_gi; ++i) {
743                 if (b->gi[i].fmt == bcf_str2int("PL", 2)) {
744                         ma->PL = (uint8_t*)b->gi[i].data;
745                         ma->PL_len = b->gi[i].len;
746                         break;
747                 }
748         }
749         if (i == b->n_gi) return -1; // no PL
750         if (b->n_alleles < 2) return -1; // FIXME: find a better solution
751         // 
752         rst->rank0 = cal_pdg(b, ma);
753         rst->f_exp = mc_cal_afs(ma, &rst->p_ref_folded, &rst->p_var_folded);
754         rst->p_ref = ma->afs1[ma->M];
755         for (k = 0, sum = 0.; k < ma->M; ++k)
756                 sum += ma->afs1[k];
757         rst->p_var = (double)sum;
758         { // compute the allele count
759                 double max = -1;
760                 rst->ac = -1;
761                 for (k = 0; k <= ma->M; ++k)
762                         if (max < ma->z[k]) max = ma->z[k], rst->ac = k;
763                 rst->ac = ma->M - rst->ac;
764         }
765         // calculate f_flat and f_em
766         for (k = 0, sum = 0.; k <= ma->M; ++k)
767                 sum += (long double)ma->z[k];
768         rst->f_flat = 0.;
769         for (k = 0; k <= ma->M; ++k) {
770                 double p = ma->z[k] / sum;
771                 rst->f_flat += k * p;
772         }
773         rst->f_flat /= ma->M;
774         { // estimate equal-tail credible interval (95% level)
775                 int l, h;
776                 double p;
777                 for (i = 0, p = 0.; i <= ma->M; ++i)
778                         if (p + ma->afs1[i] > 0.025) break;
779                         else p += ma->afs1[i];
780                 l = i;
781                 for (i = ma->M, p = 0.; i >= 0; --i)
782                         if (p + ma->afs1[i] > 0.025) break;
783                         else p += ma->afs1[i];
784                 h = i;
785                 rst->cil = (double)(ma->M - h) / ma->M; rst->cih = (double)(ma->M - l) / ma->M;
786         }
787         if (ma->n1 > 0) { // compute LRT
788                 double max0, max1, max2;
789                 for (k = 0, max0 = -1; k <= ma->M; ++k)
790                         if (max0 < ma->z[k]) max0 = ma->z[k];
791                 for (k = 0, max1 = -1; k <= ma->n1 * 2; ++k)
792                         if (max1 < ma->z1[k]) max1 = ma->z1[k];
793                 for (k = 0, max2 = -1; k <= ma->M - ma->n1 * 2; ++k)
794                         if (max2 < ma->z2[k]) max2 = ma->z2[k];
795                 rst->lrt = log(max1 * max2 / max0);
796                 rst->lrt = rst->lrt < 0? 1 : kf_gammaq(.5, rst->lrt);
797         } else rst->lrt = -1.0;
798         rst->cmp[0] = rst->cmp[1] = rst->cmp[2] = rst->p_chi2 = -1.0;
799         if (do_contrast && rst->p_var > 0.5) // skip contrast2() if the locus is a strong non-variant
800                 rst->p_chi2 = contrast2(ma, rst->cmp);
801         return 0;
802 }
803
804 void bcf_p1_dump_afs(bcf_p1aux_t *ma)
805 {
806         int k;
807         fprintf(stderr, "[afs]");
808         for (k = 0; k <= ma->M; ++k)
809                 fprintf(stderr, " %d:%.3lf", k, ma->afs[ma->M - k]);
810         fprintf(stderr, "\n");
811         memset(ma->afs, 0, sizeof(double) * (ma->M + 1));
812 }