]> git.donarmstrong.com Git - samtools.git/blob - bcftools/call1.c
Fixed brackets missing in the previous commit
[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 #ifdef _WIN32
12 #define srand48(x) srand(x)
13 #define lrand48() rand()
14 #endif
15
16 #include "kseq.h"
17 KSTREAM_INIT(gzFile, gzread, 16384)
18
19 #define VC_NO_GENO 2
20 #define VC_BCFOUT  4
21 #define VC_CALL    8
22 #define VC_VARONLY 16
23 #define VC_VCFIN   32
24 #define VC_UNCOMP  64
25 #define VC_KEEPALT 256
26 #define VC_ACGT_ONLY 512
27 #define VC_QCALL   1024
28 #define VC_CALL_GT 2048
29 #define VC_ADJLD   4096
30 #define VC_NO_INDEL 8192
31 #define VC_ANNO_MAX 16384
32 #define VC_FIX_PL   32768
33 #define VC_EM       0x10000
34 #define VC_PAIRCALL 0x20000
35 #define VC_QCNT     0x40000
36 #define VC_INDEL_ONLY 0x80000
37
38 typedef struct {
39         int flag, prior_type, n1, n_sub, *sublist, n_perm;
40         uint32_t *trio_aux;
41         char *prior_file, **subsam, *fn_dict;
42         uint8_t *ploidy;
43         double theta, pref, indel_frac, min_perm_p, min_smpl_frac, min_lrt, min_ma_lrt;
44         void *bed;
45 } viewconf_t;
46
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);
50
51 static double ttest(int n1, int n2, int a[4])
52 {
53         extern double kf_betai(double a, double b, double x);
54         double t, v, u1, u2;
55         if (n1 == 0 || n2 == 0 || n1 + n2 < 3) return 1.0;
56         u1 = (double)a[0] / n1; u2 = (double)a[2] / n2;
57         if (u1 <= u2) return 1.;
58         t = (u1 - u2) / sqrt(((a[1] - n1 * u1 * u1) + (a[3] - n2 * u2 * u2)) / (n1 + n2 - 2) * (1./n1 + 1./n2));
59         v = n1 + n2 - 2;
60 //      printf("%d,%d,%d,%d,%lf,%lf,%lf\n", a[0], a[1], a[2], a[3], t, u1, u2);
61         return t < 0.? 1. : .5 * kf_betai(.5*v, .5, v/(v+t*t));
62 }
63
64 static int test16_core(int anno[16], anno16_t *a)
65 {
66         extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
67         double left, right;
68         int i;
69         a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
70         memcpy(a->d, anno, 4 * sizeof(int));
71         a->depth = anno[0] + anno[1] + anno[2] + anno[3];
72         a->is_tested = (anno[0] + anno[1] > 0 && anno[2] + anno[3] > 0);
73         if (a->depth == 0) return -1;
74         a->mq = (int)(sqrt((anno[9] + anno[11]) / a->depth) + .499);
75         kt_fisher_exact(anno[0], anno[1], anno[2], anno[3], &left, &right, &a->p[0]);
76         for (i = 1; i < 4; ++i)
77                 a->p[i] = ttest(anno[0] + anno[1], anno[2] + anno[3], anno+4*i);
78         return 0;
79 }
80
81 int test16(bcf1_t *b, anno16_t *a)
82 {
83         char *p;
84         int i, anno[16];
85         a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
86         a->d[0] = a->d[1] = a->d[2] = a->d[3] = 0.;
87         a->mq = a->depth = a->is_tested = 0;
88         if ((p = strstr(b->info, "I16=")) == 0) return -1;
89         p += 4;
90         for (i = 0; i < 16; ++i) {
91                 errno = 0; anno[i] = strtol(p, &p, 10);
92                 if (anno[i] == 0 && (errno == EINVAL || errno == ERANGE)) return -2;
93                 ++p;
94         }
95         return test16_core(anno, a);
96 }
97
98 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)
99 {
100         kstring_t s;
101         int has_I16, is_var;
102         double fq, r;
103         anno16_t a;
104
105         has_I16 = test16(b, &a) >= 0? 1 : 0;
106         //rm_info(b, "I16="); // FIXME: probably this function has a bug. If I move it below, I16 will not be removed!
107
108         memset(&s, 0, sizeof(kstring_t));
109         kputc('\0', &s); kputs(b->ref, &s); kputc('\0', &s);
110         kputs(b->alt, &s); kputc('\0', &s); kputc('\0', &s);
111         kputs(b->info, &s);
112         if (b->info[0]) kputc(';', &s);
113         { // print EM
114                 if (em[0] >= 0) ksprintf(&s, "AF1=%.4g", 1 - em[0]);
115                 if (em[4] >= 0 && em[4] <= 0.05) ksprintf(&s, ";G3=%.4g,%.4g,%.4g;HWE=%.3g", em[3], em[2], em[1], em[4]);
116                 if (em[5] >= 0 && em[6] >= 0) ksprintf(&s, ";AF2=%.4g,%.4g", 1 - em[5], 1 - em[6]);
117                 if (em[7] >= 0) ksprintf(&s, ";LRT=%.3g", em[7]);
118                 if (em[8] >= 0) ksprintf(&s, ";LRT2=%.3g", em[8]);
119         }
120         if (cons_llr > 0) {
121                 ksprintf(&s, ";CLR=%d", cons_llr);
122                 if (cons_gt > 0)
123                         ksprintf(&s, ";UGT=%c%c%c;CGT=%c%c%c", cons_gt&0xff, cons_gt>>8&0xff, cons_gt>>16&0xff,
124                                      cons_gt>>32&0xff, cons_gt>>40&0xff, cons_gt>>48&0xff);
125         }
126         if (pr == 0) { // if pr is unset, return
127                 kputc('\0', &s); kputs(b->fmt, &s); kputc('\0', &s);
128                 free(b->str);
129                 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
130                 bcf_sync(b);
131                 return 1;
132         }
133
134         is_var = (pr->p_ref < pref);
135         r = is_var? pr->p_ref : pr->p_var;
136
137 //      ksprintf(&s, ";CI95=%.4g,%.4g", pr->cil, pr->cih); // FIXME: when EM is not used, ";" should be omitted!
138         ksprintf(&s, ";AC1=%d", pr->ac);
139         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);
140         fq = pr->p_ref_folded < 0.5? -4.343 * log(pr->p_ref_folded) : 4.343 * log(pr->p_var_folded);
141         if (fq < -999) fq = -999;
142         if (fq > 999) fq = 999;
143         ksprintf(&s, ";FQ=%.3g", fq);
144         if (pr->cmp[0] >= 0.) { // two sample groups
145                 int i, q[3];
146                 for (i = 1; i < 3; ++i) {
147                         double x = pr->cmp[i] + pr->cmp[0]/2.;
148                         q[i] = x == 0? 255 : (int)(-4.343 * log(x) + .499);
149                         if (q[i] > 255) q[i] = 255;
150                 }
151                 if (pr->perm_rank >= 0) ksprintf(&s, ";PR=%d", pr->perm_rank);
152                 // ksprintf(&s, ";LRT3=%.3g", pr->lrt);
153                 ksprintf(&s, ";PCHI2=%.3g;PC2=%d,%d", q[1], q[2], pr->p_chi2);
154         }
155         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]);
156         kputc('\0', &s);
157     rm_info(&s, "QS=");
158     rm_info(&s, "I16=");
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=AN,"))
240                 kputs("##INFO=<ID=AN,Number=1,Type=Integer,Description=\"Total number of alleles in called genotypes\">\n", &str);
241         if (!strstr(str.s, "##INFO=<ID=IS,"))
242                 kputs("##INFO=<ID=IS,Number=2,Type=Float,Description=\"Maximum number of reads supporting an indel and fraction of indel reads\">\n", &str);
243         if (!strstr(str.s, "##INFO=<ID=AC,"))
244                 kputs("##INFO=<ID=AC,Number=A,Type=Integer,Description=\"Allele count in genotypes for each ALT allele, in the same order as listed\">\n", &str);
245         if (!strstr(str.s, "##INFO=<ID=G3,"))
246                 kputs("##INFO=<ID=G3,Number=3,Type=Float,Description=\"ML estimate of genotype frequencies\">\n", &str);
247         if (!strstr(str.s, "##INFO=<ID=HWE,"))
248                 kputs("##INFO=<ID=HWE,Number=1,Type=Float,Description=\"Chi^2 based HWE test P-value based on G3\">\n", &str);
249         if (!strstr(str.s, "##INFO=<ID=CLR,"))
250                 kputs("##INFO=<ID=CLR,Number=1,Type=Integer,Description=\"Log ratio of genotype likelihoods with and without the constraint\">\n", &str);
251         if (!strstr(str.s, "##INFO=<ID=UGT,"))
252                 kputs("##INFO=<ID=UGT,Number=1,Type=String,Description=\"The most probable unconstrained genotype configuration in the trio\">\n", &str);
253         if (!strstr(str.s, "##INFO=<ID=CGT,"))
254                 kputs("##INFO=<ID=CGT,Number=1,Type=String,Description=\"The most probable constrained genotype configuration in the trio\">\n", &str);
255 //      if (!strstr(str.s, "##INFO=<ID=CI95,"))
256 //              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);
257         if (!strstr(str.s, "##INFO=<ID=PV4,"))
258                 kputs("##INFO=<ID=PV4,Number=4,Type=Float,Description=\"P-values for strand bias, baseQ bias, mapQ bias and tail distance bias\">\n", &str);
259     if (!strstr(str.s, "##INFO=<ID=INDEL,"))
260         kputs("##INFO=<ID=INDEL,Number=0,Type=Flag,Description=\"Indicates that the variant is an INDEL.\">\n", &str);
261     if (!strstr(str.s, "##INFO=<ID=PC2,"))
262         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);
263     if (!strstr(str.s, "##INFO=<ID=PCHI2,"))
264         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);
265     if (!strstr(str.s, "##INFO=<ID=QCHI2,"))
266         kputs("##INFO=<ID=QCHI2,Number=1,Type=Integer,Description=\"Phred scaled PCHI2.\">\n", &str);
267     if (!strstr(str.s, "##INFO=<ID=RP,"))
268         kputs("##INFO=<ID=PR,Number=1,Type=Integer,Description=\"# permutations yielding a smaller PCHI2.\">\n", &str);
269     if (!strstr(str.s, "##INFO=<ID=VDB,"))
270         kputs("##INFO=<ID=VDB,Number=1,Type=Float,Description=\"Variant Distance Bias\">\n", &str);
271     if (!strstr(str.s, "##FORMAT=<ID=GT,"))
272         kputs("##FORMAT=<ID=GT,Number=1,Type=String,Description=\"Genotype\">\n", &str);
273     if (!strstr(str.s, "##FORMAT=<ID=GQ,"))
274         kputs("##FORMAT=<ID=GQ,Number=1,Type=Integer,Description=\"Genotype Quality\">\n", &str);
275     if (!strstr(str.s, "##FORMAT=<ID=GL,"))
276         kputs("##FORMAT=<ID=GL,Number=3,Type=Float,Description=\"Likelihoods for RR,RA,AA genotypes (R=ref,A=alt)\">\n", &str);
277         if (!strstr(str.s, "##FORMAT=<ID=DP,"))
278                 kputs("##FORMAT=<ID=DP,Number=1,Type=Integer,Description=\"# high-quality bases\">\n", &str);
279         if (!strstr(str.s, "##FORMAT=<ID=DV,"))
280                 kputs("##FORMAT=<ID=DV,Number=1,Type=Integer,Description=\"# high-quality non-reference bases\">\n", &str);
281         if (!strstr(str.s, "##FORMAT=<ID=SP,"))
282                 kputs("##FORMAT=<ID=SP,Number=1,Type=Integer,Description=\"Phred-scaled strand bias P-value\">\n", &str);
283         if (!strstr(str.s, "##FORMAT=<ID=PL,"))
284                 kputs("##FORMAT=<ID=PL,Number=G,Type=Integer,Description=\"List of Phred-scaled genotype likelihoods\">\n", &str);
285         h->l_txt = str.l + 1; h->txt = str.s;
286 }
287
288 double bcf_pair_freq(const bcf1_t *b0, const bcf1_t *b1, double f[4]);
289
290 int bcfview(int argc, char *argv[])
291 {
292         extern int bcf_2qcall(bcf_hdr_t *h, bcf1_t *b);
293         extern void bcf_p1_indel_prior(bcf_p1aux_t *ma, double x);
294         extern int bcf_fix_gt(bcf1_t *b);
295         extern int bcf_anno_max(bcf1_t *b);
296         extern int bcf_shuffle(bcf1_t *b, int seed);
297         extern uint32_t *bcf_trio_prep(int is_x, int is_son);
298         extern int bcf_trio_call(uint32_t *prep, const bcf1_t *b, int *llr, int64_t *gt);
299         extern int bcf_pair_call(const bcf1_t *b);
300         extern int bcf_min_diff(const bcf1_t *b);
301         extern int bcf_p1_get_M(bcf_p1aux_t *b);
302
303         extern gzFile bcf_p1_fp_lk;
304
305         bcf_t *bp, *bout = 0;
306         bcf1_t *b, *blast;
307         int c, *seeds = 0;
308         uint64_t n_processed = 0, qcnt[256];
309         viewconf_t vc;
310         bcf_p1aux_t *p1 = 0;
311         bcf_hdr_t *hin, *hout;
312         int tid, begin, end;
313         char moder[4], modew[4];
314
315         tid = begin = end = -1;
316         memset(&vc, 0, sizeof(viewconf_t));
317         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; vc.min_ma_lrt = -1;
318         memset(qcnt, 0, 8 * 256);
319         while ((c = getopt(argc, argv, "FN1:l:cC:eHAGvbSuP:t:p:QgLi:IMs:D:U:X:d:T:Ywm:K:")) >= 0) {
320                 switch (c) {
321                 case '1': vc.n1 = atoi(optarg); break;
322                 case 'l': vc.bed = bed_read(optarg); if (!vc.bed) { fprintf(stderr,"Could not read \"%s\"\n", optarg); return 1; } break;
323                 case 'D': vc.fn_dict = strdup(optarg); break;
324                 case 'F': vc.flag |= VC_FIX_PL; break;
325                 case 'N': vc.flag |= VC_ACGT_ONLY; break;
326                 case 'G': vc.flag |= VC_NO_GENO; break;
327                 case 'A': vc.flag |= VC_KEEPALT; break;
328                 case 'b': vc.flag |= VC_BCFOUT; break;
329                 case 'S': vc.flag |= VC_VCFIN; break;
330                 case 'c': vc.flag |= VC_CALL; break;
331                 case 'e': vc.flag |= VC_EM; break;
332                 case 'v': vc.flag |= VC_VARONLY | VC_CALL; break;
333                 case 'u': vc.flag |= VC_UNCOMP | VC_BCFOUT; break;
334                 case 'g': vc.flag |= VC_CALL_GT | VC_CALL; break;
335                 case 'I': vc.flag |= VC_NO_INDEL; break;
336                 case 'w': vc.flag |= VC_INDEL_ONLY; break;
337                 case 'M': vc.flag |= VC_ANNO_MAX; break;
338                 case 'Y': vc.flag |= VC_QCNT; break;
339         case 'm': vc.min_ma_lrt = atof(optarg); break;
340                 case 't': vc.theta = atof(optarg); break;
341                 case 'p': vc.pref = atof(optarg); break;
342                 case 'i': vc.indel_frac = atof(optarg); break;
343                 case 'Q': vc.flag |= VC_QCALL; break;
344                 case 'L': vc.flag |= VC_ADJLD; break;
345                 case 'U': vc.n_perm = atoi(optarg); break;
346                 case 'C': vc.min_lrt = atof(optarg); break;
347                 case 'X': vc.min_perm_p = atof(optarg); break;
348                 case 'd': vc.min_smpl_frac = atof(optarg); break;
349                 case 'K': bcf_p1_fp_lk = gzopen(optarg, "w"); break;
350                 case 's': vc.subsam = read_samples(optarg, &vc.n_sub);
351                         vc.ploidy = calloc(vc.n_sub + 1, 1);
352                         for (tid = 0; tid < vc.n_sub; ++tid) vc.ploidy[tid] = vc.subsam[tid][strlen(vc.subsam[tid]) + 1];
353                         tid = -1;
354                         break;
355                 case 'T':
356                         if (strcmp(optarg, "trioauto") == 0) vc.trio_aux = bcf_trio_prep(0, 0);
357                         else if (strcmp(optarg, "trioxd") == 0) vc.trio_aux = bcf_trio_prep(1, 0);
358                         else if (strcmp(optarg, "trioxs") == 0) vc.trio_aux = bcf_trio_prep(1, 1);
359                         else if (strcmp(optarg, "pair") == 0) vc.flag |= VC_PAIRCALL;
360                         else {
361                                 fprintf(stderr, "[%s] Option '-T' can only take value trioauto, trioxd or trioxs.\n", __func__);
362                                 return 1;
363                         }
364                         break;
365                 case 'P':
366                         if (strcmp(optarg, "full") == 0) vc.prior_type = MC_PTYPE_FULL;
367                         else if (strcmp(optarg, "cond2") == 0) vc.prior_type = MC_PTYPE_COND2;
368                         else if (strcmp(optarg, "flat") == 0) vc.prior_type = MC_PTYPE_FLAT;
369                         else vc.prior_file = strdup(optarg);
370                         break;
371                 }
372         }
373         if (argc == optind) {
374                 fprintf(stderr, "\n");
375                 fprintf(stderr, "Usage: bcftools view [options] <in.bcf> [reg]\n\n");
376                 fprintf(stderr, "Input/output options:\n\n");
377                 fprintf(stderr, "       -A        keep all possible alternate alleles at variant sites\n");
378                 fprintf(stderr, "       -b        output BCF instead of VCF\n");
379                 fprintf(stderr, "       -D FILE   sequence dictionary for VCF->BCF conversion [null]\n");
380                 fprintf(stderr, "       -F        PL generated by r921 or before (which generate old ordering)\n");
381                 fprintf(stderr, "       -G        suppress all individual genotype information\n");
382                 fprintf(stderr, "       -l FILE   list of sites (chr pos) or regions (BED) to output [all sites]\n");
383                 fprintf(stderr, "       -L        calculate LD for adjacent sites\n");
384                 fprintf(stderr, "       -N        skip sites where REF is not A/C/G/T\n");
385                 fprintf(stderr, "       -Q        output the QCALL likelihood format\n");
386                 fprintf(stderr, "       -s FILE   list of samples to use [all samples]\n");
387                 fprintf(stderr, "       -S        input is VCF\n");
388                 fprintf(stderr, "       -u        uncompressed BCF output (force -b)\n");
389                 fprintf(stderr, "\nConsensus/variant calling options:\n\n");
390                 fprintf(stderr, "       -c        SNP calling (force -e)\n");
391                 fprintf(stderr, "       -d FLOAT  skip loci where less than FLOAT fraction of samples covered [0]\n");
392                 fprintf(stderr, "       -e        likelihood based analyses\n");
393                 fprintf(stderr, "       -g        call genotypes at variant sites (force -c)\n");
394                 fprintf(stderr, "       -i FLOAT  indel-to-substitution ratio [%.4g]\n", vc.indel_frac);
395                 fprintf(stderr, "       -I        skip indels\n");
396                 fprintf(stderr, "       -m FLOAT  alternative model for multiallelic and rare-variant calling, include if P(chi^2)>=FLOAT\n");
397                 fprintf(stderr, "       -p FLOAT  variant if P(ref|D)<FLOAT [%.3g]\n", vc.pref);
398                 fprintf(stderr, "       -P STR    type of prior: full, cond2, flat [full]\n");
399                 fprintf(stderr, "       -t FLOAT  scaled substitution mutation rate [%.4g]\n", vc.theta);
400                 fprintf(stderr, "       -T STR    constrained calling; STR can be: pair, trioauto, trioxd and trioxs (see manual) [null]\n");
401                 fprintf(stderr, "       -v        output potential variant sites only (force -c)\n");
402                 fprintf(stderr, "\nContrast calling and association test options:\n\n");
403                 fprintf(stderr, "       -1 INT    number of group-1 samples [0]\n");
404                 fprintf(stderr, "       -C FLOAT  posterior constrast for LRT<FLOAT and P(ref|D)<0.5 [%g]\n", vc.min_lrt);
405                 fprintf(stderr, "       -U INT    number of permutations for association testing (effective with -1) [0]\n");
406                 fprintf(stderr, "       -X FLOAT  only perform permutations for P(chi^2)<FLOAT [%g]\n", vc.min_perm_p);
407                 fprintf(stderr, "\n");
408                 return 1;
409         }
410
411         if (vc.flag & VC_CALL) vc.flag |= VC_EM;
412         if ((vc.flag & VC_VCFIN) && (vc.flag & VC_BCFOUT) && vc.fn_dict == 0) {
413                 fprintf(stderr, "[%s] For VCF->BCF conversion please specify the sequence dictionary with -D\n", __func__);
414                 return 1;
415         }
416         if (vc.n1 <= 0) vc.n_perm = 0; // TODO: give a warning here!
417         if (vc.n_perm > 0) {
418                 seeds = malloc(vc.n_perm * sizeof(int));
419                 srand48(time(0));
420                 for (c = 0; c < vc.n_perm; ++c) seeds[c] = lrand48();
421         }
422         b = calloc(1, sizeof(bcf1_t));
423         blast = calloc(1, sizeof(bcf1_t));
424         strcpy(moder, "r");
425         if (!(vc.flag & VC_VCFIN)) strcat(moder, "b");
426         strcpy(modew, "w");
427         if (vc.flag & VC_BCFOUT) strcat(modew, "b");
428         if (vc.flag & VC_UNCOMP) strcat(modew, "u");
429         bp = vcf_open(argv[optind], moder);
430         hin = hout = vcf_hdr_read(bp);
431         if (vc.fn_dict && (vc.flag & VC_VCFIN))
432                 vcf_dictread(bp, hin, vc.fn_dict);
433         bout = vcf_open("-", modew);
434         if (!(vc.flag & VC_QCALL)) {
435                 if (vc.n_sub) {
436                         vc.sublist = calloc(vc.n_sub, sizeof(int));
437                         hout = bcf_hdr_subsam(hin, vc.n_sub, vc.subsam, vc.sublist);
438                 }
439                 write_header(hout); // always print the header
440                 vcf_hdr_write(bout, hout);
441         }
442         if (vc.flag & VC_CALL) {
443                 p1 = bcf_p1_init(hout->n_smpl, vc.ploidy);
444                 if (vc.prior_file) {
445                         if (bcf_p1_read_prior(p1, vc.prior_file) < 0) {
446                                 fprintf(stderr, "[%s] fail to read the prior AFS.\n", __func__);
447                                 return 1;
448                         }
449                 } else bcf_p1_init_prior(p1, vc.prior_type, vc.theta);
450                 if (vc.n1 > 0 && vc.min_lrt > 0.) { // set n1
451                         bcf_p1_set_n1(p1, vc.n1);
452                         bcf_p1_init_subprior(p1, vc.prior_type, vc.theta);
453                 }
454                 if (vc.indel_frac > 0.) bcf_p1_indel_prior(p1, vc.indel_frac); // otherwise use the default indel_frac
455         }
456         if (optind + 1 < argc && !(vc.flag&VC_VCFIN)) {
457                 void *str2id = bcf_build_refhash(hout);
458                 if (bcf_parse_region(str2id, argv[optind+1], &tid, &begin, &end) >= 0) {
459                         bcf_idx_t *idx;
460                         idx = bcf_idx_load(argv[optind]);
461                         if (idx) {
462                                 uint64_t off;
463                                 off = bcf_idx_query(idx, tid, begin);
464                                 if (off == 0) {
465                                         fprintf(stderr, "[%s] no records in the query region.\n", __func__);
466                                         return 1; // FIXME: a lot of memory leaks...
467                                 }
468                                 bgzf_seek(bp->fp, off, SEEK_SET);
469                                 bcf_idx_destroy(idx);
470                         }
471                 }
472         }
473         if (bcf_p1_fp_lk && p1) {
474                 int32_t M = bcf_p1_get_M(p1);
475                 gzwrite(bcf_p1_fp_lk, &M, 4);
476         }
477         while (vcf_read(bp, hin, b) > 0) {
478                 int is_indel, cons_llr = -1;
479                 int64_t cons_gt = -1;
480                 double em[10];
481                 if ((vc.flag & VC_VARONLY) && strcmp(b->alt, "X") == 0) continue;
482                 if ((vc.flag & VC_VARONLY) && vc.min_smpl_frac > 0.) {
483                         extern int bcf_smpl_covered(const bcf1_t *b);
484                         int n = bcf_smpl_covered(b);
485                         if ((double)n / b->n_smpl < vc.min_smpl_frac) continue;
486                 }
487                 if (vc.n_sub) bcf_subsam(vc.n_sub, vc.sublist, b);
488                 if (vc.flag & VC_FIX_PL) bcf_fix_pl(b);
489                 is_indel = bcf_is_indel(b);
490                 if ((vc.flag & VC_NO_INDEL) && is_indel) continue;
491                 if ((vc.flag & VC_INDEL_ONLY) && !is_indel) continue;
492                 if ((vc.flag & VC_ACGT_ONLY) && !is_indel) {
493                         int x;
494                         if (b->ref[0] == 0 || b->ref[1] != 0) continue;
495                         x = toupper(b->ref[0]);
496                         if (x != 'A' && x != 'C' && x != 'G' && x != 'T') continue;
497                 }
498                 if (vc.bed && !bed_overlap(vc.bed, hin->ns[b->tid], b->pos, b->pos + strlen(b->ref))) continue;
499                 if (tid >= 0) {
500                         int l = strlen(b->ref);
501                         l = b->pos + (l > 0? l : 1);
502                         if (b->tid != tid || b->pos >= end) break;
503                         if (!(l > begin && end > b->pos)) continue;
504                 }
505                 ++n_processed;
506                 if ((vc.flag & VC_QCNT) && !is_indel) { // summarize the difference
507                         int x = bcf_min_diff(b);
508                         if (x > 255) x = 255;
509                         if (x >= 0) ++qcnt[x];
510                 }
511                 if (vc.flag & VC_QCALL) { // output QCALL format; STOP here
512                         bcf_2qcall(hout, b);
513                         continue;
514                 }
515                 if (vc.trio_aux) // do trio calling
516                         bcf_trio_call(vc.trio_aux, b, &cons_llr, &cons_gt);
517                 else if (vc.flag & VC_PAIRCALL)
518                         cons_llr = bcf_pair_call(b);
519                 if (vc.flag & (VC_CALL|VC_ADJLD|VC_EM)) bcf_gl2pl(b);
520                 if (vc.flag & VC_EM) bcf_em1(b, vc.n1, 0x1ff, em);
521                 else {
522                         int i;
523                         for (i = 0; i < 9; ++i) em[i] = -1.;
524                 }
525         if ( !(vc.flag&VC_KEEPALT) && (vc.flag&VC_CALL) && vc.min_ma_lrt>=0 )
526         {
527             bcf_p1_set_ploidy(b, p1); // could be improved: do this per site to allow pseudo-autosomal regions
528             int gts = call_multiallelic_gt(b,p1,vc.min_ma_lrt);
529             if ( gts<=1 && vc.flag & VC_VARONLY ) continue;
530         }
531                 else if (vc.flag & VC_CALL) { // call variants
532                         bcf_p1rst_t pr;
533                         int calret;
534                         gzwrite(bcf_p1_fp_lk, &b->tid, 4);
535                         gzwrite(bcf_p1_fp_lk, &b->pos, 4);
536                         gzwrite(bcf_p1_fp_lk, &em[0], sizeof(double));
537                         calret = bcf_p1_cal(b, (em[7] >= 0 && em[7] < vc.min_lrt), p1, &pr);
538                         if (n_processed % 100000 == 0) {
539                                 fprintf(stderr, "[%s] %ld sites processed.\n", __func__, (long)n_processed);
540                                 bcf_p1_dump_afs(p1);
541                         }
542                         if (pr.p_ref >= vc.pref && (vc.flag & VC_VARONLY)) continue;
543                         if (vc.n_perm && vc.n1 > 0 && pr.p_chi2 < vc.min_perm_p) { // permutation test
544                                 bcf_p1rst_t r;
545                                 int i, n = 0;
546                                 for (i = 0; i < vc.n_perm; ++i) {
547 #ifdef BCF_PERM_LRT // LRT based permutation is much faster but less robust to artifacts
548                                         double x[10];
549                                         bcf_shuffle(b, seeds[i]);
550                                         bcf_em1(b, vc.n1, 1<<7, x);
551                                         if (x[7] < em[7]) ++n;
552 #else
553                                         bcf_shuffle(b, seeds[i]);
554                                         bcf_p1_cal(b, 1, p1, &r);
555                                         if (pr.p_chi2 >= r.p_chi2) ++n;
556 #endif
557                                 }
558                                 pr.perm_rank = n;
559                         }
560                         if (calret >= 0) update_bcf1(b, p1, &pr, vc.pref, vc.flag, em, cons_llr, cons_gt);
561                 } else if (vc.flag & VC_EM) update_bcf1(b, 0, 0, 0, vc.flag, em, cons_llr, cons_gt);
562                 if (vc.flag & VC_ADJLD) { // compute LD
563                         double f[4], r2;
564                         if ((r2 = bcf_pair_freq(blast, b, f)) >= 0) {
565                                 kstring_t s;
566                                 s.m = s.l = 0; s.s = 0;
567                                 if (*b->info) kputc(';', &s);
568                                 ksprintf(&s, "NEIR=%.3f;NEIF4=%.3f,%.3f,%.3f,%.3f", r2, f[0], f[1], f[2], f[3]);
569                                 bcf_append_info(b, s.s, s.l);
570                                 free(s.s);
571                         }
572                         bcf_cpy(blast, b);
573                 }
574                 if (vc.flag & VC_ANNO_MAX) bcf_anno_max(b);
575                 if (vc.flag & VC_NO_GENO) { // do not output GENO fields
576                         b->n_gi = 0;
577                         b->fmt[0] = '\0';
578                         b->l_str = b->fmt - b->str + 1;
579                 } else bcf_fix_gt(b);
580                 vcf_write(bout, hout, b);
581         }
582
583         if (bcf_p1_fp_lk) gzclose(bcf_p1_fp_lk);
584         if (vc.prior_file) free(vc.prior_file);
585         if (vc.flag & VC_CALL) bcf_p1_dump_afs(p1);
586         if (hin != hout) bcf_hdr_destroy(hout);
587         bcf_hdr_destroy(hin);
588         bcf_destroy(b); bcf_destroy(blast);
589         vcf_close(bp); vcf_close(bout);
590         if (vc.fn_dict) free(vc.fn_dict);
591         if (vc.ploidy) free(vc.ploidy);
592         if (vc.trio_aux) free(vc.trio_aux);
593         if (vc.n_sub) {
594                 int i;
595                 for (i = 0; i < vc.n_sub; ++i) free(vc.subsam[i]);
596                 free(vc.subsam); free(vc.sublist);
597         }
598         if (vc.bed) bed_destroy(vc.bed);
599         if (vc.flag & VC_QCNT)
600                 for (c = 0; c < 256; ++c)
601                         fprintf(stderr, "QT\t%d\t%lld\n", c, (long long)qcnt[c]);
602         if (seeds) free(seeds);
603         if (p1) bcf_p1_destroy(p1);
604         return 0;
605 }