]> git.donarmstrong.com Git - samtools.git/blob - bcftools/bcf.c
use float for QUAL
[samtools.git] / bcftools / bcf.c
1 #include <string.h>
2 #include <ctype.h>
3 #include <stdio.h>
4 #include "kstring.h"
5 #include "bcf.h"
6
7 void bcf_hdr_clear(bcf_hdr_t *b)
8 {
9         free(b->name); free(b->sname); free(b->txt); free(b->ns); free(b->sns);
10         memset(b, 0, sizeof(bcf_hdr_t));
11 }
12
13 bcf_t *bcf_open(const char *fn, const char *mode)
14 {
15         bcf_t *b;
16         b = calloc(1, sizeof(bcf_t));
17         if (strchr(mode, 'w')) {
18                 b->fp = strcmp(fn, "-")? bgzf_open(fn, mode) : bgzf_fdopen(fileno(stdout), mode);
19         } else {
20                 b->fp = strcmp(fn, "-")? bgzf_open(fn, mode) : bgzf_fdopen(fileno(stdin), mode);
21         }
22         b->fp->owned_file = 1;
23         return b;
24 }
25
26 int bcf_close(bcf_t *b)
27 {
28         int ret;
29         if (b == 0) return 0;
30         ret = bgzf_close(b->fp);
31         free(b);
32         return ret;
33 }
34
35 int bcf_hdr_write(bcf_t *b, const bcf_hdr_t *h)
36 {
37         if (b == 0 || h == 0) return -1;
38         bgzf_write(b->fp, "BCF\4", 4);
39         bgzf_write(b->fp, &h->l_nm, 4);
40         bgzf_write(b->fp, h->name, h->l_nm);
41         bgzf_write(b->fp, &h->l_smpl, 4);
42         bgzf_write(b->fp, h->sname, h->l_smpl);
43         bgzf_write(b->fp, &h->l_txt, 4);
44         bgzf_write(b->fp, h->txt, h->l_txt);
45         bgzf_flush(b->fp);
46         return 16 + h->l_nm + h->l_smpl + h->l_txt;
47 }
48
49 bcf_hdr_t *bcf_hdr_read(bcf_t *b)
50 {
51         uint8_t magic[4];
52         bcf_hdr_t *h;
53         if (b == 0) return 0;
54         h = calloc(1, sizeof(bcf_hdr_t));
55         bgzf_read(b->fp, magic, 4);
56         bgzf_read(b->fp, &h->l_nm, 4);
57         h->name = malloc(h->l_nm);
58         bgzf_read(b->fp, h->name, h->l_nm);
59         bgzf_read(b->fp, &h->l_smpl, 4);
60         h->sname = malloc(h->l_smpl);
61         bgzf_read(b->fp, h->sname, h->l_smpl);
62         bgzf_read(b->fp, &h->l_txt, 4);
63         h->txt = malloc(h->l_txt);
64         bgzf_read(b->fp, h->txt, h->l_txt);
65         bcf_hdr_sync(h);
66         return h;
67 }
68
69 void bcf_hdr_destroy(bcf_hdr_t *h)
70 {
71         if (h == 0) return;
72         free(h->name); free(h->sname); free(h->txt); free(h->ns); free(h->sns);
73         free(h);
74 }
75
76 static inline char **cnt_null(int l, char *str, int *_n)
77 {
78         int n = 0;
79         char *p, **list;
80         *_n = 0;
81         if (l == 0 || str == 0) return 0;
82         for (p = str; p != str + l; ++p)
83                 if (*p == 0) ++n;
84         *_n = n;
85         list = calloc(n, sizeof(void*));
86         list[0] = str;
87         for (p = str, n = 1; p < str + l - 1; ++p)
88                 if (*p == 0) list[n++] = p + 1;
89         return list;
90 }
91
92 int bcf_hdr_sync(bcf_hdr_t *b)
93 {
94         if (b == 0) return -1;
95         if (b->ns) free(b->ns);
96         if (b->sns) free(b->sns);
97         if (b->l_nm) b->ns = cnt_null(b->l_nm, b->name, &b->n_ref);
98         else b->ns = 0, b->n_ref = 0;
99         b->sns = cnt_null(b->l_smpl, b->sname, &b->n_smpl);
100         return 0;
101 }
102
103 int bcf_sync(int n_smpl, bcf1_t *b)
104 {
105         char *p, *tmp[5];
106         int i, n;
107         ks_tokaux_t aux;
108         // set ref, alt, flt, info, fmt
109         b->ref = b->alt = b->flt = b->info = b->fmt = 0;
110         for (p = b->str, n = 0; p < b->str + b->l_str; ++p)
111                 if (*p == 0 && p+1 != b->str + b->l_str) tmp[n++] = p + 1;
112         if (n != 5) return -1;
113         b->ref = tmp[0]; b->alt = tmp[1]; b->flt = tmp[2]; b->info = tmp[3]; b->fmt = tmp[4];
114         // set n_alleles
115         for (p = b->alt, n = 1; *p; ++p)
116                 if (*p == ',') ++n;
117         b->n_alleles = n + 1;
118         // set n_gi and gi[i].fmt
119         for (p = b->fmt, n = 1; *p; ++p)
120                 if (*p == ':') ++n;
121         if (n > b->m_gi) {
122                 int old_m = b->m_gi;
123                 b->m_gi = n;
124                 kroundup32(b->m_gi);
125                 b->gi = realloc(b->gi, b->m_gi * sizeof(bcf_ginfo_t));
126                 memset(b->gi + old_m, 0, (b->m_gi - old_m) * sizeof(bcf_ginfo_t));
127         }
128         b->n_gi = n;
129         for (p = kstrtok(b->fmt, ":", &aux), n = 0; p; p = kstrtok(0, 0, &aux))
130                 b->gi[n++].fmt = bcf_str2int(p, aux.p - p);
131         // set gi[i].len
132         for (i = 0; i < b->n_gi; ++i) {
133                 if (b->gi[i].fmt == bcf_str2int("PL", 2)) {
134                         b->gi[i].len = b->n_alleles * (b->n_alleles + 1) / 2;
135                 } else if (b->gi[i].fmt == bcf_str2int("DP", 2) || b->gi[i].fmt == bcf_str2int("HQ", 2)) {
136                         b->gi[i].len = 2;
137                 } else if (b->gi[i].fmt == bcf_str2int("GQ", 2) || b->gi[i].fmt == bcf_str2int("GT", 2)) {
138                         b->gi[i].len = 1;
139                 } else if (b->gi[i].fmt == bcf_str2int("GL", 2)) {
140                         b->gi[i].len = 4;
141                 }
142                 b->gi[i].data = realloc(b->gi[i].data, n_smpl * b->gi[i].len);
143         }
144         return 0;
145 }
146
147 int bcf_write(bcf_t *bp, const bcf_hdr_t *h, const bcf1_t *b)
148 {
149         int i, l = 0;
150         if (b == 0) return -1;
151         bgzf_write(bp->fp, &b->tid, 4);
152         bgzf_write(bp->fp, &b->pos, 4);
153         bgzf_write(bp->fp, &b->qual, 4);
154         bgzf_write(bp->fp, &b->l_str, 4);
155         bgzf_write(bp->fp, b->str, b->l_str);
156         l = 12 + b->l_str;
157         for (i = 0; i < b->n_gi; ++i) {
158                 bgzf_write(bp->fp, b->gi[i].data, b->gi[i].len * h->n_smpl);
159                 l += b->gi[i].len * h->n_smpl;
160         }
161         return l;
162 }
163
164 int bcf_read(bcf_t *bp, const bcf_hdr_t *h, bcf1_t *b)
165 {
166         int i, l = 0;
167         if (b == 0) return -1;
168         if (bgzf_read(bp->fp, &b->tid, 4) == 0) return -1;
169         bgzf_read(bp->fp, &b->pos, 4);
170         bgzf_read(bp->fp, &b->qual, 4);
171         bgzf_read(bp->fp, &b->l_str, 4);
172         if (b->l_str > b->m_str) {
173                 b->m_str = b->l_str;
174                 kroundup32(b->m_str);
175                 b->str = realloc(b->str, b->m_str);
176         }
177         bgzf_read(bp->fp, b->str, b->l_str);
178         l = 12 + b->l_str;
179         bcf_sync(h->n_smpl, b);
180         for (i = 0; i < b->n_gi; ++i) {
181                 bgzf_read(bp->fp, b->gi[i].data, b->gi[i].len * h->n_smpl);
182                 l += b->gi[i].len * h->n_smpl;
183         }
184         return l;
185 }
186
187 int bcf_destroy(bcf1_t *b)
188 {
189         int i;
190         if (b == 0) return -1;
191         free(b->str);
192         for (i = 0; i < b->n_gi; ++i)
193                 free(b->gi[i].data);
194         free(b->gi);
195         free(b);
196         return 0;
197 }
198
199 static inline void fmt_str(const char *p, kstring_t *s)
200 {
201         if (*p == 0) kputc('.', s);
202         else kputs(p, s);
203 }
204
205 void bcf_fmt_core(const bcf_hdr_t *h, bcf1_t *b, kstring_t *s)
206 {
207         int i, j, x;
208         s->l = 0;
209         if (h->n_ref) kputs(h->ns[b->tid], s);
210         else kputw(b->tid, s);
211         kputc('\t', s);
212         kputw(b->pos + 1, s); kputc('\t', s);
213         fmt_str(b->str, s); kputc('\t', s);
214         fmt_str(b->ref, s); kputc('\t', s);
215         fmt_str(b->alt, s); kputc('\t', s);
216         kputw(b->qual, s); kputc('\t', s);
217         fmt_str(b->flt, s); kputc('\t', s);
218         fmt_str(b->info, s);
219         if (b->fmt[0]) {
220                 kputc('\t', s);
221                 fmt_str(b->fmt, s);
222         }
223         x = b->n_alleles * (b->n_alleles + 1) / 2;
224         if (b->n_gi == 0) return;
225         for (j = 0; j < h->n_smpl; ++j) {
226                 kputc('\t', s);
227                 for (i = 0; i < b->n_gi; ++i) {
228                         if (i) kputc(':', s);
229                         if (b->gi[i].fmt == bcf_str2int("PL", 2)) {
230                                 uint8_t *d = (uint8_t*)b->gi[i].data + j * x;
231                                 int k;
232                                 for (k = 0; k < x; ++k) {
233                                         if (k > 0) kputc(',', s);
234                                         kputw(d[k], s);
235                                 }
236                         } else if (b->gi[i].fmt == bcf_str2int("DP", 2)) {
237                                 kputw(((uint16_t*)b->gi[i].data)[j], s);
238                         } else if (b->gi[i].fmt == bcf_str2int("GQ", 2)) {
239                                 kputw(((uint8_t*)b->gi[i].data)[j], s);
240                         } else if (b->gi[i].fmt == bcf_str2int("GT", 2)) {
241                                 int y = ((uint8_t*)b->gi[i].data)[j];
242                                 kputc('0' + (y>>3&7), s);
243                                 kputc("/|"[y>>6&1], s);
244                                 kputc('0' + (y&7), s);
245                         }
246                 }
247         }
248 }
249
250 char *bcf_fmt(const bcf_hdr_t *h, bcf1_t *b)
251 {
252         kstring_t s;
253         s.l = s.m = 0; s.s = 0;
254         bcf_fmt_core(h, b, &s);
255         return s.s;
256 }