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