Imported Upstream version 0.1.18
[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;
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                         for (i = 0; i < n_plp[s]; ++i) {
152                                 const bam_pileup1_t *p = plp[s] + i;
153                                 if (rghash == 0 || p->aux == 0) {
154                                         ++n_tot;
155                                         if (p->indel != 0) {
156                                                 ++n_alt;
157                                                 aux[m++] = MINUS_CONST + p->indel;
158                                         }
159                                 }
160                                 j = bam_cigar2qlen(&p->b->core, bam1_cigar(p->b));
161                                 if (j > max_rd_len) max_rd_len = j;
162                         }
163                 }
164                 ks_introsort(uint32_t, m, aux);
165                 // squeeze out identical types
166                 for (i = 1, n_types = 1; i < m; ++i)
167                         if (aux[i] != aux[i-1]) ++n_types;
168                 if (n_types == 1 || (double)n_alt / n_tot < bca->min_frac || n_alt < bca->min_support) { // then skip
169                         free(aux); return -1;
170                 }
171                 if (n_types >= 64) {
172                         free(aux);
173                         if (bam_verbose >= 2) 
174                                 fprintf(stderr, "[%s] excessive INDEL alleles at position %d. Skip the position.\n", __func__, pos + 1);
175                         return -1;
176                 }
177                 types = (int*)calloc(n_types, sizeof(int));
178                 t = 0;
179                 types[t++] = aux[0] - MINUS_CONST; 
180                 for (i = 1; i < m; ++i)
181                         if (aux[i] != aux[i-1])
182                                 types[t++] = aux[i] - MINUS_CONST;
183                 free(aux);
184                 for (t = 0; t < n_types; ++t)
185                         if (types[t] == 0) break;
186                 ref_type = t; // the index of the reference type (0)
187         }
188         { // calculate left and right boundary
189                 left = pos > INDEL_WINDOW_SIZE? pos - INDEL_WINDOW_SIZE : 0;
190                 right = pos + INDEL_WINDOW_SIZE;
191                 if (types[0] < 0) right -= types[0];
192                 // in case the alignments stand out the reference
193                 for (i = pos; i < right; ++i)
194                         if (ref[i] == 0) break;
195                 right = i;
196         }
197         /* The following block fixes a long-existing flaw in the INDEL
198          * calling model: the interference of nearby SNPs. However, it also
199          * reduces the power because sometimes, substitutions caused by
200          * indels are not distinguishable from true mutations. Multiple
201          * sequence realignment helps to increase the power.
202          */
203         { // construct per-sample consensus
204                 int L = right - left + 1, max_i, max2_i;
205                 uint32_t *cns, max, max2;
206                 char *ref0, *r;
207                 ref_sample = calloc(n, sizeof(void*));
208                 cns = calloc(L, 4);
209                 ref0 = calloc(L, 1);
210                 for (i = 0; i < right - left; ++i)
211                         ref0[i] = bam_nt16_table[(int)ref[i+left]];
212                 for (s = 0; s < n; ++s) {
213                         r = ref_sample[s] = calloc(L, 1);
214                         memset(cns, 0, sizeof(int) * L);
215                         // collect ref and non-ref counts
216                         for (i = 0; i < n_plp[s]; ++i) {
217                                 bam_pileup1_t *p = plp[s] + i;
218                                 bam1_t *b = p->b;
219                                 uint32_t *cigar = bam1_cigar(b);
220                                 uint8_t *seq = bam1_seq(b);
221                                 int x = b->core.pos, y = 0;
222                                 for (k = 0; k < b->core.n_cigar; ++k) {
223                                         int op = cigar[k]&0xf;
224                                         int j, l = cigar[k]>>4;
225                                         if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
226                                                 for (j = 0; j < l; ++j)
227                                                         if (x + j >= left && x + j < right)
228                                                                 cns[x+j-left] += (bam1_seqi(seq, y+j) == ref0[x+j-left])? 1 : 0x10000;
229                                                 x += l; y += l;
230                                         } else if (op == BAM_CDEL || op == BAM_CREF_SKIP) x += l;
231                                         else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += l;
232                                 }
233                         }
234                         // determine the consensus
235                         for (i = 0; i < right - left; ++i) r[i] = ref0[i];
236                         max = max2 = 0; max_i = max2_i = -1;
237                         for (i = 0; i < right - left; ++i) {
238                                 if (cns[i]>>16 >= max>>16) max2 = max, max2_i = max_i, max = cns[i], max_i = i;
239                                 else if (cns[i]>>16 >= max2>>16) max2 = cns[i], max2_i = i;
240                         }
241                         if ((double)(max&0xffff) / ((max&0xffff) + (max>>16)) >= 0.7) max_i = -1;
242                         if ((double)(max2&0xffff) / ((max2&0xffff) + (max2>>16)) >= 0.7) max2_i = -1;
243                         if (max_i >= 0) r[max_i] = 15;
244                         if (max2_i >= 0) r[max2_i] = 15;
245 //                      for (i = 0; i < right - left; ++i) fputc("=ACMGRSVTWYHKDBN"[(int)r[i]], stderr); fputc('\n', stderr);
246                 }
247                 free(ref0); free(cns);
248         }
249         { // the length of the homopolymer run around the current position
250                 int c = bam_nt16_table[(int)ref[pos + 1]];
251                 if (c == 15) l_run = 1;
252                 else {
253                         for (i = pos + 2; ref[i]; ++i)
254                                 if (bam_nt16_table[(int)ref[i]] != c) break;
255                         l_run = i;
256                         for (i = pos; i >= 0; --i)
257                                 if (bam_nt16_table[(int)ref[i]] != c) break;
258                         l_run -= i + 1;
259                 }
260         }
261         // construct the consensus sequence
262         max_ins = types[n_types - 1]; // max_ins is at least 0
263         if (max_ins > 0) {
264                 int *inscns_aux = calloc(4 * n_types * max_ins, sizeof(int));
265                 // count the number of occurrences of each base at each position for each type of insertion
266                 for (t = 0; t < n_types; ++t) {
267                         if (types[t] > 0) {
268                                 for (s = 0; s < n; ++s) {
269                                         for (i = 0; i < n_plp[s]; ++i) {
270                                                 bam_pileup1_t *p = plp[s] + i;
271                                                 if (p->indel == types[t]) {
272                                                         uint8_t *seq = bam1_seq(p->b);
273                                                         for (k = 1; k <= p->indel; ++k) {
274                                                                 int c = bam_nt16_nt4_table[bam1_seqi(seq, p->qpos + k)];
275                                                                 if (c < 4) ++inscns_aux[(t*max_ins+(k-1))*4 + c];
276                                                         }
277                                                 }
278                                         }
279                                 }
280                         }
281                 }
282                 // use the majority rule to construct the consensus
283                 inscns = calloc(n_types * max_ins, 1);
284                 for (t = 0; t < n_types; ++t) {
285                         for (j = 0; j < types[t]; ++j) {
286                                 int max = 0, max_k = -1, *ia = &inscns_aux[(t*max_ins+j)*4];
287                                 for (k = 0; k < 4; ++k)
288                                         if (ia[k] > max)
289                                                 max = ia[k], max_k = k;
290                                 inscns[t*max_ins + j] = max? max_k : 4;
291                         }
292                 }
293                 free(inscns_aux);
294         }
295         // compute the likelihood given each type of indel for each read
296         max_ref2 = right - left + 2 + 2 * (max_ins > -types[0]? max_ins : -types[0]);
297         ref2  = calloc(max_ref2, 1);
298         query = calloc(right - left + max_rd_len + max_ins + 2, 1);
299         score1 = calloc(N * n_types, sizeof(int));
300         score2 = calloc(N * n_types, sizeof(int));
301         bca->indelreg = 0;
302         for (t = 0; t < n_types; ++t) {
303                 int l, ir;
304                 kpa_par_t apf1 = { 1e-4, 1e-2, 10 }, apf2 = { 1e-6, 1e-3, 10 };
305                 apf1.bw = apf2.bw = abs(types[t]) + 3;
306                 // compute indelreg
307                 if (types[t] == 0) ir = 0;
308                 else if (types[t] > 0) ir = est_indelreg(pos, ref, types[t], &inscns[t*max_ins]);
309                 else ir = est_indelreg(pos, ref, -types[t], 0);
310                 if (ir > bca->indelreg) bca->indelreg = ir;
311 //              fprintf(stderr, "%d, %d, %d\n", pos, types[t], ir);
312                 // realignment
313                 for (s = K = 0; s < n; ++s) {
314                         // write ref2
315                         for (k = 0, j = left; j <= pos; ++j)
316                                 ref2[k++] = bam_nt16_nt4_table[(int)ref_sample[s][j-left]];
317                         if (types[t] <= 0) j += -types[t];
318                         else for (l = 0; l < types[t]; ++l)
319                                          ref2[k++] = inscns[t*max_ins + l];
320                         for (; j < right && ref[j]; ++j)
321                                 ref2[k++] = bam_nt16_nt4_table[(int)ref_sample[s][j-left]];
322                         for (; k < max_ref2; ++k) ref2[k] = 4;
323                         if (j < right) right = j;
324                         // align each read to ref2
325                         for (i = 0; i < n_plp[s]; ++i, ++K) {
326                                 bam_pileup1_t *p = plp[s] + i;
327                                 int qbeg, qend, tbeg, tend, sc, kk;
328                                 uint8_t *seq = bam1_seq(p->b);
329                                 uint32_t *cigar = bam1_cigar(p->b);
330                                 if (p->b->core.flag&4) continue; // unmapped reads
331                                 // FIXME: the following loop should be better moved outside; nonetheless, realignment should be much slower anyway.
332                                 for (kk = 0; kk < p->b->core.n_cigar; ++kk)
333                                         if ((cigar[kk]&BAM_CIGAR_MASK) == BAM_CREF_SKIP) break;
334                                 if (kk < p->b->core.n_cigar) continue;
335                                 // FIXME: the following skips soft clips, but using them may be more sensitive.
336                                 // determine the start and end of sequences for alignment
337                                 qbeg = tpos2qpos(&p->b->core, bam1_cigar(p->b), left,  0, &tbeg);
338                                 qend = tpos2qpos(&p->b->core, bam1_cigar(p->b), right, 1, &tend);
339                                 if (types[t] < 0) {
340                                         int l = -types[t];
341                                         tbeg = tbeg - l > left?  tbeg - l : left;
342                                 }
343                                 // write the query sequence
344                                 for (l = qbeg; l < qend; ++l)
345                                         query[l - qbeg] = bam_nt16_nt4_table[bam1_seqi(seq, l)];
346                                 { // do realignment; this is the bottleneck
347                                         const uint8_t *qual = bam1_qual(p->b), *bq;
348                                         uint8_t *qq;
349                                         qq = calloc(qend - qbeg, 1);
350                                         bq = (uint8_t*)bam_aux_get(p->b, "ZQ");
351                                         if (bq) ++bq; // skip type
352                                         for (l = qbeg; l < qend; ++l) {
353                                                 qq[l - qbeg] = bq? qual[l] + (bq[l] - 64) : qual[l];
354                                                 if (qq[l - qbeg] > 30) qq[l - qbeg] = 30;
355                                                 if (qq[l - qbeg] < 7) qq[l - qbeg] = 7;
356                                         }
357                                         sc = kpa_glocal((uint8_t*)ref2 + tbeg - left, tend - tbeg + abs(types[t]),
358                                                                         (uint8_t*)query, qend - qbeg, qq, &apf1, 0, 0);
359                                         l = (int)(100. * sc / (qend - qbeg) + .499); // used for adjusting indelQ below
360                                         if (l > 255) l = 255;
361                                         score1[K*n_types + t] = score2[K*n_types + t] = sc<<8 | l;
362                                         if (sc > 5) {
363                                                 sc = kpa_glocal((uint8_t*)ref2 + tbeg - left, tend - tbeg + abs(types[t]),
364                                                                                 (uint8_t*)query, qend - qbeg, qq, &apf2, 0, 0);
365                                                 l = (int)(100. * sc / (qend - qbeg) + .499);
366                                                 if (l > 255) l = 255;
367                                                 score2[K*n_types + t] = sc<<8 | l;
368                                         }
369                                         free(qq);
370                                 }
371 /*
372                                 for (l = 0; l < tend - tbeg + abs(types[t]); ++l)
373                                         fputc("ACGTN"[(int)ref2[tbeg-left+l]], stderr);
374                                 fputc('\n', stderr);
375                                 for (l = 0; l < qend - qbeg; ++l) fputc("ACGTN"[(int)query[l]], stderr);
376                                 fputc('\n', stderr);
377                                 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);
378 */
379                         }
380                 }
381         }
382         free(ref2); free(query);
383         { // compute indelQ
384                 int *sc, tmp, *sumq;
385                 sc   = alloca(n_types * sizeof(int));
386                 sumq = alloca(n_types * sizeof(int));
387                 memset(sumq, 0, sizeof(int) * n_types);
388                 for (s = K = 0; s < n; ++s) {
389                         for (i = 0; i < n_plp[s]; ++i, ++K) {
390                                 bam_pileup1_t *p = plp[s] + i;
391                                 int *sct = &score1[K*n_types], indelQ1, indelQ2, seqQ, indelQ;
392                                 for (t = 0; t < n_types; ++t) sc[t] = sct[t]<<6 | t;
393                                 for (t = 1; t < n_types; ++t) // insertion sort
394                                         for (j = t; j > 0 && sc[j] < sc[j-1]; --j)
395                                                 tmp = sc[j], sc[j] = sc[j-1], sc[j-1] = tmp;
396                                 /* errmod_cal() assumes that if the call is wrong, the
397                                  * likelihoods of other events are equal. This is about
398                                  * right for substitutions, but is not desired for
399                                  * indels. To reuse errmod_cal(), I have to make
400                                  * compromise for multi-allelic indels.
401                                  */
402                                 if ((sc[0]&0x3f) == ref_type) {
403                                         indelQ1 = (sc[1]>>14) - (sc[0]>>14);
404                                         seqQ = est_seqQ(bca, types[sc[1]&0x3f], l_run);
405                                 } else {
406                                         for (t = 0; t < n_types; ++t) // look for the reference type
407                                                 if ((sc[t]&0x3f) == ref_type) break;
408                                         indelQ1 = (sc[t]>>14) - (sc[0]>>14);
409                                         seqQ = est_seqQ(bca, types[sc[0]&0x3f], l_run);
410                                 }
411                                 tmp = sc[0]>>6 & 0xff;
412                                 indelQ1 = tmp > 111? 0 : (int)((1. - tmp/111.) * indelQ1 + .499); // reduce indelQ
413                                 sct = &score2[K*n_types];
414                                 for (t = 0; t < n_types; ++t) sc[t] = sct[t]<<6 | t;
415                                 for (t = 1; t < n_types; ++t) // insertion sort
416                                         for (j = t; j > 0 && sc[j] < sc[j-1]; --j)
417                                                 tmp = sc[j], sc[j] = sc[j-1], sc[j-1] = tmp;
418                                 if ((sc[0]&0x3f) == ref_type) {
419                                         indelQ2 = (sc[1]>>14) - (sc[0]>>14);
420                                 } else {
421                                         for (t = 0; t < n_types; ++t) // look for the reference type
422                                                 if ((sc[t]&0x3f) == ref_type) break;
423                                         indelQ2 = (sc[t]>>14) - (sc[0]>>14);
424                                 }
425                                 tmp = sc[0]>>6 & 0xff;
426                                 indelQ2 = tmp > 111? 0 : (int)((1. - tmp/111.) * indelQ2 + .499);
427                                 // pick the smaller between indelQ1 and indelQ2
428                                 indelQ = indelQ1 < indelQ2? indelQ1 : indelQ2;
429                                 if (indelQ > 255) indelQ = 255;
430                                 if (seqQ > 255) seqQ = 255;
431                                 p->aux = (sc[0]&0x3f)<<16 | seqQ<<8 | indelQ; // use 22 bits in total
432                                 sumq[sc[0]&0x3f] += indelQ < seqQ? indelQ : seqQ;
433 //                              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);
434                         }
435                 }
436                 // determine bca->indel_types[] and bca->inscns
437                 bca->maxins = max_ins;
438                 bca->inscns = realloc(bca->inscns, bca->maxins * 4);
439                 for (t = 0; t < n_types; ++t)
440                         sumq[t] = sumq[t]<<6 | t;
441                 for (t = 1; t < n_types; ++t) // insertion sort
442                         for (j = t; j > 0 && sumq[j] > sumq[j-1]; --j)
443                                 tmp = sumq[j], sumq[j] = sumq[j-1], sumq[j-1] = tmp;
444                 for (t = 0; t < n_types; ++t) // look for the reference type
445                         if ((sumq[t]&0x3f) == ref_type) break;
446                 if (t) { // then move the reference type to the first
447                         tmp = sumq[t];
448                         for (; t > 0; --t) sumq[t] = sumq[t-1];
449                         sumq[0] = tmp;
450                 }
451                 for (t = 0; t < 4; ++t) bca->indel_types[t] = B2B_INDEL_NULL;
452                 for (t = 0; t < 4 && t < n_types; ++t) {
453                         bca->indel_types[t] = types[sumq[t]&0x3f];
454                         memcpy(&bca->inscns[t * bca->maxins], &inscns[(sumq[t]&0x3f) * max_ins], bca->maxins);
455                 }
456                 // update p->aux
457                 for (s = n_alt = 0; s < n; ++s) {
458                         for (i = 0; i < n_plp[s]; ++i) {
459                                 bam_pileup1_t *p = plp[s] + i;
460                                 int x = types[p->aux>>16&0x3f];
461                                 for (j = 0; j < 4; ++j)
462                                         if (x == bca->indel_types[j]) break;
463                                 p->aux = j<<16 | (j == 4? 0 : (p->aux&0xffff));
464                                 if ((p->aux>>16&0x3f) > 0) ++n_alt;
465 //                              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);
466                         }
467                 }               
468         }
469         free(score1); free(score2);
470         // free
471         for (i = 0; i < n; ++i) free(ref_sample[i]);
472         free(ref_sample);
473         free(types); free(inscns);
474         return n_alt > 0? 0 : -1;
475 }