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