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