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