]> git.donarmstrong.com Git - samtools.git/blobdiff - bam_aux.c
Release samtools-0.1.5a (for compatibility with Bio::DB::Sam)
[samtools.git] / bam_aux.c
index f9be398f0e820d9b48b2002073039af1ebeb95d5..7482500ef08a53f3471858793112d9c72263e33a 100644 (file)
--- a/bam_aux.c
+++ b/bam_aux.c
@@ -1,9 +1,31 @@
 #include <ctype.h>
 #include "bam.h"
 #include "khash.h"
+typedef char *str_p;
 KHASH_MAP_INIT_STR(s, int)
+KHASH_MAP_INIT_STR(r2l, str_p)
 
-uint8_t *bam_aux_get(bam1_t *b, const char tag[2])
+void bam_aux_append(bam1_t *b, const char tag[2], char type, int len, uint8_t *data)
+{
+       int ori_len = b->data_len;
+       b->data_len += 3 + len;
+       b->l_aux += 3 + len;
+       if (b->m_data < b->data_len) {
+               b->m_data = b->data_len;
+               kroundup32(b->m_data);
+               b->data = (uint8_t*)realloc(b->data, b->m_data);
+       }
+       b->data[ori_len] = tag[0]; b->data[ori_len + 1] = tag[1];
+       b->data[ori_len + 2] = type;
+       memcpy(b->data + ori_len + 3, data, len);
+}
+
+uint8_t *bam_aux_get_core(bam1_t *b, const char tag[2])
+{
+       return bam_aux_get(b, tag);
+}
+
+uint8_t *bam_aux_get(const bam1_t *b, const char tag[2])
 {
        uint8_t *s;
        int y = tag[0]<<8 | tag[1];
@@ -17,7 +39,7 @@ uint8_t *bam_aux_get(bam1_t *b, const char tag[2])
                else if (type == 'S') s += 2;
                else if (type == 'I' || type == 'F') s += 4;
                else if (type == 'D') s += 8;
-               else if (type == 'Z' || type == 'H') { while (*s) putchar(*s++); ++s; }
+               else if (type == 'Z' || type == 'H') { while (*s) ++s; ++s; }
        }
        return 0;
 }
@@ -50,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;
@@ -59,7 +81,7 @@ void bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *b
 
        bam_init_header_hash(header);
        h = (khash_t(s)*)header->hash;
-       
+
        l = strlen(str);
        p = s = (char*)malloc(l+1);
        /* squeeze out "," */
@@ -71,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);
@@ -85,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)
@@ -137,3 +163,70 @@ char *bam_aux2Z(const uint8_t *s)
        if (type == 'Z' || type == 'H') return (char*)s;
        else return 0;
 }
+
+/******************
+ * rg2lib related *
+ ******************/
+
+int bam_strmap_put(void *rg2lib, const char *rg, const char *lib)
+{
+       int ret;
+       khint_t k;
+       khash_t(r2l) *h = (khash_t(r2l)*)rg2lib;
+       char *key;
+       if (h == 0) return 1;
+       key = strdup(rg);
+       k = kh_put(r2l, h, key, &ret);
+       if (ret) kh_val(h, k) = strdup(lib);
+       else {
+               fprintf(stderr, "[bam_rg2lib_put] duplicated @RG ID: %s\n", rg);
+               free(key);
+       }
+       return 0;
+}
+
+const char *bam_strmap_get(const void *rg2lib, const char *rg)
+{
+       const khash_t(r2l) *h = (const khash_t(r2l)*)rg2lib;
+       khint_t k;
+       if (h == 0) return 0;
+       k = kh_get(r2l, h, rg);
+       if (k != kh_end(h)) return (const char*)kh_val(h, k);
+       else return 0;
+}
+
+void *bam_strmap_dup(const void *rg2lib)
+{
+       const khash_t(r2l) *h = (const khash_t(r2l)*)rg2lib;
+       khash_t(r2l) *g;
+       khint_t k, l;
+       int ret;
+       if (h == 0) return 0;
+       g = kh_init(r2l);
+       for (k = kh_begin(h); k < kh_end(h); ++k) {
+               if (kh_exist(h, k)) {
+                       char *key = strdup(kh_key(h, k));
+                       l = kh_put(r2l, g, key, &ret);
+                       kh_val(g, l) = strdup(kh_val(h, k));
+               }
+       }
+       return g;
+}
+
+void *bam_strmap_init()
+{
+       return (void*)kh_init(r2l);
+}
+
+void bam_strmap_destroy(void *rg2lib)
+{
+       khash_t(r2l) *h = (khash_t(r2l)*)rg2lib;
+       khint_t k;
+       if (h == 0) return;
+       for (k = kh_begin(h); k < kh_end(h); ++k) {
+               if (kh_exist(h, k)) {
+                       free((char*)kh_key(h, k)); free(kh_val(h, k));
+               }
+       }
+       kh_destroy(r2l, h);
+}