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