]> git.donarmstrong.com Git - samtools.git/blob - bam.c
more meaningful BAM truncation message
[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;
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)
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         }
164 }
165
166 int bam_read1(bamFile fp, bam1_t *b)
167 {
168         bam1_core_t *c = &b->core;
169         int32_t block_len, ret, i;
170         uint32_t x[8];
171
172         assert(BAM_CORE_SIZE == 32);
173         if ((ret = bam_read(fp, &block_len, 4)) != 4) {
174                 if (ret == 0) return -1; // normal end-of-file
175                 else return -2; // truncated
176         }
177         if (bam_read(fp, x, BAM_CORE_SIZE) != BAM_CORE_SIZE) return -3;
178         if (bam_is_be) {
179                 bam_swap_endian_4p(&block_len);
180                 for (i = 0; i < 8; ++i) bam_swap_endian_4p(x + i);
181         }
182         c->tid = x[0]; c->pos = x[1];
183         c->bin = x[2]>>16; c->qual = x[2]>>8&0xff; c->l_qname = x[2]&0xff;
184         c->flag = x[3]>>16; c->n_cigar = x[3]&0xffff;
185         c->l_qseq = x[4];
186         c->mtid = x[5]; c->mpos = x[6]; c->isize = x[7];
187         b->data_len = block_len - BAM_CORE_SIZE;
188         if (b->m_data < b->data_len) {
189                 b->m_data = b->data_len;
190                 kroundup32(b->m_data);
191                 b->data = (uint8_t*)realloc(b->data, b->m_data);
192         }
193         if (bam_read(fp, b->data, b->data_len) != b->data_len) return -4;
194         b->l_aux = b->data_len - c->n_cigar * 4 - c->l_qname - c->l_qseq - (c->l_qseq+1)/2;
195         if (bam_is_be) swap_endian_data(c, b->data_len, b->data);
196         return 4 + block_len;
197 }
198
199 inline int bam_write1_core(bamFile fp, const bam1_core_t *c, int data_len, uint8_t *data)
200 {
201         uint32_t x[8], block_len = data_len + BAM_CORE_SIZE, y;
202         int i;
203         assert(BAM_CORE_SIZE == 32);
204         x[0] = c->tid;
205         x[1] = c->pos;
206         x[2] = (uint32_t)c->bin<<16 | c->qual<<8 | c->l_qname;
207         x[3] = (uint32_t)c->flag<<16 | c->n_cigar;
208         x[4] = c->l_qseq;
209         x[5] = c->mtid;
210         x[6] = c->mpos;
211         x[7] = c->isize;
212         bgzf_flush_try(fp, 4 + block_len);
213         if (bam_is_be) {
214                 for (i = 0; i < 8; ++i) bam_swap_endian_4p(x + i);
215                 y = block_len;
216                 bam_write(fp, bam_swap_endian_4p(&y), 4);
217                 swap_endian_data(c, data_len, data);
218         } else bam_write(fp, &block_len, 4);
219         bam_write(fp, x, BAM_CORE_SIZE);
220         bam_write(fp, data, data_len);
221         if (bam_is_be) swap_endian_data(c, data_len, data);
222         return 4 + block_len;
223 }
224
225 int bam_write1(bamFile fp, const bam1_t *b)
226 {
227         return bam_write1_core(fp, &b->core, b->data_len, b->data);
228 }
229
230 char *bam_format1_core(const bam_header_t *header, const bam1_t *b, int of)
231 {
232         uint8_t *s = bam1_seq(b), *t = bam1_qual(b);
233         int i;
234         const bam1_core_t *c = &b->core;
235         kstring_t str;
236         str.l = str.m = 0; str.s = 0;
237
238         kputsn(bam1_qname(b), c->l_qname-1, &str); kputc('\t', &str);
239         if (of == BAM_OFDEC) { kputw(c->flag, &str); kputc('\t', &str); }
240         else if (of == BAM_OFHEX) ksprintf(&str, "0x%x\t", c->flag);
241         else { // BAM_OFSTR
242                 for (i = 0; i < 16; ++i)
243                         if ((c->flag & 1<<i) && bam_flag2char_table[i])
244                                 kputc(bam_flag2char_table[i], &str);
245                 kputc('\t', &str);
246         }
247         if (c->tid < 0) kputsn("*\t", 2, &str);
248         else {
249                 if (header) kputs(header->target_name[c->tid] , &str);
250                 else kputw(c->tid, &str);
251                 kputc('\t', &str);
252         }
253         kputw(c->pos + 1, &str); kputc('\t', &str); kputw(c->qual, &str); kputc('\t', &str);
254         if (c->n_cigar == 0) kputc('*', &str);
255         else {
256                 for (i = 0; i < c->n_cigar; ++i) {
257                         kputw(bam1_cigar(b)[i]>>BAM_CIGAR_SHIFT, &str);
258                         kputc("MIDNSHP"[bam1_cigar(b)[i]&BAM_CIGAR_MASK], &str);
259                 }
260         }
261         kputc('\t', &str);
262         if (c->mtid < 0) kputsn("*\t", 2, &str);
263         else if (c->mtid == c->tid) kputsn("=\t", 2, &str);
264         else {
265                 if (header) kputs(header->target_name[c->mtid], &str);
266                 else kputw(c->mtid, &str);
267                 kputc('\t', &str);
268         }
269         kputw(c->mpos + 1, &str); kputc('\t', &str); kputw(c->isize, &str); kputc('\t', &str);
270         if (c->l_qseq) {
271                 for (i = 0; i < c->l_qseq; ++i) kputc(bam_nt16_rev_table[bam1_seqi(s, i)], &str);
272                 kputc('\t', &str);
273                 if (t[0] == 0xff) kputc('*', &str);
274                 else for (i = 0; i < c->l_qseq; ++i) kputc(t[i] + 33, &str);
275         } else kputsn("*\t*", 3, &str);
276         s = bam1_aux(b);
277         while (s < b->data + b->data_len) {
278                 uint8_t type, key[2];
279                 key[0] = s[0]; key[1] = s[1];
280                 s += 2; type = *s; ++s;
281                 kputc('\t', &str); kputsn((char*)key, 2, &str); kputc(':', &str);
282                 if (type == 'A') { kputsn("A:", 2, &str); kputc(*s, &str); ++s; }
283                 else if (type == 'C') { kputsn("i:", 2, &str); kputw(*s, &str); ++s; }
284                 else if (type == 'c') { kputsn("i:", 2, &str); kputw(*(int8_t*)s, &str); ++s; }
285                 else if (type == 'S') { kputsn("i:", 2, &str); kputw(*(uint16_t*)s, &str); s += 2; }
286                 else if (type == 's') { kputsn("i:", 2, &str); kputw(*(int16_t*)s, &str); s += 2; }
287                 else if (type == 'I') { kputsn("i:", 2, &str); kputuw(*(uint32_t*)s, &str); s += 4; }
288                 else if (type == 'i') { kputsn("i:", 2, &str); kputw(*(int32_t*)s, &str); s += 4; }
289                 else if (type == 'f') { ksprintf(&str, "f:%g", *(float*)s); s += 4; }
290                 else if (type == 'd') { ksprintf(&str, "d:%lg", *(double*)s); s += 8; }
291                 else if (type == 'Z' || type == 'H') { kputc(type, &str); kputc(':', &str); while (*s) kputc(*s++, &str); ++s; }
292         }
293         return str.s;
294 }
295
296 char *bam_format1(const bam_header_t *header, const bam1_t *b)
297 {
298         return bam_format1_core(header, b, BAM_OFDEC);
299 }
300
301 void bam_view1(const bam_header_t *header, const bam1_t *b)
302 {
303         char *s = bam_format1(header, b);
304         puts(s);
305         free(s);
306 }
307
308 int bam_validate1(const bam_header_t *header, const bam1_t *b)
309 {
310         char *s;
311
312         if (b->core.tid < -1 || b->core.mtid < -1) return 0;
313         if (header && (b->core.tid >= header->n_targets || b->core.mtid >= header->n_targets)) return 0;
314
315         if (b->data_len < b->core.l_qname) return 0;
316         s = memchr(bam1_qname(b), '\0', b->core.l_qname);
317         if (s != &bam1_qname(b)[b->core.l_qname-1]) return 0;
318
319         // FIXME: Other fields could also be checked, especially the auxiliary data
320
321         return 1;
322 }
323
324 // FIXME: we should also check the LB tag associated with each alignment
325 const char *bam_get_library(bam_header_t *h, const bam1_t *b)
326 {
327         const uint8_t *rg;
328         if (h->dict == 0) h->dict = sam_header_parse2(h->text);
329         if (h->rg2lib == 0) h->rg2lib = sam_header2tbl(h->dict, "RG", "ID", "LB");
330         rg = bam_aux_get(b, "RG");
331         return (rg == 0)? 0 : sam_tbl_get(h->rg2lib, (const char*)(rg + 1));
332 }