]> git.donarmstrong.com Git - samtools.git/blob - bam2bcf_indel.c
* samtools-0.1.12-4 (r884)
[samtools.git] / bam2bcf_indel.c
1 #include <assert.h>
2 #include <ctype.h>
3 #include <string.h>
4 #include "bam.h"
5 #include "bam2bcf.h"
6 #include "ksort.h"
7 #include "kaln.h"
8 #include "kprobaln.h"
9 #include "khash.h"
10 KHASH_SET_INIT_STR(rg)
11
12 #define MINUS_CONST 0x10000000
13 #define INDEL_WINDOW_SIZE 50
14
15 void *bcf_call_add_rg(void *_hash, const char *hdtext, const char *list)
16 {
17         const char *s, *p, *q, *r, *t;
18         khash_t(rg) *hash;
19         if (list == 0 || hdtext == 0) return _hash;
20         if (_hash == 0) _hash = kh_init(rg);
21         hash = (khash_t(rg)*)_hash;
22         if ((s = strstr(hdtext, "@RG\t")) == 0) return hash;
23         do {
24                 t = strstr(s + 4, "@RG\t"); // the next @RG
25                 if ((p = strstr(s, "\tID:")) != 0) p += 4;
26                 if ((q = strstr(s, "\tPL:")) != 0) q += 4;
27                 if (p && q && (t == 0 || (p < t && q < t))) { // ID and PL are both present
28                         int lp, lq;
29                         char *x;
30                         for (r = p; *r && *r != '\t' && *r != '\n'; ++r); lp = r - p;
31                         for (r = q; *r && *r != '\t' && *r != '\n'; ++r); lq = r - q;
32                         x = calloc((lp > lq? lp : lq) + 1, 1);
33                         for (r = q; *r && *r != '\t' && *r != '\n'; ++r) x[r-q] = *r;
34                         if (strstr(list, x)) { // insert ID to the hash table
35                                 khint_t k;
36                                 int ret;
37                                 for (r = p; *r && *r != '\t' && *r != '\n'; ++r) x[r-p] = *r;
38                                 x[r-p] = 0;
39                                 k = kh_get(rg, hash, x);
40                                 if (k == kh_end(hash)) k = kh_put(rg, hash, x, &ret);
41                                 else free(x);
42                         } else free(x);
43                 }
44                 s = t;
45         } while (s);
46         return hash;
47 }
48
49 void bcf_call_del_rghash(void *_hash)
50 {
51         khint_t k;
52         khash_t(rg) *hash = (khash_t(rg)*)_hash;
53         if (hash == 0) return;
54         for (k = kh_begin(hash); k < kh_end(hash); ++k)
55                 if (kh_exist(hash, k))
56                         free((char*)kh_key(hash, k));
57         kh_destroy(rg, hash);
58 }
59
60 static int tpos2qpos(const bam1_core_t *c, const uint32_t *cigar, int32_t tpos, int is_left, int32_t *_tpos)
61 {
62         int k, x = c->pos, y = 0, last_y = 0;
63         *_tpos = c->pos;
64         for (k = 0; k < c->n_cigar; ++k) {
65                 int op = cigar[k] & BAM_CIGAR_MASK;
66                 int l = cigar[k] >> BAM_CIGAR_SHIFT;
67                 if (op == BAM_CMATCH) {
68                         if (c->pos > tpos) return y;
69                         if (x + l > tpos) {
70                                 *_tpos = tpos;
71                                 return y + (tpos - x);
72                         }
73                         x += l; y += l;
74                         last_y = y;
75                 } else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += l;
76                 else if (op == BAM_CDEL || op == BAM_CREF_SKIP) {
77                         if (x + l > tpos) {
78                                 *_tpos = is_left? x : x + l;
79                                 return y;
80                         }
81                         x += l;
82                 }
83         }
84         *_tpos = x;
85         return last_y;
86 }
87 // FIXME: check if the inserted sequence is consistent with the homopolymer run
88 // l is the relative gap length and l_run is the length of the homopolymer on the reference
89 static inline int est_seqQ(const bcf_callaux_t *bca, int l, int l_run)
90 {
91         int q, qh;
92         q = bca->openQ + bca->extQ * (abs(l) - 1);
93         qh = l_run >= 3? (int)(bca->tandemQ * (double)abs(l) / l_run + .499) : 1000;
94         return q < qh? q : qh;
95 }
96
97 static inline int est_indelreg(int pos, const char *ref, int l, char *ins4)
98 {
99         int i, j, max = 0, max_i = pos, score = 0;
100         l = abs(l);
101         for (i = pos + 1, j = 0; ref[i]; ++i, ++j) {
102                 if (ins4) score += (toupper(ref[i]) != "ACGTN"[(int)ins4[j%l]])? -10 : 1;
103                 else score += (toupper(ref[i]) != toupper(ref[pos+1+j%l]))? -10 : 1;
104                 if (score < 0) break;
105                 if (max < score) max = score, max_i = i;
106         }
107         return max_i - pos;
108 }
109
110 int bcf_call_gap_prep(int n, int *n_plp, bam_pileup1_t **plp, int pos, bcf_callaux_t *bca, const char *ref,
111                                           const void *rghash)
112 {
113         extern void ks_introsort_uint32_t(int, uint32_t*);
114         int i, s, j, k, t, n_types, *types, max_rd_len, left, right, max_ins, *score1, *score2, max_ref2;
115         int N, K, l_run, ref_type, n_alt;
116         char *inscns = 0, *ref2, *query, **ref_sample;
117         khash_t(rg) *hash = (khash_t(rg)*)rghash;
118         if (ref == 0 || bca == 0) return -1;
119         // mark filtered reads
120         if (rghash) {
121                 N = 0;
122                 for (s = N = 0; s < n; ++s) {
123                         for (i = 0; i < n_plp[s]; ++i) {
124                                 bam_pileup1_t *p = plp[s] + i;
125                                 const uint8_t *rg = bam_aux_get(p->b, "RG");
126                                 p->aux = 1; // filtered by default
127                                 if (rg) {
128                                         khint_t k = kh_get(rg, hash, (const char*)(rg + 1));
129                                         if (k != kh_end(hash)) p->aux = 0, ++N; // not filtered
130                                 }
131                         }
132                 }
133                 if (N == 0) return -1; // no reads left
134         }
135         // determine if there is a gap
136         for (s = N = 0; s < n; ++s) {
137                 for (i = 0; i < n_plp[s]; ++i)
138                         if (plp[s][i].indel != 0) break;
139                 if (i < n_plp[s]) break;
140         }
141         if (s == n) return -1; // there is no indel at this position.
142         for (s = N = 0; s < n; ++s) N += n_plp[s]; // N is the total number of reads
143         { // find out how many types of indels are present
144                 int m, n_alt = 0, n_tot = 0;
145                 uint32_t *aux;
146                 aux = calloc(N + 1, 4);
147                 m = max_rd_len = 0;
148                 aux[m++] = MINUS_CONST; // zero indel is always a type
149                 for (s = 0; s < n; ++s) {
150                         for (i = 0; i < n_plp[s]; ++i) {
151                                 const bam_pileup1_t *p = plp[s] + i;
152                                 if (rghash == 0 || p->aux == 0) {
153                                         ++n_tot;
154                                         if (p->indel != 0) {
155                                                 ++n_alt;
156                                                 aux[m++] = MINUS_CONST + p->indel;
157                                         }
158                                 }
159                                 j = bam_cigar2qlen(&p->b->core, bam1_cigar(p->b));
160                                 if (j > max_rd_len) max_rd_len = j;
161                         }
162                 }
163                 ks_introsort(uint32_t, m, aux);
164                 // squeeze out identical types
165                 for (i = 1, n_types = 1; i < m; ++i)
166                         if (aux[i] != aux[i-1]) ++n_types;
167                 if (n_types == 1 || (double)n_alt / n_tot < bca->min_frac || n_alt < bca->min_support) { // then skip
168                         free(aux); return -1;
169                 }
170                 types = (int*)calloc(n_types, sizeof(int));
171                 t = 0;
172                 types[t++] = aux[0] - MINUS_CONST; 
173                 for (i = 1; i < m; ++i)
174                         if (aux[i] != aux[i-1])
175                                 types[t++] = aux[i] - MINUS_CONST;
176                 free(aux);
177                 for (t = 0; t < n_types; ++t)
178                         if (types[t] == 0) break;
179                 ref_type = t; // the index of the reference type (0)
180                 assert(n_types < 64);
181         }
182         { // calculate left and right boundary
183                 left = pos > INDEL_WINDOW_SIZE? pos - INDEL_WINDOW_SIZE : 0;
184                 right = pos + INDEL_WINDOW_SIZE;
185                 if (types[0] < 0) right -= types[0];
186                 // in case the alignments stand out the reference
187                 for (i = pos; i < right; ++i)
188                         if (ref[i] == 0) break;
189                 right = i;
190         }
191         /* The following block fixes a long-existing flaw in the INDEL
192          * calling model: the interference of nearby SNPs. However, it also
193          * reduces the power because sometimes, substitutions caused by
194          * indels are not distinguishable from true mutations. Multiple
195          * sequence realignment helps to increase the power.
196          */
197         { // construct per-sample consensus
198                 int L = right - left + 1;
199                 uint32_t *cns;
200                 char *ref0, *r;
201                 ref_sample = calloc(n, sizeof(void*));
202                 cns = calloc(L, 4);
203                 ref0 = calloc(L, 1);
204                 for (i = 0; i < right - left; ++i)
205                         ref0[i] = bam_nt16_table[(int)ref[i+left]];
206                 for (s = 0; s < n; ++s) {
207                         r = ref_sample[s] = calloc(L, 1);
208                         memset(cns, 0, sizeof(int) * L);
209                         // collect ref and non-ref counts
210                         for (i = 0; i < n_plp[s]; ++i) {
211                                 bam_pileup1_t *p = plp[s] + i;
212                                 bam1_t *b = p->b;
213                                 uint32_t *cigar = bam1_cigar(b);
214                                 uint8_t *seq = bam1_seq(b);
215                                 int x = b->core.pos, y = 0;
216                                 for (k = 0; k < b->core.n_cigar; ++k) {
217                                         int op = cigar[k]&0xf;
218                                         int j, l = cigar[k]>>4;
219                                         if (op == BAM_CMATCH) {
220                                                 for (j = 0; j < l; ++j)
221                                                         if (x + j >= left && x + j < right)
222                                                                 cns[x+j-left] += (bam1_seqi(seq, y+j) == ref0[x+j-left])? 1 : 0x10000;
223                                                 x += l; y += l;
224                                         } else if (op == BAM_CDEL || op == BAM_CREF_SKIP) x += l;
225                                         else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += l;
226                                 }
227                         }
228                         // determine the consensus
229                         for (i = 0; i < right - left; ++i)
230                                 r[i] = (cns[i] && (double)(cns[i]&0xffff) / ((cns[i]&0xffff)+(cns[i]>>16&0xffff)) < 0.7)? 15 : ref0[i];
231 //                      for (i = 0; i < right - left; ++i) fputc("=ACMGRSVTWYHKDBN"[(int)r[i]], stderr); fputc('\n', stderr);
232                 }
233                 free(ref0); free(cns);
234         }
235         { // the length of the homopolymer run around the current position
236                 int c = bam_nt16_table[(int)ref[pos + 1]];
237                 if (c == 15) l_run = 1;
238                 else {
239                         for (i = pos + 2; ref[i]; ++i)
240                                 if (bam_nt16_table[(int)ref[i]] != c) break;
241                         l_run = i;
242                         for (i = pos; i >= 0; --i)
243                                 if (bam_nt16_table[(int)ref[i]] != c) break;
244                         l_run -= i + 1;
245                 }
246         }
247         // construct the consensus sequence
248         max_ins = types[n_types - 1]; // max_ins is at least 0
249         if (max_ins > 0) {
250                 int *inscns_aux = calloc(4 * n_types * max_ins, sizeof(int));
251                 // count the number of occurrences of each base at each position for each type of insertion
252                 for (t = 0; t < n_types; ++t) {
253                         if (types[t] > 0) {
254                                 for (s = 0; s < n; ++s) {
255                                         for (i = 0; i < n_plp[s]; ++i) {
256                                                 bam_pileup1_t *p = plp[s] + i;
257                                                 if (p->indel == types[t]) {
258                                                         uint8_t *seq = bam1_seq(p->b);
259                                                         for (k = 1; k <= p->indel; ++k) {
260                                                                 int c = bam_nt16_nt4_table[bam1_seqi(seq, p->qpos + k)];
261                                                                 if (c < 4) ++inscns_aux[(t*max_ins+(k-1))*4 + c];
262                                                         }
263                                                 }
264                                         }
265                                 }
266                         }
267                 }
268                 // use the majority rule to construct the consensus
269                 inscns = calloc(n_types * max_ins, 1);
270                 for (t = 0; t < n_types; ++t) {
271                         for (j = 0; j < types[t]; ++j) {
272                                 int max = 0, max_k = -1, *ia = &inscns_aux[(t*max_ins+j)*4];
273                                 for (k = 0; k < 4; ++k)
274                                         if (ia[k] > max)
275                                                 max = ia[k], max_k = k;
276                                 inscns[t*max_ins + j] = max? max_k : 4;
277                         }
278                 }
279                 free(inscns_aux);
280         }
281         // compute the likelihood given each type of indel for each read
282         max_ref2 = right - left + 2 + 2 * (max_ins > -types[0]? max_ins : -types[0]);
283         ref2  = calloc(max_ref2, 1);
284         query = calloc(right - left + max_rd_len + max_ins + 2, 1);
285         score1 = calloc(N * n_types, sizeof(int));
286         score2 = calloc(N * n_types, sizeof(int));
287         bca->indelreg = 0;
288         for (t = 0; t < n_types; ++t) {
289                 int l, ir;
290                 kpa_par_t apf1 = { 1e-4, 1e-2, 10 }, apf2 = { 1e-6, 1e-3, 10 };
291                 apf1.bw = apf2.bw = abs(types[t]) + 3;
292                 // compute indelreg
293                 if (types[t] == 0) ir = 0;
294                 else if (types[t] > 0) ir = est_indelreg(pos, ref, types[t], &inscns[t*max_ins]);
295                 else ir = est_indelreg(pos, ref, -types[t], 0);
296                 if (ir > bca->indelreg) bca->indelreg = ir;
297 //              fprintf(stderr, "%d, %d, %d\n", pos, types[t], ir);
298                 // realignment
299                 for (s = K = 0; s < n; ++s) {
300                         // write ref2
301                         for (k = 0, j = left; j <= pos; ++j)
302                                 ref2[k++] = bam_nt16_nt4_table[(int)ref_sample[s][j-left]];
303                         if (types[t] <= 0) j += -types[t];
304                         else for (l = 0; l < types[t]; ++l)
305                                          ref2[k++] = inscns[t*max_ins + l];
306                         for (; j < right && ref[j]; ++j)
307                                 ref2[k++] = bam_nt16_nt4_table[(int)ref_sample[s][j-left]];
308                         for (; k < max_ref2; ++k) ref2[k] = 4;
309                         if (j < right) right = j;
310                         // align each read to ref2
311                         for (i = 0; i < n_plp[s]; ++i, ++K) {
312                                 bam_pileup1_t *p = plp[s] + i;
313                                 int qbeg, qend, tbeg, tend, sc;
314                                 uint8_t *seq = bam1_seq(p->b);
315                                 // FIXME: the following skips soft clips, but using them may be more sensitive.
316                                 // determine the start and end of sequences for alignment
317                                 qbeg = tpos2qpos(&p->b->core, bam1_cigar(p->b), left,  0, &tbeg);
318                                 qend = tpos2qpos(&p->b->core, bam1_cigar(p->b), right, 1, &tend);
319                                 if (types[t] < 0) {
320                                         int l = -types[t];
321                                         tbeg = tbeg - l > left?  tbeg - l : left;
322                                 }
323                                 // write the query sequence
324                                 for (l = qbeg; l < qend; ++l)
325                                         query[l - qbeg] = bam_nt16_nt4_table[bam1_seqi(seq, l)];
326                                 { // do realignment; this is the bottleneck
327                                         const uint8_t *qual = bam1_qual(p->b), *bq;
328                                         uint8_t *qq;
329                                         qq = calloc(qend - qbeg, 1);
330                                         bq = (uint8_t*)bam_aux_get(p->b, "ZQ");
331                                         if (bq) ++bq; // skip type
332                                         for (l = qbeg; l < qend; ++l) {
333                                                 qq[l - qbeg] = bq? qual[l] + (bq[l] - 64) : qual[l];
334                                                 if (qq[l - qbeg] > 30) qq[l - qbeg] = 30;
335                                                 if (qq[l - qbeg] < 7) qq[l - qbeg] = 7;
336                                         }
337                                         sc = kpa_glocal((uint8_t*)ref2 + tbeg - left, tend - tbeg + abs(types[t]),
338                                                                         (uint8_t*)query, qend - qbeg, qq, &apf1, 0, 0);
339                                         l = (int)(100. * sc / (qend - qbeg) + .499); // used for adjusting indelQ below
340                                         if (l > 255) l = 255;
341                                         score1[K*n_types + t] = score2[K*n_types + t] = sc<<8 | l;
342                                         if (sc > 5) {
343                                                 sc = kpa_glocal((uint8_t*)ref2 + tbeg - left, tend - tbeg + abs(types[t]),
344                                                                                 (uint8_t*)query, qend - qbeg, qq, &apf2, 0, 0);
345                                                 l = (int)(100. * sc / (qend - qbeg) + .499);
346                                                 if (l > 255) l = 255;
347                                                 score2[K*n_types + t] = sc<<8 | l;
348                                         }
349                                         free(qq);
350                                 }
351 /*
352                                 for (l = 0; l < tend - tbeg + abs(types[t]); ++l)
353                                         fputc("ACGTN"[(int)ref2[tbeg-left+l]], stderr);
354                                 fputc('\n', stderr);
355                                 for (l = 0; l < qend - qbeg; ++l) fputc("ACGTN"[(int)query[l]], stderr);
356                                 fputc('\n', stderr);
357                                 fprintf(stderr, "pos=%d type=%d read=%d:%d name=%s qbeg=%d tbeg=%d score=%d\n", pos, types[t], s, i, bam1_qname(p->b), qbeg, tbeg, sc);
358 */
359                         }
360                 }
361         }
362         free(ref2); free(query);
363         { // compute indelQ
364                 int *sc, tmp, *sumq;
365                 sc   = alloca(n_types * sizeof(int));
366                 sumq = alloca(n_types * sizeof(int));
367                 memset(sumq, 0, sizeof(int) * n_types);
368                 for (s = K = 0; s < n; ++s) {
369                         for (i = 0; i < n_plp[s]; ++i, ++K) {
370                                 bam_pileup1_t *p = plp[s] + i;
371                                 int *sct = &score1[K*n_types], indelQ1, indelQ2, seqQ, indelQ;
372                                 for (t = 0; t < n_types; ++t) sc[t] = sct[t]<<6 | t;
373                                 for (t = 1; t < n_types; ++t) // insertion sort
374                                         for (j = t; j > 0 && sc[j] < sc[j-1]; --j)
375                                                 tmp = sc[j], sc[j] = sc[j-1], sc[j-1] = tmp;
376                                 /* errmod_cal() assumes that if the call is wrong, the
377                                  * likelihoods of other events are equal. This is about
378                                  * right for substitutions, but is not desired for
379                                  * indels. To reuse errmod_cal(), I have to make
380                                  * compromise for multi-allelic indels.
381                                  */
382                                 if ((sc[0]&0x3f) == ref_type) {
383                                         indelQ1 = (sc[1]>>14) - (sc[0]>>14);
384                                         seqQ = est_seqQ(bca, types[sc[1]&0x3f], l_run);
385                                 } else {
386                                         for (t = 0; t < n_types; ++t) // look for the reference type
387                                                 if ((sc[t]&0x3f) == ref_type) break;
388                                         indelQ1 = (sc[t]>>14) - (sc[0]>>14);
389                                         seqQ = est_seqQ(bca, types[sc[0]&0x3f], l_run);
390                                 }
391                                 tmp = sc[0]>>6 & 0xff;
392                                 indelQ1 = tmp > 111? 0 : (int)((1. - tmp/111.) * indelQ1 + .499); // reduce indelQ
393                                 sct = &score2[K*n_types];
394                                 for (t = 0; t < n_types; ++t) sc[t] = sct[t]<<6 | t;
395                                 for (t = 1; t < n_types; ++t) // insertion sort
396                                         for (j = t; j > 0 && sc[j] < sc[j-1]; --j)
397                                                 tmp = sc[j], sc[j] = sc[j-1], sc[j-1] = tmp;
398                                 if ((sc[0]&0x3f) == ref_type) {
399                                         indelQ2 = (sc[1]>>14) - (sc[0]>>14);
400                                 } else {
401                                         for (t = 0; t < n_types; ++t) // look for the reference type
402                                                 if ((sc[t]&0x3f) == ref_type) break;
403                                         indelQ2 = (sc[t]>>14) - (sc[0]>>14);
404                                 }
405                                 tmp = sc[0]>>6 & 0xff;
406                                 indelQ2 = tmp > 111? 0 : (int)((1. - tmp/111.) * indelQ2 + .499);
407                                 // pick the smaller between indelQ1 and indelQ2
408                                 indelQ = indelQ1 < indelQ2? indelQ1 : indelQ2;
409                                 p->aux = (sc[0]&0x3f)<<16 | seqQ<<8 | indelQ;
410                                 sumq[sc[0]&0x3f] += indelQ < seqQ? indelQ : seqQ;
411 //                              fprintf(stderr, "pos=%d read=%d:%d name=%s call=%d q=%d\n", pos, s, i, bam1_qname(p->b), types[sc[0]&0x3f], indelQ);
412                         }
413                 }
414                 // determine bca->indel_types[] and bca->inscns
415                 bca->maxins = max_ins;
416                 bca->inscns = realloc(bca->inscns, bca->maxins * 4);
417                 for (t = 0; t < n_types; ++t)
418                         sumq[t] = sumq[t]<<6 | t;
419                 for (t = 1; t < n_types; ++t) // insertion sort
420                         for (j = t; j > 0 && sumq[j] > sumq[j-1]; --j)
421                                 tmp = sumq[j], sumq[j] = sumq[j-1], sumq[j-1] = tmp;
422                 for (t = 0; t < n_types; ++t) // look for the reference type
423                         if ((sumq[t]&0x3f) == ref_type) break;
424                 if (t) { // then move the reference type to the first
425                         tmp = sumq[t];
426                         for (; t > 0; --t) sumq[t] = sumq[t-1];
427                         sumq[0] = tmp;
428                 }
429                 for (t = 0; t < 4; ++t) bca->indel_types[t] = B2B_INDEL_NULL;
430                 for (t = 0; t < 4 && t < n_types; ++t) {
431                         bca->indel_types[t] = types[sumq[t]&0x3f];
432                         memcpy(&bca->inscns[t * bca->maxins], &inscns[(sumq[t]&0x3f) * max_ins], bca->maxins);
433                 }
434                 // update p->aux
435                 for (s = n_alt = 0; s < n; ++s) {
436                         for (i = 0; i < n_plp[s]; ++i) {
437                                 bam_pileup1_t *p = plp[s] + i;
438                                 int x = types[p->aux>>16&0x3f];
439                                 for (j = 0; j < 4; ++j)
440                                         if (x == bca->indel_types[j]) break;
441                                 p->aux = j<<16 | (j == 4? 0 : (p->aux&0xffff));
442                                 if ((p->aux>>16&0x3f) > 0) ++n_alt;
443 //                              fprintf(stderr, "X pos=%d read=%d:%d name=%s call=%d type=%d q=%d seqQ=%d\n", pos, s, i, bam1_qname(p->b), p->aux>>16&63, bca->indel_types[p->aux>>16&63], p->aux&0xff, p->aux>>8&0xff);
444                         }
445                 }               
446         }
447         free(score1); free(score2);
448         // free
449         for (i = 0; i < n; ++i) free(ref_sample[i]);
450         free(ref_sample);
451         free(types); free(inscns);
452         return n_alt > 0? 0 : -1;
453 }