knetfile.c
[samtools.git] / bam_index.c
index 328f011c16ccf4a9a558702317ecdecfbcc5212a..c7f33959bd99ddef3a9a628d885e450fc2e07acf 100644 (file)
@@ -80,11 +80,17 @@ static inline void insert_offset(khash_t(i) *h, int bin, uint64_t beg, uint64_t
        l = &kh_value(h, k);
        if (ret) { // not present
                l->m = 1; l->n = 0;
-               l->list = (pair64_t*)calloc(l->m, 16);
+               if ((l->list = (pair64_t*)calloc(l->m, 16)) == NULL) {
+                       fprintf(stderr, "[%s] failed to allocate bam_binlist node.\n", __func__);
+                       abort();
+               }
        }
        if (l->n == l->m) {
                l->m <<= 1;
-               l->list = (pair64_t*)realloc(l->list, l->m * 16);
+               if ((l->list = (pair64_t*)realloc(l->list, l->m * 16)) == NULL) {
+                       fprintf(stderr, "[%s] failed to resize bam_binlist node.\n", __func__);
+                       abort();
+               }
        }
        l->list[l->n].u = beg; l->list[l->n++].v = end;
 }
@@ -98,7 +104,10 @@ static inline void insert_offset2(bam_lidx_t *index2, bam1_t *b, uint64_t offset
                int old_m = index2->m;
                index2->m = end + 1;
                kroundup32(index2->m);
-               index2->offset = (uint64_t*)realloc(index2->offset, index2->m * 8);
+               if ((index2->offset = (uint64_t*)realloc(index2->offset, index2->m * 8)) == NULL) {
+                       fprintf(stderr, "[%s] failed to resize offset array.\n", __func__);
+                       abort();
+               }
                memset(index2->offset + old_m, 0, 8 * (index2->m - old_m));
        }
        if (beg == end) {
@@ -159,32 +168,54 @@ bam_index_t *bam_index_core(bamFile fp)
        bam1_core_t *c;
        uint64_t save_off, last_off, n_mapped, n_unmapped, off_beg, off_end, n_no_coor;
 
-       idx = (bam_index_t*)calloc(1, sizeof(bam_index_t));
-       b = (bam1_t*)calloc(1, sizeof(bam1_t));
+       if ((idx = (bam_index_t*)calloc(1, sizeof(bam_index_t))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam_index structure.\n", __func__);
+               return NULL;
+       }
+       if ((b = (bam1_t*)calloc(1, sizeof(bam1_t))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam1_t buffer.\n", __func__);
+               free(idx);
+               return NULL;
+       }
        h = bam_header_read(fp);
        c = &b->core;
 
        idx->n = h->n_targets;
        bam_header_destroy(h);
-       idx->index = (khash_t(i)**)calloc(idx->n, sizeof(void*));
+       if ((idx->index = (khash_t(i)**)calloc(idx->n, sizeof(void*))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam_index bin hash.\n", __func__);
+               free(b);
+               free(idx);
+               return NULL;
+       }
        for (i = 0; i < idx->n; ++i) idx->index[i] = kh_init(i);
-       idx->index2 = (bam_lidx_t*)calloc(idx->n, sizeof(bam_lidx_t));
+       if ((idx->index2 = (bam_lidx_t*)calloc(idx->n, sizeof(bam_lidx_t))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam_index linear index.\n", __func__);
+               free(idx->index);
+               free(b);
+               free(idx);
+               return NULL;
+       }
 
        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;
+       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 (c->tid < 0) ++n_no_coor;
-               if (last_tid != c->tid) { // change of chromosomes
+               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;
                }
-               if (c->tid >= 0) 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);
@@ -203,7 +234,7 @@ 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;
@@ -217,8 +248,15 @@ bam_index_t *bam_index_core(bamFile fp)
        }
        merge_chunks(idx);
        fill_missing(idx);
-       if (ret >= 0)
-               while ((ret = bam_read1(fp, b)) >= 0) ++n_no_coor;
+       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;
@@ -248,11 +286,23 @@ void bam_index_save(const bam_index_t *idx, FILE *fp)
 {
        int32_t i, size;
        khint_t k;
-       fwrite("BAI\1", 1, 4, fp);
+       if (fwrite("BAI\1", 1, 4, fp) < 4) {
+               fprintf(stderr, "[%s] failed to write magic number.\n", __func__);
+               return;
+       }
        if (bam_is_be) {
                uint32_t x = idx->n;
-               fwrite(bam_swap_endian_4p(&x), 4, 1, fp);
-       } else fwrite(&idx->n, 4, 1, fp);
+               if (fwrite(bam_swap_endian_4p(&x), 4, 1, fp) < 1) {
+                       fprintf(stderr, "[%s] failed to write n_ref.\n", __func__);
+                       return;
+               }
+       }
+       else {
+               if (fwrite(&idx->n, 4, 1, fp) < 1) {
+                       fprintf(stderr, "[%s] failed to write n_ref.\n", __func__);
+                       return;
+               }
+       }
        for (i = 0; i < idx->n; ++i) {
                khash_t(i) *index = idx->index[i];
                bam_lidx_t *index2 = idx->index2 + i;
@@ -260,49 +310,100 @@ void bam_index_save(const bam_index_t *idx, FILE *fp)
                size = kh_size(index);
                if (bam_is_be) { // big endian
                        uint32_t x = size;
-                       fwrite(bam_swap_endian_4p(&x), 4, 1, fp);
-               } else fwrite(&size, 4, 1, fp);
+                       if (fwrite(bam_swap_endian_4p(&x), 4, 1, fp) < 1) {
+                               fprintf(stderr, "[%s] failed to write n_bin.\n", __func__);
+                               return;
+                       }
+               }
+               else {
+                       if (fwrite(&size, 4, 1, fp) < 1) {
+                               fprintf(stderr, "[%s] failed to write n_bin.\n", __func__);
+                               return;
+                       }
+               }
                for (k = kh_begin(index); k != kh_end(index); ++k) {
                        if (kh_exist(index, k)) {
                                bam_binlist_t *p = &kh_value(index, k);
                                if (bam_is_be) { // big endian
                                        uint32_t x;
-                                       x = kh_key(index, k); fwrite(bam_swap_endian_4p(&x), 4, 1, fp);
-                                       x = p->n; fwrite(bam_swap_endian_4p(&x), 4, 1, fp);
+                                       x = kh_key(index, k);
+                                       if (fwrite(bam_swap_endian_4p(&x), 4, 1, fp) < 1) {
+                                               fprintf(stderr, "[%s] failed to write bin.\n", __func__);
+                                               return;
+                                       }
+                                       x = p->n;
+                                       if (fwrite(bam_swap_endian_4p(&x), 4, 1, fp) < 1) {
+                                               fprintf(stderr, "[%s] failed to write n_chunk.\n", __func__);
+                                               return;
+                                       }
                                        for (x = 0; (int)x < p->n; ++x) {
                                                bam_swap_endian_8p(&p->list[x].u);
                                                bam_swap_endian_8p(&p->list[x].v);
                                        }
-                                       fwrite(p->list, 16, p->n, fp);
+                                       if (fwrite(p->list, 16, p->n, fp) < p->n) {
+                                               fprintf(stderr, "[%s] failed to write %d chunks.\n", __func__, p->n);
+                                               return;
+                                       }
+                                       /*
                                        for (x = 0; (int)x < p->n; ++x) {
                                                bam_swap_endian_8p(&p->list[x].u);
                                                bam_swap_endian_8p(&p->list[x].v);
+                                       } */
+                               }
+                               else {
+                                       if (fwrite(&kh_key(index, k), 4, 1, fp) < 1) {
+                                               fprintf(stderr, "[%s] failed to write bin.\n", __func__);
+                                               return;
+                                       }
+                                       if (fwrite(&p->n, 4, 1, fp) < 1) {
+                                               fprintf(stderr, "[%s] failed to write n_chunk.\n", __func__);
+                                               return;
+                                       }
+                                       if (fwrite(p->list, 16, p->n, fp) < p->n) {
+                                               fprintf(stderr, "[%s] failed to write %d chunks.\n", __func__, p->n);
+                                               return;
                                        }
-                               } else {
-                                       fwrite(&kh_key(index, k), 4, 1, fp);
-                                       fwrite(&p->n, 4, 1, fp);
-                                       fwrite(p->list, 16, p->n, fp);
                                }
                        }
                }
                // write linear index (index2)
                if (bam_is_be) {
                        int x = index2->n;
-                       fwrite(bam_swap_endian_4p(&x), 4, 1, fp);
-               } else fwrite(&index2->n, 4, 1, fp);
+                       if (fwrite(bam_swap_endian_4p(&x), 4, 1, fp) < 1) {
+                               fprintf(stderr, "[%s] failed to write n_intv.\n", __func__);
+                               return;
+                       }
+               }
+               else {
+                       if (fwrite(&index2->n, 4, 1, fp) < 1) {
+                               fprintf(stderr, "[%s] failed to write n_intv.\n", __func__);
+                               return;
+                       }
+               }
                if (bam_is_be) { // big endian
                        int x;
                        for (x = 0; (int)x < index2->n; ++x)
                                bam_swap_endian_8p(&index2->offset[x]);
-                       fwrite(index2->offset, 8, index2->n, fp);
+                       if (fwrite(index2->offset, 8, index2->n, fp) < index2->n) {
+                               fprintf(stderr, "[%s] failed to write ioffset.\n", __func__);
+                               return;
+                       }
                        for (x = 0; (int)x < index2->n; ++x)
                                bam_swap_endian_8p(&index2->offset[x]);
-               } else fwrite(index2->offset, 8, index2->n, fp);
+               }
+               else {
+                       if (fwrite(index2->offset, 8, index2->n, fp) < index2->n) {
+                               fprintf(stderr, "[%s] failed to write ioffset.\n", __func__);
+                               return;
+                       }
+               }
        }
        { // 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);
+               if (fwrite(&x, 8, 1, fp) < 1) {
+                       fprintf(stderr, "[%s] failed to write n_no_coor.\n", __func__);
+               }
        }
        fflush(fp);
 }
@@ -316,17 +417,42 @@ static bam_index_t *bam_index_load_core(FILE *fp)
                fprintf(stderr, "[bam_index_load_core] fail to load index.\n");
                return 0;
        }
-       fread(magic, 1, 4, fp);
+       if (fread(magic, 1, 4, fp) < 4) {
+               fprintf(stderr, "[%s] failed to read magic number.\n", __func__);
+               fclose(fp);
+               return 0;
+       }
        if (strncmp(magic, "BAI\1", 4)) {
                fprintf(stderr, "[bam_index_load] wrong magic number.\n");
                fclose(fp);
                return 0;
        }
-       idx = (bam_index_t*)calloc(1, sizeof(bam_index_t));     
-       fread(&idx->n, 4, 1, fp);
+       if ((idx = (bam_index_t*)calloc(1, sizeof(bam_index_t))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam_index structure.\n", __func__);
+               fclose(fp);
+               return NULL;
+       }
+       if (fread(&idx->n, 4, 1, fp) < 1) {
+               fprintf(stderr, "[%s] failed to read n_ref.\n", __func__);
+               free(idx);
+               fclose(fp);
+               return 0;
+       }
        if (bam_is_be) bam_swap_endian_4p(&idx->n);
-       idx->index = (khash_t(i)**)calloc(idx->n, sizeof(void*));
-       idx->index2 = (bam_lidx_t*)calloc(idx->n, sizeof(bam_lidx_t));
+       if ((idx->index = (khash_t(i)**)calloc(idx->n, sizeof(void*))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam_index bin hash.\n", __func__);
+               free(idx);
+               fclose(fp);
+               return NULL;
+       }
+       if ((idx->index2 = (bam_lidx_t*)calloc(idx->n, sizeof(bam_lidx_t))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam_index linear index.\n", __func__);
+               free(idx->index);
+               free(idx);
+               fclose(fp);
+               return NULL;
+       }
+
        for (i = 0; i < idx->n; ++i) {
                khash_t(i) *index;
                bam_lidx_t *index2 = idx->index2 + i;
@@ -336,18 +462,48 @@ static bam_index_t *bam_index_load_core(FILE *fp)
                bam_binlist_t *p;
                index = idx->index[i] = kh_init(i);
                // load binning index
-               fread(&size, 4, 1, fp);
+               if (fread(&size, 4, 1, fp) < 1) {
+                       fprintf(stderr, "[%s] failed to read n_bin.\n", __func__);
+                       free(idx->index);
+                       free(idx);
+                       fclose(fp);
+                       return 0;
+               }
                if (bam_is_be) bam_swap_endian_4p(&size);
                for (j = 0; j < (int)size; ++j) {
-                       fread(&key, 4, 1, fp);
+                       if (fread(&key, 4, 1, fp) < 1) {
+                               fprintf(stderr, "[%s] failed to read bin.\n", __func__);
+                               free(idx->index);
+                               free(idx);
+                               fclose(fp);
+                               return 0;
+                       }
                        if (bam_is_be) bam_swap_endian_4p(&key);
                        k = kh_put(i, index, key, &ret);
                        p = &kh_value(index, k);
-                       fread(&p->n, 4, 1, fp);
+                       if (fread(&p->n, 4, 1, fp) < 1) {
+                               fprintf(stderr, "[%s] failed to read n_chunk.\n", __func__);
+                               free(idx->index);
+                               free(idx);
+                               fclose(fp);
+                               return 0;
+                       }
                        if (bam_is_be) bam_swap_endian_4p(&p->n);
                        p->m = p->n;
-                       p->list = (pair64_t*)malloc(p->m * 16);
-                       fread(p->list, 16, p->n, fp);
+                       if ((p->list = (pair64_t*)malloc(p->m * 16)) == NULL) {
+                               fprintf(stderr, "[%s] failed to allocate chunk array.\n", __func__);
+                               free(idx->index);
+                               free(idx);
+                               fclose(fp);
+                               return 0;
+                       }
+                       if (fread(p->list, 16, p->n, fp) < p->n) {
+                               fprintf(stderr, "[%s] failed to read %d chunks.\n", __func__, p->n);
+                               free(idx->index);
+                               free(idx);
+                               fclose(fp);
+                               return 0;
+                       }
                        if (bam_is_be) {
                                int x;
                                for (x = 0; x < p->n; ++x) {
@@ -357,13 +513,37 @@ static bam_index_t *bam_index_load_core(FILE *fp)
                        }
                }
                // load linear index
-               fread(&index2->n, 4, 1, fp);
+               if (fread(&index2->n, 4, 1, fp) < 1) {
+                       fprintf(stderr, "[%s] failed to read n_intv", __func__);
+                       free(idx->index);
+                       free(idx);
+                       fclose(fp);
+                       return 0;
+               }
                if (bam_is_be) bam_swap_endian_4p(&index2->n);
-               index2->m = index2->n;
-               index2->offset = (uint64_t*)calloc(index2->m, 8);
-               fread(index2->offset, index2->n, 8, fp);
-               if (bam_is_be)
-                       for (j = 0; j < index2->n; ++j) bam_swap_endian_8p(&index2->offset[j]);
+               if (index2->n > 0) {
+                       index2->m = index2->n;
+                       if ((index2->offset = (uint64_t*)calloc(index2->m, 8)) == NULL) {
+                               fprintf(stderr, "[%s] falied to allocate interval array.\n", __func__);
+                               free(idx->index);
+                               free(idx);
+                               fclose(fp);
+                               return 0;
+                       }
+                       if (fread(index2->offset, index2->n, 8, fp) < index2->n) {
+                               fprintf(stderr, "[%s] failed to read %d intervals.", __func__, index2->n);
+                               free(index2->offset);
+                               free(idx->index);
+                               free(idx);
+                               fclose(fp);
+                               return 0;
+                       }
+                       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);
@@ -382,7 +562,11 @@ bam_index_t *bam_index_load_local(const char *_fn)
                        if (*p == '/') break;
                fn = strdup(p + 1);
        } else fn = strdup(_fn);
-       fnidx = (char*)calloc(strlen(fn) + 5, 1);
+       if ((fnidx = (char*)calloc(strlen(fn) + 5, 1)) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate space for index file name.\n", __func__);
+               free(fn);
+               return NULL;
+       }
        strcpy(fnidx, fn); strcat(fnidx, ".bai");
        fp = fopen(fnidx, "rb");
        if (fp == 0) { // try "{base}.bai"
@@ -398,7 +582,11 @@ bam_index_t *bam_index_load_local(const char *_fn)
                bam_index_t *idx = bam_index_load_core(fp);
                fclose(fp);
                return idx;
-       } else return 0;
+       }
+       else {
+               fprintf(stderr, "[%s] failed to open index file.\n", __func__);
+               return 0;
+       }
 }
 
 #ifdef _USE_KNETFILE
@@ -417,17 +605,25 @@ static void download_from_remote(const char *url)
        ++fn; // fn now points to the file name
        fp_remote = knet_open(url, "r");
        if (fp_remote == 0) {
-               fprintf(stderr, "[download_from_remote] fail to open remote file.\n");
+               fprintf(stderr, "[%s] failed to open remote file.\n", __func__);
                return;
        }
        if ((fp = fopen(fn, "wb")) == 0) {
-               fprintf(stderr, "[download_from_remote] fail to create file in the working directory.\n");
+               fprintf(stderr, "[%s] failed to create file in the working directory.\n", __func__);
                knet_close(fp_remote);
                return;
        }
-       buf = (uint8_t*)calloc(buf_size, 1);
-       while ((l = knet_read(fp_remote, buf, buf_size)) != 0)
-               fwrite(buf, 1, l, fp);
+       if ((buf = (uint8_t*)calloc(buf_size, 1)) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate buffer.\n", __func__);
+               fclose(fp);
+               return;
+       }
+       while ((l = knet_read(fp_remote, buf, buf_size)) != 0) {
+               if (fwrite(buf, 1, l, fp) < l) {
+                       fprintf(stderr, "[%s] failed to write to destination file.\n", __func__);
+                       break;
+               }
+       }
        free(buf);
        fclose(fp);
        knet_close(fp_remote);
@@ -445,12 +641,16 @@ bam_index_t *bam_index_load(const char *fn)
        idx = bam_index_load_local(fn);
        if (idx == 0 && (strstr(fn, "ftp://") == fn || strstr(fn, "http://") == fn)) {
                char *fnidx = calloc(strlen(fn) + 5, 1);
+               if (fnidx == NULL) {
+                       fprintf(stderr, "[%s] failed to allocate space for index filename.\n", __func__);
+                       return NULL;
+               }
                strcat(strcpy(fnidx, fn), ".bai");
-               fprintf(stderr, "[bam_index_load] attempting to download the remote index file.\n");
+               fprintf(stderr, "[%s] attempting to download the remote index file.\n", __func__);
                download_from_remote(fnidx);
                idx = bam_index_load_local(fn);
        }
-       if (idx == 0) fprintf(stderr, "[bam_index_load] fail to load BAM index.\n");
+       if (idx == 0) fprintf(stderr, "[%s] fail to load BAM index.\n", __func__);
        return idx;
 }
 
@@ -461,18 +661,22 @@ int bam_index_build2(const char *fn, const char *_fnidx)
        bamFile fp;
        bam_index_t *idx;
        if ((fp = bam_open(fn, "r")) == 0) {
-               fprintf(stderr, "[bam_index_build2] fail to open the BAM file.\n");
+               fprintf(stderr, "[%s] failure to open the BAM file.\n", __func__);
                return -1;
        }
        idx = bam_index_core(fp);
        bam_close(fp);
+       if(idx == 0) {
+               fprintf(stderr, "[%s] fail to index the BAM file.\n", __func__);
+               return -1;
+       }
        if (_fnidx == 0) {
                fnidx = (char*)calloc(strlen(fn) + 5, 1);
                strcpy(fnidx, fn); strcat(fnidx, ".bai");
        } else fnidx = strdup(_fnidx);
        fpidx = fopen(fnidx, "wb");
        if (fpidx == 0) {
-               fprintf(stderr, "[bam_index_build2] fail to create the index file.\n");
+               fprintf(stderr, "[%s] fail to create the index file.\n", __func__);
                free(fnidx);
                return -1;
        }
@@ -574,17 +778,23 @@ bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end)
        if (beg < 0) beg = 0;
        if (end < beg) return 0;
        // initialize iter
-       iter = calloc(1, sizeof(struct __bam_iter_t));
+       if ((iter = calloc(1, sizeof(struct __bam_iter_t))) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bam iterator.\n", __func__);
+               abort();
+       }
        iter->tid = tid, iter->beg = beg, iter->end = end; iter->i = -1;
        //
-       bins = (uint16_t*)calloc(BAM_MAX_BIN, 2);
+       if ((bins = (uint16_t*)calloc(BAM_MAX_BIN, 2)) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate bin array.\n", __func__);
+               abort();
+       }
        n_bins = reg2bins(beg, end, bins);
        index = idx->index[tid];
        if (idx->index2[tid].n > 0) {
                min_off = (beg>>BAM_LIDX_SHIFT >= idx->index2[tid].n)? idx->index2[tid].offset[idx->index2[tid].n-1]
                        : idx->index2[tid].offset[beg>>BAM_LIDX_SHIFT];
                if (min_off == 0) { // improvement for index files built by tabix prior to 0.1.4
-                       int n = beg>>BAM_LIDX_SHIFT;
+                       int n = beg >> BAM_LIDX_SHIFT;
                        if (n > idx->index2[tid].n) n = idx->index2[tid].n;
                        for (i = n - 1; i >= 0; --i)
                                if (idx->index2[tid].offset[i] != 0) break;
@@ -598,7 +808,10 @@ bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end)
        if (n_off == 0) {
                free(bins); return iter;
        }
-       off = (pair64_t*)calloc(n_off, 16);
+       if ((off = (pair64_t*)calloc(n_off, 16)) == NULL) {
+               fprintf(stderr, "[%s] failed to allocate offset pair array.\n", __func__);
+               abort();
+       }
        for (i = n_off = 0; i < n_bins; ++i) {
                if ((k = kh_get(i, index, bins[i])) != kh_end(index)) {
                        int j;
@@ -613,6 +826,10 @@ bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end)
        }
        {
                bam1_t *b = (bam1_t*)calloc(1, sizeof(bam1_t));
+               if (b == NULL) {
+                       fprintf(stderr, "[%s] failure to allocate bam1_t buffer.\n", __func__);
+                       abort();
+               }
                int l;
                ks_introsort(off, n_off, off);
                // resolve completely contained adjacent blocks