]> git.donarmstrong.com Git - samtools.git/blobdiff - bam_index.c
works
[samtools.git] / bam_index.c
index 6fffd75e1ed5c58c363847abd42fed04e638c5b0..d6b94e2413e1dab2d2dc4ce2c91e471854ddd35a 100644 (file)
@@ -42,6 +42,8 @@
 // 1<<14 is the size of minimum bin.
 #define BAM_LIDX_SHIFT    14
 
+#define BAM_MAX_BIN 37450 // =(8^6-1)/7+1
+
 typedef struct {
        uint64_t u, v;
 } pair64_t;
@@ -63,6 +65,7 @@ KHASH_MAP_INIT_INT(i, bam_binlist_t)
 
 struct __bam_index_t {
        int32_t n;
+       uint64_t n_no_coor; // unmapped reads without coordinate
        khash_t(i) **index;
        bam_lidx_t *index2;
 };
@@ -117,7 +120,7 @@ static void merge_chunks(bam_index_t *idx)
                index = idx->index[i];
                for (k = kh_begin(index); k != kh_end(index); ++k) {
                        bam_binlist_t *p;
-                       if (!kh_exist(index, k)) continue;
+                       if (!kh_exist(index, k) || kh_key(index, k) == BAM_MAX_BIN) continue;
                        p = &kh_value(index, k);
                        m = 0;
                        for (l = 1; l < p->n; ++l) {
@@ -154,11 +157,16 @@ bam_index_t *bam_index_core(bamFile fp)
        uint32_t last_bin, save_bin;
        int32_t last_coor, last_tid, save_tid;
        bam1_core_t *c;
-       uint64_t save_off, last_off;
+       uint64_t save_off, last_off, n_mapped, n_unmapped, off_beg, off_end, n_no_coor;
+
+       h = bam_header_read(fp);
+       if(h == 0) {
+           fprintf(stderr, "[bam_index_core] Invalid BAM header.");
+           return NULL;
+       }
 
        idx = (bam_index_t*)calloc(1, sizeof(bam_index_t));
        b = (bam1_t*)calloc(1, sizeof(bam1_t));
-       h = bam_header_read(fp);
        c = &b->core;
 
        idx->n = h->n_targets;
@@ -169,19 +177,33 @@ bam_index_t *bam_index_core(bamFile fp)
 
        save_bin = save_tid = last_tid = last_bin = 0xffffffffu;
        save_off = last_off = bam_tell(fp); last_coor = 0xffffffffu;
+       n_mapped = n_unmapped = n_no_coor = off_end = 0;
+       off_beg = off_end = bam_tell(fp);
        while ((ret = bam_read1(fp, b)) >= 0) {
-               if (last_tid != c->tid) { // change of chromosomes
+               if (c->tid < 0) ++n_no_coor;
+               if (last_tid < c->tid || (last_tid >= 0 && c->tid < 0)) { // change of chromosomes
                        last_tid = c->tid;
                        last_bin = 0xffffffffu;
-               } else if (last_coor > c->pos) {
+               } else if ((uint32_t)last_tid > (uint32_t)c->tid) {
+                       fprintf(stderr, "[bam_index_core] the alignment is not sorted (%s): %d-th chr > %d-th chr\n",
+                                       bam1_qname(b), last_tid+1, c->tid+1);
+                       return NULL;
+               } else if ((int32_t)c->tid >= 0 && last_coor > c->pos) {
                        fprintf(stderr, "[bam_index_core] the alignment is not sorted (%s): %u > %u in %d-th chr\n",
                                        bam1_qname(b), last_coor, c->pos, c->tid+1);
-                       exit(1);
+                       return NULL;
                }
-               insert_offset2(&idx->index2[b->core.tid], b, last_off);
+               if (c->tid >= 0 && !(c->flag & BAM_FUNMAP)) insert_offset2(&idx->index2[b->core.tid], b, last_off);
                if (c->bin != last_bin) { // then possibly write the binning index
                        if (save_bin != 0xffffffffu) // save_bin==0xffffffffu only happens to the first record
                                insert_offset(idx->index[save_tid], save_bin, save_off, last_off);
+                       if (last_bin == 0xffffffffu && save_tid != 0xffffffffu) { // write the meta element
+                               off_end = last_off;
+                               insert_offset(idx->index[save_tid], BAM_MAX_BIN, off_beg, off_end);
+                               insert_offset(idx->index[save_tid], BAM_MAX_BIN, n_mapped, n_unmapped);
+                               n_mapped = n_unmapped = 0;
+                               off_beg = off_end;
+                       }
                        save_off = last_off;
                        save_bin = last_bin = c->bin;
                        save_tid = c->tid;
@@ -190,16 +212,32 @@ bam_index_t *bam_index_core(bamFile fp)
                if (bam_tell(fp) <= last_off) {
                        fprintf(stderr, "[bam_index_core] bug in BGZF/RAZF: %llx < %llx\n",
                                        (unsigned long long)bam_tell(fp), (unsigned long long)last_off);
-                       exit(1);
+                       return NULL;
                }
+               if (c->flag & BAM_FUNMAP) ++n_unmapped;
+               else ++n_mapped;
                last_off = bam_tell(fp);
                last_coor = b->core.pos;
        }
-       if (save_tid >= 0) insert_offset(idx->index[save_tid], save_bin, save_off, bam_tell(fp));
+       if (save_tid >= 0) {
+               insert_offset(idx->index[save_tid], save_bin, save_off, bam_tell(fp));
+               insert_offset(idx->index[save_tid], BAM_MAX_BIN, off_beg, bam_tell(fp));
+               insert_offset(idx->index[save_tid], BAM_MAX_BIN, n_mapped, n_unmapped);
+       }
        merge_chunks(idx);
        fill_missing(idx);
+       if (ret >= 0) {
+               while ((ret = bam_read1(fp, b)) >= 0) {
+                       ++n_no_coor;
+                       if (c->tid >= 0 && n_no_coor) {
+                               fprintf(stderr, "[bam_index_core] the alignment is not sorted: reads without coordinates prior to reads with coordinates.\n");
+                               return NULL;
+                       }
+               }
+       }
        if (ret < -1) fprintf(stderr, "[bam_index_core] truncated file? Continue anyway. (%d)\n", ret);
        free(b->data); free(b);
+       idx->n_no_coor = n_no_coor;
        return idx;
 }
 
@@ -277,6 +315,11 @@ void bam_index_save(const bam_index_t *idx, FILE *fp)
                                bam_swap_endian_8p(&index2->offset[x]);
                } else fwrite(index2->offset, 8, index2->n, fp);
        }
+       { // write the number of reads coor-less records.
+               uint64_t x = idx->n_no_coor;
+               if (bam_is_be) bam_swap_endian_8p(&x);
+               fwrite(&x, 8, 1, fp);
+       }
        fflush(fp);
 }
 
@@ -338,6 +381,8 @@ static bam_index_t *bam_index_load_core(FILE *fp)
                if (bam_is_be)
                        for (j = 0; j < index2->n; ++j) bam_swap_endian_8p(&index2->offset[j]);
        }
+       if (fread(&idx->n_no_coor, 8, 1, fp) == 0) idx->n_no_coor = 0;
+       if (bam_is_be) bam_swap_endian_8p(&idx->n_no_coor);
        return idx;
 }
 
@@ -437,6 +482,10 @@ int bam_index_build2(const char *fn, const char *_fnidx)
        }
        idx = bam_index_core(fp);
        bam_close(fp);
+       if(idx == 0) {
+               fprintf(stderr, "[bam_index_build2] fail to index the BAM file.\n");
+               return -1;
+       }
        if (_fnidx == 0) {
                fnidx = (char*)calloc(strlen(fn) + 5, 1);
                strcpy(fnidx, fn); strcat(fnidx, ".bai");
@@ -462,7 +511,7 @@ int bam_index_build(const char *fn)
 int bam_index(int argc, char *argv[])
 {
        if (argc < 2) {
-               fprintf(stderr, "Usage: samtools index <in.bam> [<out.index>]\n");
+               fprintf(stderr, "Usage: samtools index <in.bam> [out.index]\n");
                return 1;
        }
        if (argc >= 3) bam_index_build2(argv[1], argv[2]);
@@ -470,9 +519,39 @@ int bam_index(int argc, char *argv[])
        return 0;
 }
 
-#define MAX_BIN 37450 // =(8^6-1)/7+1
+int bam_idxstats(int argc, char *argv[])
+{
+       bam_index_t *idx;
+       bam_header_t *header;
+       bamFile fp;
+       int i;
+       if (argc < 2) {
+               fprintf(stderr, "Usage: samtools idxstats <in.bam>\n");
+               return 1;
+       }
+       fp = bam_open(argv[1], "r");
+       if (fp == 0) { fprintf(stderr, "[%s] fail to open BAM.\n", __func__); return 1; }
+       header = bam_header_read(fp);
+       bam_close(fp);
+       idx = bam_index_load(argv[1]);
+       if (idx == 0) { fprintf(stderr, "[%s] fail to load the index.\n", __func__); return 1; }
+       for (i = 0; i < idx->n; ++i) {
+               khint_t k;
+               khash_t(i) *h = idx->index[i];
+               printf("%s\t%d", header->target_name[i], header->target_len[i]);
+               k = kh_get(i, h, BAM_MAX_BIN);
+               if (k != kh_end(h))
+                       printf("\t%llu\t%llu", (long long)kh_val(h, k).list[1].u, (long long)kh_val(h, k).list[1].v);
+               else printf("\t0\t0");
+               putchar('\n');
+       }
+       printf("*\t0\t0\t%llu\n", (long long)idx->n_no_coor);
+       bam_header_destroy(header);
+       bam_index_destroy(idx);
+       return 0;
+}
 
-static inline int reg2bins(uint32_t beg, uint32_t end, uint16_t list[MAX_BIN])
+static inline int reg2bins(uint32_t beg, uint32_t end, uint16_t list[BAM_MAX_BIN])
 {
        int i = 0, k;
        if (beg >= end) return 0;
@@ -518,7 +597,7 @@ bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end)
        iter = calloc(1, sizeof(struct __bam_iter_t));
        iter->tid = tid, iter->beg = beg, iter->end = end; iter->i = -1;
        //
-       bins = (uint16_t*)calloc(MAX_BIN, 2);
+       bins = (uint16_t*)calloc(BAM_MAX_BIN, 2);
        n_bins = reg2bins(beg, end, bins);
        index = idx->index[tid];
        if (idx->index2[tid].n > 0) {
@@ -549,6 +628,9 @@ bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end)
                }
        }
        free(bins);
+       if (n_off == 0) {
+               free(off); return iter;
+       }
        {
                bam1_t *b = (bam1_t*)calloc(1, sizeof(bam1_t));
                int l;
@@ -597,17 +679,17 @@ void bam_iter_destroy(bam_iter_t iter)
 
 int bam_iter_read(bamFile fp, bam_iter_t iter, bam1_t *b)
 {
-       if (iter->finished) return -1;
-       if (iter->from_first) {
-               int ret = bam_read1(fp, b);
-               if (ret < 0) iter->finished = 1;
+       int ret;
+       if (iter && iter->finished) return -1;
+       if (iter == 0 || iter->from_first) {
+               ret = bam_read1(fp, b);
+               if (ret < 0 && iter) iter->finished = 1;
                return ret;
        }
        if (iter->off == 0) return -1;
        for (;;) {
-               int ret;
                if (iter->curr_off == 0 || iter->curr_off >= iter->off[iter->i].v) { // then jump to the next chunk
-                       if (iter->i == iter->n_off - 1) break; // no more chunks
+                       if (iter->i == iter->n_off - 1) { ret = -1; break; } // no more chunks
                        if (iter->i >= 0) assert(iter->curr_off == iter->off[iter->i].v); // otherwise bug
                        if (iter->i < 0 || iter->off[iter->i].v != iter->off[iter->i+1].u) { // not adjacent chunks; then seek
                                bam_seek(fp, iter->off[iter->i+1].u, SEEK_SET);
@@ -615,23 +697,28 @@ int bam_iter_read(bamFile fp, bam_iter_t iter, bam1_t *b)
                        }
                        ++iter->i;
                }
-               if ((ret = bam_read1(fp, b)) > 0) {
+               if ((ret = bam_read1(fp, b)) >= 0) {
                        iter->curr_off = bam_tell(fp);
-                       if (b->core.tid != iter->tid || b->core.pos >= iter->end) break; // no need to proceed
+                       if (b->core.tid != iter->tid || b->core.pos >= iter->end) { // no need to proceed
+                               ret = bam_validate1(NULL, b)? -1 : -5; // determine whether end of region or error
+                               break;
+                       }
                        else if (is_overlap(iter->beg, iter->end, b)) return ret;
-               } else break; // end of file
+               } else break; // end of file or error
        }
        iter->finished = 1;
-       return -1;
+       return ret;
 }
 
 int bam_fetch(bamFile fp, const bam_index_t *idx, int tid, int beg, int end, void *data, bam_fetch_f func)
 {
+       int ret;
        bam_iter_t iter;
        bam1_t *b;
        b = bam_init1();
        iter = bam_iter_query(idx, tid, beg, end);
-       while (bam_iter_read(fp, iter, b) >= 0) func(b, data);
+       while ((ret = bam_iter_read(fp, iter, b)) >= 0) func(b, data);
+       bam_iter_destroy(iter);
        bam_destroy1(b);
-       return 0;
+       return (ret == -1)? 0 : ret;
 }