]> git.donarmstrong.com Git - samtools.git/blobdiff - bam2bcf_indel.c
* samtools-0.1.9-8 (r809)
[samtools.git] / bam2bcf_indel.c
index 2d0b5c6f8a66a4c5c006dbec560718ab1f7a0377..c28b6859dd0101dc5f9ded34a42ded1add58f781 100644 (file)
@@ -110,7 +110,7 @@ int bcf_call_gap_prep(int n, int *n_plp, bam_pileup1_t **plp, int pos, bcf_calla
                                          const void *rghash)
 {
        extern void ks_introsort_uint32_t(int, uint32_t*);
-       int i, s, j, k, t, n_types, *types, max_rd_len, left, right, max_ins, *score, N, K, l_run, ref_type;
+       int i, s, j, k, t, n_types, *types, max_rd_len, left, right, max_ins, *score, N, K, l_run, ref_type, n_alt;
        char *inscns = 0, *ref2, *query;
        khash_t(rg) *hash = (khash_t(rg)*)rghash;
        if (ref == 0 || bca == 0) return -1;
@@ -342,13 +342,14 @@ int bcf_call_gap_prep(int n, int *n_plp, bam_pileup1_t **plp, int pos, bcf_calla
                        memcpy(&bca->inscns[t * bca->maxins], &inscns[(sumq[t]&0x3f) * max_ins], bca->maxins);
                }
                // update p->aux
-               for (s = K = 0; s < n; ++s) {
-                       for (i = 0; i < n_plp[s]; ++i, ++K) {
+               for (s = n_alt = 0; s < n; ++s) {
+                       for (i = 0; i < n_plp[s]; ++i) {
                                bam_pileup1_t *p = plp[s] + i;
                                int x = types[p->aux>>16&0x3f];
                                for (j = 0; j < 4; ++j)
                                        if (x == bca->indel_types[j]) break;
                                p->aux = j<<16 | (j == 4? 0 : (p->aux&0xffff));
+                               if ((p->aux>>16&0x3f) > 0) ++n_alt;
 //                             fprintf(stderr, "pos=%d read=%d:%d name=%s call=%d q=%d\n", pos, s, i, bam1_qname(p->b), p->aux>>16&63, p->aux&0xff);
                        }
                }               
@@ -356,5 +357,5 @@ int bcf_call_gap_prep(int n, int *n_plp, bam_pileup1_t **plp, int pos, bcf_calla
        free(score);
        // free
        free(types); free(inscns);
-       return 0;
+       return n_alt > 0? 0 : -1;
 }