]> git.donarmstrong.com Git - samtools.git/blob - bam.c
* Updated samtools to the latest
[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 #include "sam_header.h"
9
10 int bam_is_be = 0, bam_verbose = 2;
11 char *bam_flag2char_table = "pPuUrR12sfd\0\0\0\0\0";
12
13 /**************************
14  * CIGAR related routines *
15  **************************/
16
17 uint32_t bam_calend(const bam1_core_t *c, const uint32_t *cigar)
18 {
19         uint32_t k, end;
20         end = c->pos;
21         for (k = 0; k < c->n_cigar; ++k) {
22                 int op = cigar[k] & BAM_CIGAR_MASK;
23                 if (op == BAM_CMATCH || op == BAM_CDEL || op == BAM_CREF_SKIP)
24                         end += cigar[k] >> BAM_CIGAR_SHIFT;
25         }
26         return end;
27 }
28
29 int32_t bam_cigar2qlen(const bam1_core_t *c, const uint32_t *cigar)
30 {
31         uint32_t k;
32         int32_t l = 0;
33         for (k = 0; k < c->n_cigar; ++k) {
34                 int op = cigar[k] & BAM_CIGAR_MASK;
35                 if (op == BAM_CMATCH || op == BAM_CINS || op == BAM_CSOFT_CLIP || op == BAM_CEQUAL || op == BAM_CDIFF)
36                         l += cigar[k] >> BAM_CIGAR_SHIFT;
37         }
38         return l;
39 }
40
41 /********************
42  * BAM I/O routines *
43  ********************/
44
45 bam_header_t *bam_header_init()
46 {
47         bam_is_be = bam_is_big_endian();
48         return (bam_header_t*)calloc(1, sizeof(bam_header_t));
49 }
50
51 void bam_header_destroy(bam_header_t *header)
52 {
53         int32_t i;
54         extern void bam_destroy_header_hash(bam_header_t *header);
55         if (header == 0) return;
56         if (header->target_name) {
57                 for (i = 0; i < header->n_targets; ++i)
58                         free(header->target_name[i]);
59                 free(header->target_name);
60                 free(header->target_len);
61         }
62         free(header->text);
63         if (header->dict) sam_header_free(header->dict);
64         if (header->rg2lib) sam_tbl_destroy(header->rg2lib);
65         bam_destroy_header_hash(header);
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         int magic_len;
74         int32_t i = 1, name_len;
75         // check EOF
76         i = bgzf_check_EOF(fp);
77         if (i < 0) {
78                 // If the file is a pipe, checking the EOF marker will *always* fail
79                 // with ESPIPE.  Suppress the error message in this case.
80                 if (errno != ESPIPE) perror("[bam_header_read] bgzf_check_EOF");
81         }
82         else if (i == 0) fprintf(stderr, "[bam_header_read] EOF marker is absent. The input is probably truncated.\n");
83         // read "BAM1"
84         magic_len = bam_read(fp, buf, 4);
85         if (magic_len != 4 || strncmp(buf, "BAM\001", 4) != 0) {
86                 fprintf(stderr, "[bam_header_read] invalid BAM binary header (this is not a BAM file).\n");
87                 return 0;
88         }
89         header = bam_header_init();
90         // read plain text and the number of reference sequences
91         bam_read(fp, &header->l_text, 4);
92         if (bam_is_be) bam_swap_endian_4p(&header->l_text);
93         header->text = (char*)calloc(header->l_text + 1, 1);
94         bam_read(fp, header->text, header->l_text);
95         bam_read(fp, &header->n_targets, 4);
96         if (bam_is_be) bam_swap_endian_4p(&header->n_targets);
97         // read reference sequence names and lengths
98         header->target_name = (char**)calloc(header->n_targets, sizeof(char*));
99         header->target_len = (uint32_t*)calloc(header->n_targets, 4);
100         for (i = 0; i != header->n_targets; ++i) {
101                 bam_read(fp, &name_len, 4);
102                 if (bam_is_be) bam_swap_endian_4p(&name_len);
103                 header->target_name[i] = (char*)calloc(name_len, 1);
104                 bam_read(fp, header->target_name[i], name_len);
105                 bam_read(fp, &header->target_len[i], 4);
106                 if (bam_is_be) bam_swap_endian_4p(&header->target_len[i]);
107         }
108         return header;
109 }
110
111 int bam_header_write(bamFile fp, const bam_header_t *header)
112 {
113         char buf[4];
114         int32_t i, name_len, x;
115         // write "BAM1"
116         strncpy(buf, "BAM\001", 4);
117         bam_write(fp, buf, 4);
118         // write plain text and the number of reference sequences
119         if (bam_is_be) {
120                 x = bam_swap_endian_4(header->l_text);
121                 bam_write(fp, &x, 4);
122                 if (header->l_text) bam_write(fp, header->text, header->l_text);
123                 x = bam_swap_endian_4(header->n_targets);
124                 bam_write(fp, &x, 4);
125         } else {
126                 bam_write(fp, &header->l_text, 4);
127                 if (header->l_text) bam_write(fp, header->text, header->l_text);
128                 bam_write(fp, &header->n_targets, 4);
129         }
130         // write sequence names and lengths
131         for (i = 0; i != header->n_targets; ++i) {
132                 char *p = header->target_name[i];
133                 name_len = strlen(p) + 1;
134                 if (bam_is_be) {
135                         x = bam_swap_endian_4(name_len);
136                         bam_write(fp, &x, 4);
137                 } else bam_write(fp, &name_len, 4);
138                 bam_write(fp, p, name_len);
139                 if (bam_is_be) {
140                         x = bam_swap_endian_4(header->target_len[i]);
141                         bam_write(fp, &x, 4);
142                 } else bam_write(fp, &header->target_len[i], 4);
143         }
144         bgzf_flush(fp);
145         return 0;
146 }
147
148 static void swap_endian_data(const bam1_core_t *c, int data_len, uint8_t *data)
149 {
150         uint8_t *s;
151         uint32_t i, *cigar = (uint32_t*)(data + c->l_qname);
152         s = data + c->n_cigar*4 + c->l_qname + c->l_qseq + (c->l_qseq + 1)/2;
153         for (i = 0; i < c->n_cigar; ++i) bam_swap_endian_4p(&cigar[i]);
154         while (s < data + data_len) {
155                 uint8_t type;
156                 s += 2; // skip key
157                 type = toupper(*s); ++s; // skip type
158                 if (type == 'C' || type == 'A') ++s;
159                 else if (type == 'S') { bam_swap_endian_2p(s); s += 2; }
160                 else if (type == 'I' || type == 'F') { bam_swap_endian_4p(s); s += 4; }
161                 else if (type == 'D') { bam_swap_endian_8p(s); s += 8; }
162                 else if (type == 'Z' || type == 'H') { while (*s) ++s; ++s; }
163                 else if (type == 'B') {
164                         int32_t n, Bsize = bam_aux_type2size(*s);
165                         memcpy(&n, s + 1, 4);
166                         if (1 == Bsize) {
167                         } else if (2 == Bsize) {
168                                 for (i = 0; i < n; i += 2)
169                                         bam_swap_endian_2p(s + 5 + i);
170                         } else if (4 == Bsize) {
171                                 for (i = 0; i < n; i += 4)
172                                         bam_swap_endian_4p(s + 5 + i);
173                         }
174                         bam_swap_endian_4p(s+1); 
175                 }
176         }
177 }
178
179 int bam_read1(bamFile fp, bam1_t *b)
180 {
181         bam1_core_t *c = &b->core;
182         int32_t block_len, ret, i;
183         uint32_t x[8];
184
185         assert(BAM_CORE_SIZE == 32);
186         if ((ret = bam_read(fp, &block_len, 4)) != 4) {
187                 if (ret == 0) return -1; // normal end-of-file
188                 else return -2; // truncated
189         }
190         if (bam_read(fp, x, BAM_CORE_SIZE) != BAM_CORE_SIZE) return -3;
191         if (bam_is_be) {
192                 bam_swap_endian_4p(&block_len);
193                 for (i = 0; i < 8; ++i) bam_swap_endian_4p(x + i);
194         }
195         c->tid = x[0]; c->pos = x[1];
196         c->bin = x[2]>>16; c->qual = x[2]>>8&0xff; c->l_qname = x[2]&0xff;
197         c->flag = x[3]>>16; c->n_cigar = x[3]&0xffff;
198         c->l_qseq = x[4];
199         c->mtid = x[5]; c->mpos = x[6]; c->isize = x[7];
200         b->data_len = block_len - BAM_CORE_SIZE;
201         if (b->m_data < b->data_len) {
202                 b->m_data = b->data_len;
203                 kroundup32(b->m_data);
204                 b->data = (uint8_t*)realloc(b->data, b->m_data);
205         }
206         if (bam_read(fp, b->data, b->data_len) != b->data_len) return -4;
207         b->l_aux = b->data_len - c->n_cigar * 4 - c->l_qname - c->l_qseq - (c->l_qseq+1)/2;
208         if (bam_is_be) swap_endian_data(c, b->data_len, b->data);
209         return 4 + block_len;
210 }
211
212 inline int bam_write1_core(bamFile fp, const bam1_core_t *c, int data_len, uint8_t *data)
213 {
214         uint32_t x[8], block_len = data_len + BAM_CORE_SIZE, y;
215         int i;
216         assert(BAM_CORE_SIZE == 32);
217         x[0] = c->tid;
218         x[1] = c->pos;
219         x[2] = (uint32_t)c->bin<<16 | c->qual<<8 | c->l_qname;
220         x[3] = (uint32_t)c->flag<<16 | c->n_cigar;
221         x[4] = c->l_qseq;
222         x[5] = c->mtid;
223         x[6] = c->mpos;
224         x[7] = c->isize;
225         bgzf_flush_try(fp, 4 + block_len);
226         if (bam_is_be) {
227                 for (i = 0; i < 8; ++i) bam_swap_endian_4p(x + i);
228                 y = block_len;
229                 bam_write(fp, bam_swap_endian_4p(&y), 4);
230                 swap_endian_data(c, data_len, data);
231         } else bam_write(fp, &block_len, 4);
232         bam_write(fp, x, BAM_CORE_SIZE);
233         bam_write(fp, data, data_len);
234         if (bam_is_be) swap_endian_data(c, data_len, data);
235         return 4 + block_len;
236 }
237
238 int bam_write1(bamFile fp, const bam1_t *b)
239 {
240         return bam_write1_core(fp, &b->core, b->data_len, b->data);
241 }
242
243 char *bam_format1_core(const bam_header_t *header, const bam1_t *b, int of)
244 {
245         uint8_t *s = bam1_seq(b), *t = bam1_qual(b);
246         int i;
247         const bam1_core_t *c = &b->core;
248         kstring_t str;
249         str.l = str.m = 0; str.s = 0;
250
251         kputsn(bam1_qname(b), c->l_qname-1, &str); kputc('\t', &str);
252         if (of == BAM_OFDEC) { kputw(c->flag, &str); kputc('\t', &str); }
253         else if (of == BAM_OFHEX) ksprintf(&str, "0x%x\t", c->flag);
254         else { // BAM_OFSTR
255                 for (i = 0; i < 16; ++i)
256                         if ((c->flag & 1<<i) && bam_flag2char_table[i])
257                                 kputc(bam_flag2char_table[i], &str);
258                 kputc('\t', &str);
259         }
260         if (c->tid < 0) kputsn("*\t", 2, &str);
261         else {
262                 if (header) kputs(header->target_name[c->tid] , &str);
263                 else kputw(c->tid, &str);
264                 kputc('\t', &str);
265         }
266         kputw(c->pos + 1, &str); kputc('\t', &str); kputw(c->qual, &str); kputc('\t', &str);
267         if (c->n_cigar == 0) kputc('*', &str);
268         else {
269                 for (i = 0; i < c->n_cigar; ++i) {
270                         kputw(bam1_cigar(b)[i]>>BAM_CIGAR_SHIFT, &str);
271                         kputc("MIDNSHP=X"[bam1_cigar(b)[i]&BAM_CIGAR_MASK], &str);
272                 }
273         }
274         kputc('\t', &str);
275         if (c->mtid < 0) kputsn("*\t", 2, &str);
276         else if (c->mtid == c->tid) kputsn("=\t", 2, &str);
277         else {
278                 if (header) kputs(header->target_name[c->mtid], &str);
279                 else kputw(c->mtid, &str);
280                 kputc('\t', &str);
281         }
282         kputw(c->mpos + 1, &str); kputc('\t', &str); kputw(c->isize, &str); kputc('\t', &str);
283         if (c->l_qseq) {
284                 for (i = 0; i < c->l_qseq; ++i) kputc(bam_nt16_rev_table[bam1_seqi(s, i)], &str);
285                 kputc('\t', &str);
286                 if (t[0] == 0xff) kputc('*', &str);
287                 else for (i = 0; i < c->l_qseq; ++i) kputc(t[i] + 33, &str);
288         } else kputsn("*\t*", 3, &str);
289         s = bam1_aux(b);
290         while (s < b->data + b->data_len) {
291                 uint8_t type, key[2];
292                 key[0] = s[0]; key[1] = s[1];
293                 s += 2; type = *s; ++s;
294                 kputc('\t', &str); kputsn((char*)key, 2, &str); kputc(':', &str);
295                 if (type == 'A') { kputsn("A:", 2, &str); kputc(*s, &str); ++s; }
296                 else if (type == 'C') { kputsn("i:", 2, &str); kputw(*s, &str); ++s; }
297                 else if (type == 'c') { kputsn("i:", 2, &str); kputw(*(int8_t*)s, &str); ++s; }
298                 else if (type == 'S') { kputsn("i:", 2, &str); kputw(*(uint16_t*)s, &str); s += 2; }
299                 else if (type == 's') { kputsn("i:", 2, &str); kputw(*(int16_t*)s, &str); s += 2; }
300                 else if (type == 'I') { kputsn("i:", 2, &str); kputuw(*(uint32_t*)s, &str); s += 4; }
301                 else if (type == 'i') { kputsn("i:", 2, &str); kputw(*(int32_t*)s, &str); s += 4; }
302                 else if (type == 'f') { ksprintf(&str, "f:%g", *(float*)s); s += 4; }
303                 else if (type == 'd') { ksprintf(&str, "d:%lg", *(double*)s); s += 8; }
304                 else if (type == 'Z' || type == 'H') { kputc(type, &str); kputc(':', &str); while (*s) kputc(*s++, &str); ++s; }
305                 else if (type == 'B') {
306                         uint8_t sub_type = *(s++);
307                         int32_t n;
308                         memcpy(&n, s, 4);
309                         s += 4; // no point to the start of the array
310                         kputc(type, &str); kputc(':', &str); kputc(sub_type, &str); // write the typing
311                         for (i = 0; i < n; ++i) {
312                                 kputc(',', &str);
313                                 if ('c' == sub_type || 'c' == sub_type) { kputw(*(int8_t*)s, &str); ++s; }
314                                 else if ('C' == sub_type) { kputw(*(uint8_t*)s, &str); ++s; }
315                                 else if ('s' == sub_type) { kputw(*(int16_t*)s, &str); s += 2; }
316                                 else if ('S' == sub_type) { kputw(*(uint16_t*)s, &str); s += 2; }
317                                 else if ('i' == sub_type) { kputw(*(int32_t*)s, &str); s += 4; }
318                                 else if ('I' == sub_type) { kputuw(*(uint32_t*)s, &str); s += 4; }
319                                 else if ('f' == sub_type) { ksprintf(&str, "%g", *(float*)s); s += 4; }
320                         }
321                 }
322         }
323         return str.s;
324 }
325
326 char *bam_format1(const bam_header_t *header, const bam1_t *b)
327 {
328         return bam_format1_core(header, b, BAM_OFDEC);
329 }
330
331 void bam_view1(const bam_header_t *header, const bam1_t *b)
332 {
333         char *s = bam_format1(header, b);
334         puts(s);
335         free(s);
336 }
337
338 int bam_validate1(const bam_header_t *header, const bam1_t *b)
339 {
340         char *s;
341
342         if (b->core.tid < -1 || b->core.mtid < -1) return 0;
343         if (header && (b->core.tid >= header->n_targets || b->core.mtid >= header->n_targets)) return 0;
344
345         if (b->data_len < b->core.l_qname) return 0;
346         s = memchr(bam1_qname(b), '\0', b->core.l_qname);
347         if (s != &bam1_qname(b)[b->core.l_qname-1]) return 0;
348
349         // FIXME: Other fields could also be checked, especially the auxiliary data
350
351         return 1;
352 }
353
354 // FIXME: we should also check the LB tag associated with each alignment
355 const char *bam_get_library(bam_header_t *h, const bam1_t *b)
356 {
357         const uint8_t *rg;
358         if (h->dict == 0) h->dict = sam_header_parse2(h->text);
359         if (h->rg2lib == 0) h->rg2lib = sam_header2tbl(h->dict, "RG", "ID", "LB");
360         rg = bam_aux_get(b, "RG");
361         return (rg == 0)? 0 : sam_tbl_get(h->rg2lib, (const char*)(rg + 1));
362 }