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