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