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