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