]> git.donarmstrong.com Git - samtools.git/blob - bam.c
* samtools-0.1.6-12 (r478)
[samtools.git] / bam.c
1 #include <stdio.h>
2 #include <ctype.h>
3 #include <errno.h>
4 #include <assert.h>
5 #include "bam.h"
6 #include "bam_endian.h"
7 #include "kstring.h"
8
9 int bam_is_be = 0;
10 char *bam_flag2char_table = "pPuUrR12sfd\0\0\0\0\0";
11
12 /**************************
13  * CIGAR related routines *
14  **************************/
15
16 uint32_t bam_calend(const bam1_core_t *c, const uint32_t *cigar)
17 {
18         uint32_t k, end;
19         end = c->pos;
20         for (k = 0; k < c->n_cigar; ++k) {
21                 int op = cigar[k] & BAM_CIGAR_MASK;
22                 if (op == BAM_CMATCH || op == BAM_CDEL || op == BAM_CREF_SKIP)
23                         end += cigar[k] >> BAM_CIGAR_SHIFT;
24         }
25         return end;
26 }
27
28 int32_t bam_cigar2qlen(const bam1_core_t *c, const uint32_t *cigar)
29 {
30         uint32_t k;
31         int32_t l = 0;
32         for (k = 0; k < c->n_cigar; ++k) {
33                 int op = cigar[k] & BAM_CIGAR_MASK;
34                 if (op == BAM_CMATCH || op == BAM_CINS || op == BAM_CSOFT_CLIP)
35                         l += cigar[k] >> BAM_CIGAR_SHIFT;
36         }
37         return l;
38 }
39
40 /********************
41  * BAM I/O routines *
42  ********************/
43
44 bam_header_t *bam_header_init()
45 {
46         bam_is_be = bam_is_big_endian();
47         return (bam_header_t*)calloc(1, sizeof(bam_header_t));
48 }
49
50 void bam_header_destroy(bam_header_t *header)
51 {
52         int32_t i;
53         extern void bam_destroy_header_hash(bam_header_t *header);
54         if (header == 0) return;
55         if (header->target_name) {
56                 for (i = 0; i < header->n_targets; ++i)
57                         free(header->target_name[i]);
58                 free(header->target_name);
59                 free(header->target_len);
60         }
61         free(header->text);
62 #ifndef BAM_NO_HASH
63         if (header->rg2lib) bam_strmap_destroy(header->rg2lib);
64         bam_destroy_header_hash(header);
65 #endif
66         free(header);
67 }
68
69 bam_header_t *bam_header_read(bamFile fp)
70 {
71         bam_header_t *header;
72         char buf[4];
73         int32_t i = 1, name_len;
74         // check EOF
75         i = bgzf_check_EOF(fp);
76         if (i < 0) {
77                 // If the file is a pipe, checking the EOF marker will *always* fail
78                 // with ESPIPE.  Suppress the error message in this case.
79                 if (errno != ESPIPE) perror("[bam_header_read] bgzf_check_EOF");
80         }
81         else if (i == 0) fprintf(stderr, "[bam_header_read] EOF marker is absent.\n");
82         // read "BAM1"
83         if (bam_read(fp, buf, 4) != 4) return 0;
84         if (strncmp(buf, "BAM\001", 4)) {
85                 fprintf(stderr, "[bam_header_read] wrong header\n");
86                 return 0;
87         }
88         header = bam_header_init();
89         // read plain text and the number of reference sequences
90         bam_read(fp, &header->l_text, 4);
91         if (bam_is_be) bam_swap_endian_4p(&header->l_text);
92         header->text = (char*)calloc(header->l_text + 1, 1);
93         bam_read(fp, header->text, header->l_text);
94         bam_read(fp, &header->n_targets, 4);
95         if (bam_is_be) bam_swap_endian_4p(&header->n_targets);
96         // read reference sequence names and lengths
97         header->target_name = (char**)calloc(header->n_targets, sizeof(char*));
98         header->target_len = (uint32_t*)calloc(header->n_targets, 4);
99         for (i = 0; i != header->n_targets; ++i) {
100                 bam_read(fp, &name_len, 4);
101                 if (bam_is_be) bam_swap_endian_4p(&name_len);
102                 header->target_name[i] = (char*)calloc(name_len, 1);
103                 bam_read(fp, header->target_name[i], name_len);
104                 bam_read(fp, &header->target_len[i], 4);
105                 if (bam_is_be) bam_swap_endian_4p(&header->target_len[i]);
106         }
107         return header;
108 }
109
110 int bam_header_write(bamFile fp, const bam_header_t *header)
111 {
112         char buf[4];
113         int32_t i, name_len, x;
114         // write "BAM1"
115         strncpy(buf, "BAM\001", 4);
116         bam_write(fp, buf, 4);
117         // write plain text and the number of reference sequences
118         if (bam_is_be) {
119                 x = bam_swap_endian_4(header->l_text);
120                 bam_write(fp, &x, 4);
121                 if (header->l_text) bam_write(fp, header->text, header->l_text);
122                 x = bam_swap_endian_4(header->n_targets);
123                 bam_write(fp, &x, 4);
124         } else {
125                 bam_write(fp, &header->l_text, 4);
126                 if (header->l_text) bam_write(fp, header->text, header->l_text);
127                 bam_write(fp, &header->n_targets, 4);
128         }
129         // write sequence names and lengths
130         for (i = 0; i != header->n_targets; ++i) {
131                 char *p = header->target_name[i];
132                 name_len = strlen(p) + 1;
133                 if (bam_is_be) {
134                         x = bam_swap_endian_4(name_len);
135                         bam_write(fp, &x, 4);
136                 } else bam_write(fp, &name_len, 4);
137                 bam_write(fp, p, name_len);
138                 if (bam_is_be) {
139                         x = bam_swap_endian_4(header->target_len[i]);
140                         bam_write(fp, &x, 4);
141                 } else bam_write(fp, &header->target_len[i], 4);
142         }
143         return 0;
144 }
145
146 static void swap_endian_data(const bam1_core_t *c, int data_len, uint8_t *data)
147 {
148         uint8_t *s;
149         uint32_t i, *cigar = (uint32_t*)(data + c->l_qname);
150         s = data + c->n_cigar*4 + c->l_qname + c->l_qseq + (c->l_qseq + 1)/2;
151         for (i = 0; i < c->n_cigar; ++i) bam_swap_endian_4p(&cigar[i]);
152         while (s < data + data_len) {
153                 uint8_t type;
154                 s += 2; // skip key
155                 type = toupper(*s); ++s; // skip type
156                 if (type == 'C' || type == 'A') ++s;
157                 else if (type == 'S') { bam_swap_endian_2p(s); s += 2; }
158                 else if (type == 'I' || type == 'F') { bam_swap_endian_4p(s); s += 4; }
159                 else if (type == 'D') { bam_swap_endian_8p(s); s += 8; }
160                 else if (type == 'Z' || type == 'H') { while (*s) ++s; ++s; }
161         }
162 }
163
164 int bam_read1(bamFile fp, bam1_t *b)
165 {
166         bam1_core_t *c = &b->core;
167         int32_t block_len, ret, i;
168         uint32_t x[8];
169
170         assert(BAM_CORE_SIZE == 32);
171         if ((ret = bam_read(fp, &block_len, 4)) != 4) {
172                 if (ret == 0) return -1; // normal end-of-file
173                 else return -2; // truncated
174         }
175         if (bam_read(fp, x, BAM_CORE_SIZE) != BAM_CORE_SIZE) return -3;
176         if (bam_is_be) {
177                 bam_swap_endian_4p(&block_len);
178                 for (i = 0; i < 8; ++i) bam_swap_endian_4p(x + i);
179         }
180         c->tid = x[0]; c->pos = x[1];
181         c->bin = x[2]>>16; c->qual = x[2]>>8&0xff; c->l_qname = x[2]&0xff;
182         c->flag = x[3]>>16; c->n_cigar = x[3]&0xffff;
183         c->l_qseq = x[4];
184         c->mtid = x[5]; c->mpos = x[6]; c->isize = x[7];
185         b->data_len = block_len - BAM_CORE_SIZE;
186         if (b->m_data < b->data_len) {
187                 b->m_data = b->data_len;
188                 kroundup32(b->m_data);
189                 b->data = (uint8_t*)realloc(b->data, b->m_data);
190         }
191         if (bam_read(fp, b->data, b->data_len) != b->data_len) return -4;
192         b->l_aux = b->data_len - c->n_cigar * 4 - c->l_qname - c->l_qseq - (c->l_qseq+1)/2;
193         if (bam_is_be) swap_endian_data(c, b->data_len, b->data);
194         return 4 + block_len;
195 }
196
197 inline int bam_write1_core(bamFile fp, const bam1_core_t *c, int data_len, uint8_t *data)
198 {
199         uint32_t x[8], block_len = data_len + BAM_CORE_SIZE, y;
200         int i;
201         assert(BAM_CORE_SIZE == 32);
202         x[0] = c->tid;
203         x[1] = c->pos;
204         x[2] = (uint32_t)c->bin<<16 | c->qual<<8 | c->l_qname;
205         x[3] = (uint32_t)c->flag<<16 | c->n_cigar;
206         x[4] = c->l_qseq;
207         x[5] = c->mtid;
208         x[6] = c->mpos;
209         x[7] = c->isize;
210         if (bam_is_be) {
211                 for (i = 0; i < 8; ++i) bam_swap_endian_4p(x + i);
212                 y = block_len;
213                 bam_write(fp, bam_swap_endian_4p(&y), 4);
214                 swap_endian_data(c, data_len, data);
215         } else bam_write(fp, &block_len, 4);
216         bam_write(fp, x, BAM_CORE_SIZE);
217         bam_write(fp, data, data_len);
218         if (bam_is_be) swap_endian_data(c, data_len, data);
219         return 4 + block_len;
220 }
221
222 int bam_write1(bamFile fp, const bam1_t *b)
223 {
224         return bam_write1_core(fp, &b->core, b->data_len, b->data);
225 }
226
227 char *bam_format1_core(const bam_header_t *header, const bam1_t *b, int of)
228 {
229         uint8_t *s = bam1_seq(b), *t = bam1_qual(b);
230         int i;
231         const bam1_core_t *c = &b->core;
232         kstring_t str;
233         str.l = str.m = 0; str.s = 0;
234
235         ksprintf(&str, "%s\t", bam1_qname(b));
236         if (of == BAM_OFDEC) ksprintf(&str, "%d\t", c->flag);
237         else if (of == BAM_OFHEX) ksprintf(&str, "0x%x\t", c->flag);
238         else { // BAM_OFSTR
239                 for (i = 0; i < 16; ++i)
240                         if ((c->flag & 1<<i) && bam_flag2char_table[i])
241                                 kputc(bam_flag2char_table[i], &str);
242                 kputc('\t', &str);
243         }
244         if (c->tid < 0) kputs("*\t", &str);
245         else ksprintf(&str, "%s\t", header->target_name[c->tid]);
246         ksprintf(&str, "%d\t%d\t", c->pos + 1, c->qual);
247         if (c->n_cigar == 0) kputc('*', &str);
248         else {
249                 for (i = 0; i < c->n_cigar; ++i)
250                         ksprintf(&str, "%d%c", bam1_cigar(b)[i]>>BAM_CIGAR_SHIFT, "MIDNSHP"[bam1_cigar(b)[i]&BAM_CIGAR_MASK]);
251         }
252         kputc('\t', &str);
253         if (c->mtid < 0) kputs("*\t", &str);
254         else if (c->mtid == c->tid) kputs("=\t", &str);
255         else ksprintf(&str, "%s\t", header->target_name[c->mtid]);
256         ksprintf(&str, "%d\t%d\t", c->mpos + 1, c->isize);
257         if (c->l_qseq) {
258                 for (i = 0; i < c->l_qseq; ++i) kputc(bam_nt16_rev_table[bam1_seqi(s, i)], &str);
259                 kputc('\t', &str);
260                 if (t[0] == 0xff) kputc('*', &str);
261                 else for (i = 0; i < c->l_qseq; ++i) kputc(t[i] + 33, &str);
262         } else ksprintf(&str, "*\t*");
263         s = bam1_aux(b);
264         while (s < b->data + b->data_len) {
265                 uint8_t type, key[2];
266                 key[0] = s[0]; key[1] = s[1];
267                 s += 2; type = *s; ++s;
268                 ksprintf(&str, "\t%c%c:", key[0], key[1]);
269                 if (type == 'A') { ksprintf(&str, "A:%c", *s); ++s; }
270                 else if (type == 'C') { ksprintf(&str, "i:%u", *s); ++s; }
271                 else if (type == 'c') { ksprintf(&str, "i:%d", *s); ++s; }
272                 else if (type == 'S') { ksprintf(&str, "i:%u", *(uint16_t*)s); s += 2; }
273                 else if (type == 's') { ksprintf(&str, "i:%d", *(int16_t*)s); s += 2; }
274                 else if (type == 'I') { ksprintf(&str, "i:%u", *(uint32_t*)s); s += 4; }
275                 else if (type == 'i') { ksprintf(&str, "i:%d", *(int32_t*)s); s += 4; }
276                 else if (type == 'f') { ksprintf(&str, "f:%g", *(float*)s); s += 4; }
277                 else if (type == 'd') { ksprintf(&str, "d:%lg", *(double*)s); s += 8; }
278                 else if (type == 'Z' || type == 'H') { ksprintf(&str, "%c:", type); while (*s) kputc(*s++, &str); ++s; }
279         }
280         return str.s;
281 }
282
283 char *bam_format1(const bam_header_t *header, const bam1_t *b)
284 {
285         return bam_format1_core(header, b, BAM_OFDEC);
286 }
287
288 void bam_view1(const bam_header_t *header, const bam1_t *b)
289 {
290         char *s = bam_format1(header, b);
291         printf("%s\n", s);
292         free(s);
293 }