]> git.donarmstrong.com Git - samtools.git/blob - bam2bcf.c
for backup. no effective changes
[samtools.git] / bam2bcf.c
1 #include <math.h>
2 #include <stdint.h>
3 #include "bam.h"
4 #include "kstring.h"
5 #include "bam2bcf.h"
6 #include "bcftools/bcf.h"
7
8 extern  void ks_introsort_uint32_t(size_t n, uint32_t a[]);
9
10 #define CALL_ETA 0.03f
11 #define CALL_MAX 256
12 #define CALL_DEFTHETA 0.83f
13
14 #define CAP_DIST 25
15
16 bcf_callaux_t *bcf_call_init(double theta, int min_baseQ)
17 {
18         bcf_callaux_t *bca;
19         if (theta <= 0.) theta = CALL_DEFTHETA;
20         bca = calloc(1, sizeof(bcf_callaux_t));
21         bca->capQ = 60;
22         bca->openQ = 40;
23         bca->extQ = 20;
24         bca->tandemQ = 100;
25         bca->min_baseQ = min_baseQ;
26         bca->e = errmod_init(1. - theta);
27         return bca;
28 }
29
30 void bcf_call_destroy(bcf_callaux_t *bca)
31 {
32         if (bca == 0) return;
33         errmod_destroy(bca->e);
34         free(bca->bases); free(bca);
35 }
36
37 int bcf_call_glfgen(int _n, const bam_pileup1_t *pl, int ref_base /*4-bit*/, bcf_callaux_t *bca, bcf_callret1_t *r)
38 {
39         int i, n, ref4;
40         memset(r, 0, sizeof(bcf_callret1_t));
41         ref4 = bam_nt16_nt4_table[ref_base];
42         if (_n == 0) return -1;
43
44         // enlarge the bases array if necessary
45         if (bca->max_bases < _n) {
46                 bca->max_bases = _n;
47                 kroundup32(bca->max_bases);
48                 bca->bases = (uint16_t*)realloc(bca->bases, 2 * bca->max_bases);
49         }
50         // fill the bases array
51         memset(r, 0, sizeof(bcf_callret1_t));
52         for (i = n = 0; i < _n; ++i) {
53                 const bam_pileup1_t *p = pl + i;
54                 int q, b, mapQ, baseQ, is_diff, min_dist;
55                 // set base
56                 if (p->is_del || p->is_refskip || (p->b->core.flag&BAM_FUNMAP)) continue;
57                 baseQ = q = (int)bam1_qual(p->b)[p->qpos]; // base quality
58                 if (q < bca->min_baseQ) continue;
59                 mapQ = p->b->core.qual < bca->capQ? p->b->core.qual : bca->capQ;
60                 if (q > mapQ) q = mapQ;
61                 if (q > 63) q = 63;
62                 if (q < 4) q = 4;
63                 b = bam1_seqi(bam1_seq(p->b), p->qpos); // base
64                 b = bam_nt16_nt4_table[b? b : ref_base]; // b is the 2-bit base
65                 bca->bases[n++] = q<<5 | (int)bam1_strand(p->b)<<4 | b;
66                 // collect annotations
67                 r->qsum[b] += q;
68                 is_diff = (ref4 < 4 && b == ref4)? 0 : 1;
69                 ++r->anno[0<<2|is_diff<<1|bam1_strand(p->b)];
70                 min_dist = p->b->core.l_qseq - 1 - p->qpos;
71                 if (min_dist > p->qpos) min_dist = p->qpos;
72                 if (min_dist > CAP_DIST) min_dist = CAP_DIST;
73                 r->anno[1<<2|is_diff<<1|0] += baseQ;
74                 r->anno[1<<2|is_diff<<1|1] += baseQ * baseQ;
75                 r->anno[2<<2|is_diff<<1|0] += mapQ;
76                 r->anno[2<<2|is_diff<<1|1] += mapQ * mapQ;
77                 r->anno[3<<2|is_diff<<1|0] += min_dist;
78                 r->anno[3<<2|is_diff<<1|1] += min_dist * min_dist;
79         }
80         r->depth = n;
81         // glfgen
82         errmod_cal(bca->e, n, 5, bca->bases, r->p);
83         return r->depth;
84 }
85
86 int bcf_call_glfgen_gap(int pos, int _n, const bam_pileup1_t *pl, bcf_callaux_t *bca, bcf_callret1_t *r)
87 {
88         int i, n, n_ins, n_del;
89         memset(r, 0, sizeof(bcf_callret1_t));
90         if (_n == 0) return -1;
91
92         // enlarge the bases array if necessary
93         if (bca->max_bases < _n) {
94                 bca->max_bases = _n;
95                 kroundup32(bca->max_bases);
96                 bca->bases = (uint16_t*)realloc(bca->bases, 2 * bca->max_bases);
97         }
98         // fill the bases array
99         memset(r, 0, sizeof(bcf_callret1_t));
100         r->indelreg = 10000;
101         for (i = n = 0; i < _n; ++i) {
102                 const bam_pileup1_t *p = pl + i;
103                 int q, b, mapQ, indelQ, is_diff, min_dist;
104                 if (p->is_del || p->is_refskip || (p->b->core.flag&BAM_FUNMAP)) continue;
105                 { // compute indel (base) quality
106                         // this can be made more efficient, but realignment is the bottleneck anyway
107                         int j, k, x, y, op, len = 0, max_left, max_rght, seqQ, indelreg;
108                         bam1_core_t *c = &p->b->core;
109                         uint32_t *cigar = bam1_cigar(p->b);
110                         uint8_t *qual = bam1_qual(p->b);
111                         for (k = y = 0, x = c->pos; k < c->n_cigar && y <= p->qpos; ++k) {
112                                 op = cigar[k]&0xf;
113                                 len = cigar[k]>>4;
114                                 if (op == BAM_CMATCH) {
115                                         if (pos > x && pos < x + len) break;
116                                         x += len; y += len;
117                                 } else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += len;
118                                 else if (op == BAM_CDEL || op == BAM_CREF_SKIP) x += len;
119                         }
120                         if (k == c->n_cigar) continue; // this actually should not happen
121                         max_left = max_rght = 0; indelreg = 0;
122                         if (pos == x + len - 1 && k+2 < c->n_cigar && ((cigar[k+1]&0xf) == BAM_CINS || (cigar[k+1]&0xf) == BAM_CDEL)
123                                 && (cigar[k+2]&0xf) == BAM_CMATCH)
124                         {
125                                 for (j = y; j < y + len; ++j)
126                                         if (max_left < qual[j]) max_left = qual[j];
127                                 if ((cigar[k+1]&0xf) == BAM_CINS) y += cigar[k+1]>>4;
128                                 else x += cigar[k+1]>>4;
129                                 op = cigar[k+2]&0xf; len = cigar[k+2]>>4;
130                                 for (j = y; j < y + len; ++j) {
131                                         if (max_rght < qual[j]) max_rght = qual[j];
132                                         if (qual[j] > BAM2BCF_INDELREG_THRES && indelreg == 0)
133                                                 indelreg = j - y + 1;
134                                 }
135                         } else {
136                                 for (j = y; j <= p->qpos; ++j)
137                                         if (max_left < qual[j]) max_left = qual[j];
138                                 for (j = p->qpos + 1; j < y + len; ++j)
139                                         if (max_rght < qual[j]) max_rght = qual[j];
140                                         
141                         }
142                         indelQ = max_left < max_rght? max_left : max_rght;
143                         // estimate the sequencing error rate
144                         seqQ = bca->openQ;
145                         if (p->indel != 0) seqQ += bca->extQ * (abs(p->indel) - 1); // FIXME: better to model homopolymer
146                         if (p->indel != 0) { // a different model for tandem repeats
147                                 uint8_t *seq = bam1_seq(p->b);
148                                 int tandemQ, qb = bam1_seqi(seq, p->qpos), l;
149                                 for (j = p->qpos + 1; j < c->l_qseq; ++j)
150                                         if (qb != bam1_seqi(seq, j)) break;
151                                 l = j;
152                                 for (j = (int)p->qpos - 1; j >= 0; --j)
153                                         if (qb != bam1_seqi(seq, j)) break;
154                                 l = l - (j + 1);
155                                 tandemQ = (int)((double)(abs(p->indel)) / l * bca->tandemQ + .499);
156                                 if (seqQ > tandemQ) seqQ = tandemQ;
157                         }
158 //                      fprintf(stderr, "%s\t%d\t%d\t%d\t%d\t%d\t%d\n", bam1_qname(p->b), pos+1, p->indel, indelQ, seqQ, max_left, max_rght);
159                         if (indelQ > seqQ) indelQ = seqQ;
160                         q = indelQ;
161                 }
162                 if (q < bca->min_baseQ) continue;
163                 mapQ = p->b->core.qual < bca->capQ? p->b->core.qual : bca->capQ;
164                 if (q > mapQ) q = mapQ;
165                 if (q > 63) q = 63;
166                 if (q < 4) q = 4;
167                 b = p->indel? 1 : 0;
168                 bca->bases[n++] = q<<5 | (int)bam1_strand(p->b)<<4 | b;
169                 // collect annotations
170                 r->qsum[b] += q;
171                 is_diff = b;
172                 ++r->anno[0<<2|is_diff<<1|bam1_strand(p->b)];
173                 min_dist = p->b->core.l_qseq - 1 - p->qpos;
174                 if (min_dist > p->qpos) min_dist = p->qpos;
175                 if (min_dist > CAP_DIST) min_dist = CAP_DIST;
176                 r->anno[1<<2|is_diff<<1|0] += indelQ;
177                 r->anno[1<<2|is_diff<<1|1] += indelQ * indelQ;
178                 r->anno[2<<2|is_diff<<1|0] += mapQ;
179                 r->anno[2<<2|is_diff<<1|1] += mapQ * mapQ;
180                 r->anno[3<<2|is_diff<<1|0] += min_dist;
181                 r->anno[3<<2|is_diff<<1|1] += min_dist * min_dist;
182         }
183         r->depth = n;
184         // glfgen
185         errmod_cal(bca->e, n, 2, bca->bases, r->p);
186         return r->depth;        
187 }
188
189 int bcf_call_combine(int n, const bcf_callret1_t *calls, int ref_base /*4-bit*/, bcf_call_t *call)
190 {
191         int ref4, i, j, qsum[4];
192         int64_t tmp;
193         call->ori_ref = ref4 = bam_nt16_nt4_table[ref_base];
194         if (ref4 > 4) ref4 = 4;
195         // calculate qsum
196         memset(qsum, 0, 4 * sizeof(int));
197         for (i = 0; i < n; ++i)
198                 for (j = 0; j < 4; ++j)
199                         qsum[j] += calls[i].qsum[j];
200         for (j = 0; j < 4; ++j) qsum[j] = qsum[j] << 2 | j;
201         // find the top 2 alleles
202         for (i = 1; i < 4; ++i) // insertion sort
203                 for (j = i; j > 0 && qsum[j] < qsum[j-1]; --j)
204                         tmp = qsum[j], qsum[j] = qsum[j-1], qsum[j-1] = tmp;
205         // set the reference allele and alternative allele(s)
206         for (i = 0; i < 5; ++i) call->a[i] = -1;
207         call->unseen = -1;
208         call->a[0] = ref4;
209         for (i = 3, j = 1; i >= 0; --i) {
210                 if ((qsum[i]&3) != ref4) {
211                         if (qsum[i]>>2 != 0) call->a[j++] = qsum[i]&3;
212                         else break;
213                 }
214         }
215         if (((ref4 < 4 && j < 4) || (ref4 == 4 && j < 5)) && i >= 0)
216                 call->unseen = j, call->a[j++] = qsum[i]&3;
217         call->n_alleles = j;
218         // set the PL array
219         if (call->n < n) {
220                 call->n = n;
221                 call->PL = realloc(call->PL, 15 * n);
222         }
223         {
224                 int x, g[15], z;
225                 double sum_min = 0.;
226                 x = call->n_alleles * (call->n_alleles + 1) / 2;
227                 // get the possible genotypes
228                 for (i = z = 0; i < call->n_alleles; ++i)
229                         for (j = i; j < call->n_alleles; ++j)
230                                 g[z++] = call->a[i] * 5 + call->a[j];
231                 for (i = 0; i < n; ++i) {
232                         uint8_t *PL = call->PL + x * i;
233                         const bcf_callret1_t *r = calls + i;
234                         float min = 1e37;
235                         for (j = 0; j < x; ++j)
236                                 if (min > r->p[g[j]]) min = r->p[g[j]];
237                         sum_min += min;
238                         for (j = 0; j < x; ++j) {
239                                 int y;
240                                 y = (int)(r->p[g[j]] - min + .499);
241                                 if (y > 255) y = 255;
242                                 PL[j] = y;
243                         }
244                 }
245                 call->shift = (int)(sum_min + .499);
246         }
247         // combine annotations
248         memset(call->anno, 0, 16 * sizeof(int));
249         for (i = call->depth = 0, tmp = 0; i < n; ++i) {
250                 call->depth += calls[i].depth;
251                 for (j = 0; j < 16; ++j) call->anno[j] += calls[i].anno[j];
252         }
253         return 0;
254 }
255
256 int bcf_call2bcf(int tid, int pos, bcf_call_t *bc, bcf1_t *b, bcf_callret1_t *bcr, int is_SP)
257 {
258         extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
259         kstring_t s;
260         int i;
261         b->n_smpl = bc->n;
262         b->tid = tid; b->pos = pos; b->qual = 0;
263         s.s = b->str; s.m = b->m_str; s.l = 0;
264         kputc('\0', &s);
265         kputc("ACGTN"[bc->ori_ref], &s); kputc('\0', &s);
266         for (i = 1; i < 5; ++i) {
267                 if (bc->a[i] < 0) break;
268                 if (i > 1) kputc(',', &s);
269                 kputc(bc->unseen == i? 'X' : "ACGT"[bc->a[i]], &s);
270         }
271         kputc('\0', &s);
272         kputc('\0', &s);
273         // INFO
274         kputs("I16=", &s);
275         for (i = 0; i < 16; ++i) {
276                 if (i) kputc(',', &s);
277                 kputw(bc->anno[i], &s);
278         }
279         kputc('\0', &s);
280         // FMT
281         kputs("PL", &s);
282         if (bcr) {
283                 kputs(":DP", &s);
284                 if (is_SP) kputs(":SP", &s);
285         }
286         kputc('\0', &s);
287         b->m_str = s.m; b->str = s.s; b->l_str = s.l;
288         bcf_sync(b);
289         memcpy(b->gi[0].data, bc->PL, b->gi[0].len * bc->n);
290         if (bcr) {
291                 uint16_t *dp = (uint16_t*)b->gi[1].data;
292                 uint8_t *sp = is_SP? b->gi[2].data : 0;
293                 for (i = 0; i < bc->n; ++i) {
294                         bcf_callret1_t *p = bcr + i;
295                         dp[i] = p->depth < 0xffff? p->depth : 0xffff;
296                         if (is_SP) {
297                                 if (p->anno[0] + p->anno[1] < 2 || p->anno[2] + p->anno[3] < 2
298                                         || p->anno[0] + p->anno[2] < 2 || p->anno[1] + p->anno[3] < 2)
299                                 {
300                                         sp[i] = 0;
301                                 } else {
302                                         double left, right, two;
303                                         int x;
304                                         kt_fisher_exact(p->anno[0], p->anno[1], p->anno[2], p->anno[3], &left, &right, &two);
305                                         x = (int)(-4.343 * log(two) + .499);
306                                         if (x > 255) x = 255;
307                                         sp[i] = x;
308                                 }
309                         }
310                 }
311         }
312         return 0;
313 }