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