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