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