]> git.donarmstrong.com Git - samtools.git/blob - bam_mcns.c
supposedly this is THE correct implementation, but more testing is needed
[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;
18         double *z, *zswap; // aux for afs
19         double *CMk; // \binom{M}{k}
20         double *afs, *afs1; // afs: accumulative AFS; afs1: site posterior distribution
21         int *qsum, *bcnt;
22 };
23
24 void mc_init_prior(mc_aux_t *ma, int type, double theta)
25 {
26         int i;
27         if (type == MC_PTYPE_COND2) {
28                 for (i = 0; i <= 2 * ma->n; ++i)
29                         ma->phi[i] = 2. * (i + 1) / (2 * ma->n + 1) / (2 * ma->n + 2);
30         } else if (type == MC_PTYPE_FLAT) {
31                 for (i = 0; i <= ma->M; ++i)
32                         ma->phi[i] = 1. / (ma->M + 1);
33         } else {
34                 double sum;
35                 for (i = 0, sum = 0.; i < 2 * ma->n; ++i)
36                         sum += (ma->phi[i] = theta / (2 * ma->n - i));
37                 ma->phi[2 * ma->n] = 1. - sum;
38         }
39 }
40
41 mc_aux_t *mc_init(int n) // FIXME: assuming diploid
42 {
43         mc_aux_t *ma;
44         int i;
45         ma = calloc(1, sizeof(mc_aux_t));
46         ma->n = n; ma->M = 2 * n;
47         ma->q2p = calloc(MC_MAX_SUMQ + 1, sizeof(double));
48         ma->qsum = calloc(4 * ma->n, sizeof(int));
49         ma->bcnt = calloc(4 * ma->n, sizeof(int));
50         ma->pdg = calloc(3 * ma->n, sizeof(double));
51         ma->phi = calloc(2 * ma->n + 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         ma->CMk = calloc(ma->M + 1, sizeof(double));
57         for (i = 0; i <= MC_MAX_SUMQ; ++i)
58                 ma->q2p[i] = pow(10., -i / 10.);
59         for (i = 0; i <= ma->M; ++i)
60                 ma->CMk[i] = exp(lgamma(ma->M+1) - lgamma(ma->M-i+1) - lgamma(i+1));
61         mc_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior
62         return ma;
63 }
64
65 void mc_destroy(mc_aux_t *ma)
66 {
67         if (ma) {
68                 free(ma->qsum); free(ma->bcnt);
69                 free(ma->q2p); free(ma->pdg);
70                 free(ma->phi);
71                 free(ma->z); free(ma->zswap);
72                 free(ma->afs); free(ma->afs1);
73                 free(ma->CMk);
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; ma->alt2 = -1;
115         if (ma->ref != ref) { // the best base is not ref
116                 if (ref >= 0 && ref <= 3) { // ref is not N
117                         if (ma->alt == ref) tmp = ma->ref, ma->ref = ma->alt, ma->alt = tmp; // then switch alt and ref
118                         else ma->alt2 = ma->alt, ma->alt = ma->ref, ma->ref = ref; // then set ref as ref
119                 } else ma->alt2 = ma->alt, ma->alt = ma->ref, ma->ref = sum[0]&3; // then set the weakest as ref
120         }
121 }
122
123 static void cal_pdg(mc_aux_t *ma)
124 {
125         int i, j;
126         for (j = 0; j < ma->n; ++j) {
127                 int pi[3], *qsum, *bcnt;
128                 double *pdg = ma->pdg + j * 3;
129                 qsum = ma->qsum + j * 4;
130                 bcnt = ma->bcnt + j * 4;
131                 pi[1] = 3 * (bcnt[ma->ref] + bcnt[ma->alt]);
132                 pi[0] = qsum[ma->ref];
133                 pi[2] = qsum[ma->alt];
134                 for (i = 0; i < 3; ++i)
135                         pdg[i] = pi[i] > MC_MAX_SUMQ? MC_MAX_SUMQP : ma->q2p[pi[i]];
136         }
137 }
138 // this calculates the naive allele frequency and Nielsen's frequency
139 static double mc_freq0(const mc_aux_t *ma, double *_f)
140 {
141         int i, cnt;
142         double f, f_nielsen, w_sum;
143         *_f = -1.;
144         for (i = cnt = 0, f = f_nielsen = w_sum = 0.; i < ma->n; ++i) {
145                 int *bcnt = ma->bcnt + i * 4;
146                 int x = bcnt[ma->ref] + bcnt[ma->alt];
147                 if (x) {
148                         double w, p;
149                         ++cnt;
150                         f += (double)bcnt[ma->ref] / x;
151                         p = (bcnt[ma->ref] - MC_AVG_ERR * x) / (1. - 2. * MC_AVG_ERR) / x;
152                         w = 2. * x / (1. + x);
153                         w_sum += w;
154                         f_nielsen += p * w;
155                 }
156         }
157         if (cnt) {
158                 f_nielsen /= w_sum;
159                 if (f_nielsen < 0.) f_nielsen = 0.;
160                 if (f_nielsen > 1.) f_nielsen = 1.;
161                 *_f = f_nielsen;
162                 return f / cnt;
163         } else return -1.;
164 }
165 // f0 is the reference allele frequency
166 static double mc_freq_iter(double f0, const mc_aux_t *ma)
167 {
168         double f, f3[3];
169         int i;
170         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
171         for (i = 0, f = 0.; i < ma->n; ++i) {
172                 double *pdg;
173                 pdg = ma->pdg + i * 3;
174                 f += (pdg[1] * f3[1] + 2. * pdg[2] * f3[2])
175                         / (pdg[0] * f3[0] + pdg[1] * f3[1] + pdg[2] * f3[2]);
176         }
177         f /= ma->n * 2.;
178         return f;
179 }
180
181 int mc_call_gt(const mc_aux_t *ma, double f0, int k)
182 {
183         double sum, g[3];
184         double max, f3[3], *pdg = ma->pdg + k * 3;
185         int q, i, max_i;
186         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
187         for (i = 0, sum = 0.; i < 3; ++i)
188                 sum += (g[i] = pdg[i] * f3[i]);
189         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
190                 g[i] /= sum;
191                 if (g[i] > max) max = g[i], max_i = i;
192         }
193         max = 1. - max;
194         if (max < 1e-308) max = 1e-308;
195         q = (int)(-3.434 * log(max) + .499);
196         if (q > 99) q = 99;
197         return q<<2|max_i;
198 }
199
200 static void mc_cal_z(mc_aux_t *ma)
201 {
202         double *z[2], *tmp, *pdg;
203         int i, j;
204         z[0] = ma->z;
205         z[1] = ma->zswap;
206         pdg = ma->pdg;
207         z[0][0] = 1.; z[0][1] = z[0][2] = 0.;
208         for (j = 0; j < ma->n; ++j) {
209                 int max = (j + 1) * 2;
210                 double p[3];
211                 pdg = ma->pdg + j * 3;
212                 p[0] = pdg[0]; p[1] = 2. * pdg[1]; p[2] = pdg[2];
213                 z[1][0] = p[0] * z[0][0];
214                 z[1][1] = p[0] * z[0][1] + p[1] * z[0][0];
215                 for (i = 2; i <= max; ++i)
216                         z[1][i] = p[0] * z[0][i] + p[1] * z[0][i-1] + p[2] * z[0][i-2];
217                 if (j < ma->n - 1) z[1][max+1] = z[1][max+2] = 0.;
218                 tmp = z[0]; z[0] = z[1]; z[1] = tmp;
219         }
220         if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (2 * ma->n + 1));
221 }
222
223 static double mc_add_afs(mc_aux_t *ma)
224 {
225         int k, l;
226         double sum = 0., avg = 0.;
227         memset(ma->afs1, 0, sizeof(double) * (ma->M + 1));
228         mc_cal_z(ma);
229         for (k = 0; k <= ma->M; ++k) {
230                 for (l = 0, sum = 0.; l <= ma->M; ++l)
231                         sum += ma->phi[l] * ma->z[l];
232                 ma->afs1[k] = ma->phi[k] * ma->z[k] / sum;
233         }
234         for (k = 0; k <= ma->M; ++k)
235                 if (isnan(ma->afs1[k]) || isinf(ma->afs1[k])) return -1.;
236         for (k = 0, sum = avg = 0.; k <= ma->M; ++k) {
237                 ma->afs[k] += ma->afs1[k];
238                 sum += ma->afs1[k];
239                 avg += k * ma->afs1[k];
240         }
241 //      for (k = 0; k <= ma->M; ++k) printf("^%lg:%lg:%lg ", ma->z[k], ma->phi[k], ma->afs1[k]);
242         return avg / ma->M;
243 }
244
245 int mc_cal(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, mc_rst_t *rst, int level)
246 {
247         int i, tot;
248         memset(rst, 0, sizeof(mc_rst_t));
249         rst->f_em = rst->f_exp = -1.; rst->ref = rst->alt = -1;
250         // precalculation
251         tot = sum_err(n, plp, ma);
252         if (tot == 0) return 0; // no good bases
253         set_allele(ref, ma);
254         cal_pdg(ma);
255         // set ref/major allele
256         rst->ref = ma->ref; rst->alt = ma->alt; rst->alt2 = ma->alt2;
257         // calculate naive and Nielsen's freq
258         rst->f_naive = mc_freq0(ma, &rst->f_nielsen);
259         { // calculate f_em
260                 double flast = rst->f_naive;
261                 for (i = 0; i < MC_MAX_EM_ITER; ++i) {
262                         rst->f_em = mc_freq_iter(flast, ma);
263                         if (fabs(rst->f_em - flast) < MC_EM_EPS) break;
264                         flast = rst->f_em;
265                 }
266         }
267         if (level >= 2) {
268                 rst->f_exp = mc_add_afs(ma);
269                 rst->p_ref = ma->afs1[ma->M];
270         }
271         return tot;
272 }
273
274 void mc_dump_afs(mc_aux_t *ma)
275 {
276         int k;
277         fprintf(stderr, "[afs]");
278         for (k = 0; k <= ma->M; ++k)
279                 fprintf(stderr, " %d:%.3lf", k, ma->afs[ma->M - k]);
280         fprintf(stderr, "\n");
281         memset(ma->afs, 0, sizeof(double) * (ma->M + 1));
282 }