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