]> git.donarmstrong.com Git - samtools.git/blob - bam_mcns.c
posterior expectation FINALLY working. I am so tired...
[samtools.git] / bam_mcns.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <stdio.h>
4 #include "bam_mcns.h"
5
6 #define MC_MIN_QUAL 13
7 #define MC_AVG_ERR 0.007
8 #define MC_MAX_SUMQ 3000
9 #define MC_MAX_SUMQP 1e-300
10 #define MC_MAX_EM_ITER 16
11 #define MC_EM_EPS 1e-4
12
13 struct __mc_aux_t {
14         int n, M;
15         int ref, alt;
16         double *q2p, *pdg; // pdg -> P(D|g)
17         double *alpha, *beta;
18         double *z, *zswap; // aux for afs
19         double *afs, *afs1; // afs: accumulative AFS; afs1: site posterior distribution
20         int *qsum, *bcnt;
21 };
22
23 void mc_init_prior(mc_aux_t *ma, int type, double theta)
24 {
25         int i;
26         if (type == MC_PTYPE_COND2) {
27                 for (i = 0; i <= 2 * ma->n; ++i)
28                         ma->alpha[i] = 2. * (i + 1) / (2 * ma->n + 1) / (2 * ma->n + 2);
29         } else {
30                 double sum;
31                 for (i = 0, sum = 0.; i < 2 * ma->n; ++i)
32                         sum += (ma->alpha[i] = theta / (2 * ma->n - i));
33                 ma->alpha[2 * ma->n] = 1. - sum;
34         }
35 }
36
37 mc_aux_t *mc_init(int n) // FIXME: assuming diploid
38 {
39         mc_aux_t *ma;
40         int i;
41         ma = calloc(1, sizeof(mc_aux_t));
42         ma->n = n; ma->M = 2 * n;
43         ma->q2p = calloc(MC_MAX_SUMQ + 1, sizeof(double));
44         ma->qsum = calloc(4 * ma->n, sizeof(int));
45         ma->bcnt = calloc(4 * ma->n, sizeof(int));
46         ma->pdg = calloc(3 * ma->n, sizeof(double));
47         ma->alpha = calloc(2 * ma->n + 1, sizeof(double));
48         ma->beta = calloc((2 * ma->n + 1) * 3, sizeof(double));
49         ma->z = calloc(2 * ma->n + 1, sizeof(double));
50         ma->zswap = calloc(2 * ma->n + 1, sizeof(double));
51         ma->afs = calloc(2 * ma->n + 1, sizeof(double));
52         ma->afs1 = calloc(2 * ma->n + 1, sizeof(double));
53         for (i = 0; i <= MC_MAX_SUMQ; ++i)
54                 ma->q2p[i] = pow(10., -i / 10.);
55         for (i = 0; i <= ma->M; ++i) { // beta[k][g]=P(g|k/M)
56                 double *bi = ma->beta + 3 * i;
57                 double f = (double)i / ma->M;
58                 bi[0] = (1. - f) * (1. - f);
59                 bi[1] = 2 * f * (1. - f);
60                 bi[2] = f * f;
61         }
62         mc_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior
63         return ma;
64 }
65
66 void mc_destroy(mc_aux_t *ma)
67 {
68         if (ma) {
69                 free(ma->qsum); free(ma->bcnt);
70                 free(ma->q2p); free(ma->pdg);
71                 free(ma->alpha); free(ma->beta);
72                 free(ma->z); free(ma->zswap);
73                 free(ma->afs); free(ma->afs1);
74                 free(ma);
75         }
76 }
77
78 static int sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
79 {
80         int i, j, tot = 0;
81         memset(ma->qsum, 0, sizeof(int) * 4 * ma->n);
82         memset(ma->bcnt, 0, sizeof(int) * 4 * ma->n);
83         for (j = 0; j < ma->n; ++j) {
84                 int *qsum = ma->qsum + j * 4;
85                 int *bcnt = ma->bcnt + j * 4;
86                 for (i = 0; i < n[j]; ++i) {
87                         const bam_pileup1_t *p = plp[j] + i;
88                         int q, b;
89                         if (p->is_del || (p->b->core.flag&BAM_FUNMAP)) continue;
90                         q = bam1_qual(p->b)[p->qpos];
91                         if (p->b->core.qual < q) q = p->b->core.qual;
92                         if (q < MC_MIN_QUAL) continue; // small qual
93                         b = bam_nt16_nt4_table[(int)bam1_seqi(bam1_seq(p->b), p->qpos)];
94                         if (b > 3) continue; // N
95                         qsum[b] += q;
96                         ++bcnt[b];
97                         ++tot;
98                 }
99         }
100         return tot;
101 }
102
103 static void set_allele(int ref, mc_aux_t *ma)
104 {
105         int i, j, sum[4], tmp;
106         sum[0] = sum[1] = sum[2] = sum[3] = 0;
107         for (i = 0; i < ma->n; ++i)
108                 for (j = 0; j < 4; ++j)
109                         sum[j] += ma->qsum[i * 4 + j];
110         for (j = 0; j < 4; ++j) sum[j] = sum[j]<<2 | j;
111         for (i = 1; i < 4; ++i) // insertion sort
112                 for (j = i; j > 0 && sum[j] < sum[j-1]; --j)
113                         tmp = sum[j], sum[j] = sum[j-1], sum[j-1] = tmp;
114         ma->ref = sum[3]&3; ma->alt = sum[2]&3;
115         if (ref == ma->alt) tmp = ma->ref, ma->ref = ma->alt, ma->alt = tmp;
116         // note that ma->ref might not be ref in case of triallele
117 }
118
119 static void cal_pdg(mc_aux_t *ma)
120 {
121         int i, j;
122         for (j = 0; j < ma->n; ++j) {
123                 int pi[3], *qsum, *bcnt;
124                 double *pdg = ma->pdg + j * 3;
125                 qsum = ma->qsum + j * 4;
126                 bcnt = ma->bcnt + j * 4;
127                 pi[1] = 3 * (bcnt[ma->ref] + bcnt[ma->alt]);
128                 pi[0] = qsum[ma->ref];
129                 pi[2] = qsum[ma->alt];
130                 for (i = 0; i < 3; ++i)
131                         pdg[i] = pi[i] > MC_MAX_SUMQ? MC_MAX_SUMQP : ma->q2p[pi[i]];
132         }
133 }
134 // this calculates the naive allele frequency and Nielsen's frequency
135 static double mc_freq0(const mc_aux_t *ma, double *_f)
136 {
137         int i, cnt;
138         double f, f_nielsen, w_sum;
139         *_f = -1.;
140         for (i = cnt = 0, f = f_nielsen = w_sum = 0.; i < ma->n; ++i) {
141                 int *bcnt = ma->bcnt + i * 4;
142                 int x = bcnt[ma->ref] + bcnt[ma->alt];
143                 if (x) {
144                         double w, p;
145                         ++cnt;
146                         f += (double)bcnt[ma->ref] / x;
147                         p = (bcnt[ma->ref] - MC_AVG_ERR * x) / (1. - 2. * MC_AVG_ERR) / x;
148                         w = 2. * x / (1. + x);
149                         w_sum += w;
150                         f_nielsen += p * w;
151                 }
152         }
153         if (cnt) {
154                 f_nielsen /= w_sum;
155                 if (f_nielsen < 0.) f_nielsen = 0.;
156                 if (f_nielsen > 1.) f_nielsen = 1.;
157                 *_f = f_nielsen;
158                 return f / cnt;
159         } else return -1.;
160 }
161 // f0 is the reference allele frequency
162 static double mc_freq_iter(double f0, const mc_aux_t *ma)
163 {
164         double f, f3[3];
165         int i;
166         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
167         for (i = 0, f = 0.; i < ma->n; ++i) {
168                 double *pdg;
169                 pdg = ma->pdg + i * 3;
170                 f += (pdg[1] * f3[1] + 2. * pdg[2] * f3[2])
171                         / (pdg[0] * f3[0] + pdg[1] * f3[1] + pdg[2] * f3[2]);
172         }
173         f /= ma->n * 2.;
174         return f;
175 }
176
177 static double mc_ref_prob(const mc_aux_t *ma, double *_PD, double *f_exp)
178 {
179         int k, i;
180         long double PD = 0., Pref = 0., Ef = 0.;
181         for (k = 0; k <= ma->M; ++k) {
182                 long double x = 1., y = 0.;
183                 double *bk = ma->beta + k * 3;
184                 for (i = 0; i < ma->n; ++i) {
185                         double *pdg = ma->pdg + i * 3;
186                         double z = pdg[0] * bk[0] + pdg[1] * bk[1] + pdg[2] * bk[2];
187                         x *= z;
188                         y += (pdg[1] * bk[1] + 2. * pdg[2] * bk[2]) / z;
189                 }
190                 PD += x * ma->alpha[k];
191                 Ef += x * y * ma->alpha[k];
192         }
193         for (k = 0; k <= ma->n * 2; ++k) {
194                 long double x = 1.0;
195                 for (i = 0; i < ma->n; ++i)
196                         x *= ma->pdg[i * 3 + 2] * ma->beta[k * 3 + 2];
197                 Pref += x * ma->alpha[k];
198         }
199         *f_exp = (double)(Ef / PD / ma->M);
200         *_PD = PD;
201         return Pref / PD;
202 }
203
204 int mc_call_gt(const mc_aux_t *ma, double f0, int k)
205 {
206         double sum, g[3];
207         double max, f3[3], *pdg = ma->pdg + k * 3;
208         int q, i, max_i;
209         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
210         for (i = 0, sum = 0.; i < 3; ++i)
211                 sum += (g[i] = pdg[i] * f3[i]);
212         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
213                 g[i] /= sum;
214                 if (g[i] > max) max = g[i], max_i = i;
215         }
216         max = 1. - max;
217         if (max < 1e-308) max = 1e-308;
218         q = (int)(-3.434 * log(max) + .499);
219         if (q > 99) q = 99;
220         return q<<2|max_i;
221 }
222 // calculate z_{nr}^{(k)}
223 static void mc_cal_z(mc_aux_t *ma, int k)
224 {
225         double *z[2], *tmp, *bk, *pdg;
226         int i, j;
227         z[0] = ma->z;
228         z[1] = ma->zswap;
229         bk = ma->beta + k * 3; pdg = ma->pdg;
230         z[0][0] = 1.; z[0][1] = z[0][2] = 0.;
231         for (j = 0; j < ma->n; ++j) {
232                 int max = (j + 1) * 2;
233                 double p[3];
234                 pdg = ma->pdg + j * 3;
235                 p[0] = bk[0] * pdg[0]; p[1] = bk[1] * pdg[1]; p[2] = bk[2] * pdg[2];
236                 z[1][0] = p[0] * z[0][0];
237                 z[1][1] = p[0] * z[0][1] + p[1] * z[0][0];
238                 for (i = 2; i <= max; ++i)
239                         z[1][i] = p[0] * z[0][i] + p[1] * z[0][i-1] + p[2] * z[0][i-2];
240                 if (j < ma->n - 1) z[1][max+1] = z[1][max+2] = 0.;
241                 tmp = z[0]; z[0] = z[1]; z[1] = tmp;
242         }
243         if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (2 * ma->n + 1));
244 }
245 // Warning: this is cubic in ma->n, very sloooooow
246 static void mc_add_afs(mc_aux_t *ma, double PD, double *f_map, double *p_map)
247 {
248         int k, l;
249         double sum = 0.;
250         memset(ma->afs1, 0, sizeof(double) * (2 * ma->n + 1));
251         for (k = 0; k <= 2 * ma->n; ++k) {
252                 mc_cal_z(ma, k);
253                 for (l = 0; l <= 2 * ma->n; ++l)
254                         ma->afs1[l] += ma->alpha[k] * ma->z[l] / PD;
255         }
256         for (k = 0; k <= 2 * ma->n; ++k) {
257                 ma->afs[k] += ma->afs1[k];
258                 sum += ma->afs1[k];
259         }
260         {
261                 int max_k = 0;
262                 double max = -1., e = 0.;
263                 for (k = 0; k <= 2 * ma->n; ++k) {
264                         if (ma->afs1[k] > max) max = ma->afs1[k], max_k = k;
265                         e += k * ma->afs1[k];
266                 }
267                 *f_map = .5 * max_k / ma->n; *p_map = max;
268                 printf(" * %.3lg:%.3lg:%.3lg:%.3lg * ", sum, 1.-.5*max_k/ma->n, max, 1.-.5*e/ma->n);
269         }
270 }
271
272 int mc_cal(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, mc_rst_t *rst, int level)
273 {
274         int i, tot;
275         memset(rst, 0, sizeof(mc_rst_t));
276         rst->f_em = rst->f_exp = -1.; rst->ref = rst->alt = -1;
277         // precalculation
278         tot = sum_err(n, plp, ma);
279         if (tot == 0) return 0; // no good bases
280         set_allele(ref, ma);
281         cal_pdg(ma);
282         // set ref/major allele
283         rst->ref = ma->ref; rst->alt = ma->alt;
284         // calculate naive and Nielsen's freq
285         rst->f_naive = mc_freq0(ma, &rst->f_nielsen);
286         { // calculate f_em
287                 double flast = rst->f_naive;
288                 for (i = 0; i < MC_MAX_EM_ITER; ++i) {
289                         rst->f_em = mc_freq_iter(flast, ma);
290                         if (fabs(rst->f_em - flast) < MC_EM_EPS) break;
291                         flast = rst->f_em;
292                 }
293         }
294         if (level >= 2) // quadratic-time calculations; necessary for genotyping
295                 rst->p_ref = mc_ref_prob(ma, &rst->PD, &rst->f_exp);
296         if (level >= 3)
297                 mc_add_afs(ma, rst->PD, &rst->f_map, &rst->p_map);
298         return tot;
299 }