]> git.donarmstrong.com Git - rsem.git/blob - sam/bcftools/call1.c
Updated samtools to 0.1.19
[rsem.git] / sam / 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) 
194     {
195         // interpret as sample names, not as a file name
196         const char *t = fn, *p = t;
197         while (*t)
198         {
199             t++;
200             if ( *t==',' || !*t )
201             { 
202                 sam = realloc(sam, sizeof(void*)*(n+1));
203                 sam[n] = (char*) malloc(sizeof(char)*(t-p+2));
204                 memcpy(sam[n], p, t-p);
205                 sam[n][t-p]   = 0;
206                 sam[n][t-p+1] = 2;    // assume diploid
207                 p = t+1;
208                 n++; 
209             }
210         }
211         *_n = n;
212         return sam; // fail to open file
213     }
214         ks = ks_init(fp);
215         while (ks_getuntil(ks, 0, &s, &dret) >= 0) {
216                 int l;
217                 if (max == n) {
218                         max = max? max<<1 : 4;
219                         sam = realloc(sam, sizeof(void*)*max);
220                 }
221                 l = s.l;
222                 sam[n] = malloc(s.l + 2);
223                 strcpy(sam[n], s.s);
224                 sam[n][l+1] = 2; // by default, diploid
225                 if (dret != '\n') {
226                         if (ks_getuntil(ks, 0, &s, &dret) >= 0) { // read ploidy, 1 or 2
227                                 int x = (int)s.s[0] - '0';
228                                 if (x == 1 || x == 2) sam[n][l+1] = x;
229                                 else fprintf(stderr, "(%s) ploidy can only be 1 or 2; assume diploid\n", __func__);
230                         }
231                         if (dret != '\n') ks_getuntil(ks, '\n', &s, &dret);
232                 }
233                 ++n;
234         }
235         ks_destroy(ks);
236         gzclose(fp);
237         free(s.s);
238         *_n = n;
239         return sam;
240 }
241
242 static void write_header(bcf_hdr_t *h)
243 {
244         kstring_t str;
245         str.l = h->l_txt? h->l_txt - 1 : 0;
246         str.m = str.l + 1; str.s = h->txt;
247         if (!strstr(str.s, "##INFO=<ID=DP,"))
248                 kputs("##INFO=<ID=DP,Number=1,Type=Integer,Description=\"Raw read depth\">\n", &str);
249         if (!strstr(str.s, "##INFO=<ID=DP4,"))
250                 kputs("##INFO=<ID=DP4,Number=4,Type=Integer,Description=\"# high-quality ref-forward bases, ref-reverse, alt-forward and alt-reverse bases\">\n", &str);
251         if (!strstr(str.s, "##INFO=<ID=MQ,"))
252                 kputs("##INFO=<ID=MQ,Number=1,Type=Integer,Description=\"Root-mean-square mapping quality of covering reads\">\n", &str);
253         if (!strstr(str.s, "##INFO=<ID=FQ,"))
254                 kputs("##INFO=<ID=FQ,Number=1,Type=Float,Description=\"Phred probability of all samples being the same\">\n", &str);
255         if (!strstr(str.s, "##INFO=<ID=AF1,"))
256                 kputs("##INFO=<ID=AF1,Number=1,Type=Float,Description=\"Max-likelihood estimate of the first ALT allele frequency (assuming HWE)\">\n", &str);
257         if (!strstr(str.s, "##INFO=<ID=AC1,"))
258                 kputs("##INFO=<ID=AC1,Number=1,Type=Float,Description=\"Max-likelihood estimate of the first ALT allele count (no HWE assumption)\">\n", &str);
259         if (!strstr(str.s, "##INFO=<ID=AN,"))
260                 kputs("##INFO=<ID=AN,Number=1,Type=Integer,Description=\"Total number of alleles in called genotypes\">\n", &str);
261         if (!strstr(str.s, "##INFO=<ID=IS,"))
262                 kputs("##INFO=<ID=IS,Number=2,Type=Float,Description=\"Maximum number of reads supporting an indel and fraction of indel reads\">\n", &str);
263         if (!strstr(str.s, "##INFO=<ID=AC,"))
264                 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);
265         if (!strstr(str.s, "##INFO=<ID=G3,"))
266                 kputs("##INFO=<ID=G3,Number=3,Type=Float,Description=\"ML estimate of genotype frequencies\">\n", &str);
267         if (!strstr(str.s, "##INFO=<ID=HWE,"))
268                 kputs("##INFO=<ID=HWE,Number=1,Type=Float,Description=\"Chi^2 based HWE test P-value based on G3\">\n", &str);
269         if (!strstr(str.s, "##INFO=<ID=CLR,"))
270                 kputs("##INFO=<ID=CLR,Number=1,Type=Integer,Description=\"Log ratio of genotype likelihoods with and without the constraint\">\n", &str);
271         if (!strstr(str.s, "##INFO=<ID=UGT,"))
272                 kputs("##INFO=<ID=UGT,Number=1,Type=String,Description=\"The most probable unconstrained genotype configuration in the trio\">\n", &str);
273         if (!strstr(str.s, "##INFO=<ID=CGT,"))
274                 kputs("##INFO=<ID=CGT,Number=1,Type=String,Description=\"The most probable constrained genotype configuration in the trio\">\n", &str);
275 //      if (!strstr(str.s, "##INFO=<ID=CI95,"))
276 //              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);
277         if (!strstr(str.s, "##INFO=<ID=PV4,"))
278                 kputs("##INFO=<ID=PV4,Number=4,Type=Float,Description=\"P-values for strand bias, baseQ bias, mapQ bias and tail distance bias\">\n", &str);
279     if (!strstr(str.s, "##INFO=<ID=INDEL,"))
280         kputs("##INFO=<ID=INDEL,Number=0,Type=Flag,Description=\"Indicates that the variant is an INDEL.\">\n", &str);
281     if (!strstr(str.s, "##INFO=<ID=PC2,"))
282         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);
283     if (!strstr(str.s, "##INFO=<ID=PCHI2,"))
284         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);
285     if (!strstr(str.s, "##INFO=<ID=QCHI2,"))
286         kputs("##INFO=<ID=QCHI2,Number=1,Type=Integer,Description=\"Phred scaled PCHI2.\">\n", &str);
287     if (!strstr(str.s, "##INFO=<ID=RP,"))
288         kputs("##INFO=<ID=PR,Number=1,Type=Integer,Description=\"# permutations yielding a smaller PCHI2.\">\n", &str);
289     if (!strstr(str.s, "##INFO=<ID=QBD,"))
290         kputs("##INFO=<ID=QBD,Number=1,Type=Float,Description=\"Quality by Depth: QUAL/#reads\">\n", &str);
291     //if (!strstr(str.s, "##INFO=<ID=RPS,"))
292     //    kputs("##INFO=<ID=RPS,Number=3,Type=Float,Description=\"Read Position Stats: depth, average, stddev\">\n", &str);
293     if (!strstr(str.s, "##INFO=<ID=RPB,"))
294         kputs("##INFO=<ID=RPB,Number=1,Type=Float,Description=\"Read Position Bias\">\n", &str);
295     if (!strstr(str.s, "##INFO=<ID=MDV,"))
296         kputs("##INFO=<ID=MDV,Number=1,Type=Integer,Description=\"Maximum number of high-quality nonRef reads in samples\">\n", &str);
297     if (!strstr(str.s, "##INFO=<ID=VDB,"))
298         kputs("##INFO=<ID=VDB,Number=1,Type=Float,Description=\"Variant Distance Bias (v2) for filtering splice-site artefacts in RNA-seq data. Note: this version may be broken.\">\n", &str);
299     if (!strstr(str.s, "##FORMAT=<ID=GT,"))
300         kputs("##FORMAT=<ID=GT,Number=1,Type=String,Description=\"Genotype\">\n", &str);
301     if (!strstr(str.s, "##FORMAT=<ID=GQ,"))
302         kputs("##FORMAT=<ID=GQ,Number=1,Type=Integer,Description=\"Genotype Quality\">\n", &str);
303     if (!strstr(str.s, "##FORMAT=<ID=GL,"))
304         kputs("##FORMAT=<ID=GL,Number=3,Type=Float,Description=\"Likelihoods for RR,RA,AA genotypes (R=ref,A=alt)\">\n", &str);
305         if (!strstr(str.s, "##FORMAT=<ID=DP,"))
306                 kputs("##FORMAT=<ID=DP,Number=1,Type=Integer,Description=\"# high-quality bases\">\n", &str);
307         if (!strstr(str.s, "##FORMAT=<ID=DV,"))
308                 kputs("##FORMAT=<ID=DV,Number=1,Type=Integer,Description=\"# high-quality non-reference bases\">\n", &str);
309         if (!strstr(str.s, "##FORMAT=<ID=SP,"))
310                 kputs("##FORMAT=<ID=SP,Number=1,Type=Integer,Description=\"Phred-scaled strand bias P-value\">\n", &str);
311         if (!strstr(str.s, "##FORMAT=<ID=PL,"))
312                 kputs("##FORMAT=<ID=PL,Number=G,Type=Integer,Description=\"List of Phred-scaled genotype likelihoods\">\n", &str);
313         h->l_txt = str.l + 1; h->txt = str.s;
314 }
315
316 double bcf_pair_freq(const bcf1_t *b0, const bcf1_t *b1, double f[4]);
317
318 int bcfview(int argc, char *argv[])
319 {
320         extern int bcf_2qcall(bcf_hdr_t *h, bcf1_t *b);
321         extern void bcf_p1_indel_prior(bcf_p1aux_t *ma, double x);
322         extern int bcf_fix_gt(bcf1_t *b);
323         extern int bcf_anno_max(bcf1_t *b);
324         extern int bcf_shuffle(bcf1_t *b, int seed);
325         extern uint32_t *bcf_trio_prep(int is_x, int is_son);
326         extern int bcf_trio_call(uint32_t *prep, const bcf1_t *b, int *llr, int64_t *gt);
327         extern int bcf_pair_call(const bcf1_t *b);
328         extern int bcf_min_diff(const bcf1_t *b);
329         extern int bcf_p1_get_M(bcf_p1aux_t *b);
330
331         extern gzFile bcf_p1_fp_lk;
332
333         bcf_t *bp, *bout = 0;
334         bcf1_t *b, *blast;
335         int c, *seeds = 0;
336         uint64_t n_processed = 0, qcnt[256];
337         viewconf_t vc;
338         bcf_p1aux_t *p1 = 0;
339         bcf_hdr_t *hin, *hout;
340         int tid, begin, end;
341         char moder[4], modew[4];
342
343         tid = begin = end = -1;
344         memset(&vc, 0, sizeof(viewconf_t));
345         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;
346         memset(qcnt, 0, 8 * 256);
347         while ((c = getopt(argc, argv, "FN1:l:cC:eHAGvbSuP:t:p:QgLi:IMs:D:U:X:d:T:Ywm:K:")) >= 0) {
348                 switch (c) {
349                 case '1': vc.n1 = atoi(optarg); break;
350                 case 'l': vc.bed = bed_read(optarg); if (!vc.bed) { fprintf(stderr,"Could not read \"%s\"\n", optarg); return 1; } break;
351                 case 'D': vc.fn_dict = strdup(optarg); break;
352                 case 'F': vc.flag |= VC_FIX_PL; break;
353                 case 'N': vc.flag |= VC_ACGT_ONLY; break;
354                 case 'G': vc.flag |= VC_NO_GENO; break;
355                 case 'A': vc.flag |= VC_KEEPALT; break;
356                 case 'b': vc.flag |= VC_BCFOUT; break;
357                 case 'S': vc.flag |= VC_VCFIN; break;
358                 case 'c': vc.flag |= VC_CALL; break;
359                 case 'e': vc.flag |= VC_EM; break;
360                 case 'v': vc.flag |= VC_VARONLY | VC_CALL; break;
361                 case 'u': vc.flag |= VC_UNCOMP | VC_BCFOUT; break;
362                 case 'g': vc.flag |= VC_CALL_GT | VC_CALL; break;
363                 case 'I': vc.flag |= VC_NO_INDEL; break;
364                 case 'w': vc.flag |= VC_INDEL_ONLY; break;
365                 case 'M': vc.flag |= VC_ANNO_MAX; break;
366                 case 'Y': vc.flag |= VC_QCNT; break;
367         case 'm': vc.min_ma_lrt = atof(optarg); break;
368                 case 't': vc.theta = atof(optarg); break;
369                 case 'p': vc.pref = atof(optarg); break;
370                 case 'i': vc.indel_frac = atof(optarg); break;
371                 case 'Q': vc.flag |= VC_QCALL; break;
372                 case 'L': vc.flag |= VC_ADJLD; break;
373                 case 'U': vc.n_perm = atoi(optarg); break;
374                 case 'C': vc.min_lrt = atof(optarg); break;
375                 case 'X': vc.min_perm_p = atof(optarg); break;
376                 case 'd': vc.min_smpl_frac = atof(optarg); break;
377                 case 'K': bcf_p1_fp_lk = gzopen(optarg, "w"); break;
378                 case 's': vc.subsam = read_samples(optarg, &vc.n_sub);
379                         vc.ploidy = calloc(vc.n_sub + 1, 1);
380                         for (tid = 0; tid < vc.n_sub; ++tid) vc.ploidy[tid] = vc.subsam[tid][strlen(vc.subsam[tid]) + 1];
381                         tid = -1;
382                         break;
383                 case 'T':
384                         if (strcmp(optarg, "trioauto") == 0) vc.trio_aux = bcf_trio_prep(0, 0);
385                         else if (strcmp(optarg, "trioxd") == 0) vc.trio_aux = bcf_trio_prep(1, 0);
386                         else if (strcmp(optarg, "trioxs") == 0) vc.trio_aux = bcf_trio_prep(1, 1);
387                         else if (strcmp(optarg, "pair") == 0) vc.flag |= VC_PAIRCALL;
388                         else {
389                                 fprintf(stderr, "[%s] Option '-T' can only take value trioauto, trioxd or trioxs.\n", __func__);
390                                 return 1;
391                         }
392                         break;
393                 case 'P':
394                         if (strcmp(optarg, "full") == 0) vc.prior_type = MC_PTYPE_FULL;
395                         else if (strcmp(optarg, "cond2") == 0) vc.prior_type = MC_PTYPE_COND2;
396                         else if (strcmp(optarg, "flat") == 0) vc.prior_type = MC_PTYPE_FLAT;
397                         else vc.prior_file = strdup(optarg);
398                         break;
399                 }
400         }
401         if (argc == optind) {
402                 fprintf(stderr, "\n");
403                 fprintf(stderr, "Usage: bcftools view [options] <in.bcf> [reg]\n\n");
404                 fprintf(stderr, "Input/output options:\n\n");
405                 fprintf(stderr, "       -A        keep all possible alternate alleles at variant sites\n");
406                 fprintf(stderr, "       -b        output BCF instead of VCF\n");
407                 fprintf(stderr, "       -D FILE   sequence dictionary for VCF->BCF conversion [null]\n");
408                 fprintf(stderr, "       -F        PL generated by r921 or before (which generate old ordering)\n");
409                 fprintf(stderr, "       -G        suppress all individual genotype information\n");
410                 fprintf(stderr, "       -l FILE   list of sites (chr pos) or regions (BED) to output [all sites]\n");
411                 fprintf(stderr, "       -L        calculate LD for adjacent sites\n");
412                 fprintf(stderr, "       -N        skip sites where REF is not A/C/G/T\n");
413                 fprintf(stderr, "       -Q        output the QCALL likelihood format\n");
414                 fprintf(stderr, "       -s FILE   list of samples to use [all samples]\n");
415                 fprintf(stderr, "       -S        input is VCF\n");
416                 fprintf(stderr, "       -u        uncompressed BCF output (force -b)\n");
417                 fprintf(stderr, "\nConsensus/variant calling options:\n\n");
418                 fprintf(stderr, "       -c        SNP calling (force -e)\n");
419                 fprintf(stderr, "       -d FLOAT  skip loci where less than FLOAT fraction of samples covered [0]\n");
420                 fprintf(stderr, "       -e        likelihood based analyses\n");
421                 fprintf(stderr, "       -g        call genotypes at variant sites (force -c)\n");
422                 fprintf(stderr, "       -i FLOAT  indel-to-substitution ratio [%.4g]\n", vc.indel_frac);
423                 fprintf(stderr, "       -I        skip indels\n");
424                 fprintf(stderr, "       -m FLOAT  alternative model for multiallelic and rare-variant calling, include if P(chi^2)>=FLOAT\n");
425                 fprintf(stderr, "       -p FLOAT  variant if P(ref|D)<FLOAT [%.3g]\n", vc.pref);
426                 fprintf(stderr, "       -P STR    type of prior: full, cond2, flat [full]\n");
427                 fprintf(stderr, "       -t FLOAT  scaled substitution mutation rate [%.4g]\n", vc.theta);
428                 fprintf(stderr, "       -T STR    constrained calling; STR can be: pair, trioauto, trioxd and trioxs (see manual) [null]\n");
429                 fprintf(stderr, "       -v        output potential variant sites only (force -c)\n");
430                 fprintf(stderr, "\nContrast calling and association test options:\n\n");
431                 fprintf(stderr, "       -1 INT    number of group-1 samples [0]\n");
432                 fprintf(stderr, "       -C FLOAT  posterior constrast for LRT<FLOAT and P(ref|D)<0.5 [%g]\n", vc.min_lrt);
433                 fprintf(stderr, "       -U INT    number of permutations for association testing (effective with -1) [0]\n");
434                 fprintf(stderr, "       -X FLOAT  only perform permutations for P(chi^2)<FLOAT [%g]\n", vc.min_perm_p);
435                 fprintf(stderr, "\n");
436                 return 1;
437         }
438
439         if (vc.flag & VC_CALL) vc.flag |= VC_EM;
440         if ((vc.flag & VC_VCFIN) && (vc.flag & VC_BCFOUT) && vc.fn_dict == 0) {
441                 fprintf(stderr, "[%s] For VCF->BCF conversion please specify the sequence dictionary with -D\n", __func__);
442                 return 1;
443         }
444         if (vc.n1 <= 0) vc.n_perm = 0; // TODO: give a warning here!
445         if (vc.n_perm > 0) {
446                 seeds = malloc(vc.n_perm * sizeof(int));
447                 srand48(time(0));
448                 for (c = 0; c < vc.n_perm; ++c) seeds[c] = lrand48();
449         }
450         b = calloc(1, sizeof(bcf1_t));
451         blast = calloc(1, sizeof(bcf1_t));
452         strcpy(moder, "r");
453         if (!(vc.flag & VC_VCFIN)) strcat(moder, "b");
454         strcpy(modew, "w");
455         if (vc.flag & VC_BCFOUT) strcat(modew, "b");
456         if (vc.flag & VC_UNCOMP) strcat(modew, "u");
457         bp = vcf_open(argv[optind], moder);
458         hin = hout = vcf_hdr_read(bp);
459         if (vc.fn_dict && (vc.flag & VC_VCFIN))
460                 vcf_dictread(bp, hin, vc.fn_dict);
461         bout = vcf_open("-", modew);
462         if (!(vc.flag & VC_QCALL)) {
463                 if (vc.n_sub) {
464                         vc.sublist = calloc(vc.n_sub, sizeof(int));
465                         hout = bcf_hdr_subsam(hin, vc.n_sub, vc.subsam, vc.sublist);
466                 }
467                 write_header(hout); // always print the header
468                 vcf_hdr_write(bout, hout);
469         }
470         if (vc.flag & VC_CALL) {
471                 p1 = bcf_p1_init(hout->n_smpl, vc.ploidy);
472                 if (vc.prior_file) {
473                         if (bcf_p1_read_prior(p1, vc.prior_file) < 0) {
474                                 fprintf(stderr, "[%s] fail to read the prior AFS.\n", __func__);
475                                 return 1;
476                         }
477                 } else bcf_p1_init_prior(p1, vc.prior_type, vc.theta);
478                 if (vc.n1 > 0 && vc.min_lrt > 0.) { // set n1
479                         bcf_p1_set_n1(p1, vc.n1);
480                         bcf_p1_init_subprior(p1, vc.prior_type, vc.theta);
481                 }
482                 if (vc.indel_frac > 0.) bcf_p1_indel_prior(p1, vc.indel_frac); // otherwise use the default indel_frac
483         }
484         if (optind + 1 < argc && !(vc.flag&VC_VCFIN)) {
485                 void *str2id = bcf_build_refhash(hout);
486                 if (bcf_parse_region(str2id, argv[optind+1], &tid, &begin, &end) >= 0) {
487                         bcf_idx_t *idx;
488                         idx = bcf_idx_load(argv[optind]);
489                         if (idx) {
490                                 uint64_t off;
491                                 off = bcf_idx_query(idx, tid, begin);
492                                 if (off == 0) {
493                                         fprintf(stderr, "[%s] no records in the query region.\n", __func__);
494                                         return 1; // FIXME: a lot of memory leaks...
495                                 }
496                                 bgzf_seek(bp->fp, off, SEEK_SET);
497                                 bcf_idx_destroy(idx);
498                         }
499                 }
500         }
501         if (bcf_p1_fp_lk && p1) {
502                 int32_t M = bcf_p1_get_M(p1);
503                 gzwrite(bcf_p1_fp_lk, &M, 4);
504         }
505         while (vcf_read(bp, hin, b) > 0) {
506                 int is_indel, cons_llr = -1;
507                 int64_t cons_gt = -1;
508                 double em[10];
509                 if ((vc.flag & VC_VARONLY) && strcmp(b->alt, "X") == 0) continue;
510                 if ((vc.flag & VC_VARONLY) && vc.min_smpl_frac > 0.) {
511                         extern int bcf_smpl_covered(const bcf1_t *b);
512                         int n = bcf_smpl_covered(b);
513                         if ((double)n / b->n_smpl < vc.min_smpl_frac) continue;
514                 }
515                 if (vc.n_sub) bcf_subsam(vc.n_sub, vc.sublist, b);
516                 if (vc.flag & VC_FIX_PL) bcf_fix_pl(b);
517                 is_indel = bcf_is_indel(b);
518                 if ((vc.flag & VC_NO_INDEL) && is_indel) continue;
519                 if ((vc.flag & VC_INDEL_ONLY) && !is_indel) continue;
520                 if ((vc.flag & VC_ACGT_ONLY) && !is_indel) {
521                         int x;
522                         if (b->ref[0] == 0 || b->ref[1] != 0) continue;
523                         x = toupper(b->ref[0]);
524                         if (x != 'A' && x != 'C' && x != 'G' && x != 'T') continue;
525                 }
526                 if (vc.bed && !bed_overlap(vc.bed, hin->ns[b->tid], b->pos, b->pos + strlen(b->ref))) continue;
527                 if (tid >= 0) {
528                         int l = strlen(b->ref);
529                         l = b->pos + (l > 0? l : 1);
530                         if (b->tid != tid || b->pos >= end) break;
531                         if (!(l > begin && end > b->pos)) continue;
532                 }
533                 ++n_processed;
534                 if ((vc.flag & VC_QCNT) && !is_indel) { // summarize the difference
535                         int x = bcf_min_diff(b);
536                         if (x > 255) x = 255;
537                         if (x >= 0) ++qcnt[x];
538                 }
539                 if (vc.flag & VC_QCALL) { // output QCALL format; STOP here
540                         bcf_2qcall(hout, b);
541                         continue;
542                 }
543                 if (vc.trio_aux) // do trio calling
544                         bcf_trio_call(vc.trio_aux, b, &cons_llr, &cons_gt);
545                 else if (vc.flag & VC_PAIRCALL)
546                         cons_llr = bcf_pair_call(b);
547                 if (vc.flag & (VC_CALL|VC_ADJLD|VC_EM)) bcf_gl2pl(b);
548                 if (vc.flag & VC_EM) bcf_em1(b, vc.n1, 0x1ff, em);
549                 else {
550                         int i;
551                         for (i = 0; i < 9; ++i) em[i] = -1.;
552                 }
553         if ( !(vc.flag&VC_KEEPALT) && (vc.flag&VC_CALL) && vc.min_ma_lrt>=0 )
554         {
555             bcf_p1_set_ploidy(b, p1); // could be improved: do this per site to allow pseudo-autosomal regions
556             int gts = call_multiallelic_gt(b, p1, vc.min_ma_lrt, vc.flag&VC_VARONLY);
557             if ( gts<=1 && vc.flag & VC_VARONLY ) continue;
558         }
559                 else if (vc.flag & VC_CALL) { // call variants
560                         bcf_p1rst_t pr;
561                         int calret;
562                         gzwrite(bcf_p1_fp_lk, &b->tid, 4);
563                         gzwrite(bcf_p1_fp_lk, &b->pos, 4);
564                         gzwrite(bcf_p1_fp_lk, &em[0], sizeof(double));
565                         calret = bcf_p1_cal(b, (em[7] >= 0 && em[7] < vc.min_lrt), p1, &pr);
566                         if (n_processed % 100000 == 0) {
567                                 fprintf(stderr, "[%s] %ld sites processed.\n", __func__, (long)n_processed);
568                                 bcf_p1_dump_afs(p1);
569                         }
570                         if (pr.p_ref >= vc.pref && (vc.flag & VC_VARONLY)) continue;
571                         if (vc.n_perm && vc.n1 > 0 && pr.p_chi2 < vc.min_perm_p) { // permutation test
572                                 bcf_p1rst_t r;
573                                 int i, n = 0;
574                                 for (i = 0; i < vc.n_perm; ++i) {
575 #ifdef BCF_PERM_LRT // LRT based permutation is much faster but less robust to artifacts
576                                         double x[10];
577                                         bcf_shuffle(b, seeds[i]);
578                                         bcf_em1(b, vc.n1, 1<<7, x);
579                                         if (x[7] < em[7]) ++n;
580 #else
581                                         bcf_shuffle(b, seeds[i]);
582                                         bcf_p1_cal(b, 1, p1, &r);
583                                         if (pr.p_chi2 >= r.p_chi2) ++n;
584 #endif
585                                 }
586                                 pr.perm_rank = n;
587                         }
588                         if (calret >= 0) update_bcf1(b, p1, &pr, vc.pref, vc.flag, em, cons_llr, cons_gt);
589                 } else if (vc.flag & VC_EM) update_bcf1(b, 0, 0, 0, vc.flag, em, cons_llr, cons_gt);
590                 if (vc.flag & VC_ADJLD) { // compute LD
591                         double f[4], r2;
592                         if ((r2 = bcf_pair_freq(blast, b, f)) >= 0) {
593                                 kstring_t s;
594                                 s.m = s.l = 0; s.s = 0;
595                                 if (*b->info) kputc(';', &s);
596                                 ksprintf(&s, "NEIR=%.3f;NEIF4=%.3f,%.3f,%.3f,%.3f", r2, f[0], f[1], f[2], f[3]);
597                                 bcf_append_info(b, s.s, s.l);
598                                 free(s.s);
599                         }
600                         bcf_cpy(blast, b);
601                 }
602                 if (vc.flag & VC_ANNO_MAX) bcf_anno_max(b);
603                 if (vc.flag & VC_NO_GENO) { // do not output GENO fields
604                         b->n_gi = 0;
605                         b->fmt[0] = '\0';
606                         b->l_str = b->fmt - b->str + 1;
607                 } else bcf_fix_gt(b);
608                 vcf_write(bout, hout, b);
609         }
610
611         if (bcf_p1_fp_lk) gzclose(bcf_p1_fp_lk);
612         if (vc.prior_file) free(vc.prior_file);
613         if (vc.flag & VC_CALL) bcf_p1_dump_afs(p1);
614         if (hin != hout) bcf_hdr_destroy(hout);
615         bcf_hdr_destroy(hin);
616         bcf_destroy(b); bcf_destroy(blast);
617         vcf_close(bp); vcf_close(bout);
618         if (vc.fn_dict) free(vc.fn_dict);
619         if (vc.ploidy) free(vc.ploidy);
620         if (vc.trio_aux) free(vc.trio_aux);
621         if (vc.n_sub) {
622                 int i;
623                 for (i = 0; i < vc.n_sub; ++i) free(vc.subsam[i]);
624                 free(vc.subsam); free(vc.sublist);
625         }
626         if (vc.bed) bed_destroy(vc.bed);
627         if (vc.flag & VC_QCNT)
628                 for (c = 0; c < 256; ++c)
629                         fprintf(stderr, "QT\t%d\t%lld\n", c, (long long)qcnt[c]);
630         if (seeds) free(seeds);
631         if (p1) bcf_p1_destroy(p1);
632         return 0;
633 }