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