]> git.donarmstrong.com Git - samtools.git/blob - bcftools/prob1.c
* changed 3.434 to 4.343 (typo!)
[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 "prob1.h"
7
8 #include "kseq.h"
9 KSTREAM_INIT(gzFile, gzread, 16384)
10
11 #define MC_AVG_ERR 0.007
12 #define MC_MAX_EM_ITER 16
13 #define MC_EM_EPS 1e-4
14
15 //#define _BCF_QUAD
16
17 unsigned char seq_nt4_table[256] = {
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, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
22         4, 0, 4, 1,  4, 4, 4, 2,  4, 4, 4, 4,  4, 4, 4, 4, 
23         4, 4, 4, 4,  3, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4, 
24         4, 0, 4, 1,  4, 4, 4, 2,  4, 4, 4, 4,  4, 4, 4, 4, 
25         4, 4, 4, 4,  3, 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         4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4,  4, 4, 4, 4
34 };
35
36 struct __bcf_p1aux_t {
37         int n, M, n1;
38         double *q2p, *pdg; // pdg -> P(D|g)
39         double *phi;
40         double *z, *zswap; // aux for afs
41         double *z1, *z2; // only calculated when n1 is set
42         double t, t1, t2;
43         double *afs, *afs1; // afs: accumulative AFS; afs1: site posterior distribution
44         double *k1k2;
45         const uint8_t *PL; // point to PL
46         int PL_len;
47 };
48
49 void bcf_p1_init_prior(bcf_p1aux_t *ma, int type, double theta)
50 {
51         int i;
52         if (type == MC_PTYPE_COND2) {
53                 for (i = 0; i <= ma->M; ++i)
54                         ma->phi[i] = 2. * (i + 1) / (ma->M + 1) / (ma->M + 2);
55         } else if (type == MC_PTYPE_FLAT) {
56                 for (i = 0; i <= ma->M; ++i)
57                         ma->phi[i] = 1. / (ma->M + 1);
58         } else {
59                 double sum;
60                 for (i = 0, sum = 0.; i < ma->M; ++i)
61                         sum += (ma->phi[i] = theta / (ma->M - i));
62                 ma->phi[ma->M] = 1. - sum;
63         }
64 }
65
66 int bcf_p1_read_prior(bcf_p1aux_t *ma, const char *fn)
67 {
68         gzFile fp;
69         kstring_t s;
70         kstream_t *ks;
71         long double sum;
72         int dret, k;
73         memset(&s, 0, sizeof(kstring_t));
74         fp = strcmp(fn, "-")? gzopen(fn, "r") : gzdopen(fileno(stdin), "r");
75         ks = ks_init(fp);
76         memset(ma->phi, 0, sizeof(double) * (ma->M + 1));
77         while (ks_getuntil(ks, '\n', &s, &dret) >= 0) {
78                 if (strstr(s.s, "[afs] ") == s.s) {
79                         char *p = s.s + 6;
80                         for (k = 0; k <= ma->M; ++k) {
81                                 int x;
82                                 double y;
83                                 x = strtol(p, &p, 10);
84                                 if (x != k && (errno == EINVAL || errno == ERANGE)) return -1;
85                                 ++p;
86                                 y = strtod(p, &p);
87                                 if (y == 0. && (errno == EINVAL || errno == ERANGE)) return -1;
88                                 ma->phi[ma->M - k] += y;
89                         }
90                 }
91         }
92         ks_destroy(ks);
93         gzclose(fp);
94         free(s.s);
95         for (sum = 0., k = 0; k <= ma->M; ++k) sum += ma->phi[k];
96         fprintf(stderr, "[prior]");
97         for (k = 0; k <= ma->M; ++k) ma->phi[k] /= sum;
98         for (k = 0; k <= ma->M; ++k) fprintf(stderr, " %d:%.3lg", k, ma->phi[ma->M - k]);
99         fputc('\n', stderr);
100         for (sum = 0., k = 1; k <= ma->M; ++k) sum += k * ma->phi[ma->M - k];
101         fprintf(stderr, "[heterozygosity] %lf\n", (double)sum / ma->M);
102         return 0;
103 }
104
105 bcf_p1aux_t *bcf_p1_init(int n)
106 {
107         bcf_p1aux_t *ma;
108         int i;
109         ma = calloc(1, sizeof(bcf_p1aux_t));
110         ma->n1 = -1;
111         ma->n = n; ma->M = 2 * n;
112         ma->q2p = calloc(256, sizeof(double));
113         ma->pdg = calloc(3 * ma->n, sizeof(double));
114         ma->phi = calloc(ma->M + 1, sizeof(double));
115         ma->z = calloc(2 * ma->n + 1, sizeof(double));
116         ma->zswap = calloc(2 * ma->n + 1, sizeof(double));
117         ma->z1 = calloc(ma->M + 1, sizeof(double)); // actually we do not need this large
118         ma->z2 = calloc(ma->M + 1, sizeof(double));
119         ma->afs = calloc(2 * ma->n + 1, sizeof(double));
120         ma->afs1 = calloc(2 * ma->n + 1, sizeof(double));
121         for (i = 0; i < 256; ++i)
122                 ma->q2p[i] = pow(10., -i / 10.);
123         bcf_p1_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior
124         return ma;
125 }
126
127 #ifdef _BCF_QUAD
128 static double lbinom(int n, int k)
129 {
130         return lgamma(n+1) - lgamma(k+1) - lgamma(n-k+1);
131 }
132 #endif
133
134 int bcf_p1_set_n1(bcf_p1aux_t *b, int n1)
135 {
136         if (n1 == 0 || n1 >= b->n) return -1;
137         b->n1 = n1;
138 #ifdef _BCF_QUAD
139         {
140                 int k1, k2, n2 = b->n - b->n1;
141                 b->k1k2 = calloc((2*n1+1) * (2*n2+1), sizeof(double));
142                 for (k1 = 0; k1 <= 2*n1; ++k1)
143                         for (k2 = 0; k2 <= 2*n2; ++k2)
144                                 b->k1k2[k1*(2*n2+1)+k2] = exp(lbinom(2*n1,k1) + lbinom(2*n2,k2) - lbinom(b->M,k1+k2));
145         }
146 #endif
147         return 0;
148 }
149
150 void bcf_p1_destroy(bcf_p1aux_t *ma)
151 {
152         if (ma) {
153                 free(ma->q2p); free(ma->pdg);
154                 free(ma->phi);
155                 free(ma->z); free(ma->zswap); free(ma->z1); free(ma->z2);
156                 free(ma->afs); free(ma->afs1);
157                 free(ma->k1k2);
158                 free(ma);
159         }
160 }
161
162 #define char2int(s) (((int)s[0])<<8|s[1])
163
164 static int cal_pdg(const bcf1_t *b, bcf_p1aux_t *ma)
165 {
166         int i, j, k;
167         long *p, tmp;
168         p = alloca(b->n_alleles * sizeof(long));
169         memset(p, 0, sizeof(long) * b->n_alleles);
170         for (j = 0; j < ma->n; ++j) {
171                 const uint8_t *pi = ma->PL + j * ma->PL_len;
172                 double *pdg = ma->pdg + j * 3;
173                 pdg[0] = ma->q2p[pi[b->n_alleles]]; pdg[1] = ma->q2p[pi[1]]; pdg[2] = ma->q2p[pi[0]];
174                 for (i = k = 0; i < b->n_alleles; ++i) {
175                         p[i] += (int)pi[k];
176                         k += b->n_alleles - i;
177                 }
178         }
179         for (i = 0; i < b->n_alleles; ++i) p[i] = p[i]<<4 | i;
180         for (i = 1; i < b->n_alleles; ++i) // insertion sort
181                 for (j = i; j > 0 && p[j] < p[j-1]; --j)
182                         tmp = p[j], p[j] = p[j-1], p[j-1] = tmp;
183         for (i = b->n_alleles - 1; i >= 0; --i)
184                 if ((p[i]&0xf) == 0) break;
185         return i;
186 }
187 // f0 is the reference allele frequency
188 static double mc_freq_iter(double f0, const bcf_p1aux_t *ma)
189 {
190         double f, f3[3];
191         int i;
192         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
193         for (i = 0, f = 0.; i < ma->n; ++i) {
194                 double *pdg;
195                 pdg = ma->pdg + i * 3;
196                 f += (pdg[1] * f3[1] + 2. * pdg[2] * f3[2])
197                         / (pdg[0] * f3[0] + pdg[1] * f3[1] + pdg[2] * f3[2]);
198         }
199         f /= ma->n * 2.;
200         return f;
201 }
202
203 int bcf_p1_call_gt(const bcf_p1aux_t *ma, double f0, int k)
204 {
205         double sum, g[3];
206         double max, f3[3], *pdg = ma->pdg + k * 3;
207         int q, i, max_i;
208         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
209         for (i = 0, sum = 0.; i < 3; ++i)
210                 sum += (g[i] = pdg[i] * f3[i]);
211         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
212                 g[i] /= sum;
213                 if (g[i] > max) max = g[i], max_i = i;
214         }
215         max = 1. - max;
216         if (max < 1e-308) max = 1e-308;
217         q = (int)(-4.343 * log(max) + .499);
218         if (q > 99) q = 99;
219         return q<<2|max_i;
220 }
221
222 #define TINY 1e-20
223
224 static void mc_cal_y_core(bcf_p1aux_t *ma, int beg)
225 {
226         double *z[2], *tmp, *pdg;
227         int _j, last_min, last_max;
228         z[0] = ma->z;
229         z[1] = ma->zswap;
230         pdg = ma->pdg;
231         memset(z[0], 0, sizeof(double) * (ma->M + 1));
232         memset(z[1], 0, sizeof(double) * (ma->M + 1));
233         z[0][0] = 1.;
234         last_min = last_max = 0;
235         ma->t = 0.;
236         for (_j = beg; _j < ma->n; ++_j) {
237                 int k, j = _j - beg, _min = last_min, _max = last_max;
238                 double p[3], sum;
239                 pdg = ma->pdg + _j * 3;
240                 p[0] = pdg[0]; p[1] = 2. * pdg[1]; p[2] = pdg[2];
241                 for (; _min < _max && z[0][_min] < TINY; ++_min) z[0][_min] = z[1][_min] = 0.;
242                 for (; _max > _min && z[0][_max] < TINY; --_max) z[0][_max] = z[1][_max] = 0.;
243                 _max += 2;
244                 if (_min == 0) 
245                         k = 0, z[1][k] = (2*j+2-k)*(2*j-k+1) * p[0] * z[0][k];
246                 if (_min <= 1)
247                         k = 1, z[1][k] = (2*j+2-k)*(2*j-k+1) * p[0] * z[0][k] + k*(2*j+2-k) * p[1] * z[0][k-1];
248                 for (k = _min < 2? 2 : _min; k <= _max; ++k)
249                         z[1][k] = (2*j+2-k)*(2*j-k+1) * p[0] * z[0][k]
250                                 + k*(2*j+2-k) * p[1] * z[0][k-1]
251                                 + k*(k-1)* p[2] * z[0][k-2];
252                 for (k = _min, sum = 0.; k <= _max; ++k) sum += z[1][k];
253                 ma->t += log(sum / ((2. * j + 2) * (2. * j + 1)));
254                 for (k = _min; k <= _max; ++k) z[1][k] /= sum;
255                 if (_min >= 1) z[1][_min-1] = 0.;
256                 if (_min >= 2) z[1][_min-2] = 0.;
257                 if (j < ma->n - 1) z[1][_max+1] = z[1][_max+2] = 0.;
258                 if (_j == ma->n1 - 1) { // set pop1
259                         ma->t1 = ma->t;
260                         memcpy(ma->z1, z[1], sizeof(double) * (ma->n1 * 2 + 1));
261                 }
262                 tmp = z[0]; z[0] = z[1]; z[1] = tmp;
263                 last_min = _min; last_max = _max;
264         }
265         if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (ma->M + 1));
266 }
267
268 static void mc_cal_y(bcf_p1aux_t *ma)
269 {
270         if (ma->n1 > 0 && ma->n1 < ma->n) {
271                 int k;
272                 long double x;
273                 memset(ma->z1, 0, sizeof(double) * (2 * ma->n1 + 1));
274                 memset(ma->z2, 0, sizeof(double) * (2 * (ma->n - ma->n1) + 1));
275                 ma->t1 = ma->t2 = 0.;
276                 mc_cal_y_core(ma, ma->n1);
277                 ma->t2 = ma->t;
278                 memcpy(ma->z2, ma->z, sizeof(double) * (2 * (ma->n - ma->n1) + 1));
279                 mc_cal_y_core(ma, 0);
280                 // rescale z
281                 x = expl(ma->t - (ma->t1 + ma->t2));
282                 for (k = 0; k <= ma->M; ++k) ma->z[k] *= x;
283         } else mc_cal_y_core(ma, 0);
284 #ifdef _BCF_QUAD
285 /*
286         if (ma->n1 > 0 && ma->n1 < ma->n) { // DEBUG: consistency check; z[i] should equal y[i]
287                 int i, k1, k2, n1 = ma->n1, n2 = ma->n - n1;
288                 double *y;
289                 printf("*** ");
290                 y = calloc(ma->M + 1, sizeof(double));
291                 for (k1 = 0; k1 <= 2*n1; ++k1)
292                         for (k2 = 0; k2 <= 2*n2; ++k2)
293                                 y[k1+k2] += ma->k1k2[k1*(2*n2+1)+k2] * ma->z1[k1] * ma->z2[k2];
294                 for (i = 0; i <= ma->M; ++i) printf("(%lf,%lf) ", ma->z[i], y[i]);
295                 printf("\n");
296                 free(y);
297         }
298 */
299 #endif
300 }
301
302 static void contrast(bcf_p1aux_t *ma, double pc[4]) // mc_cal_y() must be called before hand
303 {
304         int k, n1 = ma->n1, n2 = ma->n - ma->n1;
305         long double sum = -1., x, sum_alt;
306         double y;
307         pc[0] = pc[1] = pc[2] = pc[3] = -1.;
308         if (n1 <= 0 || n2 <= 0) return;
309 #ifdef _BCF_QUAD
310         { // FIXME: can be improved by skipping zero cells
311                 int k1, k2;
312                 long double z[3];
313                 z[0] = z[1] = z[2] = 0.;
314                 for (k1 = 0; k1 <= 2*n1; ++k1)
315                         for (k2 = 0; k2 <= 2*n2; ++k2) {
316                                 double zz = ma->phi[k1+k2] * ma->z1[k1] * ma->z2[k2] * ma->k1k2[k1*(2*n2+1)+k2];
317                                 if ((double)k1/n1 < (double)k2/n2) z[0] += zz;
318                                 else if ((double)k1/n1 > (double)k2/n2) z[1] += zz;
319                                 else z[2] += zz;
320                         }
321                 sum = z[0] + z[1] + z[2];
322                 pc[2] = z[0] / sum; pc[3] = z[1] / sum;
323         }
324 #else
325         pc[2] = pc[3] = 0.;
326 #endif
327         for (k = 0, sum_alt = 0.; k <= ma->M; ++k)
328                 sum_alt += (long double)ma->phi[k] * ma->z[k];
329 //      printf("* %lg, %lg *\n", (double)sum, (double)sum_alt); // DEBUG: sum should equal sum_alt
330         sum = sum_alt;
331         // the variant is specific to group2
332 //      printf("%lg %lg %lg %lg\n", ma->z[2*(n1+n2)]/exp(ma->t - (ma->t1 + ma->t2)), ma->z1[2*n1], ma->z2[2*n2], (double)sum);
333         y = lgamma(2*n2 + 1) - lgamma(ma->M + 1);
334         for (k = 0, x = 0.; k < 2 * n2; ++k)
335                 x += ma->phi[2*n1+k] * ma->z2[k] * expl(lgamma(2*n1 + k + 1) - lgamma(k + 1) + y);
336         pc[1] = ma->z1[2*n1] * x / sum;
337         for (k = 1, x = 0.; k <= 2 * n2; ++k)
338                 x += ma->phi[k] * ma->z2[k] * expl(lgamma(ma->M - k + 1) - lgamma(2*n2 - k + 1) + y);
339         pc[1] += ma->z1[0] * x / sum;
340         // the variant is specific to group1
341         y = lgamma(2*n1 + 1) - lgamma(ma->M + 1);
342         for (k = 0, x = 0.; k < 2 * n1; ++k)
343                 x += ma->phi[2*n2+k] * ma->z1[k] * expl(lgamma(2*n2 + k + 1) - lgamma(k + 1) + y);
344         pc[0] = ma->z2[2*n2] * x / sum;
345         for (k = 1, x = 0.; k <= 2 * n1; ++k)
346                 x += ma->phi[k] * ma->z1[k] * expl(lgamma(ma->M - k + 1) - lgamma(2*n1 - k + 1) + y);
347         pc[0] += ma->z2[0] * x / sum;
348         // rescale
349         for (k = 2; k < 4; ++k) {
350                 y = 1. - pc[k];
351                 if (y <= 0.) y = 1e-100;
352                 pc[k] = (int)(-4.343 * log(y) + .499);
353                 if (pc[k] > 99.) pc[k] = 99.;
354         }
355 }
356
357 static double mc_cal_afs(bcf_p1aux_t *ma)
358 {
359         int k;
360         long double sum = 0.;
361         memset(ma->afs1, 0, sizeof(double) * (ma->M + 1));
362         mc_cal_y(ma);
363         for (k = 0, sum = 0.; k <= ma->M; ++k)
364                 sum += (long double)ma->phi[k] * ma->z[k];
365         for (k = 0; k <= ma->M; ++k) {
366                 ma->afs1[k] = ma->phi[k] * ma->z[k] / sum;
367                 if (isnan(ma->afs1[k]) || isinf(ma->afs1[k])) return -1.;
368         }
369         for (k = 0, sum = 0.; k <= ma->M; ++k) {
370                 ma->afs[k] += ma->afs1[k];
371                 sum += k * ma->afs1[k];
372         }
373         return sum / ma->M;
374 }
375
376 long double bcf_p1_cal_g3(bcf_p1aux_t *p1a, double g[3])
377 {
378         long double pd = 0., g2[3];
379         int i, k;
380         memset(g2, 0, sizeof(long double) * 3);
381         for (k = 0; k < p1a->M; ++k) {
382                 double f = (double)k / p1a->M, f3[3], g1[3];
383                 long double z = 1.;
384                 g1[0] = g1[1] = g1[2] = 0.;
385                 f3[0] = (1. - f) * (1. - f); f3[1] = 2. * f * (1. - f); f3[2] = f * f;
386                 for (i = 0; i < p1a->n; ++i) {
387                         double *pdg = p1a->pdg + i * 3;
388                         double x = pdg[0] * f3[0] + pdg[1] * f3[1] + pdg[2] * f3[2];
389                         z *= x;
390                         g1[0] += pdg[0] * f3[0] / x;
391                         g1[1] += pdg[1] * f3[1] / x;
392                         g1[2] += pdg[2] * f3[2] / x;
393                 }
394                 pd += p1a->phi[k] * z;
395                 for (i = 0; i < 3; ++i)
396                         g2[i] += p1a->phi[k] * z * g1[i];
397         }
398         for (i = 0; i < 3; ++i) g[i] = g2[i] / pd;
399         return pd;
400 }
401
402 int bcf_p1_cal(bcf1_t *b, bcf_p1aux_t *ma, bcf_p1rst_t *rst)
403 {
404         int i, k;
405         long double sum = 0.;
406         // set PL and PL_len
407         for (i = 0; i < b->n_gi; ++i) {
408                 if (b->gi[i].fmt == char2int("PL")) {
409                         ma->PL = (uint8_t*)b->gi[i].data;
410                         ma->PL_len = b->gi[i].len;
411                         break;
412                 }
413         }
414         if (b->n_alleles < 2) return -1; // FIXME: find a better solution
415         // 
416         rst->rank0 = cal_pdg(b, ma);
417         rst->f_exp = mc_cal_afs(ma);
418         rst->p_ref = ma->afs1[ma->M];
419         // calculate f_flat and f_em
420         for (k = 0, sum = 0.; k <= ma->M; ++k)
421                 sum += (long double)ma->z[k];
422         rst->f_flat = 0.;
423         for (k = 0; k <= ma->M; ++k) {
424                 double p = ma->z[k] / sum;
425                 rst->f_flat += k * p;
426         }
427         rst->f_flat /= ma->M;
428         { // calculate f_em
429                 double flast = rst->f_flat;
430                 for (i = 0; i < MC_MAX_EM_ITER; ++i) {
431                         rst->f_em = mc_freq_iter(flast, ma);
432                         if (fabs(rst->f_em - flast) < MC_EM_EPS) break;
433                         flast = rst->f_em;
434                 }
435         }
436         rst->g[0] = rst->g[1] = rst->g[2] = -1.;
437         contrast(ma, rst->pc);
438 //      bcf_p1_cal_g3(ma, rst->g);
439         return 0;
440 }
441
442 void bcf_p1_dump_afs(bcf_p1aux_t *ma)
443 {
444         int k;
445         fprintf(stderr, "[afs]");
446         for (k = 0; k <= ma->M; ++k)
447                 fprintf(stderr, " %d:%.3lf", k, ma->afs[ma->M - k]);
448         fprintf(stderr, "\n");
449         memset(ma->afs, 0, sizeof(double) * (ma->M + 1));
450 }