11 KHASH_SET_INIT_INT64(set64)
14 KSTREAM_INIT(gzFile, gzread, 16384)
23 #define VC_KEEPALT 256
24 #define VC_ACGT_ONLY 512
26 #define VC_CALL_GT 2048
28 #define VC_NO_INDEL 8192
29 #define VC_ANNO_MAX 16384
32 int flag, prior_type, n1, n_sub, *sublist;
33 char *fn_list, *prior_file, **subsam, *fn_dict;
34 double theta, pref, indel_frac;
37 khash_t(set64) *bcf_load_pos(const char *fn, bcf_hdr_t *_h)
42 int ret, dret, lineno = 1;
44 khash_t(set64) *hash = 0;
46 hash = kh_init(set64);
47 str2id = bcf_build_refhash(_h);
48 str = calloc(1, sizeof(kstring_t));
49 fp = strcmp(fn, "-")? gzopen(fn, "r") : gzdopen(fileno(stdin), "r");
51 while (ks_getuntil(ks, 0, str, &dret) >= 0) {
52 int tid = bcf_str2id(str2id, str->s);
53 if (tid >= 0 && dret != '\n') {
54 if (ks_getuntil(ks, 0, str, &dret) >= 0) {
55 uint64_t x = (uint64_t)tid<<32 | (atoi(str->s) - 1);
56 kh_put(set64, hash, x, &ret);
58 } else fprintf(stderr, "[%s] %s is not a reference name (line %d).\n", __func__, str->s, lineno);
59 if (dret != '\n') while ((dret = ks_getc(ks)) > 0 && dret != '\n');
63 bcf_str2id_destroy(str2id);
66 free(str->s); free(str);
70 static double test_hwe(const double g[3])
72 extern double kf_gammaq(double p, double x);
73 double fexp, chi2, f[3], n;
75 n = g[0] + g[1] + g[2];
76 fexp = (2. * g[2] + g[1]) / (2. * n);
77 if (fexp > 1. - 1e-10) fexp = 1. - 1e-10;
78 if (fexp < 1e-10) fexp = 1e-10;
79 f[0] = n * (1. - fexp) * (1. - fexp);
80 f[1] = n * 2. * fexp * (1. - fexp);
81 f[2] = n * fexp * fexp;
82 for (i = 0, chi2 = 0.; i < 3; ++i)
83 chi2 += (g[i] - f[i]) * (g[i] - f[i]) / f[i];
84 return kf_gammaq(.5, chi2 / 2.);
89 int mq, depth, is_tested, d[4];
92 static double ttest(int n1, int n2, int a[4])
94 extern double kf_betai(double a, double b, double x);
96 if (n1 == 0 || n2 == 0 || n1 + n2 < 3) return 1.0;
97 u1 = (double)a[0] / n1; u2 = (double)a[2] / n2;
98 if (u1 <= u2) return 1.;
99 t = (u1 - u2) / sqrt(((a[1] - n1 * u1 * u1) + (a[3] - n2 * u2 * u2)) / (n1 + n2 - 2) * (1./n1 + 1./n2));
101 // printf("%d,%d,%d,%d,%lf,%lf,%lf\n", a[0], a[1], a[2], a[3], t, u1, u2);
102 return t < 0.? 1. : .5 * kf_betai(.5*v, .5, v/(v+t*t));
105 static int test16_core(int anno[16], anno16_t *a)
107 extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
110 a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
111 memcpy(a->d, anno, 4 * sizeof(int));
112 a->depth = anno[0] + anno[1] + anno[2] + anno[3];
113 a->is_tested = (anno[0] + anno[1] > 0 && anno[2] + anno[3] > 0);
114 if (a->depth == 0) return -1;
115 a->mq = (int)(sqrt((anno[9] + anno[11]) / a->depth) + .499);
116 kt_fisher_exact(anno[0], anno[1], anno[2], anno[3], &left, &right, &a->p[0]);
117 for (i = 1; i < 4; ++i)
118 a->p[i] = ttest(anno[0] + anno[1], anno[2] + anno[3], anno+4*i);
122 static int test16(bcf1_t *b, anno16_t *a)
126 a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
127 a->d[0] = a->d[1] = a->d[2] = a->d[3] = 0.;
128 a->mq = a->depth = a->is_tested = 0;
129 if ((p = strstr(b->info, "I16=")) == 0) return -1;
131 for (i = 0; i < 16; ++i) {
132 errno = 0; anno[i] = strtol(p, &p, 10);
133 if (anno[i] == 0 && (errno == EINVAL || errno == ERANGE)) return -2;
136 return test16_core(anno, a);
139 static void rm_info(bcf1_t *b, const char *key)
142 if ((p = strstr(b->info, key)) == 0) return;
143 for (q = p; *q && *q != ';'; ++q);
144 if (p > b->info && *(p-1) == ';') --p;
145 memmove(p, q, b->l_str - (q - b->str));
150 static int update_bcf1(int n_smpl, bcf1_t *b, const bcf_p1aux_t *pa, const bcf_p1rst_t *pr, double pref, int flag)
153 int is_var = (pr->p_ref < pref);
154 double p_hwe, r = is_var? pr->p_ref : pr->p_var, fq;
157 p_hwe = pr->g[0] >= 0.? test_hwe(pr->g) : 1.0; // only do HWE g[] is calculated
161 memset(&s, 0, sizeof(kstring_t));
162 kputc('\0', &s); kputs(b->ref, &s); kputc('\0', &s);
163 kputs(b->alt, &s); kputc('\0', &s); kputc('\0', &s);
165 if (b->info[0]) kputc(';', &s);
166 // ksprintf(&s, "AF1=%.4lg;AFE=%.4lg;CI95=%.4lg,%.4lg", 1.-pr->f_em, 1.-pr->f_exp, pr->cil, pr->cih);
167 ksprintf(&s, "AF1=%.4g;CI95=%.4g,%.4g", 1.-pr->f_em, pr->cil, pr->cih);
168 ksprintf(&s, ";DP4=%d,%d,%d,%d;MQ=%d", a.d[0], a.d[1], a.d[2], a.d[3], a.mq);
169 fq = pr->p_ref_folded < 0.5? -4.343 * log(pr->p_ref_folded) : 4.343 * log(pr->p_var_folded);
170 if (fq < -999) fq = -999;
171 if (fq > 999) fq = 999;
172 ksprintf(&s, ";FQ=%.3g", fq);
174 if (pr->pc[0] >= 0.) ksprintf(&s, ";PC4=%g,%g,%g,%g", pr->pc[0], pr->pc[1], pr->pc[2], pr->pc[3]);
175 ksprintf(&s, ";PV4=%.2g,%.2g,%.2g,%.2g", a.p[0], a.p[1], a.p[2], a.p[3]);
177 if (pr->g[0] >= 0. && p_hwe <= .2)
178 ksprintf(&s, ";GC=%.2f,%.2f,%.2f;HWE=%.3f", pr->g[2], pr->g[1], pr->g[0], p_hwe);
180 kputs(b->fmt, &s); kputc('\0', &s);
182 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
183 b->qual = r < 1e-100? 999 : -4.343 * log(r);
184 if (b->qual > 999) b->qual = 999;
186 if (!is_var) bcf_shrink_alt(b, 1);
187 else if (!(flag&VC_KEEPALT))
188 bcf_shrink_alt(b, pr->rank0 < 2? 2 : pr->rank0+1);
189 if (is_var && (flag&VC_CALL_GT)) { // call individual genotype
190 int i, x, old_n_gi = b->n_gi;
191 s.m = b->m_str; s.l = b->l_str - 1; s.s = b->str;
192 kputs(":GT:GQ", &s); kputc('\0', &s);
193 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
195 for (i = 0; i < b->n_smpl; ++i) {
196 x = bcf_p1_call_gt(pa, pr->f_em, i);
197 ((uint8_t*)b->gi[old_n_gi].data)[i] = (x&3) == 0? 1<<3|1 : (x&3) == 1? 1 : 0;
198 ((uint8_t*)b->gi[old_n_gi+1].data)[i] = x>>2;
204 static char **read_samples(const char *fn, int *_n)
209 int dret, n = 0, max = 0;
212 s.l = s.m = 0; s.s = 0;
213 fp = gzopen(fn, "r");
214 if (fp == 0) return 0; // fail to open file
216 while (ks_getuntil(ks, 0, &s, &dret) >= 0) {
218 max = max? max<<1 : 4;
219 sam = realloc(sam, sizeof(void*)*max);
221 sam[n++] = strdup(s.s);
230 static void write_header(bcf_hdr_t *h)
233 str.l = h->l_txt? h->l_txt - 1 : 0;
234 str.m = str.l + 1; str.s = h->txt;
235 if (!strstr(str.s, "##INFO=<ID=DP,"))
236 kputs("##INFO=<ID=DP,Number=1,Type=Integer,Description=\"Raw read depth\">\n", &str);
237 if (!strstr(str.s, "##INFO=<ID=DP4,"))
238 kputs("##INFO=<ID=DP4,Number=4,Type=Integer,Description=\"# high-quality ref-forward bases, ref-reverse, alt-forward and alt-reverse bases\">\n", &str);
239 if (!strstr(str.s, "##INFO=<ID=MQ,"))
240 kputs("##INFO=<ID=MQ,Number=1,Type=Integer,Description=\"Root-mean-square mapping quality of covering reads\">\n", &str);
241 if (!strstr(str.s, "##INFO=<ID=FQ,"))
242 kputs("##INFO=<ID=FQ,Number=1,Type=Float,Description=\"Phred probability that sample chromosomes are not all the same\">\n", &str);
243 if (!strstr(str.s, "##INFO=<ID=AF1,"))
244 kputs("##INFO=<ID=AF1,Number=1,Type=Float,Description=\"Max-likelihood estimate of the site allele frequency of the first ALT allele\">\n", &str);
245 if (!strstr(str.s, "##INFO=<ID=CI95,"))
246 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);
247 if (!strstr(str.s, "##INFO=<ID=PV4,"))
248 kputs("##INFO=<ID=PV4,Number=4,Type=Float,Description=\"P-values for strand bias, baseQ bias, mapQ bias and tail distance bias\">\n", &str);
249 if (!strstr(str.s, "##INFO=<ID=INDEL,"))
250 kputs("##INFO=<ID=INDEL,Number=0,Type=Flag,Description=\"Indicates that the variant is an INDEL.\">\n", &str);
251 if (!strstr(str.s, "##INFO=<ID=GT,"))
252 kputs("##FORMAT=<ID=GT,Number=1,Type=String,Description=\"Genotype\">\n", &str);
253 if (!strstr(str.s, "##FORMAT=<ID=GQ,"))
254 kputs("##FORMAT=<ID=GQ,Number=1,Type=Integer,Description=\"Genotype Quality\">\n", &str);
255 if (!strstr(str.s, "##FORMAT=<ID=GL,"))
256 kputs("##FORMAT=<ID=GL,Number=3,Type=Float,Description=\"Likelihoods for RR,RA,AA genotypes (R=ref,A=alt)\">\n", &str);
257 if (!strstr(str.s, "##FORMAT=<ID=DP,"))
258 kputs("##FORMAT=<ID=DP,Number=1,Type=Integer,Description=\"# high-quality bases\">\n", &str);
259 if (!strstr(str.s, "##FORMAT=<ID=SP,"))
260 kputs("##FORMAT=<ID=SP,Number=1,Type=Integer,Description=\"Phred-scaled strand bias P-value\">\n", &str);
261 if (!strstr(str.s, "##FORMAT=<ID=PL,"))
262 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);
263 h->l_txt = str.l + 1; h->txt = str.s;
266 double bcf_ld_freq(const bcf1_t *b0, const bcf1_t *b1, double f[4]);
268 int bcfview(int argc, char *argv[])
270 extern int bcf_2qcall(bcf_hdr_t *h, bcf1_t *b);
271 extern void bcf_p1_indel_prior(bcf_p1aux_t *ma, double x);
272 extern int bcf_fix_gt(bcf1_t *b);
273 extern int bcf_anno_max(bcf1_t *b);
274 bcf_t *bp, *bout = 0;
277 uint64_t n_processed = 0;
280 bcf_hdr_t *hin, *hout;
282 char moder[4], modew[4];
283 khash_t(set64) *hash = 0;
285 tid = begin = end = -1;
286 memset(&vc, 0, sizeof(viewconf_t));
287 vc.prior_type = vc.n1 = -1; vc.theta = 1e-3; vc.pref = 0.5; vc.indel_frac = -1.;
288 vc.n_sub = 0; vc.subsam = 0; vc.sublist = 0;
289 while ((c = getopt(argc, argv, "N1:l:cHAGvbSuP:t:p:QgLi:IMs:D:")) >= 0) {
291 case '1': vc.n1 = atoi(optarg); break;
292 case 'l': vc.fn_list = strdup(optarg); break;
293 case 'D': vc.fn_dict = strdup(optarg); break;
294 case 'N': vc.flag |= VC_ACGT_ONLY; break;
295 case 'G': vc.flag |= VC_NO_GENO; break;
296 case 'A': vc.flag |= VC_KEEPALT; break;
297 case 'b': vc.flag |= VC_BCFOUT; break;
298 case 'S': vc.flag |= VC_VCFIN; break;
299 case 'c': vc.flag |= VC_CALL; break;
300 case 'v': vc.flag |= VC_VARONLY | VC_CALL; break;
301 case 'u': vc.flag |= VC_UNCOMP | VC_BCFOUT; break;
302 case 'H': vc.flag |= VC_HWE; break;
303 case 'g': vc.flag |= VC_CALL_GT | VC_CALL; break;
304 case 'I': vc.flag |= VC_NO_INDEL; break;
305 case 'M': vc.flag |= VC_ANNO_MAX; break;
306 case 't': vc.theta = atof(optarg); break;
307 case 'p': vc.pref = atof(optarg); break;
308 case 'i': vc.indel_frac = atof(optarg); break;
309 case 'Q': vc.flag |= VC_QCALL; break;
310 case 'L': vc.flag |= VC_ADJLD; break;
311 case 's': vc.subsam = read_samples(optarg, &vc.n_sub); break;
313 if (strcmp(optarg, "full") == 0) vc.prior_type = MC_PTYPE_FULL;
314 else if (strcmp(optarg, "cond2") == 0) vc.prior_type = MC_PTYPE_COND2;
315 else if (strcmp(optarg, "flat") == 0) vc.prior_type = MC_PTYPE_FLAT;
316 else vc.prior_file = strdup(optarg);
320 if (argc == optind) {
321 fprintf(stderr, "\n");
322 fprintf(stderr, "Usage: bcftools view [options] <in.bcf> [reg]\n\n");
323 fprintf(stderr, "Options: -c SNP calling\n");
324 fprintf(stderr, " -v output potential variant sites only (force -c)\n");
325 fprintf(stderr, " -g call genotypes at variant sites (force -c)\n");
326 fprintf(stderr, " -b output BCF instead of VCF\n");
327 fprintf(stderr, " -u uncompressed BCF output (force -b)\n");
328 fprintf(stderr, " -S input is VCF\n");
329 fprintf(stderr, " -A keep all possible alternate alleles at variant sites\n");
330 fprintf(stderr, " -G suppress all individual genotype information\n");
331 fprintf(stderr, " -H perform Hardy-Weinberg test (slower)\n");
332 fprintf(stderr, " -N skip sites where REF is not A/C/G/T\n");
333 fprintf(stderr, " -Q output the QCALL likelihood format\n");
334 fprintf(stderr, " -L calculate LD for adjacent sites\n");
335 fprintf(stderr, " -I skip indels\n");
336 fprintf(stderr, " -D FILE sequence dictionary for VCF->BCF conversion [null]\n");
337 fprintf(stderr, " -1 INT number of group-1 samples [0]\n");
338 fprintf(stderr, " -l FILE list of sites to output [all sites]\n");
339 fprintf(stderr, " -t FLOAT scaled substitution mutation rate [%.4g]\n", vc.theta);
340 fprintf(stderr, " -i FLOAT indel-to-substitution ratio [%.4g]\n", vc.indel_frac);
341 fprintf(stderr, " -p FLOAT variant if P(ref|D)<FLOAT [%.3g]\n", vc.pref);
342 fprintf(stderr, " -P STR type of prior: full, cond2, flat [full]\n");
343 fprintf(stderr, " -s FILE list of samples to use [all samples]\n");
344 fprintf(stderr, "\n");
348 if ((vc.flag & VC_VCFIN) && (vc.flag & VC_BCFOUT) && vc.fn_dict == 0) {
349 fprintf(stderr, "[%s] For VCF->BCF conversion please specify the sequence dictionary with -D\n", __func__);
352 b = calloc(1, sizeof(bcf1_t));
353 blast = calloc(1, sizeof(bcf1_t));
355 if (!(vc.flag & VC_VCFIN)) strcat(moder, "b");
357 if (vc.flag & VC_BCFOUT) strcat(modew, "b");
358 if (vc.flag & VC_UNCOMP) strcat(modew, "u");
359 bp = vcf_open(argv[optind], moder);
360 hin = hout = vcf_hdr_read(bp);
361 if (vc.fn_dict && (vc.flag & VC_VCFIN))
362 vcf_dictread(bp, hin, vc.fn_dict);
363 bout = vcf_open("-", modew);
364 if (!(vc.flag & VC_QCALL)) {
366 vc.sublist = calloc(vc.n_sub, sizeof(int));
367 hout = bcf_hdr_subsam(hin, vc.n_sub, vc.subsam, vc.sublist);
369 if (vc.flag & VC_CALL) write_header(hout);
370 vcf_hdr_write(bout, hout);
372 if (vc.flag & VC_CALL) {
373 p1 = bcf_p1_init(hout->n_smpl);
375 if (bcf_p1_read_prior(p1, vc.prior_file) < 0) {
376 fprintf(stderr, "[%s] fail to read the prior AFS.\n", __func__);
379 } else bcf_p1_init_prior(p1, vc.prior_type, vc.theta);
381 bcf_p1_set_n1(p1, vc.n1);
382 bcf_p1_init_subprior(p1, vc.prior_type, vc.theta);
384 if (vc.indel_frac > 0.) bcf_p1_indel_prior(p1, vc.indel_frac); // otherwise use the default indel_frac
386 if (vc.fn_list) hash = bcf_load_pos(vc.fn_list, hin);
387 if (optind + 1 < argc && !(vc.flag&VC_VCFIN)) {
388 void *str2id = bcf_build_refhash(hout);
389 if (bcf_parse_region(str2id, argv[optind+1], &tid, &begin, &end) >= 0) {
391 idx = bcf_idx_load(argv[optind]);
394 off = bcf_idx_query(idx, tid, begin);
396 fprintf(stderr, "[%s] no records in the query region.\n", __func__);
397 return 1; // FIXME: a lot of memory leaks...
399 bgzf_seek(bp->fp, off, SEEK_SET);
400 bcf_idx_destroy(idx);
404 while (vcf_read(bp, hin, b) > 0) {
406 if (vc.n_sub) bcf_subsam(vc.n_sub, vc.sublist, b);
407 is_indel = bcf_is_indel(b);
408 if ((vc.flag & VC_NO_INDEL) && is_indel) continue;
409 if ((vc.flag & VC_ACGT_ONLY) && !is_indel) {
411 if (b->ref[0] == 0 || b->ref[1] != 0) continue;
412 x = toupper(b->ref[0]);
413 if (x != 'A' && x != 'C' && x != 'G' && x != 'T') continue;
416 uint64_t x = (uint64_t)b->tid<<32 | b->pos;
417 khint_t k = kh_get(set64, hash, x);
418 if (kh_size(hash) == 0) break;
419 if (k == kh_end(hash)) continue;
420 kh_del(set64, hash, k);
423 int l = strlen(b->ref);
424 l = b->pos + (l > 0? l : 1);
425 if (b->tid != tid || b->pos >= end) break;
426 if (!(l > begin && end > b->pos)) continue;
429 if (vc.flag & VC_QCALL) { // output QCALL format; STOP here
433 if (vc.flag & (VC_CALL|VC_ADJLD)) bcf_gl2pl(b);
434 if (vc.flag & VC_CALL) { // call variants
436 bcf_p1_cal(b, p1, &pr); // pr.g[3] is not calculated here
437 if (vc.flag&VC_HWE) bcf_p1_cal_g3(p1, pr.g);
438 if (n_processed % 100000 == 0) {
439 fprintf(stderr, "[%s] %ld sites processed.\n", __func__, (long)n_processed);
442 if (pr.p_ref >= vc.pref && (vc.flag & VC_VARONLY)) continue;
443 update_bcf1(hout->n_smpl, b, p1, &pr, vc.pref, vc.flag);
445 if (vc.flag & VC_ADJLD) { // compute LD
447 if ((r2 = bcf_ld_freq(blast, b, f)) >= 0) {
449 s.m = s.l = 0; s.s = 0;
450 if (*b->info) kputc(';', &s);
451 ksprintf(&s, "NEIR=%.3f;NEIF=%.3f,%.3f", r2, f[0]+f[2], f[0]+f[1]);
452 bcf_append_info(b, s.s, s.l);
457 if (vc.flag & VC_ANNO_MAX) bcf_anno_max(b);
458 if (vc.flag & VC_NO_GENO) { // do not output GENO fields
461 b->l_str = b->fmt - b->str + 1;
462 } else bcf_fix_gt(b);
463 vcf_write(bout, hout, b);
465 if (vc.prior_file) free(vc.prior_file);
466 if (vc.flag & VC_CALL) bcf_p1_dump_afs(p1);
467 if (hin != hout) bcf_hdr_destroy(hout);
468 bcf_hdr_destroy(hin);
469 bcf_destroy(b); bcf_destroy(blast);
470 vcf_close(bp); vcf_close(bout);
471 if (hash) kh_destroy(set64, hash);
472 if (vc.fn_list) free(vc.fn_list);
473 if (vc.fn_dict) free(vc.fn_dict);
476 for (i = 0; i < vc.n_sub; ++i) free(vc.subsam[i]);
477 free(vc.subsam); free(vc.sublist);
479 if (p1) bcf_p1_destroy(p1);