]> git.donarmstrong.com Git - samtools.git/blobdiff - bam_import.c
* samtools-0.1.4-9 (r340)
[samtools.git] / bam_import.c
index 746dc5bbbde0415911794b28d666a4a53ef86be8..c7f866701edd87f37e18fbc9768ed63f68be7f9e 100644 (file)
@@ -5,6 +5,7 @@
 #include <stdlib.h>
 #include <unistd.h>
 #include <assert.h>
+#include "kstring.h"
 #include "bam.h"
 #include "kseq.h"
 #include "khash.h"
@@ -42,6 +43,7 @@ struct __tamFile_t {
        kstream_t *ks;
        kstring_t *str;
        uint64_t n_lines;
+       int is_first;
 };
 
 char **__bam_get_lines(const char *fn, int *_n) // for bam_plcmd.c only
@@ -144,13 +146,114 @@ static inline void append_text(bam_header_t *header, kstring_t *str)
        header->l_text += str->l + 1;
        header->text[header->l_text] = 0;
 }
-int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
+
+int sam_header_parse_rg(bam_header_t *h)
 {
-       int ret, doff, doff0, dret;
-       bam1_core_t *c = &b->core;
-       kstring_t *str = fp->str;
-       kstream_t *ks = fp->ks;
+       kstring_t *rgid, *rglib;
+       char *p, *q, *s, *r;
+       int n = 0;
+
+       // free
+       bam_strmap_destroy(h->rg2lib); h->rg2lib = 0;
+       if (h->l_text < 3) return 0;
+       // parse @RG lines
+       h->rg2lib = bam_strmap_init();
+       rgid = calloc(1, sizeof(kstring_t));
+       rglib = calloc(1, sizeof(kstring_t));
+       s = h->text;
+       while ((s = strstr(s, "@RG")) != 0) {
+               if (rgid->l && rglib->l) {
+                       bam_strmap_put(h->rg2lib, rgid->s, rglib->s);
+                       ++n;
+               }
+               rgid->l = rglib->l = 0;
+               s += 3;
+               r = s;
+               if ((p = strstr(s, "ID:")) != 0) {
+                       q = p + 3;
+                       for (p = q; *p && *p != '\t' && *p != '\r' && *p != '\n'; ++p);
+                       kputsn(q, p - q, rgid);
+               } else {
+                       fprintf(stderr, "[bam_header_parse] missing ID tag in @RG lines.\n");
+                       break;
+               }
+               if (r < p) r = p;
+               if ((p = strstr(s, "LB:")) != 0) {
+                       q = p + 3;
+                       for (p = q; *p && *p != '\t' && *p != '\r' && *p != '\n'; ++p);
+                       kputsn(q, p - q, rglib);
+               } else {
+                       fprintf(stderr, "[bam_header_parse] missing LB tag in @RG lines.\n");
+                       break;
+               }
+               if (r < p) r = p;
+               s = r + 3;
+       }
+       if (rgid->l && rglib->l) {
+               bam_strmap_put(h->rg2lib, rgid->s, rglib->s);
+               ++n;
+       }
+       free(rgid->s); free(rgid);
+       free(rglib->s); free(rglib);
+       if (n == 0) {
+               bam_strmap_destroy(h->rg2lib);
+               h->rg2lib = 0;
+       }
+       return n;
+}
 
+int sam_header_parse(bam_header_t *h)
+{
+       int i;
+       char *s, *p, *q, *r;
+
+       // free
+       free(h->target_len); free(h->target_name);
+       h->n_targets = 0; h->target_len = 0; h->target_name = 0;
+       if (h->l_text < 3) return 0;
+       // count number of @SQ
+       s = h->text;
+       while ((s = strstr(s, "@SQ")) != 0) {
+               ++h->n_targets;
+               s += 3;
+       }
+       if (h->n_targets == 0) return 0;
+       h->target_len = (uint32_t*)calloc(h->n_targets, 4);
+       h->target_name = (char**)calloc(h->n_targets, sizeof(void*));
+       // parse @SQ lines
+       i = 0;
+       s = h->text;
+       while ((s = strstr(s, "@SQ")) != 0) {
+               s += 3;
+               r = s;
+               if ((p = strstr(s, "SN:")) != 0) {
+                       q = p + 3;
+                       for (p = q; *p && *p != '\t' && *p != '\r' && *p != '\n'; ++p);
+                       h->target_name[i] = (char*)calloc(p - q + 1, 1);
+                       strncpy(h->target_name[i], q, p - q);
+               } else goto header_err_ret;
+               if (r < p) r = p;
+               if ((p = strstr(s, "LN:")) != 0) h->target_len[i] = strtol(p + 3, 0, 10);
+               else goto header_err_ret;
+               if (r < p) r = p;
+               s = r + 3;
+               ++i;
+       }
+       sam_header_parse_rg(h);
+       return h->n_targets;
+
+header_err_ret:
+       fprintf(stderr, "[bam_header_parse] missing SN or LN tag in @SQ lines.\n");
+       free(h->target_len); free(h->target_name);
+       h->n_targets = 0; h->target_len = 0; h->target_name = 0;
+       return 0;
+}
+
+bam_header_t *sam_header_read(tamFile fp)
+{
+       int ret, dret;
+       bam_header_t *header = bam_header_init();
+       kstring_t *str = fp->str;
        while ((ret = ks_getuntil(fp->ks, KS_SEP_TAB, str, &dret)) >= 0 && str->s[0] == '@') { // skip header
                str->s[str->l] = dret; // note that str->s is NOT null terminated!!
                append_text(header, str);
@@ -161,7 +264,28 @@ int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
                }
                ++fp->n_lines;
        }
-       while (ret == 0) ret = ks_getuntil(fp->ks, KS_SEP_TAB, str, &dret); // special consideration for "\r\n"
+       sam_header_parse(header);
+       bam_init_header_hash(header);
+       fp->is_first = 1;
+       return header;
+}
+
+int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
+{
+       int ret, doff, doff0, dret, z = 0;
+       bam1_core_t *c = &b->core;
+       kstring_t *str = fp->str;
+       kstream_t *ks = fp->ks;
+
+       if (fp->is_first) {
+               fp->is_first = 0;
+               ret = str->l;
+       } else {
+               do { // special consideration for empty lines
+                       ret = ks_getuntil(fp->ks, KS_SEP_TAB, str, &dret);
+                       if (ret >= 0) z += str->l + 1;
+               } while (ret == 0);
+       }
        if (ret < 0) return -1;
        ++fp->n_lines;
        doff = 0;
@@ -172,10 +296,16 @@ int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
                doff += c->l_qname;
        }
        { // flag, tid, pos, qual
-               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); c->flag = atoi(str->s);
-               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); c->tid = bam_get_tid(header, str->s);
-               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); c->pos = isdigit(str->s[0])? atoi(str->s) - 1 : -1;
-               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); c->qual = isdigit(str->s[0])? atoi(str->s) : 0;
+               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); z += str->l + 1; c->flag = atoi(str->s);
+               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); z += str->l + 1; c->tid = bam_get_tid(header, str->s);
+               if (c->tid < 0 && strcmp(str->s, "*")) {
+                       if (header->n_targets == 0) {
+                               fprintf(stderr, "[sam_read1] missing header? Abort!\n");
+                               exit(1);
+                       } else fprintf(stderr, "[sam_read1] reference '%s' is recognized as '*'.\n", str->s);
+               }
+               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); z += str->l + 1; c->pos = isdigit(str->s[0])? atoi(str->s) - 1 : -1;
+               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); z += str->l + 1; c->qual = isdigit(str->s[0])? atoi(str->s) : 0;
                if (ret < 0) return -2;
        }
        { // cigar
@@ -184,6 +314,7 @@ int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
                long x;
                c->n_cigar = 0;
                if (ks_getuntil(ks, KS_SEP_TAB, str, &dret) < 0) return -3;
+               z += str->l + 1;
                if (str->s[0] != '*') {
                        for (s = str->s; *s; ++s) {
                                if (isalpha(*s)) ++c->n_cigar;
@@ -210,15 +341,19 @@ int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
                } else c->bin = bam_reg2bin(c->pos, c->pos + 1);
        }
        { // mtid, mpos, isize
-               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); c->mtid = strcmp(str->s, "=")? bam_get_tid(header, str->s) : c->tid;
-               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); c->mpos = isdigit(str->s[0])? atoi(str->s) - 1 : -1;
-               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); c->isize = (str->s[0] == '-' || isdigit(str->s[0]))? atoi(str->s) : 0;
+               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); z += str->l + 1;
+               c->mtid = strcmp(str->s, "=")? bam_get_tid(header, str->s) : c->tid;
+               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); z += str->l + 1;
+               c->mpos = isdigit(str->s[0])? atoi(str->s) - 1 : -1;
+               ret = ks_getuntil(ks, KS_SEP_TAB, str, &dret); z += str->l + 1;
+               c->isize = (str->s[0] == '-' || isdigit(str->s[0]))? atoi(str->s) : 0;
                if (ret < 0) return -4;
        }
        { // seq and qual
                int i;
                uint8_t *p;
                if (ks_getuntil(ks, KS_SEP_TAB, str, &dret) < 0) return -5; // seq
+               z += str->l + 1;
                c->l_qseq = strlen(str->s);
                if (c->n_cigar && c->l_qseq != (int32_t)bam_cigar2qlen(c, bam1_cigar(b)))
                        parse_error(fp->n_lines, "CIGAR and sequence length are inconsistent");
@@ -227,6 +362,7 @@ int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
                for (i = 0; i < c->l_qseq; ++i)
                        p[i/2] |= bam_nt16_table[(int)str->s[i]] << 4*(1-i%2);
                if (ks_getuntil(ks, KS_SEP_TAB, str, &dret) < 0) return -6; // qual
+               z += str->l + 1;
                if (strcmp(str->s, "*") && c->l_qseq != strlen(str->s))
                        parse_error(fp->n_lines, "sequence and quality are inconsistent");
                p += (c->l_qseq+1)/2;
@@ -238,6 +374,7 @@ int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
        if (dret != '\n' && dret != '\r') { // aux
                while (ks_getuntil(ks, KS_SEP_TAB, str, &dret) >= 0) {
                        uint8_t *s, type, key[2];
+                       z += str->l + 1;
                        if (str->l < 6 || str->s[2] != ':' || str->s[4] != ':')
                                parse_error(fp->n_lines, "missing colon in auxiliary data");
                        key[0] = str->s[0]; key[1] = str->s[1];
@@ -313,7 +450,7 @@ int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b)
        }
        b->l_aux = doff - doff0;
        b->data_len = doff;
-       return 0;
+       return z;
 }
 
 tamFile sam_open(const char *fn)
@@ -323,7 +460,6 @@ tamFile sam_open(const char *fn)
        fp->str = (kstring_t*)calloc(1, sizeof(kstring_t));
        fp->fp = (strcmp(fn, "-") == 0)? gzdopen(fileno(stdin), "r") : gzopen(fn, "r");
        fp->ks = ks_init(fp->fp);
-       fp->n_lines = 0;
        return fp;
 }
 
@@ -336,41 +472,3 @@ void sam_close(tamFile fp)
                free(fp);
        }
 }
-
-static void taf2baf_core(const char *fntaf, const char *fnbaf, bam_header_t *header)
-{
-       bamFile fpbaf;
-       bam1_t *b;
-       tamFile fp;
-       int ret;
-
-       b = (bam1_t*)calloc(1, sizeof(bam1_t));
-       fpbaf = (strcmp(fnbaf, "-") == 0)? bam_dopen(fileno(stdout), "w") : bam_open(fnbaf, "w");
-       fp = sam_open(fntaf);
-       ret = sam_read1(fp, header, b);
-       bam_header_write(fpbaf, header);
-       if (ret >= 0) {
-               bam_write1(fpbaf, b);
-               while (sam_read1(fp, header, b) >= 0) bam_write1(fpbaf, b);
-       }
-       bam_close(fpbaf);
-       free(b->data); free(b);
-       sam_close(fp);
-}
-
-int bam_taf2baf(int argc, char *argv[])
-{
-       int c;
-       bam_header_t *header;
-
-       while ((c = getopt(argc, argv, "")) >= 0) {
-       }
-       if (optind + 3 > argc) {
-               fprintf(stderr, "Usage: bamtk import <in.ref_list> <in.sam> <out.bam>\n");
-               return 1;
-       }
-       header = sam_header_read2(argv[optind]);
-       taf2baf_core(argv[optind+1], argv[optind+2], header);
-       bam_header_destroy(header);
-       return 0;
-}