]> git.donarmstrong.com Git - samtools.git/commitdiff
* samtools-0.1.4-21 (r368)
authorHeng Li <lh3@live.co.uk>
Thu, 2 Jul 2009 10:32:26 +0000 (10:32 +0000)
committerHeng Li <lh3@live.co.uk>
Thu, 2 Jul 2009 10:32:26 +0000 (10:32 +0000)
 * propagate errors rather than exit or complain assertion failure. Assertion
   should be only used for checking internal bugs, but not for external input
   inconsistency. I was just a bit lazy.
 * small memory leak may be present on failure, though

14 files changed:
bam.c
bam.h
bam_aux.c
bam_index.c
bam_lpileup.c
bam_md.c
bam_pileup.c
bam_rmdup.c
bam_stat.c
bam_tview.c
bamtk.c
faidx.c
faidx.h
glf.c

diff --git a/bam.c b/bam.c
index 7f1ecfd26217cd5f2e3b6ea4b0e1c58cb6d4d06b..1ff4a5aba25028f4fccdd0735a9dcae594b6e9a8 100644 (file)
--- a/bam.c
+++ b/bam.c
@@ -1,5 +1,6 @@
 #include <stdio.h>
 #include <ctype.h>
+#include <assert.h>
 #include "bam.h"
 #include "bam_endian.h"
 #include "kstring.h"
diff --git a/bam.h b/bam.h
index 6316c9b438320dc5df8751f5c1a73ed6c185b6e4..a2ceb170e18e6ce135d126ef6111f7eccb54620e 100644 (file)
--- a/bam.h
+++ b/bam.h
@@ -41,7 +41,6 @@
  */
 
 #include <stdint.h>
-#include <assert.h>
 #include <stdlib.h>
 #include <string.h>
 #include <stdio.h>
@@ -600,8 +599,9 @@ extern "C" {
          @param  ref_id  the returned chromosome ID
          @param  begin   the returned start coordinate
          @param  end     the returned end coordinate
+         @return         0 on success; -1 on failure
         */
-       void bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *begin, int *end);
+       int bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *begin, int *end);
 
        /*!
          @abstract       Retrieve data of a tag
index a63e2aeadccf35e95b786097cd33394cb31a7a06..9b918424eeaf046719105956ea7baaef29b6cd03 100644 (file)
--- a/bam_aux.c
+++ b/bam_aux.c
@@ -72,7 +72,7 @@ int32_t bam_get_tid(const bam_header_t *header, const char *seq_name)
        return k == kh_end(h)? -1 : kh_value(h, k);
 }
 
-void bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *begin, int *end)
+int bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *begin, int *end)
 {
        char *s, *p;
        int i, l, k;
@@ -93,12 +93,12 @@ void bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *b
        iter = kh_get(s, h, s); /* get the ref_id */
        if (iter == kh_end(h)) { // name not found
                *ref_id = -1; free(s);
-               return;
+               return -1;
        }
        *ref_id = kh_value(h, iter);
        if (i == k) { /* dump the whole sequence */
                *begin = 0; *end = 1<<29; free(s);
-               return;
+               return -1;
        }
        for (p = s + i + 1; i != k; ++i) if (s[i] == '-') break;
        *begin = atoi(p);
@@ -107,8 +107,12 @@ void bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *b
                *end = atoi(p);
        } else *end = 1<<29;
        if (*begin > 0) --*begin;
-       assert(*begin <= *end);
        free(s);
+       if (*begin > *end) {
+               fprintf(stderr, "[bam_parse_region] invalid region.\n");
+               return -1;
+       }
+       return 0;
 }
 
 int32_t bam_aux2i(const uint8_t *s)
index 329df323851f7bd61ecbbc345790f122b4303b6f..72ef270f9c7afd143759e1f98d681ce69147bb13 100644 (file)
@@ -1,4 +1,5 @@
 #include <ctype.h>
+#include <assert.h>
 #include "bam.h"
 #include "khash.h"
 #include "ksort.h"
@@ -324,7 +325,6 @@ static bam_index_t *bam_index_load_core(FILE *fp)
 
 bam_index_t *bam_index_load_local(const char *_fn)
 {
-       bam_index_t *idx;
        FILE *fp;
        char *fnidx, *fn;
 
@@ -347,9 +347,11 @@ bam_index_t *bam_index_load_local(const char *_fn)
                }
        }
        free(fnidx); free(fn);
-       idx = bam_index_load_core(fp);
-       fclose(fp);
-       return idx;
+       if (fp) {
+               bam_index_t *idx = bam_index_load_core(fp);
+               fclose(fp);
+               return idx;
+       } else return 0;
 }
 
 static void download_from_remote(const char *url)
@@ -394,6 +396,7 @@ bam_index_t *bam_index_load(const char *fn)
                download_from_remote(fnidx);
                idx = bam_index_load_local(fn);
        }
+       if (idx == 0) fprintf(stderr, "[bam_index_load] fail to load BAM index.\n");
        return idx;
 }
 
@@ -403,7 +406,10 @@ int bam_index_build2(const char *fn, const char *_fnidx)
        FILE *fpidx;
        bamFile fp;
        bam_index_t *idx;
-       assert(fp = bam_open(fn, "r"));
+       if ((fp = bam_open(fn, "r")) == 0) {
+               fprintf(stderr, "[bam_index_build2] fail to open the BAM file.\n");
+               return -1;
+       }
        idx = bam_index_core(fp);
        bam_close(fp);
        if (_fnidx == 0) {
@@ -414,7 +420,7 @@ int bam_index_build2(const char *fn, const char *_fnidx)
        if (fpidx == 0) {
                fprintf(stderr, "[bam_index_build2] fail to create the index file.\n");
                free(fnidx);
-               return 1;
+               return -1;
        }
        bam_index_save(idx, fpidx);
        bam_index_destroy(idx);
index 1562170e9159e8b53f2d65d7c87e4a59e8acecac..425290e26f11d29cd9b3cc4633287ca5e8ea41e2 100644 (file)
@@ -1,5 +1,6 @@
 #include <stdlib.h>
 #include <stdio.h>
+#include <assert.h>
 #include "bam.h"
 #include "ksort.h"
 
index 9abff46290f183fa4e6f8fdd3323b7345513a3a8..a20f9b354d477c1f11a8d98854ff9ec7682d01c1 100644 (file)
--- a/bam_md.c
+++ b/bam_md.c
@@ -1,4 +1,5 @@
 #include <unistd.h>
+#include <assert.h>
 #include <string.h>
 #include <ctype.h>
 #include "faidx.h"
index 46906dba429032c49920f61a4149c29ec43eae11..3ffd52851e3be64bed628081b2b57d8ca6378b8b 100644 (file)
@@ -1,6 +1,7 @@
 #include <stdio.h>
 #include <stdlib.h>
 #include <ctype.h>
+#include <assert.h>
 #include "sam.h"
 
 typedef struct __linkbuf_t {
@@ -61,7 +62,7 @@ static inline int resolve_cigar(bam_pileup1_t *p, uint32_t pos)
        int ret = 1, is_restart = 1;
 
        if (c->flag&BAM_FUNMAP) return 0; // unmapped read
-       assert(x <= pos);
+       assert(x <= pos); // otherwise a bug
        p->qpos = -1; p->indel = 0; p->is_del = p->is_head = p->is_tail = 0;
        for (k = 0; k < c->n_cigar; ++k) {
                int op = bam1_cigar(b)[k] & BAM_CIGAR_MASK; // operation
@@ -97,7 +98,7 @@ static inline int resolve_cigar(bam_pileup1_t *p, uint32_t pos)
                        break;
                }
        }
-       assert(x > pos);
+       assert(x > pos); // otherwise a bug
        return ret;
 }
 
@@ -196,7 +197,12 @@ int bam_plbuf_push(const bam1_t *b, bam_plbuf_t *buf)
                        buf->func(buf->tid, buf->pos, n_pu, buf->pu, buf->func_data);
                }
                // update tid and pos
-               if (buf->head->next) assert(buf->tid <= buf->head->b.core.tid); // otherwise, not sorted
+               if (buf->head->next) {
+                       if (buf->tid > buf->head->b.core.tid) {
+                               fprintf(stderr, "[bam_plbuf_push] unsorted input. Pileup aborts.\n");
+                               return 1;
+                       }
+               }
                if (buf->tid < buf->head->b.core.tid) { // come to a new reference sequence
                        buf->tid = buf->head->b.core.tid; buf->pos = buf->head->beg; // jump to the next reference
                } else if (buf->pos < buf->head->beg) { // here: tid == head->b.core.tid
index 50269010d7e6760ddfa8f908df3b7349ab4d83cd..1fa6cad5d7582b9b7a1d148c47df4e00aa5b638e 100644 (file)
@@ -1,7 +1,6 @@
 #include <stdlib.h>
 #include <string.h>
 #include <stdio.h>
-#include <assert.h>
 #include <zlib.h>
 #include "bam.h"
 
index 81b7b6356ed2f32c255ddcf6129c38d84ca4ab15..c1c4a432ce8167448f31577d2ce0f650610ff626 100644 (file)
@@ -1,4 +1,5 @@
 #include <unistd.h>
+#include <assert.h>
 #include "bam.h"
 
 typedef struct {
index e5b83a11c17e6ce393e119945cf43e048ad07e8a..be2579ca0ee28f59792e863bd102b2f7937566d0 100644 (file)
@@ -161,6 +161,7 @@ tview_t *tv_init(const char *fn, const char *fn_fa)
        tview_t *tv = (tview_t*)calloc(1, sizeof(tview_t));
        tv->is_dot = 1;
        tv->idx = bam_index_load(fn);
+       if (tv->idx == 0) exit(1);
        tv->fp = bam_open(fn, "r");
        bgzf_set_cache_size(tv->fp, 8 * 1024 *1024);
        assert(tv->fp);
diff --git a/bamtk.c b/bamtk.c
index 61ead1058c818b5ee6b00df27fad2864db3f664f..9876279da9e23dc7d1da8ed3cd6861ae5ec89eb6 100644 (file)
--- a/bamtk.c
+++ b/bamtk.c
@@ -1,9 +1,10 @@
 #include <stdio.h>
 #include <unistd.h>
+#include <assert.h>
 #include "bam.h"
 
 #ifndef PACKAGE_VERSION
-#define PACKAGE_VERSION "0.1.4-20 (r365)"
+#define PACKAGE_VERSION "0.1.4-21 (r368)"
 #endif
 
 int bam_taf2baf(int argc, char *argv[]);
diff --git a/faidx.c b/faidx.c
index 6b4514ffa2f24946bdfed2228c8044cd8268e692..36366c20db595ee9ffeaac4d87628dc40c3690d2 100644 (file)
--- a/faidx.c
+++ b/faidx.c
@@ -1,5 +1,4 @@
 #include <ctype.h>
-#include <assert.h>
 #include <string.h>
 #include <stdlib.h>
 #include <stdio.h>
@@ -85,14 +84,19 @@ faidx_t *fai_build_core(RAZF *rz)
                                name[l_name++] = c;
                        }
                        name[l_name] = '\0';
-                       assert(ret);
+                       if (ret == 0) {
+                               fprintf(stderr, "[fai_build_core] the last entry has no sequence\n");
+                               free(name); fai_destroy(idx);
+                               return 0;
+                       }
                        if (c != '\n') while (razf_read(rz, &c, 1) && c != '\n');
                        state = 1; len = 0;
                        offset = razf_tell(rz);
                } else {
                        if (state == 3) {
-                               fprintf(stderr, "[fai_build_core] inlined empty line is not allowed in sequence '%s'. Abort!\n", name);
-                               exit(1);
+                               fprintf(stderr, "[fai_build_core] inlined empty line is not allowed in sequence '%s'.\n", name);
+                               free(name); fai_destroy(idx);
+                               return 0;
                        }
                        if (state == 2) state = 3;
                        l1 = l2 = 0;
@@ -101,13 +105,15 @@ faidx_t *fai_build_core(RAZF *rz)
                                if (isgraph(c)) ++l2;
                        } while ((ret = razf_read(rz, &c, 1)) && c != '\n');
                        if (state == 3 && l2) {
-                               fprintf(stderr, "[fai_build_core] different line length in sequence '%s'. Abort!\n", name);
-                               exit(1);
+                               fprintf(stderr, "[fai_build_core] different line length in sequence '%s'.\n", name);
+                               free(name); fai_destroy(idx);
+                               return 0;
                        }
                        ++l1; len += l2;
                        if (l2 >= 0x10000) {
-                               fprintf(stderr, "[fai_build_core] line length exceeds 65535 in sequence '%s'. Abort!\n", name);
-                               exit(1);
+                               fprintf(stderr, "[fai_build_core] line length exceeds 65535 in sequence '%s'.\n", name);
+                               free(name); fai_destroy(idx);
+                               return 0;
                        }
                        if (state == 1) line_len = l1, line_blen = l2, state = 0;
                        else if (state == 0) {
@@ -161,7 +167,7 @@ void fai_destroy(faidx_t *fai)
        free(fai);
 }
 
-void fai_build(const char *fn)
+int fai_build(const char *fn)
 {
        char *str;
        RAZF *rz;
@@ -170,15 +176,24 @@ void fai_build(const char *fn)
        str = (char*)calloc(strlen(fn) + 5, 1);
        sprintf(str, "%s.fai", fn);
        rz = razf_open(fn, "r");
-       assert(rz);
+       if (rz == 0) {
+               fprintf(stderr, "[fai_build] fail to open the FASTA file.\n");
+               free(str);
+               return -1;
+       }
        fai = fai_build_core(rz);
        razf_close(rz);
        fp = fopen(str, "w");
-       assert(fp);
+       if (fp == 0) {
+               fprintf(stderr, "[fai_build] fail to write FASTA index.\n");
+               fai_destroy(fai); free(str);
+               return -1;
+       }
        fai_save(fai, fp);
        fclose(fp);
        free(str);
        fai_destroy(fai);
+       return 0;
 }
 
 faidx_t *fai_load(const char *fn)
@@ -194,6 +209,7 @@ faidx_t *fai_load(const char *fn)
                fai_build(fn);
                fp = fopen(str, "r");
                if (fp == 0) {
+                       fprintf(stderr, "[fai_load] fail to open FASTA index.\n");
                        free(str);
                        return 0;
                }
@@ -201,9 +217,11 @@ faidx_t *fai_load(const char *fn)
        fai = fai_read(fp);
        fclose(fp);
        fai->rz = razf_open(fn, "r");
-       if (fai->rz == 0) return 0;
-       assert(fai->rz);
        free(str);
+       if (fai->rz == 0) {
+               fprintf(stderr, "[fai_load] fail to open FASTA file.\n");
+               return 0;
+       }
        return fai;
 }
 
@@ -271,7 +289,7 @@ int faidx_main(int argc, char *argv[])
                        char *s;
                        faidx_t *fai;
                        fai = fai_load(argv[1]);
-                       assert(fai);
+                       if (fai == 0) return 1;
                        for (i = 2; i != argc; ++i) {
                                printf(">%s\n", argv[i]);
                                s = fai_fetch(fai, argv[i], &l);
diff --git a/faidx.h b/faidx.h
index 98c60e4cf1dfbc69e9d6014d96174ca2adaca2e8..1a52fb7eb73afc36529c263061c93c46cc525217 100644 (file)
--- a/faidx.h
+++ b/faidx.h
@@ -46,9 +46,10 @@ extern "C" {
        /*!
          @abstract   Build index for a FASTA or razip compressed FASTA file.
          @param  fn  FASTA file name
+         @return     0 on success; or -1 on failure
          @discussion File "fn.fai" will be generated.
         */
-       void fai_build(const char *fn);
+       int fai_build(const char *fn);
 
        /*!
          @abstract    Distroy a faidx_t struct.
diff --git a/glf.c b/glf.c
index 4a6c6675112668625b8f74b1159c75a09b6d6806..8d5346ae70af825afeeefaa2d50f276749d5301f 100644 (file)
--- a/glf.c
+++ b/glf.c
@@ -38,8 +38,9 @@ glf3_header_t *glf3_header_read(glfFile fp)
        h = glf3_header_init();
        bgzf_read(fp, magic, 4);
        if (strncmp(magic, "GLF\3", 4)) {
-               fprintf(stderr, "[glf3_header_read] invalid magic. Abort.\n");
-               exit(1);
+               fprintf(stderr, "[glf3_header_read] invalid magic.\n");
+               glf3_header_destroy(h);
+               return 0;
        }
        bgzf_read(fp, &h->l_text, 4);
        if (glf3_is_BE) h->l_text = bam_swap_endian_4(h->l_text);