]> git.donarmstrong.com Git - samtools.git/blob - bam_mcns.c
* samtools-0.1.8-2 (r617)
[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 *q2p, *pdg; // pdg -> P(D|g)
13         int *qsum, *bcnt, rcnt[4];
14 };
15
16 mc_aux_t *mc_init(int n) // FIXME: assuming diploid
17 {
18         mc_aux_t *ma;
19         int i;
20         ma = calloc(1, sizeof(mc_aux_t));
21         ma->n = n; ma->N = 2 * n;
22         ma->q2p = calloc(MC_MAX_SUMQ + 1, sizeof(double));
23         ma->qsum = calloc(4 * ma->n, sizeof(int));
24         ma->bcnt = calloc(4 * ma->n, sizeof(int));
25         ma->pdg = calloc(3 * ma->n, sizeof(double));
26         for (i = 0; i <= MC_MAX_SUMQ; ++i)
27                 ma->q2p[i] = pow(10., -i / 10.);
28         return ma;
29 }
30
31 void mc_destroy(mc_aux_t *ma)
32 {
33         if (ma) {
34                 free(ma->qsum); free(ma->bcnt);
35                 free(ma->q2p); free(ma->pdg);
36                 free(ma);
37         }
38 }
39
40 static void sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
41 {
42         int i, j;
43         memset(ma->qsum, 0, sizeof(int) * 4 * ma->n);
44         memset(ma->bcnt, 0, sizeof(int) * 4 * ma->n);
45         memset(ma->rcnt, 0, sizeof(int) * 4);
46         for (j = 0; j < ma->n; ++j) {
47                 int tmp, *qsum = ma->qsum + j * 4;
48                 int *bcnt = ma->bcnt + j * 4;
49                 double r = drand48(), rr;
50                 for (i = tmp = 0; i < n[j]; ++i) {
51                         const bam_pileup1_t *p = plp[j] + i;
52                         int q, b;
53                         if (p->is_del || (p->b->core.flag&BAM_FUNMAP)) continue;
54                         q = bam1_qual(p->b)[p->qpos];
55                         if (p->b->core.qual < q) q = p->b->core.qual;
56                         if (q < MC_MIN_QUAL) continue; // small qual
57                         b = bam_nt16_nt4_table[(int)bam1_seqi(bam1_seq(p->b), p->qpos)];
58                         if (b > 3) continue; // N
59                         qsum[b] += q;
60                         ++bcnt[b];
61                         ++tmp;
62                 }
63                 if (tmp) {
64                         for (i = 0, rr = 0.; i < 4; ++i) {
65                                 rr += (double)bcnt[i] / tmp;
66                                 if (rr >= r) break;
67                         }
68                         if (i < 4) ++ma->rcnt[i];
69                 }
70         }
71 }
72
73 static void set_allele(int ref, mc_aux_t *ma)
74 {
75         int i, j, sum[4], tmp;
76         sum[0] = sum[1] = sum[2] = sum[3] = 0;
77         for (i = 0; i < ma->n; ++i)
78                 for (j = 0; j < 4; ++j)
79                         sum[j] += ma->qsum[i * 4 + j];
80         for (j = 0; j < 4; ++j) sum[j] = sum[j]<<2 | j;
81         for (i = 1; i < 4; ++i) // insertion sort
82                 for (j = i; j > 0 && sum[j] < sum[j-1]; --j)
83                         tmp = sum[j], sum[j] = sum[j-1], sum[j-1] = tmp;
84         ma->ref = sum[3]&3; ma->alt = sum[2]&3;
85         if (ref == ma->alt) tmp = ma->ref, ma->ref = ma->alt, ma->alt = tmp;
86         // note that ma->ref might not be ref in case of triallele
87 }
88
89 static void cal_pdg(mc_aux_t *ma)
90 {
91         int i, j;
92         for (j = 0; j < ma->n; ++j) {
93                 int pi[3], *qsum, *bcnt;
94                 double *pdg = ma->pdg + j * 3;
95                 qsum = ma->qsum + j * 4;
96                 bcnt = ma->bcnt + j * 4;
97                 pi[1] = 3 * (bcnt[ma->ref] + bcnt[ma->alt]);
98                 pi[0] = qsum[ma->ref];
99                 pi[2] = qsum[ma->alt];
100                 for (i = 0; i < 3; ++i)
101                         pdg[i] = pi[i] > MC_MAX_SUMQ? MC_MAX_SUMQP : ma->q2p[pi[i]];
102         }
103 }
104
105 double mc_freq0(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, int *_ref, int *alt)
106 {
107         int i, acnt[4], j;
108         double fr = -1., f0 = -1.;
109         sum_err(n, plp, ma);
110         set_allele(ref, ma);
111         cal_pdg(ma);
112         acnt[0] = acnt[1] = acnt[2] = acnt[3] = 0;
113         for (i = 0; i < ma->n; ++i)
114                 for (j = 0; j < 4; ++j)
115                         acnt[j] += ma->bcnt[i * 4 + j];
116         if (ma->rcnt[ma->ref] + ma->rcnt[ma->alt]) // never used...
117                 fr = (double)ma->rcnt[ref] / (ma->rcnt[ma->ref] + ma->rcnt[ma->alt]);
118         if (acnt[ma->ref] + acnt[ma->alt])
119                 f0 = (double)acnt[ref] / (acnt[ma->ref] + acnt[ma->alt]);
120         *_ref = ma->ref; *alt = ma->alt;
121         return ma->n == 1 || fr < 0.? f0 : fr;
122 }
123
124 double mc_freq_iter(double f0, mc_aux_t *ma)
125 {
126         double f, f3[3];
127         int i, j;
128         f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
129         for (i = 0, f = 0.; i < ma->n; ++i) {
130                 double up, dn, *pdg;
131                 pdg = ma->pdg + i * 3;
132                 for (j = 1, up = 0.; j < 3; ++j)
133                         up += j * pdg[j] * f3[j];
134                 for (j = 0, dn = 0.; j < 3; ++j)
135                         dn += pdg[j] * f3[j];
136                 f += up / dn;
137         }
138         f /= ma->n * 2.;
139         return f;
140 }