]> git.donarmstrong.com Git - samtools.git/blob - bam_mcns.c
* samtools-0.1.8-6 (r638)
[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 *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 if (type == MC_PTYPE_FLAT) {
30                 for (i = 0; i <= ma->M; ++i)
31                         ma->alpha[i] = 1. / (ma->M + 1);
32         } else {
33                 double sum;
34                 for (i = 0, sum = 0.; i < 2 * ma->n; ++i)
35                         sum += (ma->alpha[i] = theta / (2 * ma->n - i));
36                 ma->alpha[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->alpha = calloc(2 * ma->n + 1, sizeof(double));
51         ma->beta = calloc((2 * ma->n + 1) * 3, 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) { // beta[k][g]=P(g|k/M)
59                 double *bi = ma->beta + 3 * i;
60                 double f = (double)i / ma->M;
61                 bi[0] = (1. - f) * (1. - f);
62                 bi[1] = 2 * f * (1. - f);
63                 bi[2] = f * f;
64         }
65         mc_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior
66         return ma;
67 }
68
69 void mc_destroy(mc_aux_t *ma)
70 {
71         if (ma) {
72                 free(ma->qsum); free(ma->bcnt);
73                 free(ma->q2p); free(ma->pdg);
74                 free(ma->alpha); free(ma->beta);
75                 free(ma->z); free(ma->zswap);
76                 free(ma->afs); free(ma->afs1);
77                 free(ma);
78         }
79 }
80
81 static int sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
82 {
83         int i, j, tot = 0;
84         memset(ma->qsum, 0, sizeof(int) * 4 * ma->n);
85         memset(ma->bcnt, 0, sizeof(int) * 4 * ma->n);
86         for (j = 0; j < ma->n; ++j) {
87                 int *qsum = ma->qsum + j * 4;
88                 int *bcnt = ma->bcnt + j * 4;
89                 for (i = 0; i < n[j]; ++i) {
90                         const bam_pileup1_t *p = plp[j] + i;
91                         int q, b;
92                         if (p->is_del || (p->b->core.flag&BAM_FUNMAP)) continue;
93                         q = bam1_qual(p->b)[p->qpos];
94                         if (p->b->core.qual < q) q = p->b->core.qual;
95                         if (q < MC_MIN_QUAL) continue; // small qual
96                         b = bam_nt16_nt4_table[(int)bam1_seqi(bam1_seq(p->b), p->qpos)];
97                         if (b > 3) continue; // N
98                         qsum[b] += q;
99                         ++bcnt[b];
100                         ++tot;
101                 }
102         }
103         return tot;
104 }
105
106 static void set_allele(int ref, mc_aux_t *ma)
107 {
108         int i, j, sum[4], tmp;
109         sum[0] = sum[1] = sum[2] = sum[3] = 0;
110         for (i = 0; i < ma->n; ++i)
111                 for (j = 0; j < 4; ++j)
112                         sum[j] += ma->qsum[i * 4 + j];
113         for (j = 0; j < 4; ++j) sum[j] = sum[j]<<2 | j;
114         for (i = 1; i < 4; ++i) // insertion sort
115                 for (j = i; j > 0 && sum[j] < sum[j-1]; --j)
116                         tmp = sum[j], sum[j] = sum[j-1], sum[j-1] = tmp;
117         ma->ref = sum[3]&3; ma->alt = sum[2]&3; ma->alt2 = -1;
118         if (ma->ref != ref) {
119                 if (ma->alt == ref) tmp = ma->ref, ma->ref = ma->alt, ma->alt = tmp;
120                 else ma->alt2 = ma->alt, ma->alt = ma->ref, ma->ref = ref;
121         }
122 }
123
124 static void cal_pdg(mc_aux_t *ma)
125 {
126         int i, j;
127         for (j = 0; j < ma->n; ++j) {
128                 int pi[3], *qsum, *bcnt;
129                 double *pdg = ma->pdg + j * 3;
130                 qsum = ma->qsum + j * 4;
131                 bcnt = ma->bcnt + j * 4;
132                 pi[1] = 3 * (bcnt[ma->ref] + bcnt[ma->alt]);
133                 pi[0] = qsum[ma->ref];
134                 pi[2] = qsum[ma->alt];
135                 for (i = 0; i < 3; ++i)
136                         pdg[i] = pi[i] > MC_MAX_SUMQ? MC_MAX_SUMQP : ma->q2p[pi[i]];
137         }
138 }
139 // this calculates the naive allele frequency and Nielsen's frequency
140 static double mc_freq0(const mc_aux_t *ma, double *_f)
141 {
142         int i, cnt;
143         double f, f_nielsen, w_sum;
144         *_f = -1.;
145         for (i = cnt = 0, f = f_nielsen = w_sum = 0.; i < ma->n; ++i) {
146                 int *bcnt = ma->bcnt + i * 4;
147                 int x = bcnt[ma->ref] + bcnt[ma->alt];
148                 if (x) {
149                         double w, p;
150                         ++cnt;
151                         f += (double)bcnt[ma->ref] / x;
152                         p = (bcnt[ma->ref] - MC_AVG_ERR * x) / (1. - 2. * MC_AVG_ERR) / x;
153                         w = 2. * x / (1. + x);
154                         w_sum += w;
155                         f_nielsen += p * w;
156                 }
157         }
158         if (cnt) {
159                 f_nielsen /= w_sum;
160                 if (f_nielsen < 0.) f_nielsen = 0.;
161                 if (f_nielsen > 1.) f_nielsen = 1.;
162                 *_f = f_nielsen;
163                 return f / cnt;
164         } else return -1.;
165 }
166 // f0 is the reference allele frequency
167 static double mc_freq_iter(double f0, const mc_aux_t *ma)
168 {
169         double f, f3[3];
170         int i;
171         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
172         for (i = 0, f = 0.; i < ma->n; ++i) {
173                 double *pdg;
174                 pdg = ma->pdg + i * 3;
175                 f += (pdg[1] * f3[1] + 2. * pdg[2] * f3[2])
176                         / (pdg[0] * f3[0] + pdg[1] * f3[1] + pdg[2] * f3[2]);
177         }
178         f /= ma->n * 2.;
179         return f;
180 }
181
182 static double mc_ref_prob(const mc_aux_t *ma, double *_PD, double *f_exp)
183 {
184         int k, i;
185         long double PD = 0., Pref = 0., Ef = 0.;
186         for (k = 0; k <= ma->M; ++k) {
187                 long double x = 1., y = 0.;
188                 double *bk = ma->beta + k * 3;
189                 for (i = 0; i < ma->n; ++i) {
190                         double *pdg = ma->pdg + i * 3;
191                         double z = pdg[0] * bk[0] + pdg[1] * bk[1] + pdg[2] * bk[2];
192                         x *= z;
193                         y += (pdg[1] * bk[1] + 2. * pdg[2] * bk[2]) / z;
194                 }
195                 PD += x * ma->alpha[k];
196                 Ef += x * y * ma->alpha[k];
197         }
198         for (k = 0; k <= ma->n * 2; ++k) {
199                 long double x = 1.0;
200                 for (i = 0; i < ma->n; ++i)
201                         x *= ma->pdg[i * 3 + 2] * ma->beta[k * 3 + 2];
202                 Pref += x * ma->alpha[k];
203         }
204         *f_exp = (double)(Ef / PD / ma->M);
205         *_PD = PD;
206         return Pref / PD;
207 }
208
209 int mc_call_gt(const mc_aux_t *ma, double f0, int k)
210 {
211         double sum, g[3];
212         double max, f3[3], *pdg = ma->pdg + k * 3;
213         int q, i, max_i;
214         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
215         for (i = 0, sum = 0.; i < 3; ++i)
216                 sum += (g[i] = pdg[i] * f3[i]);
217         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
218                 g[i] /= sum;
219                 if (g[i] > max) max = g[i], max_i = i;
220         }
221         max = 1. - max;
222         if (max < 1e-308) max = 1e-308;
223         q = (int)(-3.434 * log(max) + .499);
224         if (q > 99) q = 99;
225         return q<<2|max_i;
226 }
227 // calculate z_{nr}^{(k)}
228 static void mc_cal_z(mc_aux_t *ma, int k)
229 {
230         double *z[2], *tmp, *bk, *pdg;
231         int i, j;
232         z[0] = ma->z;
233         z[1] = ma->zswap;
234         bk = ma->beta + k * 3; pdg = ma->pdg;
235         z[0][0] = 1.; z[0][1] = z[0][2] = 0.;
236         for (j = 0; j < ma->n; ++j) {
237                 int max = (j + 1) * 2;
238                 double p[3];
239                 pdg = ma->pdg + j * 3;
240                 p[0] = bk[0] * pdg[0]; p[1] = bk[1] * pdg[1]; p[2] = bk[2] * pdg[2];
241                 z[1][0] = p[0] * z[0][0];
242                 z[1][1] = p[0] * z[0][1] + p[1] * z[0][0];
243                 for (i = 2; i <= max; ++i)
244                         z[1][i] = p[0] * z[0][i] + p[1] * z[0][i-1] + p[2] * z[0][i-2];
245                 if (j < ma->n - 1) z[1][max+1] = z[1][max+2] = 0.;
246                 tmp = z[0]; z[0] = z[1]; z[1] = tmp;
247         }
248         if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (2 * ma->n + 1));
249 }
250 // Warning: this is cubic in ma->n, very sloooooow
251 static void mc_add_afs(mc_aux_t *ma, double PD, double *f_map, double *p_map)
252 {
253         int k, l;
254         double sum = 0.;
255         memset(ma->afs1, 0, sizeof(double) * (2 * ma->n + 1));
256         *f_map = *p_map = -1.;
257         for (k = 0; k <= 2 * ma->n; ++k) {
258                 mc_cal_z(ma, k);
259                 for (l = 0; l <= 2 * ma->n; ++l)
260                         ma->afs1[l] += ma->alpha[k] * ma->z[l] / PD;
261         }
262         for (k = 0; k <= ma->M; ++k)
263                 if (isnan(ma->afs1[k]) || isinf(ma->afs1[k])) return;
264         for (k = 0; k <= 2 * ma->n; ++k) {
265                 ma->afs[k] += ma->afs1[k];
266                 sum += ma->afs1[k];
267         }
268         {
269                 int max_k = 0;
270                 double max = -1., e = 0.;
271                 for (k = 0; k <= 2 * ma->n; ++k) {
272                         if (ma->afs1[k] > max) max = ma->afs1[k], max_k = k;
273                         e += k * ma->afs1[k];
274                 }
275                 *f_map = .5 * max_k / ma->n; *p_map = max; // e should equal mc_rst_t::f_exp
276 //              printf(" * %.3lg:%.3lg:%.3lg:%.3lg * ", sum, 1.-.5*max_k/ma->n, max, 1.-.5*e/ma->n);
277         }
278 }
279
280 int mc_cal(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, mc_rst_t *rst, int level)
281 {
282         int i, tot;
283         memset(rst, 0, sizeof(mc_rst_t));
284         rst->f_em = rst->f_exp = -1.; rst->ref = rst->alt = -1;
285         // precalculation
286         tot = sum_err(n, plp, ma);
287         if (tot == 0) return 0; // no good bases
288         set_allele(ref, ma);
289         cal_pdg(ma);
290         // set ref/major allele
291         rst->ref = ma->ref; rst->alt = ma->alt; rst->alt2 = ma->alt2;
292         // calculate naive and Nielsen's freq
293         rst->f_naive = mc_freq0(ma, &rst->f_nielsen);
294         { // calculate f_em
295                 double flast = rst->f_naive;
296                 for (i = 0; i < MC_MAX_EM_ITER; ++i) {
297                         rst->f_em = mc_freq_iter(flast, ma);
298                         if (fabs(rst->f_em - flast) < MC_EM_EPS) break;
299                         flast = rst->f_em;
300                 }
301         }
302         if (level >= 2) // quadratic-time calculations; necessary for genotyping
303                 rst->p_ref = mc_ref_prob(ma, &rst->PD, &rst->f_exp);
304         if (level >= 3)
305                 mc_add_afs(ma, rst->PD, &rst->f_map, &rst->p_map);
306         return tot;
307 }
308
309 void mc_dump_afs(mc_aux_t *ma)
310 {
311         int k;
312         fprintf(stderr, "[afs]");
313         for (k = 0; k <= ma->M; ++k)
314                 fprintf(stderr, " %d:%.3lf", k, ma->afs[ma->M - k]);
315         fprintf(stderr, "\n");
316         memset(ma->afs, 0, sizeof(double) * (ma->M + 1));
317 }