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