]> git.donarmstrong.com Git - samtools.git/blob - bam_mcns.c
* use the revised MAQ error model for mpileup
[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, alt2;
16         double *q2p, *pdg; // pdg -> P(D|g)
17         double *phi, *CMk; // CMk=\binom{M}{k}
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->phi[i] = 2. * (i + 1) / (2 * ma->n + 1) / (2 * ma->n + 2);
29         } else if (type == MC_PTYPE_FLAT) {
30                 for (i = 0; i <= ma->M; ++i)
31                         ma->phi[i] = 1. / (ma->M + 1);
32         } else {
33                 double sum;
34                 for (i = 0, sum = 0.; i < 2 * ma->n; ++i)
35                         sum += (ma->phi[i] = theta / (2 * ma->n - i));
36                 ma->phi[2 * ma->n] = 1. - sum;
37         }
38 }
39
40 mc_aux_t *mc_init(int n) // FIXME: assuming diploid
41 {
42         mc_aux_t *ma;
43         int i;
44         ma = calloc(1, sizeof(mc_aux_t));
45         ma->n = n; ma->M = 2 * n;
46         ma->q2p = calloc(MC_MAX_SUMQ + 1, sizeof(double));
47         ma->qsum = calloc(4 * ma->n, sizeof(int));
48         ma->bcnt = calloc(4 * ma->n, sizeof(int));
49         ma->pdg = calloc(3 * ma->n, sizeof(double));
50         ma->phi = calloc(ma->M + 1, sizeof(double));
51         ma->CMk = calloc(ma->M + 1, sizeof(double));
52         ma->z = calloc(2 * ma->n + 1, sizeof(double));
53         ma->zswap = calloc(2 * ma->n + 1, sizeof(double));
54         ma->afs = calloc(2 * ma->n + 1, sizeof(double));
55         ma->afs1 = calloc(2 * ma->n + 1, sizeof(double));
56         for (i = 0; i <= MC_MAX_SUMQ; ++i)
57                 ma->q2p[i] = pow(10., -i / 10.);
58         for (i = 0; i <= ma->M; ++i)
59                 ma->CMk[i] = exp(lgamma(ma->M + 1) - lgamma(i + 1) - lgamma(ma->M - i + 1));
60         mc_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior
61         return ma;
62 }
63
64 void mc_destroy(mc_aux_t *ma)
65 {
66         if (ma) {
67                 free(ma->qsum); free(ma->bcnt);
68                 free(ma->q2p); free(ma->pdg);
69                 free(ma->phi); free(ma->CMk);
70                 free(ma->z); free(ma->zswap);
71                 free(ma->afs); free(ma->afs1);
72                 free(ma);
73         }
74 }
75
76 static int sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
77 {
78         int i, j, tot = 0;
79         memset(ma->qsum, 0, sizeof(int) * 4 * ma->n);
80         memset(ma->bcnt, 0, sizeof(int) * 4 * ma->n);
81         for (j = 0; j < ma->n; ++j) {
82                 int *qsum = ma->qsum + j * 4;
83                 int *bcnt = ma->bcnt + j * 4;
84                 for (i = 0; i < n[j]; ++i) {
85                         const bam_pileup1_t *p = plp[j] + i;
86                         int q, b;
87                         if (p->is_del || (p->b->core.flag&BAM_FUNMAP)) continue;
88                         q = bam1_qual(p->b)[p->qpos];
89                         if (p->b->core.qual < q) q = p->b->core.qual;
90                         if (q < MC_MIN_QUAL) continue; // small qual
91                         b = bam_nt16_nt4_table[(int)bam1_seqi(bam1_seq(p->b), p->qpos)];
92                         if (b > 3) continue; // N
93                         qsum[b] += q;
94                         ++bcnt[b];
95                         ++tot;
96                 }
97         }
98         return tot;
99 }
100
101 static void set_allele(int ref, mc_aux_t *ma)
102 {
103         int i, j, sum[4], tmp;
104         sum[0] = sum[1] = sum[2] = sum[3] = 0;
105         for (i = 0; i < ma->n; ++i)
106                 for (j = 0; j < 4; ++j)
107                         sum[j] += ma->qsum[i * 4 + j];
108         for (j = 0; j < 4; ++j) sum[j] = sum[j]<<2 | j;
109         for (i = 1; i < 4; ++i) // insertion sort
110                 for (j = i; j > 0 && sum[j] < sum[j-1]; --j)
111                         tmp = sum[j], sum[j] = sum[j-1], sum[j-1] = tmp;
112         ma->ref = sum[3]&3; ma->alt = sum[2]&3; ma->alt2 = -1;
113         if (ma->ref != ref) { // the best base is not ref
114                 if (ref >= 0 && ref <= 3) { // ref is not N
115                         if (ma->alt == ref) tmp = ma->ref, ma->ref = ma->alt, ma->alt = tmp; // then switch alt and ref
116                         else ma->alt2 = ma->alt, ma->alt = ma->ref, ma->ref = ref; // then set ref as ref
117                 } else ma->alt2 = ma->alt, ma->alt = ma->ref, ma->ref = sum[0]&3; // then set the weakest as ref
118         }
119 }
120
121 static void cal_pdg(mc_aux_t *ma)
122 {
123         int i, j;
124         for (j = 0; j < ma->n; ++j) {
125                 int pi[3], *qsum, *bcnt;
126                 double *pdg = ma->pdg + j * 3;
127                 qsum = ma->qsum + j * 4;
128                 bcnt = ma->bcnt + j * 4;
129                 pi[1] = 3 * (bcnt[ma->ref] + bcnt[ma->alt]);
130                 pi[0] = qsum[ma->ref];
131                 pi[2] = qsum[ma->alt];
132                 for (i = 0; i < 3; ++i)
133                         pdg[i] = pi[i] > MC_MAX_SUMQ? MC_MAX_SUMQP : ma->q2p[pi[i]];
134         }
135 }
136 // this calculates the naive allele frequency and Nielsen's frequency
137 static double mc_freq0(const mc_aux_t *ma, double *_f)
138 {
139         int i, cnt;
140         double f, f_nielsen, w_sum;
141         *_f = -1.;
142         for (i = cnt = 0, f = f_nielsen = w_sum = 0.; i < ma->n; ++i) {
143                 int *bcnt = ma->bcnt + i * 4;
144                 int x = bcnt[ma->ref] + bcnt[ma->alt];
145                 if (x) {
146                         double w, p;
147                         ++cnt;
148                         f += (double)bcnt[ma->ref] / x;
149                         p = (bcnt[ma->ref] - MC_AVG_ERR * x) / (1. - 2. * MC_AVG_ERR) / x;
150                         w = 2. * x / (1. + x);
151                         w_sum += w;
152                         f_nielsen += p * w;
153                 }
154         }
155         if (cnt) {
156                 f_nielsen /= w_sum;
157                 if (f_nielsen < 0.) f_nielsen = 0.;
158                 if (f_nielsen > 1.) f_nielsen = 1.;
159                 *_f = f_nielsen;
160                 return f / cnt;
161         } else return -1.;
162 }
163 // f0 is the reference allele frequency
164 static double mc_freq_iter(double f0, const mc_aux_t *ma)
165 {
166         double f, f3[3];
167         int i;
168         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
169         for (i = 0, f = 0.; i < ma->n; ++i) {
170                 double *pdg;
171                 pdg = ma->pdg + i * 3;
172                 f += (pdg[1] * f3[1] + 2. * pdg[2] * f3[2])
173                         / (pdg[0] * f3[0] + pdg[1] * f3[1] + pdg[2] * f3[2]);
174         }
175         f /= ma->n * 2.;
176         return f;
177 }
178
179 int mc_call_gt(const mc_aux_t *ma, double f0, int k)
180 {
181         double sum, g[3];
182         double max, f3[3], *pdg = ma->pdg + k * 3;
183         int q, i, max_i;
184         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
185         for (i = 0, sum = 0.; i < 3; ++i)
186                 sum += (g[i] = pdg[i] * f3[i]);
187         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
188                 g[i] /= sum;
189                 if (g[i] > max) max = g[i], max_i = i;
190         }
191         max = 1. - max;
192         if (max < 1e-308) max = 1e-308;
193         q = (int)(-3.434 * log(max) + .499);
194         if (q > 99) q = 99;
195         return q<<2|max_i;
196 }
197
198 static void mc_cal_z(mc_aux_t *ma)
199 {
200         double *z[2], *tmp, *pdg;
201         int i, j;
202         z[0] = ma->z;
203         z[1] = ma->zswap;
204         pdg = ma->pdg;
205         z[0][0] = 1.; z[0][1] = z[0][2] = 0.;
206         for (j = 0; j < ma->n; ++j) {
207                 int max = (j + 1) * 2;
208                 double p[3];
209                 pdg = ma->pdg + j * 3;
210                 p[0] = pdg[0]; p[1] = 2. * pdg[1]; p[2] = pdg[2];
211                 z[1][0] = p[0] * z[0][0];
212                 z[1][1] = p[0] * z[0][1] + p[1] * z[0][0];
213                 for (i = 2; i <= max; ++i)
214                         z[1][i] = p[0] * z[0][i] + p[1] * z[0][i-1] + p[2] * z[0][i-2];
215                 if (j < ma->n - 1) z[1][max+1] = z[1][max+2] = 0.;
216 //              int k; for (k = 0; k <= max; ++k) printf("%d:%.3lg ", k, z[1][k]); putchar('\n');
217                 tmp = z[0]; z[0] = z[1]; z[1] = tmp;
218         }
219         if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (2 * ma->n + 1));
220 }
221
222 static double mc_add_afs(mc_aux_t *ma)
223 {
224         int k;
225         long double sum = 0.;
226         memset(ma->afs1, 0, sizeof(double) * (ma->M + 1));
227         mc_cal_z(ma);
228         for (k = 0, sum = 0.; k <= ma->M; ++k)
229                 sum += (long double)ma->phi[k] * ma->z[k] / ma->CMk[k];
230         for (k = 0; k <= ma->M; ++k) {
231                 ma->afs1[k] = ma->phi[k] * ma->z[k] / ma->CMk[k] / sum;
232                 if (isnan(ma->afs1[k]) || isinf(ma->afs1[k])) return -1.;
233         }
234         for (k = 0, sum = 0.; k <= ma->M; ++k) {
235                 ma->afs[k] += ma->afs1[k];
236                 sum += k * ma->afs1[k];
237         }
238         return sum / ma->M;
239 }
240
241 int mc_cal(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, mc_rst_t *rst, int level)
242 {
243         int i, tot;
244         memset(rst, 0, sizeof(mc_rst_t));
245         rst->f_em = rst->f_exp = -1.; rst->ref = rst->alt = -1;
246         // precalculation
247         tot = sum_err(n, plp, ma);
248         if (tot == 0) return 0; // no good bases
249         set_allele(ref, ma);
250         cal_pdg(ma);
251         // set ref/major allele
252         rst->ref = ma->ref; rst->alt = ma->alt; rst->alt2 = ma->alt2;
253         // calculate naive and Nielsen's freq
254         rst->f_naive = mc_freq0(ma, &rst->f_nielsen);
255         { // calculate f_em
256                 double flast = rst->f_naive;
257                 for (i = 0; i < MC_MAX_EM_ITER; ++i) {
258                         rst->f_em = mc_freq_iter(flast, ma);
259                         if (fabs(rst->f_em - flast) < MC_EM_EPS) break;
260                         flast = rst->f_em;
261                 }
262         }
263         if (level >= 2) {
264                 rst->f_exp = mc_add_afs(ma);
265                 rst->p_ref = ma->afs1[ma->M];
266         }
267         return tot;
268 }
269
270 void mc_dump_afs(mc_aux_t *ma)
271 {
272         int k;
273         fprintf(stderr, "[afs]");
274         for (k = 0; k <= ma->M; ++k)
275                 fprintf(stderr, " %d:%.3lf", k, ma->afs[ma->M - k]);
276         fprintf(stderr, "\n");
277         memset(ma->afs, 0, sizeof(double) * (ma->M + 1));
278 }