]> git.donarmstrong.com Git - samtools.git/blob - bcftools/call1.c
9c0b498c18cb46e2dd673f0527f7c02d2e76dc3d
[samtools.git] / bcftools / call1.c
1 #include <unistd.h>
2 #include <stdlib.h>
3 #include <math.h>
4 #include <zlib.h>
5 #include <errno.h>
6 #include "bcf.h"
7 #include "prob1.h"
8 #include "kstring.h"
9
10 #include "khash.h"
11 KHASH_SET_INIT_INT64(set64)
12
13 #include "kseq.h"
14 KSTREAM_INIT(gzFile, gzread, 16384)
15
16 #define VC_NO_GENO 2
17 #define VC_BCFOUT  4
18 #define VC_CALL    8
19 #define VC_VARONLY 16
20 #define VC_VCFIN   32
21 #define VC_UNCOMP  64
22 #define VC_KEEPALT 256
23 #define VC_ACGT_ONLY 512
24 #define VC_QCALL   1024
25 #define VC_CALL_GT 2048
26 #define VC_ADJLD   4096
27 #define VC_NO_INDEL 8192
28 #define VC_ANNO_MAX 16384
29
30 typedef struct {
31         int flag, prior_type, n1, n_sub, *sublist;
32         char *fn_list, *prior_file, **subsam, *fn_dict;
33         uint8_t *ploidy;
34         double theta, pref, indel_frac;
35 } viewconf_t;
36
37 khash_t(set64) *bcf_load_pos(const char *fn, bcf_hdr_t *_h)
38 {
39         void *str2id;
40         gzFile fp;
41         kstream_t *ks;
42         int ret, dret, lineno = 1;
43         kstring_t *str;
44         khash_t(set64) *hash = 0;
45
46         hash = kh_init(set64);
47         str2id = bcf_build_refhash(_h);
48         str = calloc(1, sizeof(kstring_t));
49         fp = strcmp(fn, "-")? gzopen(fn, "r") : gzdopen(fileno(stdin), "r");
50         ks = ks_init(fp);
51         while (ks_getuntil(ks, 0, str, &dret) >= 0) {
52                 int tid = bcf_str2id(str2id, str->s);
53                 if (tid >= 0 && dret != '\n') {
54                         if (ks_getuntil(ks, 0, str, &dret) >= 0) {
55                                 uint64_t x = (uint64_t)tid<<32 | (atoi(str->s) - 1);
56                                 kh_put(set64, hash, x, &ret);
57                         } else break;
58                 } else fprintf(stderr, "[%s] %s is not a reference name (line %d).\n", __func__, str->s, lineno);
59                 if (dret != '\n') while ((dret = ks_getc(ks)) > 0 && dret != '\n');
60                 if (dret < 0) break;
61                 ++lineno;
62         }
63         bcf_str2id_destroy(str2id);
64         ks_destroy(ks);
65         gzclose(fp);
66         free(str->s); free(str);
67         return hash;
68 }
69
70 typedef struct {
71         double p[4];
72         int mq, depth, is_tested, d[4];
73 } anno16_t;
74
75 static double ttest(int n1, int n2, int a[4])
76 {
77         extern double kf_betai(double a, double b, double x);
78         double t, v, u1, u2;
79         if (n1 == 0 || n2 == 0 || n1 + n2 < 3) return 1.0;
80         u1 = (double)a[0] / n1; u2 = (double)a[2] / n2;
81         if (u1 <= u2) return 1.;
82         t = (u1 - u2) / sqrt(((a[1] - n1 * u1 * u1) + (a[3] - n2 * u2 * u2)) / (n1 + n2 - 2) * (1./n1 + 1./n2));
83         v = n1 + n2 - 2;
84 //      printf("%d,%d,%d,%d,%lf,%lf,%lf\n", a[0], a[1], a[2], a[3], t, u1, u2);
85         return t < 0.? 1. : .5 * kf_betai(.5*v, .5, v/(v+t*t));
86 }
87
88 static int test16_core(int anno[16], anno16_t *a)
89 {
90         extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
91         double left, right;
92         int i;
93         a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
94         memcpy(a->d, anno, 4 * sizeof(int));
95         a->depth = anno[0] + anno[1] + anno[2] + anno[3];
96         a->is_tested = (anno[0] + anno[1] > 0 && anno[2] + anno[3] > 0);
97         if (a->depth == 0) return -1;
98         a->mq = (int)(sqrt((anno[9] + anno[11]) / a->depth) + .499);
99         kt_fisher_exact(anno[0], anno[1], anno[2], anno[3], &left, &right, &a->p[0]);
100         for (i = 1; i < 4; ++i)
101                 a->p[i] = ttest(anno[0] + anno[1], anno[2] + anno[3], anno+4*i);
102         return 0;
103 }
104
105 static int test16(bcf1_t *b, anno16_t *a)
106 {
107         char *p;
108         int i, anno[16];
109         a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
110         a->d[0] = a->d[1] = a->d[2] = a->d[3] = 0.;
111         a->mq = a->depth = a->is_tested = 0;
112         if ((p = strstr(b->info, "I16=")) == 0) return -1;
113         p += 4;
114         for (i = 0; i < 16; ++i) {
115                 errno = 0; anno[i] = strtol(p, &p, 10);
116                 if (anno[i] == 0 && (errno == EINVAL || errno == ERANGE)) return -2;
117                 ++p;
118         }
119         return test16_core(anno, a);
120 }
121
122 static void rm_info(bcf1_t *b, const char *key)
123 {
124         char *p, *q;
125         if ((p = strstr(b->info, key)) == 0) return;
126         for (q = p; *q && *q != ';'; ++q);
127         if (p > b->info && *(p-1) == ';') --p;
128         memmove(p, q, b->l_str - (q - b->str));
129         b->l_str -= q - p;
130         bcf_sync(b);
131 }
132
133 static int update_bcf1(int n_smpl, bcf1_t *b, const bcf_p1aux_t *pa, const bcf_p1rst_t *pr, double pref, int flag)
134 {
135         kstring_t s;
136         int is_var = (pr->p_ref < pref);
137         double r = is_var? pr->p_ref : pr->p_var, fq;
138         anno16_t a;
139
140         test16(b, &a);
141         rm_info(b, "I16=");
142
143         memset(&s, 0, sizeof(kstring_t));
144         kputc('\0', &s); kputs(b->ref, &s); kputc('\0', &s);
145         kputs(b->alt, &s); kputc('\0', &s); kputc('\0', &s);
146         kputs(b->info, &s);
147         if (b->info[0]) kputc(';', &s);
148 //      ksprintf(&s, "AF1=%.4lg;AFE=%.4lg;CI95=%.4lg,%.4lg", 1.-pr->f_em, 1.-pr->f_exp, pr->cil, pr->cih);
149         ksprintf(&s, "AF1=%.4g;CI95=%.4g,%.4g", 1.-pr->f_em, pr->cil, pr->cih);
150         ksprintf(&s, ";DP4=%d,%d,%d,%d;MQ=%d", a.d[0], a.d[1], a.d[2], a.d[3], a.mq);
151         fq = pr->p_ref_folded < 0.5? -4.343 * log(pr->p_ref_folded) : 4.343 * log(pr->p_var_folded);
152         if (fq < -999) fq = -999;
153         if (fq > 999) fq = 999;
154         ksprintf(&s, ";FQ=%.3g", fq);
155         if (a.is_tested) {
156                 if (pr->pc[0] >= 0.) ksprintf(&s, ";PC4=%g,%g,%g,%g", pr->pc[0], pr->pc[1], pr->pc[2], pr->pc[3]);
157                 ksprintf(&s, ";PV4=%.2g,%.2g,%.2g,%.2g", a.p[0], a.p[1], a.p[2], a.p[3]);
158         }
159         kputc('\0', &s);
160         kputs(b->fmt, &s); kputc('\0', &s);
161         free(b->str);
162         b->m_str = s.m; b->l_str = s.l; b->str = s.s;
163         b->qual = r < 1e-100? 999 : -4.343 * log(r);
164         if (b->qual > 999) b->qual = 999;
165         bcf_sync(b);
166         if (!is_var) bcf_shrink_alt(b, 1);
167         else if (!(flag&VC_KEEPALT))
168                 bcf_shrink_alt(b, pr->rank0 < 2? 2 : pr->rank0+1);
169         if (is_var && (flag&VC_CALL_GT)) { // call individual genotype
170                 int i, x, old_n_gi = b->n_gi;
171                 s.m = b->m_str; s.l = b->l_str - 1; s.s = b->str;
172                 kputs(":GT:GQ", &s); kputc('\0', &s);
173                 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
174                 bcf_sync(b);
175                 for (i = 0; i < b->n_smpl; ++i) {
176                         x = bcf_p1_call_gt(pa, pr->f_em, i);
177                         ((uint8_t*)b->gi[old_n_gi].data)[i] = (x&3) == 0? 1<<3|1 : (x&3) == 1? 1 : 0;
178                         ((uint8_t*)b->gi[old_n_gi+1].data)[i] = x>>2;
179                 }
180         }
181         return is_var;
182 }
183
184 static char **read_samples(const char *fn, int *_n)
185 {
186         gzFile fp;
187         kstream_t *ks;
188         kstring_t s;
189         int dret, n = 0, max = 0;
190         char **sam = 0;
191         *_n = 0;
192         s.l = s.m = 0; s.s = 0;
193         fp = gzopen(fn, "r");
194         if (fp == 0) return 0; // fail to open file
195         ks = ks_init(fp);
196         while (ks_getuntil(ks, 0, &s, &dret) >= 0) {
197                 int l;
198                 if (max == n) {
199                         max = max? max<<1 : 4;
200                         sam = realloc(sam, sizeof(void*)*max);
201                 }
202                 l = s.l;
203                 sam[n] = malloc(s.l + 2);
204                 strcpy(sam[n], s.s);
205                 sam[n][l+1] = 2; // by default, diploid
206                 if (dret != '\n') {
207                         if (ks_getuntil(ks, 0, &s, &dret) >= 0) { // read ploidy, 1 or 2
208                                 int x = (int)s.s[0] - '0';
209                                 if (x == 1 || x == 2) sam[n][l+1] = x;
210                                 else fprintf(stderr, "(%s) ploidy can only be 1 or 2; assume diploid\n", __func__);
211                         }
212                         if (dret != '\n') ks_getuntil(ks, '\n', &s, &dret);
213                 }
214                 ++n;
215         }
216         ks_destroy(ks);
217         gzclose(fp);
218         free(s.s);
219         *_n = n;
220         return sam;
221 }
222
223 static void write_header(bcf_hdr_t *h)
224 {
225         kstring_t str;
226         str.l = h->l_txt? h->l_txt - 1 : 0;
227         str.m = str.l + 1; str.s = h->txt;
228         if (!strstr(str.s, "##INFO=<ID=DP,"))
229                 kputs("##INFO=<ID=DP,Number=1,Type=Integer,Description=\"Raw read depth\">\n", &str);
230         if (!strstr(str.s, "##INFO=<ID=DP4,"))
231                 kputs("##INFO=<ID=DP4,Number=4,Type=Integer,Description=\"# high-quality ref-forward bases, ref-reverse, alt-forward and alt-reverse bases\">\n", &str);
232         if (!strstr(str.s, "##INFO=<ID=MQ,"))
233                 kputs("##INFO=<ID=MQ,Number=1,Type=Integer,Description=\"Root-mean-square mapping quality of covering reads\">\n", &str);
234         if (!strstr(str.s, "##INFO=<ID=FQ,"))
235                 kputs("##INFO=<ID=FQ,Number=1,Type=Float,Description=\"Phred probability that sample chromosomes are not all the same\">\n", &str);
236         if (!strstr(str.s, "##INFO=<ID=AF1,"))
237                 kputs("##INFO=<ID=AF1,Number=1,Type=Float,Description=\"Max-likelihood estimate of the site allele frequency of the first ALT allele\">\n", &str);
238         if (!strstr(str.s, "##INFO=<ID=CI95,"))
239                 kputs("##INFO=<ID=CI95,Number=2,Type=Float,Description=\"Equal-tail Bayesian credible interval of the site allele frequency at the 95% level\">\n", &str);
240         if (!strstr(str.s, "##INFO=<ID=PV4,"))
241                 kputs("##INFO=<ID=PV4,Number=4,Type=Float,Description=\"P-values for strand bias, baseQ bias, mapQ bias and tail distance bias\">\n", &str);
242     if (!strstr(str.s, "##INFO=<ID=INDEL,"))
243         kputs("##INFO=<ID=INDEL,Number=0,Type=Flag,Description=\"Indicates that the variant is an INDEL.\">\n", &str);
244     if (!strstr(str.s, "##INFO=<ID=GT,"))
245         kputs("##FORMAT=<ID=GT,Number=1,Type=String,Description=\"Genotype\">\n", &str);
246     if (!strstr(str.s, "##FORMAT=<ID=GQ,"))
247         kputs("##FORMAT=<ID=GQ,Number=1,Type=Integer,Description=\"Genotype Quality\">\n", &str);
248     if (!strstr(str.s, "##FORMAT=<ID=GL,"))
249         kputs("##FORMAT=<ID=GL,Number=3,Type=Float,Description=\"Likelihoods for RR,RA,AA genotypes (R=ref,A=alt)\">\n", &str);
250         if (!strstr(str.s, "##FORMAT=<ID=DP,"))
251                 kputs("##FORMAT=<ID=DP,Number=1,Type=Integer,Description=\"# high-quality bases\">\n", &str);
252         if (!strstr(str.s, "##FORMAT=<ID=SP,"))
253                 kputs("##FORMAT=<ID=SP,Number=1,Type=Integer,Description=\"Phred-scaled strand bias P-value\">\n", &str);
254         if (!strstr(str.s, "##FORMAT=<ID=PL,"))
255                 kputs("##FORMAT=<ID=PL,Number=-1,Type=Integer,Description=\"List of Phred-scaled genotype likelihoods, number of values is (#ALT+1)*(#ALT+2)/2\">\n", &str);
256         h->l_txt = str.l + 1; h->txt = str.s;
257 }
258
259 double bcf_ld_freq(const bcf1_t *b0, const bcf1_t *b1, double f[4]);
260
261 int bcfview(int argc, char *argv[])
262 {
263         extern int bcf_2qcall(bcf_hdr_t *h, bcf1_t *b);
264         extern void bcf_p1_indel_prior(bcf_p1aux_t *ma, double x);
265         extern int bcf_fix_gt(bcf1_t *b);
266         extern int bcf_anno_max(bcf1_t *b);
267         bcf_t *bp, *bout = 0;
268         bcf1_t *b, *blast;
269         int c;
270         uint64_t n_processed = 0;
271         viewconf_t vc;
272         bcf_p1aux_t *p1 = 0;
273         bcf_hdr_t *hin, *hout;
274         int tid, begin, end;
275         char moder[4], modew[4];
276         khash_t(set64) *hash = 0;
277
278         tid = begin = end = -1;
279         memset(&vc, 0, sizeof(viewconf_t));
280         vc.prior_type = vc.n1 = -1; vc.theta = 1e-3; vc.pref = 0.5; vc.indel_frac = -1.;
281         while ((c = getopt(argc, argv, "N1:l:cHAGvbSuP:t:p:QgLi:IMs:D:")) >= 0) {
282                 switch (c) {
283                 case '1': vc.n1 = atoi(optarg); break;
284                 case 'l': vc.fn_list = strdup(optarg); break;
285                 case 'D': vc.fn_dict = strdup(optarg); break;
286                 case 'N': vc.flag |= VC_ACGT_ONLY; break;
287                 case 'G': vc.flag |= VC_NO_GENO; break;
288                 case 'A': vc.flag |= VC_KEEPALT; break;
289                 case 'b': vc.flag |= VC_BCFOUT; break;
290                 case 'S': vc.flag |= VC_VCFIN; break;
291                 case 'c': vc.flag |= VC_CALL; break;
292                 case 'v': vc.flag |= VC_VARONLY | VC_CALL; break;
293                 case 'u': vc.flag |= VC_UNCOMP | VC_BCFOUT; break;
294                 case 'g': vc.flag |= VC_CALL_GT | VC_CALL; break;
295                 case 'I': vc.flag |= VC_NO_INDEL; break;
296                 case 'M': vc.flag |= VC_ANNO_MAX; break;
297                 case 't': vc.theta = atof(optarg); break;
298                 case 'p': vc.pref = atof(optarg); break;
299                 case 'i': vc.indel_frac = atof(optarg); break;
300                 case 'Q': vc.flag |= VC_QCALL; break;
301                 case 'L': vc.flag |= VC_ADJLD; break;
302                 case 's': vc.subsam = read_samples(optarg, &vc.n_sub);
303                         vc.ploidy = calloc(vc.n_sub + 1, 1);
304                         for (tid = 0; tid < vc.n_sub; ++tid) vc.ploidy[tid] = vc.subsam[tid][strlen(vc.subsam[tid]) + 1];
305                         tid = -1;
306                         break;
307                 case 'P':
308                         if (strcmp(optarg, "full") == 0) vc.prior_type = MC_PTYPE_FULL;
309                         else if (strcmp(optarg, "cond2") == 0) vc.prior_type = MC_PTYPE_COND2;
310                         else if (strcmp(optarg, "flat") == 0) vc.prior_type = MC_PTYPE_FLAT;
311                         else vc.prior_file = strdup(optarg);
312                         break;
313                 }
314         }
315         if (argc == optind) {
316                 fprintf(stderr, "\n");
317                 fprintf(stderr, "Usage:   bcftools view [options] <in.bcf> [reg]\n\n");
318                 fprintf(stderr, "Options: -c        SNP calling\n");
319                 fprintf(stderr, "         -v        output potential variant sites only (force -c)\n");
320                 fprintf(stderr, "         -g        call genotypes at variant sites (force -c)\n");
321                 fprintf(stderr, "         -b        output BCF instead of VCF\n");
322                 fprintf(stderr, "         -u        uncompressed BCF output (force -b)\n");
323                 fprintf(stderr, "         -S        input is VCF\n");
324                 fprintf(stderr, "         -A        keep all possible alternate alleles at variant sites\n");
325                 fprintf(stderr, "         -G        suppress all individual genotype information\n");
326                 fprintf(stderr, "         -H        perform Hardy-Weinberg test (slower)\n");
327                 fprintf(stderr, "         -N        skip sites where REF is not A/C/G/T\n");
328                 fprintf(stderr, "         -Q        output the QCALL likelihood format\n");
329                 fprintf(stderr, "         -L        calculate LD for adjacent sites\n");
330                 fprintf(stderr, "         -I        skip indels\n");
331                 fprintf(stderr, "         -D FILE   sequence dictionary for VCF->BCF conversion [null]\n");
332                 fprintf(stderr, "         -1 INT    number of group-1 samples [0]\n");
333                 fprintf(stderr, "         -l FILE   list of sites to output [all sites]\n");
334                 fprintf(stderr, "         -t FLOAT  scaled substitution mutation rate [%.4g]\n", vc.theta);
335                 fprintf(stderr, "         -i FLOAT  indel-to-substitution ratio [%.4g]\n", vc.indel_frac);
336                 fprintf(stderr, "         -p FLOAT  variant if P(ref|D)<FLOAT [%.3g]\n", vc.pref);
337                 fprintf(stderr, "         -P STR    type of prior: full, cond2, flat [full]\n");
338                 fprintf(stderr, "         -s FILE   list of samples to use [all samples]\n");
339                 fprintf(stderr, "\n");
340                 return 1;
341         }
342
343         if ((vc.flag & VC_VCFIN) && (vc.flag & VC_BCFOUT) && vc.fn_dict == 0) {
344                 fprintf(stderr, "[%s] For VCF->BCF conversion please specify the sequence dictionary with -D\n", __func__);
345                 return 1;
346         }
347         b = calloc(1, sizeof(bcf1_t));
348         blast = calloc(1, sizeof(bcf1_t));
349         strcpy(moder, "r");
350         if (!(vc.flag & VC_VCFIN)) strcat(moder, "b");
351         strcpy(modew, "w");
352         if (vc.flag & VC_BCFOUT) strcat(modew, "b");
353         if (vc.flag & VC_UNCOMP) strcat(modew, "u");
354         bp = vcf_open(argv[optind], moder);
355         hin = hout = vcf_hdr_read(bp);
356         if (vc.fn_dict && (vc.flag & VC_VCFIN))
357                 vcf_dictread(bp, hin, vc.fn_dict);
358         bout = vcf_open("-", modew);
359         if (!(vc.flag & VC_QCALL)) {
360                 if (vc.n_sub) {
361                         vc.sublist = calloc(vc.n_sub, sizeof(int));
362                         hout = bcf_hdr_subsam(hin, vc.n_sub, vc.subsam, vc.sublist);
363                 }
364                 if (vc.flag & VC_CALL) write_header(hout);
365                 vcf_hdr_write(bout, hout);
366         }
367         if (vc.flag & VC_CALL) {
368                 p1 = bcf_p1_init(hout->n_smpl, vc.ploidy);
369                 if (vc.prior_file) {
370                         if (bcf_p1_read_prior(p1, vc.prior_file) < 0) {
371                                 fprintf(stderr, "[%s] fail to read the prior AFS.\n", __func__);
372                                 return 1;
373                         }
374                 } else bcf_p1_init_prior(p1, vc.prior_type, vc.theta);
375                 if (vc.n1 > 0) {
376                         bcf_p1_set_n1(p1, vc.n1);
377                         bcf_p1_init_subprior(p1, vc.prior_type, vc.theta);
378                 }
379                 if (vc.indel_frac > 0.) bcf_p1_indel_prior(p1, vc.indel_frac); // otherwise use the default indel_frac
380         }
381         if (vc.fn_list) hash = bcf_load_pos(vc.fn_list, hin);
382         if (optind + 1 < argc && !(vc.flag&VC_VCFIN)) {
383                 void *str2id = bcf_build_refhash(hout);
384                 if (bcf_parse_region(str2id, argv[optind+1], &tid, &begin, &end) >= 0) {
385                         bcf_idx_t *idx;
386                         idx = bcf_idx_load(argv[optind]);
387                         if (idx) {
388                                 uint64_t off;
389                                 off = bcf_idx_query(idx, tid, begin);
390                                 if (off == 0) {
391                                         fprintf(stderr, "[%s] no records in the query region.\n", __func__);
392                                         return 1; // FIXME: a lot of memory leaks...
393                                 }
394                                 bgzf_seek(bp->fp, off, SEEK_SET);
395                                 bcf_idx_destroy(idx);
396                         }
397                 }
398         }
399         while (vcf_read(bp, hin, b) > 0) {
400                 int is_indel;
401                 if (vc.n_sub) bcf_subsam(vc.n_sub, vc.sublist, b);
402                 is_indel = bcf_is_indel(b);
403                 if ((vc.flag & VC_NO_INDEL) && is_indel) continue;
404                 if ((vc.flag & VC_ACGT_ONLY) && !is_indel) {
405                         int x;
406                         if (b->ref[0] == 0 || b->ref[1] != 0) continue;
407                         x = toupper(b->ref[0]);
408                         if (x != 'A' && x != 'C' && x != 'G' && x != 'T') continue;
409                 }
410                 if (hash) {
411                         uint64_t x = (uint64_t)b->tid<<32 | b->pos;
412                         khint_t k = kh_get(set64, hash, x);
413                         if (kh_size(hash) == 0) break;
414                         if (k == kh_end(hash)) continue;
415                         kh_del(set64, hash, k);
416                 }
417                 if (tid >= 0) {
418                         int l = strlen(b->ref);
419                         l = b->pos + (l > 0? l : 1);
420                         if (b->tid != tid || b->pos >= end) break;
421                         if (!(l > begin && end > b->pos)) continue;
422                 }
423                 ++n_processed;
424                 if (vc.flag & VC_QCALL) { // output QCALL format; STOP here
425                         bcf_2qcall(hout, b);
426                         continue;
427                 }
428                 if (vc.flag & (VC_CALL|VC_ADJLD)) bcf_gl2pl(b);
429                 if (vc.flag & VC_CALL) { // call variants
430                         bcf_p1rst_t pr;
431                         bcf_p1_cal(b, p1, &pr); // pr.g[3] is not calculated here
432                         if (n_processed % 100000 == 0) {
433                                 fprintf(stderr, "[%s] %ld sites processed.\n", __func__, (long)n_processed);
434                                 bcf_p1_dump_afs(p1);
435                         }
436                         if (pr.p_ref >= vc.pref && (vc.flag & VC_VARONLY)) continue;
437                         update_bcf1(hout->n_smpl, b, p1, &pr, vc.pref, vc.flag);
438                 }
439                 if (vc.flag & VC_ADJLD) { // compute LD
440                         double f[4], r2;
441                         if ((r2 = bcf_ld_freq(blast, b, f)) >= 0) {
442                                 kstring_t s;
443                                 s.m = s.l = 0; s.s = 0;
444                                 if (*b->info) kputc(';', &s);
445                                 ksprintf(&s, "NEIR=%.3f;NEIF=%.3f,%.3f", r2, f[0]+f[2], f[0]+f[1]);
446                                 bcf_append_info(b, s.s, s.l);
447                                 free(s.s);
448                         }
449                         bcf_cpy(blast, b);
450                 }
451                 if (vc.flag & VC_ANNO_MAX) bcf_anno_max(b);
452                 if (vc.flag & VC_NO_GENO) { // do not output GENO fields
453                         b->n_gi = 0;
454                         b->fmt[0] = '\0';
455                         b->l_str = b->fmt - b->str + 1;
456                 } else bcf_fix_gt(b);
457                 vcf_write(bout, hout, b);
458         }
459         if (vc.prior_file) free(vc.prior_file);
460         if (vc.flag & VC_CALL) bcf_p1_dump_afs(p1);
461         if (hin != hout) bcf_hdr_destroy(hout);
462         bcf_hdr_destroy(hin);
463         bcf_destroy(b); bcf_destroy(blast);
464         vcf_close(bp); vcf_close(bout);
465         if (hash) kh_destroy(set64, hash);
466         if (vc.fn_list) free(vc.fn_list);
467         if (vc.fn_dict) free(vc.fn_dict);
468         if (vc.ploidy) free(vc.ploidy);
469         if (vc.n_sub) {
470                 int i;
471                 for (i = 0; i < vc.n_sub; ++i) free(vc.subsam[i]);
472                 free(vc.subsam); free(vc.sublist);
473         }
474         if (p1) bcf_p1_destroy(p1);
475         return 0;
476 }