]> git.donarmstrong.com Git - samtools.git/blob - bam_mcns.c
cd64b0ac139814706632a05d2e7b75d619a36b8a
[samtools.git] / bam_mcns.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include "bam_mcns.h"
4
5 #define MC_MIN_QUAL 13
6 #define MC_MAX_SUMQ 3000
7 #define MC_MAX_SUMQP 1e-300
8
9 struct __mc_aux_t {
10         int n, N;
11         int ref, alt;
12         double theta;
13         double *q2p, *pdg; // pdg -> P(D|g)
14         double *alpha, *beta;
15         int *qsum, *bcnt, rcnt[4];
16 };
17
18 void mc_init_prior(mc_aux_t *ma, double theta)
19 {
20         double sum;
21         int i;
22         ma->theta = theta;
23         for (i = 0, sum = 0.; i < 2 * ma->n; ++i)
24                 sum += (ma->alpha[i] = ma->theta / (2 * ma->n - i));
25         ma->alpha[2 * ma->n] = 1. - sum;
26 }
27
28 mc_aux_t *mc_init(int n) // FIXME: assuming diploid
29 {
30         mc_aux_t *ma;
31         int i;
32         ma = calloc(1, sizeof(mc_aux_t));
33         ma->n = n; ma->N = 2 * n;
34         ma->q2p = calloc(MC_MAX_SUMQ + 1, sizeof(double));
35         ma->qsum = calloc(4 * ma->n, sizeof(int));
36         ma->bcnt = calloc(4 * ma->n, sizeof(int));
37         ma->pdg = calloc(3 * ma->n, sizeof(double));
38         ma->alpha = calloc(2 * ma->n + 1, sizeof(double));
39         ma->beta = calloc((2 * ma->n + 1) * 3, sizeof(double));
40         for (i = 0; i <= MC_MAX_SUMQ; ++i)
41                 ma->q2p[i] = pow(10., -i / 10.);
42         for (i = 0; i <= 2 * ma->n; ++i) {
43                 double *bi = ma->beta + 3 * i;
44                 double f = (double)i / (2 * ma->n);
45                 bi[0] = (1. - f) * (1. - f);
46                 bi[1] = 2 * f * (1. - f);
47                 bi[2] = f * f;
48         }
49         mc_init_prior(ma, 1e-3); // the simplest prior
50         return ma;
51 }
52
53 void mc_destroy(mc_aux_t *ma)
54 {
55         if (ma) {
56                 free(ma->qsum); free(ma->bcnt);
57                 free(ma->q2p); free(ma->pdg);
58                 free(ma->alpha); free(ma->beta);
59                 free(ma);
60         }
61 }
62
63 static void sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
64 {
65         int i, j;
66         memset(ma->qsum, 0, sizeof(int) * 4 * ma->n);
67         memset(ma->bcnt, 0, sizeof(int) * 4 * ma->n);
68         memset(ma->rcnt, 0, sizeof(int) * 4);
69         for (j = 0; j < ma->n; ++j) {
70                 int tmp, *qsum = ma->qsum + j * 4;
71                 int *bcnt = ma->bcnt + j * 4;
72                 double r = drand48(), rr;
73                 for (i = tmp = 0; i < n[j]; ++i) {
74                         const bam_pileup1_t *p = plp[j] + i;
75                         int q, b;
76                         if (p->is_del || (p->b->core.flag&BAM_FUNMAP)) continue;
77                         q = bam1_qual(p->b)[p->qpos];
78                         if (p->b->core.qual < q) q = p->b->core.qual;
79                         if (q < MC_MIN_QUAL) continue; // small qual
80                         b = bam_nt16_nt4_table[(int)bam1_seqi(bam1_seq(p->b), p->qpos)];
81                         if (b > 3) continue; // N
82                         qsum[b] += q;
83                         ++bcnt[b];
84                         ++tmp;
85                 }
86                 if (tmp) {
87                         for (i = 0, rr = 0.; i < 4; ++i) {
88                                 rr += (double)bcnt[i] / tmp;
89                                 if (rr >= r) break;
90                         }
91                         if (i < 4) ++ma->rcnt[i];
92                 }
93         }
94 }
95
96 static void set_allele(int ref, mc_aux_t *ma)
97 {
98         int i, j, sum[4], tmp;
99         sum[0] = sum[1] = sum[2] = sum[3] = 0;
100         for (i = 0; i < ma->n; ++i)
101                 for (j = 0; j < 4; ++j)
102                         sum[j] += ma->qsum[i * 4 + j];
103         for (j = 0; j < 4; ++j) sum[j] = sum[j]<<2 | j;
104         for (i = 1; i < 4; ++i) // insertion sort
105                 for (j = i; j > 0 && sum[j] < sum[j-1]; --j)
106                         tmp = sum[j], sum[j] = sum[j-1], sum[j-1] = tmp;
107         ma->ref = sum[3]&3; ma->alt = sum[2]&3;
108         if (ref == ma->alt) tmp = ma->ref, ma->ref = ma->alt, ma->alt = tmp;
109         // note that ma->ref might not be ref in case of triallele
110 }
111
112 static void cal_pdg(mc_aux_t *ma)
113 {
114         int i, j;
115         for (j = 0; j < ma->n; ++j) {
116                 int pi[3], *qsum, *bcnt;
117                 double *pdg = ma->pdg + j * 3;
118                 qsum = ma->qsum + j * 4;
119                 bcnt = ma->bcnt + j * 4;
120                 pi[1] = 3 * (bcnt[ma->ref] + bcnt[ma->alt]);
121                 pi[0] = qsum[ma->ref];
122                 pi[2] = qsum[ma->alt];
123                 for (i = 0; i < 3; ++i)
124                         pdg[i] = pi[i] > MC_MAX_SUMQ? MC_MAX_SUMQP : ma->q2p[pi[i]];
125         }
126 }
127 // return the reference allele frequency
128 double mc_freq0(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, int *_ref, int *alt)
129 {
130         int i, acnt[4], j;
131         double fr = -1., f0 = -1.;
132         sum_err(n, plp, ma);
133         set_allele(ref, ma);
134         cal_pdg(ma);
135         acnt[0] = acnt[1] = acnt[2] = acnt[3] = 0;
136         for (i = 0; i < ma->n; ++i)
137                 for (j = 0; j < 4; ++j)
138                         acnt[j] += ma->bcnt[i * 4 + j];
139         if (ma->rcnt[ma->ref] + ma->rcnt[ma->alt]) // never used...
140                 fr = (double)ma->rcnt[ref] / (ma->rcnt[ma->ref] + ma->rcnt[ma->alt]);
141         if (acnt[ma->ref] + acnt[ma->alt])
142                 f0 = (double)acnt[ref] / (acnt[ma->ref] + acnt[ma->alt]);
143         *_ref = ma->ref; *alt = ma->alt;
144         return ma->n == 1 || fr < 0.? f0 : fr;
145 }
146 // f0 is the reference allele frequency
147 double mc_freq_iter(double f0, mc_aux_t *ma)
148 {
149         double f, f3[3];
150         int i, j;
151         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
152         for (i = 0, f = 0.; i < ma->n; ++i) {
153                 double up, dn, *pdg;
154                 pdg = ma->pdg + i * 3;
155                 for (j = 1, up = 0.; j < 3; ++j)
156                         up += j * pdg[j] * f3[j];
157                 for (j = 0, dn = 0.; j < 3; ++j)
158                         dn += pdg[j] * f3[j];
159                 f += up / dn;
160         }
161         f /= ma->n * 2.;
162         return f;
163 }
164
165 double mc_ref_prob(mc_aux_t *ma)
166 {
167         int k, i, g;
168         long double PD = 0., Pref = 0.;
169         for (k = 0; k <= ma->n * 2; ++k) {
170                 long double x = 1.;
171                 double *bk = ma->beta + k * 3;
172                 for (i = 0; i < ma->n; ++i) {
173                         double y = 0., *pdg = ma->pdg + i * 3;
174                         for (g = 0; g < 3; ++g)
175                                 y += pdg[g] * bk[g];
176                         x *= y;
177                 }
178                 PD += x * ma->alpha[k];
179         }
180         for (k = 0; k <= ma->n * 2; ++k) {
181                 long double x = 1.0;
182                 for (i = 0; i < ma->n; ++i)
183                         x *= ma->pdg[i * 3 + 2] * ma->beta[k * 3 + 2];
184                 Pref += x * ma->alpha[k];
185         }
186         return Pref / PD;
187 }
188
189 int mc_call_gt(const mc_aux_t *ma, double f0, int k)
190 {
191         double sum, g[3];
192         double max, f3[3], *pdg = ma->pdg + k * 3;
193         int q, i, max_i;
194         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
195         for (i = 0, sum = 0.; i < 3; ++i)
196                 sum += (g[i] = pdg[i] * f3[i]);
197         for (i = 0, max = -1., max_i = 0; i < 3; ++i) {
198                 g[i] /= sum;
199                 if (g[i] > max) max = g[i], max_i = i;
200         }
201         max = 1. - max;
202         if (max < 1e-308) max = 1e-308;
203         q = (int)(-3.434 * log(max) + .499);
204         if (q > 99) q = 99;
205         return q<<2|max_i;
206 }