12 #define srand48(x) srand(x)
13 #define lrand48() rand()
17 KSTREAM_INIT(gzFile, gzread, 16384)
25 #define VC_KEEPALT 256
26 #define VC_ACGT_ONLY 512
28 #define VC_CALL_GT 2048
30 #define VC_NO_INDEL 8192
31 #define VC_ANNO_MAX 16384
32 #define VC_FIX_PL 32768
34 #define VC_PAIRCALL 0x20000
35 #define VC_QCNT 0x40000
36 #define VC_INDEL_ONLY 0x80000
39 int flag, prior_type, n1, n_sub, *sublist, n_perm;
41 char *prior_file, **subsam, *fn_dict;
43 double theta, pref, indel_frac, min_perm_p, min_smpl_frac, min_lrt;
47 void *bed_read(const char *fn);
48 void bed_destroy(void *_h);
49 int bed_overlap(const void *_h, const char *chr, int beg, int end);
53 int mq, depth, is_tested, d[4];
56 static double ttest(int n1, int n2, int a[4])
58 extern double kf_betai(double a, double b, double x);
60 if (n1 == 0 || n2 == 0 || n1 + n2 < 3) return 1.0;
61 u1 = (double)a[0] / n1; u2 = (double)a[2] / n2;
62 if (u1 <= u2) return 1.;
63 t = (u1 - u2) / sqrt(((a[1] - n1 * u1 * u1) + (a[3] - n2 * u2 * u2)) / (n1 + n2 - 2) * (1./n1 + 1./n2));
65 // printf("%d,%d,%d,%d,%lf,%lf,%lf\n", a[0], a[1], a[2], a[3], t, u1, u2);
66 return t < 0.? 1. : .5 * kf_betai(.5*v, .5, v/(v+t*t));
69 static int test16_core(int anno[16], anno16_t *a)
71 extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
74 a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
75 memcpy(a->d, anno, 4 * sizeof(int));
76 a->depth = anno[0] + anno[1] + anno[2] + anno[3];
77 a->is_tested = (anno[0] + anno[1] > 0 && anno[2] + anno[3] > 0);
78 if (a->depth == 0) return -1;
79 a->mq = (int)(sqrt((anno[9] + anno[11]) / a->depth) + .499);
80 kt_fisher_exact(anno[0], anno[1], anno[2], anno[3], &left, &right, &a->p[0]);
81 for (i = 1; i < 4; ++i)
82 a->p[i] = ttest(anno[0] + anno[1], anno[2] + anno[3], anno+4*i);
86 static int test16(bcf1_t *b, anno16_t *a)
90 a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
91 a->d[0] = a->d[1] = a->d[2] = a->d[3] = 0.;
92 a->mq = a->depth = a->is_tested = 0;
93 if ((p = strstr(b->info, "I16=")) == 0) return -1;
95 for (i = 0; i < 16; ++i) {
96 errno = 0; anno[i] = strtol(p, &p, 10);
97 if (anno[i] == 0 && (errno == EINVAL || errno == ERANGE)) return -2;
100 return test16_core(anno, a);
103 static void rm_info(bcf1_t *b, const char *key)
106 if ((p = strstr(b->info, key)) == 0) return;
107 for (q = p; *q && *q != ';'; ++q);
108 if (p > b->info && *(p-1) == ';') --p;
109 memmove(p, q, b->l_str - (q - b->str));
114 static int update_bcf1(bcf1_t *b, const bcf_p1aux_t *pa, const bcf_p1rst_t *pr, double pref, int flag, double em[10], int cons_llr, int64_t cons_gt)
121 has_I16 = test16(b, &a) >= 0? 1 : 0;
122 rm_info(b, "I16="); // FIXME: probably this function has a bug. If I move it below, I16 will not be removed!
124 memset(&s, 0, sizeof(kstring_t));
125 kputc('\0', &s); kputs(b->ref, &s); kputc('\0', &s);
126 kputs(b->alt, &s); kputc('\0', &s); kputc('\0', &s);
128 if (b->info[0]) kputc(';', &s);
130 if (em[0] >= 0) ksprintf(&s, "AF1=%.4g", 1 - em[0]);
131 if (em[4] >= 0 && em[4] <= 0.05) ksprintf(&s, ";G3=%.4g,%.4g,%.4g;HWE=%.3g", em[3], em[2], em[1], em[4]);
132 if (em[5] >= 0 && em[6] >= 0) ksprintf(&s, ";AF2=%.4g,%.4g", 1 - em[5], 1 - em[6]);
133 if (em[7] >= 0) ksprintf(&s, ";LRT=%.3g", em[7]);
134 if (em[8] >= 0) ksprintf(&s, ";LRT2=%.3g", em[8]);
137 ksprintf(&s, ";CLR=%d", cons_llr);
139 ksprintf(&s, ";UGT=%c%c%c;CGT=%c%c%c", cons_gt&0xff, cons_gt>>8&0xff, cons_gt>>16&0xff,
140 cons_gt>>32&0xff, cons_gt>>40&0xff, cons_gt>>48&0xff);
142 if (pr == 0) { // if pr is unset, return
143 kputc('\0', &s); kputs(b->fmt, &s); kputc('\0', &s);
145 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
150 is_var = (pr->p_ref < pref);
151 r = is_var? pr->p_ref : pr->p_var;
153 // ksprintf(&s, ";CI95=%.4g,%.4g", pr->cil, pr->cih); // FIXME: when EM is not used, ";" should be omitted!
154 ksprintf(&s, ";AC1=%d", pr->ac);
155 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);
156 fq = pr->p_ref_folded < 0.5? -4.343 * log(pr->p_ref_folded) : 4.343 * log(pr->p_var_folded);
157 if (fq < -999) fq = -999;
158 if (fq > 999) fq = 999;
159 ksprintf(&s, ";FQ=%.3g", fq);
160 if (pr->cmp[0] >= 0.) { // two sample groups
162 for (i = 1; i < 3; ++i) {
163 double x = pr->cmp[i] + pr->cmp[0]/2.;
164 q[i] = x == 0? 255 : (int)(-4.343 * log(x) + .499);
165 if (q[i] > 255) q[i] = 255;
167 if (pr->perm_rank >= 0) ksprintf(&s, ";PR=%d", pr->perm_rank);
168 // ksprintf(&s, ";LRT3=%.3g", pr->lrt);
169 ksprintf(&s, ";PCHI2=%.3g;PC2=%d,%d", q[1], q[2], pr->p_chi2);
171 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]);
173 kputs(b->fmt, &s); kputc('\0', &s);
175 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
176 b->qual = r < 1e-100? 999 : -4.343 * log(r);
177 if (b->qual > 999) b->qual = 999;
179 if (!is_var) bcf_shrink_alt(b, 1);
180 else if (!(flag&VC_KEEPALT))
181 bcf_shrink_alt(b, pr->rank0 < 2? 2 : pr->rank0+1);
182 if (is_var && (flag&VC_CALL_GT)) { // call individual genotype
183 int i, x, old_n_gi = b->n_gi;
184 s.m = b->m_str; s.l = b->l_str - 1; s.s = b->str;
185 kputs(":GT:GQ", &s); kputc('\0', &s);
186 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
188 for (i = 0; i < b->n_smpl; ++i) {
189 x = bcf_p1_call_gt(pa, pr->f_exp, i);
190 ((uint8_t*)b->gi[old_n_gi].data)[i] = (x&3) == 0? 1<<3|1 : (x&3) == 1? 1 : 0;
191 ((uint8_t*)b->gi[old_n_gi+1].data)[i] = x>>2;
197 static char **read_samples(const char *fn, int *_n)
202 int dret, n = 0, max = 0;
205 s.l = s.m = 0; s.s = 0;
206 fp = gzopen(fn, "r");
207 if (fp == 0) return 0; // fail to open file
209 while (ks_getuntil(ks, 0, &s, &dret) >= 0) {
212 max = max? max<<1 : 4;
213 sam = realloc(sam, sizeof(void*)*max);
216 sam[n] = malloc(s.l + 2);
218 sam[n][l+1] = 2; // by default, diploid
220 if (ks_getuntil(ks, 0, &s, &dret) >= 0) { // read ploidy, 1 or 2
221 int x = (int)s.s[0] - '0';
222 if (x == 1 || x == 2) sam[n][l+1] = x;
223 else fprintf(stderr, "(%s) ploidy can only be 1 or 2; assume diploid\n", __func__);
225 if (dret != '\n') ks_getuntil(ks, '\n', &s, &dret);
236 static void write_header(bcf_hdr_t *h)
239 str.l = h->l_txt? h->l_txt - 1 : 0;
240 str.m = str.l + 1; str.s = h->txt;
241 if (!strstr(str.s, "##INFO=<ID=DP,"))
242 kputs("##INFO=<ID=DP,Number=1,Type=Integer,Description=\"Raw read depth\">\n", &str);
243 if (!strstr(str.s, "##INFO=<ID=DP4,"))
244 kputs("##INFO=<ID=DP4,Number=4,Type=Integer,Description=\"# high-quality ref-forward bases, ref-reverse, alt-forward and alt-reverse bases\">\n", &str);
245 if (!strstr(str.s, "##INFO=<ID=MQ,"))
246 kputs("##INFO=<ID=MQ,Number=1,Type=Integer,Description=\"Root-mean-square mapping quality of covering reads\">\n", &str);
247 if (!strstr(str.s, "##INFO=<ID=FQ,"))
248 kputs("##INFO=<ID=FQ,Number=1,Type=Float,Description=\"Phred probability of all samples being the same\">\n", &str);
249 if (!strstr(str.s, "##INFO=<ID=AF1,"))
250 kputs("##INFO=<ID=AF1,Number=1,Type=Float,Description=\"Max-likelihood estimate of the first ALT allele frequency (assuming HWE)\">\n", &str);
251 if (!strstr(str.s, "##INFO=<ID=AC1,"))
252 kputs("##INFO=<ID=AC1,Number=1,Type=Float,Description=\"Max-likelihood estimate of the first ALT allele count (no HWE assumption)\">\n", &str);
253 if (!strstr(str.s, "##INFO=<ID=G3,"))
254 kputs("##INFO=<ID=G3,Number=3,Type=Float,Description=\"ML estimate of genotype frequencies\">\n", &str);
255 if (!strstr(str.s, "##INFO=<ID=HWE,"))
256 kputs("##INFO=<ID=HWE,Number=1,Type=Float,Description=\"Chi^2 based HWE test P-value based on G3\">\n", &str);
257 if (!strstr(str.s, "##INFO=<ID=CLR,"))
258 kputs("##INFO=<ID=CLR,Number=1,Type=Integer,Description=\"Log ratio of genotype likelihoods with and without the constraint\">\n", &str);
259 if (!strstr(str.s, "##INFO=<ID=UGT,"))
260 kputs("##INFO=<ID=UGT,Number=1,Type=String,Description=\"The most probable unconstrained genotype configuration in the trio\">\n", &str);
261 if (!strstr(str.s, "##INFO=<ID=CGT,"))
262 kputs("##INFO=<ID=CGT,Number=1,Type=String,Description=\"The most probable constrained genotype configuration in the trio\">\n", &str);
263 // if (!strstr(str.s, "##INFO=<ID=CI95,"))
264 // 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);
265 if (!strstr(str.s, "##INFO=<ID=PV4,"))
266 kputs("##INFO=<ID=PV4,Number=4,Type=Float,Description=\"P-values for strand bias, baseQ bias, mapQ bias and tail distance bias\">\n", &str);
267 if (!strstr(str.s, "##INFO=<ID=INDEL,"))
268 kputs("##INFO=<ID=INDEL,Number=0,Type=Flag,Description=\"Indicates that the variant is an INDEL.\">\n", &str);
269 if (!strstr(str.s, "##INFO=<ID=PC2,"))
270 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);
271 if (!strstr(str.s, "##INFO=<ID=PCHI2,"))
272 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);
273 if (!strstr(str.s, "##INFO=<ID=QCHI2,"))
274 kputs("##INFO=<ID=QCHI2,Number=1,Type=Integer,Description=\"Phred scaled PCHI2.\">\n", &str);
275 if (!strstr(str.s, "##INFO=<ID=RP,"))
276 kputs("##INFO=<ID=PR,Number=1,Type=Integer,Description=\"# permutations yielding a smaller PCHI2.\">\n", &str);
277 if (!strstr(str.s, "##INFO=<ID=VDB,"))
278 kputs("##INFO=<ID=VDB,Number=1,Type=Float,Description=\"Variant Distance Bias\">\n", &str);
279 if (!strstr(str.s, "##FORMAT=<ID=GT,"))
280 kputs("##FORMAT=<ID=GT,Number=1,Type=String,Description=\"Genotype\">\n", &str);
281 if (!strstr(str.s, "##FORMAT=<ID=GQ,"))
282 kputs("##FORMAT=<ID=GQ,Number=1,Type=Integer,Description=\"Genotype Quality\">\n", &str);
283 if (!strstr(str.s, "##FORMAT=<ID=GL,"))
284 kputs("##FORMAT=<ID=GL,Number=3,Type=Float,Description=\"Likelihoods for RR,RA,AA genotypes (R=ref,A=alt)\">\n", &str);
285 if (!strstr(str.s, "##FORMAT=<ID=DP,"))
286 kputs("##FORMAT=<ID=DP,Number=1,Type=Integer,Description=\"# high-quality bases\">\n", &str);
287 if (!strstr(str.s, "##FORMAT=<ID=DV,"))
288 kputs("##FORMAT=<ID=DV,Number=1,Type=Integer,Description=\"# high-quality non-reference bases\">\n", &str);
289 if (!strstr(str.s, "##FORMAT=<ID=SP,"))
290 kputs("##FORMAT=<ID=SP,Number=1,Type=Integer,Description=\"Phred-scaled strand bias P-value\">\n", &str);
291 if (!strstr(str.s, "##FORMAT=<ID=PL,"))
292 kputs("##FORMAT=<ID=PL,Number=G,Type=Integer,Description=\"List of Phred-scaled genotype likelihoods\">\n", &str);
293 h->l_txt = str.l + 1; h->txt = str.s;
296 double bcf_pair_freq(const bcf1_t *b0, const bcf1_t *b1, double f[4]);
298 int bcfview(int argc, char *argv[])
300 extern int bcf_2qcall(bcf_hdr_t *h, bcf1_t *b);
301 extern void bcf_p1_indel_prior(bcf_p1aux_t *ma, double x);
302 extern int bcf_fix_gt(bcf1_t *b);
303 extern int bcf_anno_max(bcf1_t *b);
304 extern int bcf_shuffle(bcf1_t *b, int seed);
305 extern uint32_t *bcf_trio_prep(int is_x, int is_son);
306 extern int bcf_trio_call(uint32_t *prep, const bcf1_t *b, int *llr, int64_t *gt);
307 extern int bcf_pair_call(const bcf1_t *b);
308 extern int bcf_min_diff(const bcf1_t *b);
310 bcf_t *bp, *bout = 0;
313 uint64_t n_processed = 0, qcnt[256];
316 bcf_hdr_t *hin, *hout;
318 char moder[4], modew[4];
320 tid = begin = end = -1;
321 memset(&vc, 0, sizeof(viewconf_t));
322 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;
323 memset(qcnt, 0, 8 * 256);
324 while ((c = getopt(argc, argv, "FN1:l:cC:eHAGvbSuP:t:p:QgLi:IMs:D:U:X:d:T:Yw")) >= 0) {
326 case '1': vc.n1 = atoi(optarg); break;
327 case 'l': vc.bed = bed_read(optarg); break;
328 case 'D': vc.fn_dict = strdup(optarg); break;
329 case 'F': vc.flag |= VC_FIX_PL; break;
330 case 'N': vc.flag |= VC_ACGT_ONLY; break;
331 case 'G': vc.flag |= VC_NO_GENO; break;
332 case 'A': vc.flag |= VC_KEEPALT; break;
333 case 'b': vc.flag |= VC_BCFOUT; break;
334 case 'S': vc.flag |= VC_VCFIN; break;
335 case 'c': vc.flag |= VC_CALL; break;
336 case 'e': vc.flag |= VC_EM; break;
337 case 'v': vc.flag |= VC_VARONLY | VC_CALL; break;
338 case 'u': vc.flag |= VC_UNCOMP | VC_BCFOUT; break;
339 case 'g': vc.flag |= VC_CALL_GT | VC_CALL; break;
340 case 'I': vc.flag |= VC_NO_INDEL; break;
341 case 'w': vc.flag |= VC_INDEL_ONLY; break;
342 case 'M': vc.flag |= VC_ANNO_MAX; break;
343 case 'Y': vc.flag |= VC_QCNT; break;
344 case 't': vc.theta = atof(optarg); break;
345 case 'p': vc.pref = atof(optarg); break;
346 case 'i': vc.indel_frac = atof(optarg); break;
347 case 'Q': vc.flag |= VC_QCALL; break;
348 case 'L': vc.flag |= VC_ADJLD; break;
349 case 'U': vc.n_perm = atoi(optarg); break;
350 case 'C': vc.min_lrt = atof(optarg); break;
351 case 'X': vc.min_perm_p = atof(optarg); break;
352 case 'd': vc.min_smpl_frac = atof(optarg); break;
353 case 's': vc.subsam = read_samples(optarg, &vc.n_sub);
354 vc.ploidy = calloc(vc.n_sub + 1, 1);
355 for (tid = 0; tid < vc.n_sub; ++tid) vc.ploidy[tid] = vc.subsam[tid][strlen(vc.subsam[tid]) + 1];
359 if (strcmp(optarg, "trioauto") == 0) vc.trio_aux = bcf_trio_prep(0, 0);
360 else if (strcmp(optarg, "trioxd") == 0) vc.trio_aux = bcf_trio_prep(1, 0);
361 else if (strcmp(optarg, "trioxs") == 0) vc.trio_aux = bcf_trio_prep(1, 1);
362 else if (strcmp(optarg, "pair") == 0) vc.flag |= VC_PAIRCALL;
364 fprintf(stderr, "[%s] Option '-T' can only take value trioauto, trioxd or trioxs.\n", __func__);
369 if (strcmp(optarg, "full") == 0) vc.prior_type = MC_PTYPE_FULL;
370 else if (strcmp(optarg, "cond2") == 0) vc.prior_type = MC_PTYPE_COND2;
371 else if (strcmp(optarg, "flat") == 0) vc.prior_type = MC_PTYPE_FLAT;
372 else vc.prior_file = strdup(optarg);
376 if (argc == optind) {
377 fprintf(stderr, "\n");
378 fprintf(stderr, "Usage: bcftools view [options] <in.bcf> [reg]\n\n");
379 fprintf(stderr, "Input/output options:\n\n");
380 fprintf(stderr, " -A keep all possible alternate alleles at variant sites\n");
381 fprintf(stderr, " -b output BCF instead of VCF\n");
382 fprintf(stderr, " -D FILE sequence dictionary for VCF->BCF conversion [null]\n");
383 fprintf(stderr, " -F PL generated by r921 or before (which generate old ordering)\n");
384 fprintf(stderr, " -G suppress all individual genotype information\n");
385 fprintf(stderr, " -l FILE list of sites (chr pos) or regions (BED) to output [all sites]\n");
386 fprintf(stderr, " -L calculate LD for adjacent sites\n");
387 fprintf(stderr, " -N skip sites where REF is not A/C/G/T\n");
388 fprintf(stderr, " -Q output the QCALL likelihood format\n");
389 fprintf(stderr, " -s FILE list of samples to use [all samples]\n");
390 fprintf(stderr, " -S input is VCF\n");
391 fprintf(stderr, " -u uncompressed BCF output (force -b)\n");
392 fprintf(stderr, "\nConsensus/variant calling options:\n\n");
393 fprintf(stderr, " -c SNP calling (force -e)\n");
394 fprintf(stderr, " -d FLOAT skip loci where less than FLOAT fraction of samples covered [0]\n");
395 fprintf(stderr, " -e likelihood based analyses\n");
396 fprintf(stderr, " -g call genotypes at variant sites (force -c)\n");
397 fprintf(stderr, " -i FLOAT indel-to-substitution ratio [%.4g]\n", vc.indel_frac);
398 fprintf(stderr, " -I skip indels\n");
399 fprintf(stderr, " -p FLOAT variant if P(ref|D)<FLOAT [%.3g]\n", vc.pref);
400 fprintf(stderr, " -P STR type of prior: full, cond2, flat [full]\n");
401 fprintf(stderr, " -t FLOAT scaled substitution mutation rate [%.4g]\n", vc.theta);
402 fprintf(stderr, " -T STR constrained calling; STR can be: pair, trioauto, trioxd and trioxs (see manual) [null]\n");
403 fprintf(stderr, " -v output potential variant sites only (force -c)\n");
404 fprintf(stderr, "\nContrast calling and association test options:\n\n");
405 fprintf(stderr, " -1 INT number of group-1 samples [0]\n");
406 fprintf(stderr, " -C FLOAT posterior constrast for LRT<FLOAT and P(ref|D)<0.5 [%g]\n", vc.min_lrt);
407 fprintf(stderr, " -U INT number of permutations for association testing (effective with -1) [0]\n");
408 fprintf(stderr, " -X FLOAT only perform permutations for P(chi^2)<FLOAT [%g]\n", vc.min_perm_p);
409 fprintf(stderr, "\n");
413 if (vc.flag & VC_CALL) vc.flag |= VC_EM;
414 if ((vc.flag & VC_VCFIN) && (vc.flag & VC_BCFOUT) && vc.fn_dict == 0) {
415 fprintf(stderr, "[%s] For VCF->BCF conversion please specify the sequence dictionary with -D\n", __func__);
418 if (vc.n1 <= 0) vc.n_perm = 0; // TODO: give a warning here!
420 seeds = malloc(vc.n_perm * sizeof(int));
422 for (c = 0; c < vc.n_perm; ++c) seeds[c] = lrand48();
424 b = calloc(1, sizeof(bcf1_t));
425 blast = calloc(1, sizeof(bcf1_t));
427 if (!(vc.flag & VC_VCFIN)) strcat(moder, "b");
429 if (vc.flag & VC_BCFOUT) strcat(modew, "b");
430 if (vc.flag & VC_UNCOMP) strcat(modew, "u");
431 bp = vcf_open(argv[optind], moder);
432 hin = hout = vcf_hdr_read(bp);
433 if (vc.fn_dict && (vc.flag & VC_VCFIN))
434 vcf_dictread(bp, hin, vc.fn_dict);
435 bout = vcf_open("-", modew);
436 if (!(vc.flag & VC_QCALL)) {
438 vc.sublist = calloc(vc.n_sub, sizeof(int));
439 hout = bcf_hdr_subsam(hin, vc.n_sub, vc.subsam, vc.sublist);
441 if (vc.flag & VC_CALL) write_header(hout);
442 vcf_hdr_write(bout, hout);
444 if (vc.flag & VC_CALL) {
445 p1 = bcf_p1_init(hout->n_smpl, vc.ploidy);
447 if (bcf_p1_read_prior(p1, vc.prior_file) < 0) {
448 fprintf(stderr, "[%s] fail to read the prior AFS.\n", __func__);
451 } else bcf_p1_init_prior(p1, vc.prior_type, vc.theta);
452 if (vc.n1 > 0 && vc.min_lrt > 0.) { // set n1
453 bcf_p1_set_n1(p1, vc.n1);
454 bcf_p1_init_subprior(p1, vc.prior_type, vc.theta);
456 if (vc.indel_frac > 0.) bcf_p1_indel_prior(p1, vc.indel_frac); // otherwise use the default indel_frac
458 if (optind + 1 < argc && !(vc.flag&VC_VCFIN)) {
459 void *str2id = bcf_build_refhash(hout);
460 if (bcf_parse_region(str2id, argv[optind+1], &tid, &begin, &end) >= 0) {
462 idx = bcf_idx_load(argv[optind]);
465 off = bcf_idx_query(idx, tid, begin);
467 fprintf(stderr, "[%s] no records in the query region.\n", __func__);
468 return 1; // FIXME: a lot of memory leaks...
470 bgzf_seek(bp->fp, off, SEEK_SET);
471 bcf_idx_destroy(idx);
475 while (vcf_read(bp, hin, b) > 0) {
476 int is_indel, cons_llr = -1;
477 int64_t cons_gt = -1;
479 if ((vc.flag & VC_VARONLY) && strcmp(b->alt, "X") == 0) continue;
480 if ((vc.flag & VC_VARONLY) && vc.min_smpl_frac > 0.) {
481 extern int bcf_smpl_covered(const bcf1_t *b);
482 int n = bcf_smpl_covered(b);
483 if ((double)n / b->n_smpl < vc.min_smpl_frac) continue;
485 if (vc.n_sub) bcf_subsam(vc.n_sub, vc.sublist, b);
486 if (vc.flag & VC_FIX_PL) bcf_fix_pl(b);
487 is_indel = bcf_is_indel(b);
488 if ((vc.flag & VC_NO_INDEL) && is_indel) continue;
489 if ((vc.flag & VC_INDEL_ONLY) && !is_indel) continue;
490 if ((vc.flag & VC_ACGT_ONLY) && !is_indel) {
492 if (b->ref[0] == 0 || b->ref[1] != 0) continue;
493 x = toupper(b->ref[0]);
494 if (x != 'A' && x != 'C' && x != 'G' && x != 'T') continue;
496 if (vc.bed && !bed_overlap(vc.bed, hin->ns[b->tid], b->pos, b->pos + strlen(b->ref))) continue;
498 int l = strlen(b->ref);
499 l = b->pos + (l > 0? l : 1);
500 if (b->tid != tid || b->pos >= end) break;
501 if (!(l > begin && end > b->pos)) continue;
504 if ((vc.flag & VC_QCNT) && !is_indel) { // summarize the difference
505 int x = bcf_min_diff(b);
506 if (x > 255) x = 255;
507 if (x >= 0) ++qcnt[x];
509 if (vc.flag & VC_QCALL) { // output QCALL format; STOP here
513 if (vc.trio_aux) // do trio calling
514 bcf_trio_call(vc.trio_aux, b, &cons_llr, &cons_gt);
515 else if (vc.flag & VC_PAIRCALL)
516 cons_llr = bcf_pair_call(b);
517 if (vc.flag & (VC_CALL|VC_ADJLD|VC_EM)) bcf_gl2pl(b);
518 if (vc.flag & VC_EM) bcf_em1(b, vc.n1, 0x1ff, em);
521 for (i = 0; i < 9; ++i) em[i] = -1.;
523 if (vc.flag & VC_CALL) { // call variants
525 int calret = bcf_p1_cal(b, (em[7] >= 0 && em[7] < vc.min_lrt), p1, &pr);
526 if (n_processed % 100000 == 0) {
527 fprintf(stderr, "[%s] %ld sites processed.\n", __func__, (long)n_processed);
530 if (pr.p_ref >= vc.pref && (vc.flag & VC_VARONLY)) continue;
531 if (vc.n_perm && vc.n1 > 0 && pr.p_chi2 < vc.min_perm_p) { // permutation test
534 for (i = 0; i < vc.n_perm; ++i) {
535 #ifdef BCF_PERM_LRT // LRT based permutation is much faster but less robust to artifacts
537 bcf_shuffle(b, seeds[i]);
538 bcf_em1(b, vc.n1, 1<<7, x);
539 if (x[7] < em[7]) ++n;
541 bcf_shuffle(b, seeds[i]);
542 bcf_p1_cal(b, 1, p1, &r);
543 if (pr.p_chi2 >= r.p_chi2) ++n;
548 if (calret >= 0) update_bcf1(b, p1, &pr, vc.pref, vc.flag, em, cons_llr, cons_gt);
549 } else if (vc.flag & VC_EM) update_bcf1(b, 0, 0, 0, vc.flag, em, cons_llr, cons_gt);
550 if (vc.flag & VC_ADJLD) { // compute LD
552 if ((r2 = bcf_pair_freq(blast, b, f)) >= 0) {
554 s.m = s.l = 0; s.s = 0;
555 if (*b->info) kputc(';', &s);
556 ksprintf(&s, "NEIR=%.3f;NEIF4=%.3f,%.3f,%.3f,%.3f", r2, f[0], f[1], f[2], f[3]);
557 bcf_append_info(b, s.s, s.l);
562 if (vc.flag & VC_ANNO_MAX) bcf_anno_max(b);
563 if (vc.flag & VC_NO_GENO) { // do not output GENO fields
566 b->l_str = b->fmt - b->str + 1;
567 } else bcf_fix_gt(b);
568 vcf_write(bout, hout, b);
570 if (vc.prior_file) free(vc.prior_file);
571 if (vc.flag & VC_CALL) bcf_p1_dump_afs(p1);
572 if (hin != hout) bcf_hdr_destroy(hout);
573 bcf_hdr_destroy(hin);
574 bcf_destroy(b); bcf_destroy(blast);
575 vcf_close(bp); vcf_close(bout);
576 if (vc.fn_dict) free(vc.fn_dict);
577 if (vc.ploidy) free(vc.ploidy);
578 if (vc.trio_aux) free(vc.trio_aux);
581 for (i = 0; i < vc.n_sub; ++i) free(vc.subsam[i]);
582 free(vc.subsam); free(vc.sublist);
584 if (vc.bed) bed_destroy(vc.bed);
585 if (vc.flag & VC_QCNT)
586 for (c = 0; c < 256; ++c)
587 fprintf(stderr, "QT\t%d\t%lld\n", c, (long long)qcnt[c]);
588 if (seeds) free(seeds);
589 if (p1) bcf_p1_destroy(p1);