]> git.donarmstrong.com Git - samtools.git/blob - bcftools/bcf.c
faster for large sample size (in principle)
[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         uint32_t x;
150         int i, l = 0;
151         if (b == 0) return -1;
152         bgzf_write(bp->fp, &b->tid, 4);
153         bgzf_write(bp->fp, &b->pos, 4);
154         x = b->qual<<24 | b->l_str;
155         bgzf_write(bp->fp, &x, 4);
156         bgzf_write(bp->fp, b->str, b->l_str);
157         l = 12 + b->l_str;
158         for (i = 0; i < b->n_gi; ++i) {
159                 bgzf_write(bp->fp, b->gi[i].data, b->gi[i].len * h->n_smpl);
160                 l += b->gi[i].len * h->n_smpl;
161         }
162         return l;
163 }
164
165 int bcf_read(bcf_t *bp, const bcf_hdr_t *h, bcf1_t *b)
166 {
167         int i, l = 0;
168         uint32_t x;
169         if (b == 0) return -1;
170         if (bgzf_read(bp->fp, &b->tid, 4) == 0) return -1;
171         bgzf_read(bp->fp, &b->pos, 4);
172         bgzf_read(bp->fp, &x, 4);
173         b->qual = x >> 24; b->l_str = x << 8 >> 8;
174         if (b->l_str > b->m_str) {
175                 b->m_str = b->l_str;
176                 kroundup32(b->m_str);
177                 b->str = realloc(b->str, b->m_str);
178         }
179         bgzf_read(bp->fp, b->str, b->l_str);
180         l = 12 + b->l_str;
181         bcf_sync(h->n_smpl, b);
182         for (i = 0; i < b->n_gi; ++i) {
183                 bgzf_read(bp->fp, b->gi[i].data, b->gi[i].len * h->n_smpl);
184                 l += b->gi[i].len * h->n_smpl;
185         }
186         return l;
187 }
188
189 int bcf_destroy(bcf1_t *b)
190 {
191         int i;
192         if (b == 0) return -1;
193         free(b->str);
194         for (i = 0; i < b->n_gi; ++i)
195                 free(b->gi[i].data);
196         free(b->gi);
197         free(b);
198         return 0;
199 }
200
201 static inline void fmt_str(const char *p, kstring_t *s)
202 {
203         if (*p == 0) kputc('.', s);
204         else kputs(p, s);
205 }
206
207 void bcf_fmt_core(const bcf_hdr_t *h, bcf1_t *b, kstring_t *s)
208 {
209         int i, j, x;
210         s->l = 0;
211         if (h->n_ref) kputs(h->ns[b->tid], s);
212         else kputw(b->tid, s);
213         kputc('\t', s);
214         kputw(b->pos + 1, s); kputc('\t', s);
215         fmt_str(b->str, s); kputc('\t', s);
216         fmt_str(b->ref, s); kputc('\t', s);
217         fmt_str(b->alt, s); kputc('\t', s);
218         kputw(b->qual, s); kputc('\t', s);
219         fmt_str(b->flt, s); kputc('\t', s);
220         fmt_str(b->info, s);
221         if (b->fmt[0]) {
222                 kputc('\t', s);
223                 fmt_str(b->fmt, s);
224         }
225         x = b->n_alleles * (b->n_alleles + 1) / 2;
226         if (b->n_gi == 0) return;
227         for (j = 0; j < h->n_smpl; ++j) {
228                 kputc('\t', s);
229                 for (i = 0; i < b->n_gi; ++i) {
230                         if (i) kputc(':', s);
231                         if (b->gi[i].fmt == bcf_str2int("PL", 2)) {
232                                 uint8_t *d = (uint8_t*)b->gi[i].data + j * x;
233                                 int k;
234                                 for (k = 0; k < x; ++k) {
235                                         if (k > 0) kputc(',', s);
236                                         kputw(d[k], s);
237                                 }
238                         } else if (b->gi[i].fmt == bcf_str2int("DP", 2)) {
239                                 kputw(((uint16_t*)b->gi[i].data)[j], s);
240                         } else if (b->gi[i].fmt == bcf_str2int("GQ", 2)) {
241                                 kputw(((uint8_t*)b->gi[i].data)[j], s);
242                         } else if (b->gi[i].fmt == bcf_str2int("GT", 2)) {
243                                 int y = ((uint8_t*)b->gi[i].data)[j];
244                                 kputc('0' + (y>>3&7), s);
245                                 kputc("/|"[y>>6&1], s);
246                                 kputc('0' + (y&7), s);
247                         }
248                 }
249         }
250 }
251
252 char *bcf_fmt(const bcf_hdr_t *h, bcf1_t *b)
253 {
254         kstring_t s;
255         s.l = s.m = 0; s.s = 0;
256         bcf_fmt_core(h, b, &s);
257         return s.s;
258 }