]> git.donarmstrong.com Git - samtools.git/blob - bcftools/vcfout.c
e02e364ca1d2a3921e6900fa0cf4e0c63908d780
[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_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 extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
81
82 static double test_fisher(bcf1_t *b, const char *key, int d[4], int is_single)
83 {
84         double left, right, two;
85         char *p;
86         int i;
87         if ((p = strstr(b->info, key)) == 0) return -1.;
88         p += 4;
89         for (i = 0; i < 4; ++i) {
90                 d[i] = strtol(p, &p, 10);
91                 if (d[i] == 0 && (errno == EINVAL || errno == ERANGE)) return -2.;
92                 ++p;
93         }
94         kt_fisher_exact(d[0], d[1], d[2], d[3], &left, &right, &two);
95         return is_single? right : two;
96 }
97
98 static void rm_info(int n_smpl, bcf1_t *b, const char *key)
99 {
100         char *p, *q;
101         if ((p = strstr(b->info, key)) == 0) return;
102         for (q = p; *q && *q != ';'; ++q);
103         if (p > b->info && *(p-1) == ';') --p;
104         memmove(p, q, b->l_str - (q - b->str));
105         b->l_str -= q - p;
106         bcf_sync(n_smpl, b);
107 }
108
109 static int update_bcf1(int n_smpl, bcf1_t *b, const bcf_p1aux_t *pa, const bcf_p1rst_t *pr, double pref, int flag)
110 {
111         kstring_t s;
112         int d[4], is_var = (pr->p_ref < pref);
113         double p_hwe, p_dp, p_ed, r = is_var? pr->p_ref : 1. - pr->p_ref;
114
115         p_hwe = test_hwe(pr->g);
116         p_ed = test_fisher(b, "ED4=", d, 1);
117         p_dp = test_fisher(b, "DP4=", d, 0);
118         rm_info(n_smpl, b, "ED4=");
119
120         memset(&s, 0, sizeof(kstring_t));
121         kputc('\0', &s); kputs(b->ref, &s); kputc('\0', &s);
122         if (is_var) {
123                 kputs(b->alt, &s);
124         }
125         kputc('\0', &s); kputc('\0', &s);
126         kputs(b->info, &s);
127         if (b->info[0]) kputc(';', &s);
128         ksprintf(&s, "AF1=%.3lf;AFE=%.3lf", 1.-pr->f_em, 1.-pr->f_exp);
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         b->qual = r < 1e-100? 99 : -3.434 * log(r);
137         if (b->qual > 99) b->qual = 99;
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         char moder[4], modew[4];
153         khash_t(set64) *hash = 0;
154
155         tid = begin = end = -1;
156         memset(&vc, 0, sizeof(viewconf_t));
157         vc.prior_type = -1; vc.theta = 1e-3; vc.pref = 0.9;
158         while ((c = getopt(argc, argv, "l:cGvLbSuP:t:p:")) >= 0) {
159                 switch (c) {
160                 case 'l': vc.fn_list = strdup(optarg); break;
161                 case 'G': vc.flag |= VC_NO_GENO; break;
162                 case 'L': vc.flag |= VC_NO_PL; break;
163                 case 'b': vc.flag |= VC_BCFOUT; break;
164                 case 'S': vc.flag |= VC_VCFIN; break;
165                 case 'c': vc.flag |= VC_CALL; break;
166                 case 'v': vc.flag |= VC_VARONLY; break;
167                 case 'u': vc.flag |= VC_UNCOMP | VC_BCFOUT; break;
168                 case 't': vc.theta = atof(optarg); break;
169                 case 'p': vc.pref = atof(optarg); break;
170                 case 'P':
171                         if (strcmp(optarg, "full") == 0) vc.prior_type = MC_PTYPE_FULL;
172                         else if (strcmp(optarg, "cond2") == 0) vc.prior_type = MC_PTYPE_COND2;
173                         else if (strcmp(optarg, "flat") == 0) vc.prior_type = MC_PTYPE_FLAT;
174                         else vc.prior_file = strdup(optarg);
175                         break;
176                 }
177         }
178         if (argc == optind) {
179                 fprintf(stderr, "\n");
180                 fprintf(stderr, "Usage:   bcftools view [options] <in.bcf> [reg]\n\n");
181                 fprintf(stderr, "Options: -c        SNP calling\n");
182                 fprintf(stderr, "         -b        output BCF instead of VCF\n");
183                 fprintf(stderr, "         -u        uncompressed BCF output\n");
184                 fprintf(stderr, "         -S        input is VCF\n");
185                 fprintf(stderr, "         -G        suppress all individual genotype information\n");
186                 fprintf(stderr, "         -L        discard the PL genotype field\n");
187                 fprintf(stderr, "         -v        output potential variant sites only\n");
188                 fprintf(stderr, "         -l FILE   list of sites to output [all sites]\n");
189                 fprintf(stderr, "         -t FLOAT  scaled mutation rate [%.4lg]\n", vc.theta);
190                 fprintf(stderr, "         -p FLOAT  variant if P(ref|D)<FLOAT [%.3lg]\n", vc.pref);
191                 fprintf(stderr, "         -P STR    type of prior: full, cond2, flat [full]\n");
192                 fprintf(stderr, "\n");
193                 return 1;
194         }
195
196         b = calloc(1, sizeof(bcf1_t));
197         strcpy(moder, "r");
198         if (!(vc.flag & VC_VCFIN)) strcat(moder, "b");
199         strcpy(modew, "w");
200         if (vc.flag & VC_BCFOUT) strcat(modew, "b");
201         if (vc.flag & VC_UNCOMP) strcat(modew, "u");
202         bp = vcf_open(argv[optind], moder);
203         h = vcf_hdr_read(bp);
204         bout = vcf_open("-", modew);
205         vcf_hdr_write(bout, h);
206         if (vc.flag & VC_CALL) {
207                 p1 = bcf_p1_init(h->n_smpl);
208                 if (vc.prior_file) {
209                         if (bcf_p1_read_prior(p1, vc.prior_file) < 0) {
210                                 fprintf(stderr, "[%s] fail to read the prior AFS.\n", __func__);
211                                 return 1;
212                         }
213                 } else bcf_p1_init_prior(p1, vc.prior_type, vc.theta);
214         }
215         if (vc.fn_list) hash = bcf_load_pos(vc.fn_list, h);
216         if (optind + 1 < argc && !(vc.flag&VC_VCFIN)) {
217                 void *str2id = bcf_build_refhash(h);
218                 if (bcf_parse_region(str2id, argv[optind+1], &tid, &begin, &end) >= 0) {
219                         bcf_idx_t *idx;
220                         idx = bcf_idx_load(argv[optind]);
221                         if (idx) {
222                                 uint64_t off;
223                                 off = bcf_idx_query(idx, tid, begin, end);
224                                 bgzf_seek(bp->fp, off, SEEK_SET);
225                                 bcf_idx_destroy(idx);
226                         }
227                 }
228         }
229         while (vcf_read(bp, h, b) > 0) {
230                 if (hash) {
231                         uint64_t x = (uint64_t)b->tid<<32 | b->pos;
232                         khint_t k = kh_get(set64, hash, x);
233                         if (k == kh_end(hash)) continue;
234                 }
235                 if (tid >= 0) {
236                         int l = strlen(b->ref);
237                         l = b->pos + (l > 0? l : 1);
238                         if (b->tid != tid || b->pos >= end) break;
239                         if (!(l > begin && end > b->pos)) continue;
240                 }
241                 if (vc.flag & VC_CALL) {
242                         bcf_p1rst_t pr;
243                         bcf_p1_cal(b, p1, &pr);
244                         if ((n_processed + 1) % 50000 == 0) bcf_p1_dump_afs(p1);
245                         if (pr.p_ref >= vc.pref && (vc.flag & VC_VARONLY)) continue;
246                         update_bcf1(h->n_smpl, b, p1, &pr, vc.pref, vc.flag);
247                 }
248                 if (vc.flag & VC_NO_GENO) {
249                         b->n_gi = 0;
250                         b->fmt[0] = '\0';
251                 }
252                 vcf_write(bout, h, b);
253                 ++n_processed;
254         }
255         if (vc.prior_file) free(vc.prior_file);
256         if (vc.flag & VC_CALL) bcf_p1_dump_afs(p1);
257         bcf_hdr_destroy(h);
258         bcf_destroy(b);
259         vcf_close(bp); vcf_close(bout);
260         if (hash) kh_destroy(set64, hash);
261         if (vc.fn_list) free(vc.fn_list);
262         if (p1) bcf_p1_destroy(p1);
263         return 0;
264 }
265