2b284527185c0e12747ca1dfb5387aa0bd615c36
[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
10 #include "khash.h"
11 KHASH_SET_INIT_INT64(set64)
12
13 #include "kseq.h"
14 KSTREAM_INIT(gzFile, gzread, 16384)
15
16 #define VC_NO_GENO 2
17 #define VC_BCFOUT  4
18 #define VC_CALL    8
19 #define VC_VARONLY 16
20 #define VC_VCFIN   32
21 #define VC_UNCOMP  64
22 #define VC_HWE     128
23 #define VC_KEEPALT 256
24 #define VC_ACGT_ONLY 512
25 #define VC_QCALL   1024
26 #define VC_CALL_GT 2048
27 #define VC_ADJLD   4096
28
29 typedef struct {
30         int flag, prior_type, n1;
31         char *fn_list, *prior_file;
32         double theta, pref;
33 } viewconf_t;
34
35 khash_t(set64) *bcf_load_pos(const char *fn, bcf_hdr_t *_h)
36 {
37         void *str2id;
38         gzFile fp;
39         kstream_t *ks;
40         int ret, dret, lineno = 1;
41         kstring_t *str;
42         khash_t(set64) *hash = 0;
43
44         hash = kh_init(set64);
45         str2id = bcf_build_refhash(_h);
46         str = calloc(1, sizeof(kstring_t));
47         fp = strcmp(fn, "-")? gzopen(fn, "r") : gzdopen(fileno(stdin), "r");
48         ks = ks_init(fp);
49         while (ks_getuntil(ks, 0, str, &dret) >= 0) {
50                 int tid = bcf_str2id(str2id, str->s);
51                 if (tid >= 0 && dret != '\n') {
52                         if (ks_getuntil(ks, 0, str, &dret) >= 0) {
53                                 uint64_t x = (uint64_t)tid<<32 | (atoi(str->s) - 1);
54                                 kh_put(set64, hash, x, &ret);
55                         } else break;
56                 } else fprintf(stderr, "[%s] %s is not a reference name (line %d).\n", __func__, str->s, lineno);
57                 if (dret != '\n') while ((dret = ks_getc(ks)) > 0 && dret != '\n');
58                 if (dret < 0) break;
59                 ++lineno;
60         }
61         bcf_str2id_destroy(str2id);
62         ks_destroy(ks);
63         gzclose(fp);
64         free(str->s); free(str);
65         return hash;
66 }
67
68 static double test_hwe(const double g[3])
69 {
70         extern double kf_gammaq(double p, double x);
71         double fexp, chi2, f[3], n;
72         int i;
73         n = g[0] + g[1] + g[2];
74         fexp = (2. * g[2] + g[1]) / (2. * n);
75         if (fexp > 1. - 1e-10) fexp = 1. - 1e-10;
76         if (fexp < 1e-10) fexp = 1e-10;
77         f[0] = n * (1. - fexp) * (1. - fexp);
78         f[1] = n * 2. * fexp * (1. - fexp);
79         f[2] = n * fexp * fexp;
80         for (i = 0, chi2 = 0.; i < 3; ++i)
81                 chi2 += (g[i] - f[i]) * (g[i] - f[i]) / f[i];
82         return kf_gammaq(.5, chi2 / 2.);
83 }
84
85 typedef struct {
86         double p[4];
87         int mq, depth, is_tested, d[4];
88 } anno16_t;
89
90 static double ttest(int n1, int n2, int a[4])
91 {
92         extern double kf_betai(double a, double b, double x);
93         double t, v, u1, u2;
94         if (n1 == 0 || n2 == 0 || n1 + n2 < 3) return 1.0;
95         u1 = (double)a[0] / n1; u2 = (double)a[2] / n2;
96         if (u1 <= u2) return 1.;
97         t = (u1 - u2) / sqrt(((a[1] - n1 * u1 * u1) + (a[3] - n2 * u2 * u2)) / (n1 + n2 - 2) * (1./n1 + 1./n2));
98         v = n1 + n2 - 2;
99 //      printf("%d,%d,%d,%d,%lf,%lf,%lf\n", a[0], a[1], a[2], a[3], t, u1, u2);
100         return t < 0.? 1. : .5 * kf_betai(.5*v, .5, v/(v+t*t));
101 }
102
103 static int test16_core(int anno[16], anno16_t *a)
104 {
105         extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
106         double left, right;
107         int i;
108         a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
109         memcpy(a->d, anno, 4 * sizeof(int));
110         a->depth = anno[0] + anno[1] + anno[2] + anno[3];
111         a->is_tested = (anno[0] + anno[1] > 0 && anno[2] + anno[3] > 0);
112         if (a->depth == 0) return -1;
113         a->mq = (int)(sqrt((anno[9] + anno[11]) / a->depth) + .499);
114         kt_fisher_exact(anno[0], anno[1], anno[2], anno[3], &left, &right, &a->p[0]);
115         for (i = 1; i < 4; ++i)
116                 a->p[i] = ttest(anno[0] + anno[1], anno[2] + anno[3], anno+4*i);
117         return 0;
118 }
119
120 static int test16(bcf1_t *b, anno16_t *a)
121 {
122         char *p;
123         int i, anno[16];
124         a->p[0] = a->p[1] = a->p[2] = a->p[3] = 1.;
125         a->d[0] = a->d[1] = a->d[2] = a->d[3] = 0.;
126         a->mq = a->depth = a->is_tested = 0;
127         if ((p = strstr(b->info, "I16=")) == 0) return -1;
128         p += 4;
129         for (i = 0; i < 16; ++i) {
130                 anno[i] = strtol(p, &p, 10);
131                 if (anno[i] == 0 && (errno == EINVAL || errno == ERANGE)) return -2;
132                 ++p;
133         }
134         return test16_core(anno, a);
135 }
136
137 static void rm_info(bcf1_t *b, const char *key)
138 {
139         char *p, *q;
140         if ((p = strstr(b->info, key)) == 0) return;
141         for (q = p; *q && *q != ';'; ++q);
142         if (p > b->info && *(p-1) == ';') --p;
143         memmove(p, q, b->l_str - (q - b->str));
144         b->l_str -= q - p;
145         bcf_sync(b);
146 }
147
148 static int update_bcf1(int n_smpl, bcf1_t *b, const bcf_p1aux_t *pa, const bcf_p1rst_t *pr, double pref, int flag)
149 {
150         kstring_t s;
151         int is_var = (pr->p_ref < pref);
152         double p_hwe, r = is_var? pr->p_ref : 1. - pr->p_ref;
153         anno16_t a;
154
155         p_hwe = pr->g[0] >= 0.? test_hwe(pr->g) : 1.0; // only do HWE g[] is calculated
156         test16(b, &a);
157         rm_info(b, "I16=");
158
159         memset(&s, 0, sizeof(kstring_t));
160         kputc('\0', &s); kputs(b->ref, &s); kputc('\0', &s);
161         kputs(b->alt, &s); kputc('\0', &s); kputc('\0', &s);
162         kputs(b->info, &s);
163         if (b->info[0]) kputc(';', &s);
164         ksprintf(&s, "AF1=%.3lf;AFE=%.3lf", 1.-pr->f_em, 1.-pr->f_exp);
165         ksprintf(&s, ";DP4=%d,%d,%d,%d;MQ=%d", a.d[0], a.d[1], a.d[2], a.d[3], a.mq);
166         if (a.is_tested) {
167                 if (pr->pc[0] >= 0.) ksprintf(&s, ";PC4=%lg,%lg,%lg,%lg", pr->pc[0], pr->pc[1], pr->pc[2], pr->pc[3]);
168                 ksprintf(&s, ";PV4=%.2lg,%.2lg,%.2lg,%.2lg", a.p[0], a.p[1], a.p[2], a.p[3]);
169         }
170         if (pr->g[0] >= 0. && p_hwe <= .2)
171                 ksprintf(&s, ";GC=%.2lf,%.2lf,%.2lf;HWE=%.3lf", pr->g[2], pr->g[1], pr->g[0], p_hwe);
172         kputc('\0', &s);
173         kputs(b->fmt, &s); kputc('\0', &s);
174         free(b->str);
175         b->m_str = s.m; b->l_str = s.l; b->str = s.s;
176         b->qual = r < 1e-100? 99 : -4.343 * log(r);
177         if (b->qual > 99) b->qual = 99;
178         bcf_sync(b);
179         if (!is_var) bcf_shrink_alt(b, 1);
180         else if (!(flag&VC_KEEPALT))
181                 bcf_shrink_alt(b, pr->rank0 < 2? 2 : pr->rank0+1);
182         if (is_var && (flag&VC_CALL_GT)) { // call individual genotype
183                 int i, x, old_n_gi = b->n_gi;
184                 s.m = b->m_str; s.l = b->l_str - 1; s.s = b->str;
185                 kputs(":GT:GQ", &s); kputc('\0', &s);
186                 b->m_str = s.m; b->l_str = s.l; b->str = s.s;
187                 bcf_sync(b);
188                 for (i = 0; i < b->n_smpl; ++i) {
189                         x = bcf_p1_call_gt(pa, pr->f_em, i);
190                         ((uint8_t*)b->gi[old_n_gi].data)[i] = (x&3) == 0? 1<<3|1 : (x&3) == 1? 1 : 0;
191                         ((uint8_t*)b->gi[old_n_gi+1].data)[i] = x>>2;
192                 }
193         }
194         return is_var;
195 }
196
197 double bcf_ld_freq(const bcf1_t *b0, const bcf1_t *b1, double f[4]);
198
199 int bcfview(int argc, char *argv[])
200 {
201         extern int bcf_2qcall(bcf_hdr_t *h, bcf1_t *b);
202         bcf_t *bp, *bout = 0;
203         bcf1_t *b, *blast;
204         int c;
205         uint64_t n_processed = 0;
206         viewconf_t vc;
207         bcf_p1aux_t *p1 = 0;
208         bcf_hdr_t *h;
209         int tid, begin, end;
210         char moder[4], modew[4];
211         khash_t(set64) *hash = 0;
212
213         tid = begin = end = -1;
214         memset(&vc, 0, sizeof(viewconf_t));
215         vc.prior_type = vc.n1 = -1; vc.theta = 1e-3; vc.pref = 0.5;
216         while ((c = getopt(argc, argv, "N1:l:cHAGvbSuP:t:p:QgL")) >= 0) {
217                 switch (c) {
218                 case '1': vc.n1 = atoi(optarg); break;
219                 case 'l': vc.fn_list = strdup(optarg); break;
220                 case 'N': vc.flag |= VC_ACGT_ONLY; break;
221                 case 'G': vc.flag |= VC_NO_GENO; break;
222                 case 'A': vc.flag |= VC_KEEPALT; break;
223                 case 'b': vc.flag |= VC_BCFOUT; break;
224                 case 'S': vc.flag |= VC_VCFIN; break;
225                 case 'c': vc.flag |= VC_CALL; break;
226                 case 'v': vc.flag |= VC_VARONLY | VC_CALL; break;
227                 case 'u': vc.flag |= VC_UNCOMP | VC_BCFOUT; break;
228                 case 'H': vc.flag |= VC_HWE; break;
229                 case 'g': vc.flag |= VC_CALL_GT | VC_CALL; break;
230                 case 't': vc.theta = atof(optarg); break;
231                 case 'p': vc.pref = atof(optarg); break;
232                 case 'Q': vc.flag |= VC_QCALL; break;
233                 case 'L': vc.flag |= VC_ADJLD; break;
234                 case 'P':
235                         if (strcmp(optarg, "full") == 0) vc.prior_type = MC_PTYPE_FULL;
236                         else if (strcmp(optarg, "cond2") == 0) vc.prior_type = MC_PTYPE_COND2;
237                         else if (strcmp(optarg, "flat") == 0) vc.prior_type = MC_PTYPE_FLAT;
238                         else vc.prior_file = strdup(optarg);
239                         break;
240                 }
241         }
242         if (argc == optind) {
243                 fprintf(stderr, "\n");
244                 fprintf(stderr, "Usage:   bcftools view [options] <in.bcf> [reg]\n\n");
245                 fprintf(stderr, "Options: -c        SNP calling\n");
246                 fprintf(stderr, "         -v        output potential variant sites only (force -c)\n");
247                 fprintf(stderr, "         -g        call genotypes at variant sites (force -c)\n");
248                 fprintf(stderr, "         -b        output BCF instead of VCF\n");
249                 fprintf(stderr, "         -u        uncompressed BCF output (force -b)\n");
250                 fprintf(stderr, "         -S        input is VCF\n");
251                 fprintf(stderr, "         -A        keep all possible alternate alleles at variant sites\n");
252                 fprintf(stderr, "         -G        suppress all individual genotype information\n");
253                 fprintf(stderr, "         -H        perform Hardy-Weinberg test (slower)\n");
254                 fprintf(stderr, "         -N        skip sites where REF is not A/C/G/T\n");
255                 fprintf(stderr, "         -Q        output the QCALL likelihood format\n");
256                 fprintf(stderr, "         -L        calculate LD for adjacent sites\n");
257                 fprintf(stderr, "         -1 INT    number of group-1 samples [0]\n");
258                 fprintf(stderr, "         -l FILE   list of sites to output [all sites]\n");
259                 fprintf(stderr, "         -t FLOAT  scaled mutation rate [%.4lg]\n", vc.theta);
260                 fprintf(stderr, "         -p FLOAT  variant if P(ref|D)<FLOAT [%.3lg]\n", vc.pref);
261                 fprintf(stderr, "         -P STR    type of prior: full, cond2, flat [full]\n");
262                 fprintf(stderr, "\n");
263                 return 1;
264         }
265
266         b = calloc(1, sizeof(bcf1_t));
267         blast = calloc(1, sizeof(bcf1_t));
268         strcpy(moder, "r");
269         if (!(vc.flag & VC_VCFIN)) strcat(moder, "b");
270         strcpy(modew, "w");
271         if (vc.flag & VC_BCFOUT) strcat(modew, "b");
272         if (vc.flag & VC_UNCOMP) strcat(modew, "u");
273         bp = vcf_open(argv[optind], moder);
274         h = vcf_hdr_read(bp);
275         bout = vcf_open("-", modew);
276         if (!(vc.flag & VC_QCALL)) vcf_hdr_write(bout, h);
277         if (vc.flag & VC_CALL) {
278                 p1 = bcf_p1_init(h->n_smpl);
279                 if (vc.prior_file) {
280                         if (bcf_p1_read_prior(p1, vc.prior_file) < 0) {
281                                 fprintf(stderr, "[%s] fail to read the prior AFS.\n", __func__);
282                                 return 1;
283                         }
284                 } else bcf_p1_init_prior(p1, vc.prior_type, vc.theta);
285                 if (vc.n1 > 0) {
286                         bcf_p1_set_n1(p1, vc.n1);
287                         bcf_p1_init_subprior(p1, vc.prior_type, vc.theta);
288                 }
289         }
290         if (vc.fn_list) hash = bcf_load_pos(vc.fn_list, h);
291         if (optind + 1 < argc && !(vc.flag&VC_VCFIN)) {
292                 void *str2id = bcf_build_refhash(h);
293                 if (bcf_parse_region(str2id, argv[optind+1], &tid, &begin, &end) >= 0) {
294                         bcf_idx_t *idx;
295                         idx = bcf_idx_load(argv[optind]);
296                         if (idx) {
297                                 uint64_t off;
298                                 off = bcf_idx_query(idx, tid, begin);
299                                 if (off == 0) {
300                                         fprintf(stderr, "[%s] no records in the query region.\n", __func__);
301                                         return 1; // FIXME: a lot of memory leaks...
302                                 }
303                                 bgzf_seek(bp->fp, off, SEEK_SET);
304                                 bcf_idx_destroy(idx);
305                         }
306                 }
307         }
308         while (vcf_read(bp, h, b) > 0) {
309                 if (vc.flag & VC_ACGT_ONLY) {
310                         int x;
311                         if (b->ref[0] == 0 || b->ref[1] != 0) continue;
312                         x = toupper(b->ref[0]);
313                         if (x != 'A' && x != 'C' && x != 'G' && x != 'T') continue;
314                 }
315                 if (hash) {
316                         uint64_t x = (uint64_t)b->tid<<32 | b->pos;
317                         khint_t k = kh_get(set64, hash, x);
318                         if (kh_size(hash) == 0) break;
319                         if (k == kh_end(hash)) continue;
320                         kh_del(set64, hash, k);
321                 }
322                 if (tid >= 0) {
323                         int l = strlen(b->ref);
324                         l = b->pos + (l > 0? l : 1);
325                         if (b->tid != tid || b->pos >= end) break;
326                         if (!(l > begin && end > b->pos)) continue;
327                 }
328                 ++n_processed;
329                 if (vc.flag & VC_QCALL) { // output QCALL format; STOP here
330                         bcf_2qcall(h, b);
331                         continue;
332                 }
333                 if (vc.flag & (VC_CALL|VC_ADJLD)) bcf_gl2pl(b);
334                 if (vc.flag & VC_CALL) { // call variants
335                         bcf_p1rst_t pr;
336                         bcf_p1_cal(b, p1, &pr); // pr.g[3] is not calculated here
337                         if (vc.flag&VC_HWE) bcf_p1_cal_g3(p1, pr.g);
338                         if (n_processed % 100000 == 0) {
339                                 fprintf(stderr, "[%s] %ld sites processed.\n", __func__, (long)n_processed);
340                                 bcf_p1_dump_afs(p1);
341                         }
342                         if (pr.p_ref >= vc.pref && (vc.flag & VC_VARONLY)) continue;
343                         update_bcf1(h->n_smpl, b, p1, &pr, vc.pref, vc.flag);
344                 }
345                 if (vc.flag & VC_ADJLD) { // compute LD
346                         double f[4], r2;
347                         if ((r2 = bcf_ld_freq(blast, b, f)) >= 0) {
348                                 kstring_t s;
349                                 s.m = s.l = 0; s.s = 0;
350                                 if (*b->info) kputc(';', &s);
351                                 ksprintf(&s, "NEIR=%.3lf;NEIF=%.3lf,%.3lf", r2, f[0]+f[2], f[0]+f[1]);
352                                 bcf_append_info(b, s.s, s.l);
353                                 free(s.s);
354                         }
355                         bcf_cpy(blast, b);
356                 }
357                 if (vc.flag & VC_NO_GENO) { // do not output GENO fields
358                         b->n_gi = 0;
359                         b->fmt[0] = '\0';
360                 }
361                 vcf_write(bout, h, b);
362         }
363         if (vc.prior_file) free(vc.prior_file);
364         if (vc.flag & VC_CALL) bcf_p1_dump_afs(p1);
365         bcf_hdr_destroy(h);
366         bcf_destroy(b); bcf_destroy(blast);
367         vcf_close(bp); vcf_close(bout);
368         if (hash) kh_destroy(set64, hash);
369         if (vc.fn_list) free(vc.fn_list);
370         if (p1) bcf_p1_destroy(p1);
371         return 0;
372 }