]> git.donarmstrong.com Git - samtools.git/blob - bam_plcmd.c
posterior expectation FINALLY working. I am so tired...
[samtools.git] / bam_plcmd.c
1 #include <math.h>
2 #include <stdio.h>
3 #include <unistd.h>
4 #include <ctype.h>
5 #include "sam.h"
6 #include "faidx.h"
7 #include "bam_maqcns.h"
8 #include "bam_mcns.h"
9 #include "khash.h"
10 #include "glf.h"
11 #include "kstring.h"
12
13 typedef int *indel_list_t;
14 KHASH_MAP_INIT_INT64(64, indel_list_t)
15
16 #define BAM_PLF_SIMPLE     0x01
17 #define BAM_PLF_CNS        0x02
18 #define BAM_PLF_INDEL_ONLY 0x04
19 #define BAM_PLF_GLF        0x08
20 #define BAM_PLF_VAR_ONLY   0x10
21 #define BAM_PLF_2ND        0x20
22 #define BAM_PLF_RANBASE    0x40
23 #define BAM_PLF_1STBASE    0x80
24 #define BAM_PLF_ALLBASE    0x100
25 #define BAM_PLF_READPOS    0x200
26
27 typedef struct {
28         bam_header_t *h;
29         bam_maqcns_t *c;
30         bam_maqindel_opt_t *ido;
31         faidx_t *fai;
32         khash_t(64) *hash;
33         uint32_t format;
34         int tid, len, last_pos;
35         int mask;
36     int max_depth;  // for indel calling, ignore reads with the depth too high. 0 for unlimited
37         char *ref;
38         glfFile fp_glf; // for glf output only
39 } pu_data_t;
40
41 char **__bam_get_lines(const char *fn, int *_n);
42 void bam_init_header_hash(bam_header_t *header);
43 int32_t bam_get_tid(const bam_header_t *header, const char *seq_name);
44
45 static khash_t(64) *load_pos(const char *fn, bam_header_t *h)
46 {
47         char **list;
48         int i, j, n, *fields, max_fields;
49         khash_t(64) *hash;
50         bam_init_header_hash(h);
51         list = __bam_get_lines(fn, &n);
52         hash = kh_init(64);
53         max_fields = 0; fields = 0;
54         for (i = 0; i < n; ++i) {
55                 char *str = list[i];
56                 int chr, n_fields, ret;
57                 khint_t k;
58                 uint64_t x;
59                 n_fields = ksplit_core(str, 0, &max_fields, &fields);
60                 if (n_fields < 2) continue;
61                 chr = bam_get_tid(h, str + fields[0]);
62                 if (chr < 0) {
63                         fprintf(stderr, "[load_pos] unknown reference sequence name: %s\n", str + fields[0]);
64                         continue;
65                 }
66                 x = (uint64_t)chr << 32 | (atoi(str + fields[1]) - 1);
67                 k = kh_put(64, hash, x, &ret);
68                 if (ret == 0) {
69                         fprintf(stderr, "[load_pos] position %s:%s has been loaded.\n", str+fields[0], str+fields[1]);
70                         continue;
71                 }
72                 kh_val(hash, k) = 0;
73                 if (n_fields > 2) {
74                         // count
75                         for (j = 2; j < n_fields; ++j) {
76                                 char *s = str + fields[j];
77                                 if ((*s != '+' && *s != '-') || !isdigit(s[1])) break;
78                         }
79                         if (j > 2) { // update kh_val()
80                                 int *q, y, z;
81                                 q = kh_val(hash, k) = (int*)calloc(j - 1, sizeof(int));
82                                 q[0] = j - 2; z = j; y = 1;
83                                 for (j = 2; j < z; ++j)
84                                         q[y++] = atoi(str + fields[j]);
85                         }
86                 }
87                 free(str);
88         }
89         free(list); free(fields);
90         return hash;
91 }
92
93 // an analogy to pileup_func() below
94 static int glt3_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
95 {
96         pu_data_t *d = (pu_data_t*)data;
97         bam_maqindel_ret_t *r = 0;
98         int rb, *proposed_indels = 0;
99         glf1_t *g;
100         glf3_t *g3;
101
102         if (d->fai == 0) {
103                 fprintf(stderr, "[glt3_func] reference sequence is required for generating GLT. Abort!\n");
104                 exit(1);
105         }
106         if (d->hash) { // only output a list of sites
107                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
108                 if (k == kh_end(d->hash)) return 0;
109                 proposed_indels = kh_val(d->hash, k);
110         }
111         g3 = glf3_init1();
112         if (d->fai && (int)tid != d->tid) {
113                 if (d->ref) { // then write the end mark
114                         g3->rtype = GLF3_RTYPE_END;
115                         glf3_write1(d->fp_glf, g3);
116                 }
117                 glf3_ref_write(d->fp_glf, d->h->target_name[tid], d->h->target_len[tid]); // write reference
118                 free(d->ref);
119                 d->ref = fai_fetch(d->fai, d->h->target_name[tid], &d->len);
120                 d->tid = tid;
121                 d->last_pos = 0;
122         }
123         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
124         g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
125         memcpy(g3, g, sizeof(glf1_t));
126         g3->rtype = GLF3_RTYPE_SUB;
127         g3->offset = pos - d->last_pos;
128         d->last_pos = pos;
129         glf3_write1(d->fp_glf, g3);
130     if (pos < d->len) {
131         int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
132                 if (proposed_indels)
133                         r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
134                 else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
135         }
136         if (r) { // then write indel line
137                 int het = 3 * n, min;
138                 min = het;
139                 if (min > r->gl[0]) min = r->gl[0];
140                 if (min > r->gl[1]) min = r->gl[1];
141                 g3->ref_base = 0;
142                 g3->rtype = GLF3_RTYPE_INDEL;
143                 memset(g3->lk, 0, 10);
144                 g3->lk[0] = r->gl[0] - min < 255? r->gl[0] - min : 255;
145                 g3->lk[1] = r->gl[1] - min < 255? r->gl[1] - min : 255;
146                 g3->lk[2] = het - min < 255? het - min : 255;
147                 g3->offset = 0;
148                 g3->indel_len[0] = r->indel1;
149                 g3->indel_len[1] = r->indel2;
150                 g3->min_lk = min < 255? min : 255;
151                 g3->max_len = (abs(r->indel1) > abs(r->indel2)? abs(r->indel1) : abs(r->indel2)) + 1;
152                 g3->indel_seq[0] = strdup(r->s[0]+1);
153                 g3->indel_seq[1] = strdup(r->s[1]+1);
154                 glf3_write1(d->fp_glf, g3);
155                 bam_maqindel_ret_destroy(r);
156         }
157         free(g);
158         glf3_destroy1(g3);
159         return 0;
160 }
161
162 static void pileup_seq(const bam_pileup1_t *p, int pos, int ref_len, const char *ref)
163 {
164         if (p->is_head) printf("^%c", p->b->core.qual > 93? 126 : p->b->core.qual + 33);
165         if (!p->is_del) {
166                 int j, rb, c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)];
167                 rb = (ref && pos < ref_len)? ref[pos] : 'N';
168                 if (c == '=' || toupper(c) == toupper(rb)) c = bam1_strand(p->b)? ',' : '.';
169                 else c = bam1_strand(p->b)? tolower(c) : toupper(c);
170                 putchar(c);
171                 if (p->indel > 0) {
172                         printf("+%d", p->indel);
173                         for (j = 1; j <= p->indel; ++j) {
174                                 c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)];
175                                 putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
176                         }
177                 } else if (p->indel < 0) {
178                         printf("%d", p->indel);
179                         for (j = 1; j <= -p->indel; ++j) {
180                                 c = (ref && (int)pos+j < ref_len)? ref[pos+j] : 'N';
181                                 putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
182                         }
183                 }
184         } else putchar('*');
185         if (p->is_tail) putchar('$');
186 }
187
188 static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
189 {
190         pu_data_t *d = (pu_data_t*)data;
191         bam_maqindel_ret_t *r = 0;
192         int i, rb, rms_mapq = -1, *proposed_indels = 0;
193         uint64_t rms_aux;
194         uint32_t cns = 0;
195
196         // if GLF is required, suppress -c completely
197         if (d->format & BAM_PLF_GLF) return glt3_func(tid, pos, n, pu, data);
198         // if d->hash is initialized, only output the sites in the hash table
199         if (d->hash) {
200                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
201                 if (k == kh_end(d->hash)) return 0;
202                 proposed_indels = kh_val(d->hash, k);
203         }
204         // update d->ref if necessary
205         if (d->fai && (int)tid != d->tid) {
206                 free(d->ref);
207                 d->ref = faidx_fetch_seq(d->fai, d->h->target_name[tid], 0, 0x7fffffff, &d->len);
208                 d->tid = tid;
209         }
210         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
211         // when the indel-only mode is asked for, return if no reads mapped with indels
212         if (d->format & BAM_PLF_INDEL_ONLY) {
213                 for (i = 0; i < n; ++i)
214                         if (pu[i].indel != 0) break;
215                 if (i == n) return 0;
216         }
217         // call the consensus and indel
218         if (d->format & BAM_PLF_CNS) { // call consensus
219                 if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE)) { // use a random base or the 1st base as the consensus call
220                         const bam_pileup1_t *p = (d->format & BAM_PLF_1STBASE)? pu : pu + (int)(drand48() * n);
221                         int q = bam1_qual(p->b)[p->qpos];
222                         int mapQ = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
223                         uint32_t b = bam1_seqi(bam1_seq(p->b), p->qpos);
224                         cns = b<<28 | 0xf<<24 | mapQ<<16 | q<<8;
225                 } else if (d->format & BAM_PLF_ALLBASE) { // collapse all bases
226                         uint64_t rmsQ = 0;
227                         uint32_t b = 0;
228                         for (i = 0; i < n; ++i) {
229                                 const bam_pileup1_t *p = pu + i;
230                                 int q = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
231                                 b |= bam1_seqi(bam1_seq(p->b), p->qpos);
232                                 rmsQ += q * q;
233                         }
234                         rmsQ = (uint64_t)(sqrt((double)rmsQ / n) + .499);
235                         cns = b<<28 | 0xf<<24 | rmsQ<<16 | 60<<8;
236                 } else cns = bam_maqcns_call(n, pu, d->c);
237         }
238     if ((d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)) && d->ref && pos < d->len) { // call indels
239         int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
240         if (proposed_indels) // the first element gives the size of the array
241             r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
242         else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
243         }
244         // when only variant sites are asked for, test if the site is a variant
245         if ((d->format & BAM_PLF_CNS) && (d->format & BAM_PLF_VAR_ONLY)) {
246                 if (!(bam_nt16_table[rb] != 15 && cns>>28 != bam_nt16_table[rb])) { // not a SNP
247                         if (!(r && (r->gt == 2 || strcmp(r->s[r->gt], "*")))) { // not an indel
248                                 if (r) bam_maqindel_ret_destroy(r);
249                                 return 0;
250                         }
251                 }
252         }
253         // print the first 3 columns
254         printf("%s\t%d\t%c\t", d->h->target_name[tid], pos + 1, rb);
255         // print consensus information if required
256         if (d->format & BAM_PLF_CNS) {
257                 int ref_q, rb4 = bam_nt16_table[rb];
258                 ref_q = 0;
259                 if (rb4 != 15 && cns>>28 != 15 && cns>>28 != rb4) { // a SNP
260                         ref_q = ((cns>>24&0xf) == rb4)? cns>>8&0xff : (cns>>8&0xff) + (cns&0xff);
261                         if (ref_q > 255) ref_q = 255;
262                 }
263                 rms_mapq = cns>>16&0xff;
264                 printf("%c\t%d\t%d\t%d\t", bam_nt16_rev_table[cns>>28], cns>>8&0xff, ref_q, rms_mapq);
265         }
266         // print pileup sequences
267         printf("%d\t", n);
268         rms_aux = 0; // we need to recalculate rms_mapq when -c is not flagged on the command line
269         for (i = 0; i < n; ++i) {
270                 const bam_pileup1_t *p = pu + i;
271                 int tmp = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
272                 rms_aux += tmp * tmp;
273                 pileup_seq(p, pos, d->len, d->ref);
274         }
275         // finalize rms_mapq
276         rms_aux = (uint64_t)(sqrt((double)rms_aux / n) + .499);
277         if (rms_mapq < 0) rms_mapq = rms_aux;
278         putchar('\t');
279         // print quality
280         for (i = 0; i < n; ++i) {
281                 const bam_pileup1_t *p = pu + i;
282                 int c = bam1_qual(p->b)[p->qpos] + 33;
283                 if (c > 126) c = 126;
284                 putchar(c);
285         }
286         if (d->format & BAM_PLF_2ND) { // print 2nd calls and qualities
287                 const unsigned char *q;
288                 putchar('\t');
289                 for (i = 0; i < n; ++i) {
290                         const bam_pileup1_t *p = pu + i;
291                         q = bam_aux_get(p->b, "E2");
292                         putchar(q? q[p->qpos + 1] : 'N');
293                 }
294                 putchar('\t');
295                 for (i = 0; i < n; ++i) {
296                         const bam_pileup1_t *p = pu + i;
297                         q = bam_aux_get(p->b, "U2");
298                         putchar(q? q[p->qpos + 1] : '!');
299                 }
300         }
301         // print mapping quality if -s is flagged on the command line
302         if (d->format & BAM_PLF_SIMPLE) {
303                 putchar('\t');
304                 for (i = 0; i < n; ++i) {
305                         int c = pu[i].b->core.qual + 33;
306                         if (c > 126) c = 126;
307                         putchar(c);
308                 }
309         }
310         // print read position
311         if (d->format & BAM_PLF_READPOS) {
312                 putchar('\t');
313                 for (i = 0; i < n; ++i) {
314                         int x = pu[i].qpos;
315                         int l = pu[i].b->core.l_qseq;
316                         printf("%d,", x < l/2? x+1 : -((l-1)-x+1));
317                 }
318         }
319         putchar('\n');
320         // print the indel line if r has been calculated. This only happens if:
321         // a) -c or -i are flagged, AND b) the reference sequence is available
322         if (r) {
323                 printf("%s\t%d\t*\t", d->h->target_name[tid], pos + 1);
324                 if (r->gt < 2) printf("%s/%s\t", r->s[r->gt], r->s[r->gt]);
325                 else printf("%s/%s\t", r->s[0], r->s[1]);
326                 printf("%d\t%d\t", r->q_cns, r->q_ref);
327                 printf("%d\t%d\t", rms_mapq, n);
328                 printf("%s\t%s\t", r->s[0], r->s[1]);
329                 //printf("%d\t%d\t", r->gl[0], r->gl[1]);
330                 printf("%d\t%d\t%d\t", r->cnt1, r->cnt2, r->cnt_anti);
331                 printf("%d\t%d\n", r->cnt_ref, r->cnt_ambi);
332                 bam_maqindel_ret_destroy(r);
333         }
334         return 0;
335 }
336
337 int bam_pileup(int argc, char *argv[])
338 {
339         int c, is_SAM = 0;
340         char *fn_list = 0, *fn_fa = 0, *fn_pos = 0;
341         pu_data_t *d = (pu_data_t*)calloc(1, sizeof(pu_data_t));
342     d->max_depth = 0;
343         d->tid = -1; d->mask = BAM_DEF_MASK;
344         d->c = bam_maqcns_init();
345         d->c->is_soap = 1; // change the default model
346         d->ido = bam_maqindel_opt_init();
347         while ((c = getopt(argc, argv, "st:f:cT:N:r:l:d:im:gI:G:vM:S2aR:PA")) >= 0) {
348                 switch (c) {
349                 case 'a': d->c->is_soap = 1; break;
350                 case 'A': d->c->is_soap = 0; break;
351                 case 's': d->format |= BAM_PLF_SIMPLE; break;
352                 case 't': fn_list = strdup(optarg); break;
353                 case 'l': fn_pos = strdup(optarg); break;
354                 case 'f': fn_fa = strdup(optarg); break;
355                 case 'T': d->c->theta = atof(optarg); break;
356                 case 'N': d->c->n_hap = atoi(optarg); break;
357                 case 'r': d->c->het_rate = atof(optarg); d->ido->r_snp = d->c->het_rate; break;
358                 case 'M': d->c->cap_mapQ = atoi(optarg); break;
359                 case 'd': d->max_depth = atoi(optarg); break;
360                 case 'c': d->format |= BAM_PLF_CNS; break;
361                 case 'i': d->format |= BAM_PLF_INDEL_ONLY; break;
362                 case 'v': d->format |= BAM_PLF_VAR_ONLY; break;
363                 case 'm': d->mask = strtol(optarg, 0, 0); break;
364                 case 'g': d->format |= BAM_PLF_GLF; break;
365                 case '2': d->format |= BAM_PLF_2ND; break;
366                 case 'P': d->format |= BAM_PLF_READPOS; break;
367                 case 'I': d->ido->q_indel = atoi(optarg); break;
368                 case 'G': d->ido->r_indel = atof(optarg); break;
369                 case 'S': is_SAM = 1; break;
370                 case 'R':
371                         if (strcmp(optarg, "random") == 0) d->format |= BAM_PLF_RANBASE;
372                         else if (strcmp(optarg, "first") == 0) d->format |= BAM_PLF_1STBASE;
373                         else if (strcmp(optarg, "all") == 0) d->format |= BAM_PLF_ALLBASE;
374                         else fprintf(stderr, "[bam_pileup] unrecognized -R\n");
375                         break;
376                 default: fprintf(stderr, "Unrecognizd option '-%c'.\n", c); return 1;
377                 }
378         }
379         if (fn_list) is_SAM = 1;
380         if (optind == argc) {
381                 fprintf(stderr, "\n");
382                 fprintf(stderr, "Usage:  samtools pileup [options] <in.bam>|<in.sam>\n\n");
383                 fprintf(stderr, "Option: -s        simple (yet incomplete) pileup format\n");
384                 fprintf(stderr, "        -S        the input is in SAM\n");
385                 fprintf(stderr, "        -A        use the MAQ model for SNP calling\n");
386                 fprintf(stderr, "        -2        output the 2nd best call and quality\n");
387                 fprintf(stderr, "        -i        only show lines/consensus with indels\n");
388                 fprintf(stderr, "        -m INT    filtering reads with bits in INT [%d]\n", d->mask);
389                 fprintf(stderr, "        -M INT    cap mapping quality at INT [%d]\n", d->c->cap_mapQ);
390         fprintf(stderr, "        -d INT    limit maximum depth for indels [unlimited]\n");
391                 fprintf(stderr, "        -t FILE   list of reference sequences (force -S)\n");
392                 fprintf(stderr, "        -l FILE   list of sites at which pileup is output\n");
393                 fprintf(stderr, "        -f FILE   reference sequence in the FASTA format\n\n");
394                 fprintf(stderr, "        -c        output the SOAPsnp consensus sequence\n");
395                 fprintf(stderr, "        -v        print variants only (for -c)\n");
396                 fprintf(stderr, "        -g        output in the GLFv3 format (suppressing -c/-i/-s)\n");
397                 fprintf(stderr, "        -T FLOAT  theta in maq consensus calling model (for -c/-g) [%f]\n", d->c->theta);
398                 fprintf(stderr, "        -N INT    number of haplotypes in the sample (for -c/-g) [%d]\n", d->c->n_hap);
399                 fprintf(stderr, "        -r FLOAT  prior of a difference between two haplotypes (for -c/-g) [%f]\n", d->c->het_rate);
400                 fprintf(stderr, "        -G FLOAT  prior of an indel between two haplotypes (for -c/-g) [%f]\n", d->ido->r_indel);
401                 fprintf(stderr, "        -I INT    phred prob. of an indel in sequencing/prep. (for -c/-g) [%d]\n", d->ido->q_indel);
402                 fprintf(stderr, "\n");
403                 free(fn_list); free(fn_fa); free(d);
404                 return 1;
405         }
406         if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE|BAM_PLF_ALLBASE)) d->format |= BAM_PLF_CNS;
407         if (fn_fa) d->fai = fai_load(fn_fa);
408         if (d->format & (BAM_PLF_CNS|BAM_PLF_GLF)) bam_maqcns_prepare(d->c); // consensus calling
409         if (d->format & BAM_PLF_GLF) { // for glf output
410                 glf3_header_t *h;
411                 h = glf3_header_init();
412                 d->fp_glf = bgzf_fdopen(fileno(stdout), "w");
413                 glf3_header_write(d->fp_glf, h);
414                 glf3_header_destroy(h);
415         }
416         if (d->fai == 0 && (d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)))
417                 fprintf(stderr, "[bam_pileup] indels will not be called when -f is absent.\n");
418         if (fn_fa && is_SAM && fn_list == 0) fn_list = samfaipath(fn_fa);
419
420         {
421                 samfile_t *fp;
422                 fp = is_SAM? samopen(argv[optind], "r", fn_list) : samopen(argv[optind], "rb", 0);
423                 if (fp == 0 || fp->header == 0) {
424                         fprintf(stderr, "[bam_pileup] fail to read the header: non-exisiting file or wrong format.\n");
425                         return 1;
426                 }
427                 d->h = fp->header;
428                 if (fn_pos) d->hash = load_pos(fn_pos, d->h);
429                 sampileup(fp, d->mask, pileup_func, d);
430                 samclose(fp); // d->h will be destroyed here
431         }
432
433         // free
434         if (d->format & BAM_PLF_GLF) bgzf_close(d->fp_glf);
435         if (fn_pos) { // free the hash table
436                 khint_t k;
437                 for (k = kh_begin(d->hash); k < kh_end(d->hash); ++k)
438                         if (kh_exist(d->hash, k)) free(kh_val(d->hash, k));
439                 kh_destroy(64, d->hash);
440         }
441         free(fn_pos); free(fn_list); free(fn_fa);
442         if (d->fai) fai_destroy(d->fai);
443         bam_maqcns_destroy(d->c);
444         free(d->ido); free(d->ref); free(d);
445         return 0;
446 }
447
448 /***********
449  * mpileup *
450  ***********/
451
452 typedef struct {
453         int vcf, max_mq, min_mq, var_only, prior_type, afs;
454         double theta;
455         char *reg, *fn_pos;
456         faidx_t *fai;
457         kh_64_t *hash;
458 } mplp_conf_t;
459
460 typedef struct {
461         bamFile fp;
462         bam_iter_t iter;
463         int min_mq;
464 } mplp_aux_t;
465
466 static int mplp_func(void *data, bam1_t *b)
467 {
468         mplp_aux_t *ma = (mplp_aux_t*)data;
469         int ret;
470         do {
471                 ret = ma->iter? bam_iter_read(ma->fp, ma->iter, b) : bam_read1(ma->fp, b);
472         } while (b->core.qual < ma->min_mq && ret >= 0);
473         return ret;
474 }
475
476 static int mpileup(mplp_conf_t *conf, int n, char **fn)
477 {
478         mplp_aux_t **data;
479         mc_aux_t *ma = 0;
480         int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid;
481         const bam_pileup1_t **plp;
482         bam_mplp_t iter;
483         bam_header_t *h = 0;
484         char *ref;
485         khash_t(64) *hash = 0;
486         // allocate
487         data = calloc(n, sizeof(void*));
488         plp = calloc(n, sizeof(void*));
489         n_plp = calloc(n, sizeof(int*));
490         // read the header and initialize data
491         for (i = 0; i < n; ++i) {
492                 bam_header_t *h_tmp;
493                 data[i] = calloc(1, sizeof(mplp_aux_t));
494                 data[i]->min_mq = conf->min_mq;
495                 data[i]->fp = bam_open(fn[i], "r");
496                 h_tmp = bam_header_read(data[i]->fp);
497                 if (conf->reg) {
498                         int beg, end;
499                         bam_index_t *idx;
500                         idx = bam_index_load(fn[i]);
501                         if (idx == 0) {
502                                 fprintf(stderr, "[%s] fail to load index for %d-th input.\n", __func__, i+1);
503                                 exit(1);
504                         }
505                         if (bam_parse_region(h_tmp, conf->reg, &tid, &beg, &end) < 0) {
506                                 fprintf(stderr, "[%s] malformatted region or wrong seqname for %d-th input.\n", __func__, i+1);
507                                 exit(1);
508                         }
509                         if (i == 0) beg0 = beg, end0 = end;
510                         data[i]->iter = bam_iter_query(idx, tid, beg, end);
511                         bam_index_destroy(idx);
512                 }
513                 if (i == 0) h = h_tmp;
514                 else {
515                         // FIXME: to check consistency
516                         bam_header_destroy(h_tmp);
517                 }
518         }
519         if (conf->fn_pos) hash = load_pos(conf->fn_pos, h);
520         // write the VCF header
521         if (conf->vcf) {
522                 kstring_t s;
523                 s.l = s.m = 0; s.s = 0;
524                 puts("##fileformat=VCFv4.0");
525                 puts("##INFO=<ID=DP,Number=1,Type=Integer,Description=\"Total read depth\">");
526 //              puts("##INFO=<ID=AF,Number=1,Type=Float,Description=\"Posterior non-reference allele frequency\">");
527 //              puts("##INFO=<ID=AFEM,Number=1,Type=Float,Description=\"Prior-free non-reference allele frequency by EM\">");
528                 puts("##FILTER=<ID=Q13,Description=\"All min{baseQ,mapQ} below 13\">");
529                 kputs("#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO\tFORMAT", &s);
530                 for (i = 0; i < n; ++i) {
531                         const char *p;
532                         kputc('\t', &s);
533                         if ((p = strstr(fn[i], ".bam")) != 0)
534                                 kputsn(fn[i], p - fn[i], &s);
535                         else kputs(fn[i], &s);
536                 }
537                 puts(s.s);
538                 free(s.s);
539         }
540         // mpileup
541         if (conf->vcf) {
542                 ma = mc_init(n);
543                 mc_init_prior(ma, conf->prior_type, conf->theta);
544         }
545         ref_tid = -1; ref = 0;
546         iter = bam_mplp_init(n, mplp_func, (void**)data);
547         while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
548                 if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested
549                 if (hash) {
550                         khint_t k;
551                         k = kh_get(64, hash, (uint64_t)tid<<32 | pos);
552                         if (k == kh_end(hash)) continue;
553                 }
554                 if (tid != ref_tid) {
555                         free(ref);
556                         if (conf->fai) ref = fai_fetch(conf->fai, h->target_name[tid], &ref_len);
557                         ref_tid = tid;
558                 }
559                 if (conf->vcf) {
560                         mc_rst_t r;
561                         int j, _ref0, depth, rms_q, _ref0b, is_var = 0, qref = 0, level = 2, tot;
562                         uint64_t sqr_sum;
563                         if (conf->afs) level = 3;
564                         _ref0 = _ref0b = (ref && pos < ref_len)? ref[pos] : 'N';
565                         _ref0 = bam_nt16_nt4_table[bam_nt16_table[_ref0]];
566                         tot = mc_cal(_ref0, n_plp, plp, ma, &r, level);
567                         if (tot) { // has good bases
568                                 double q;
569                                 is_var = (r.p_ref < .5);
570                                 q = is_var? r.p_ref : 1. - r.p_ref;
571                                 if (q < 1e-308) q = 1e-308;
572                                 qref = (int)(-3.434 * log(q) + .499);
573                                 if (qref > 99) qref = 99;
574                         }
575                         if (conf->var_only && !is_var) continue;
576                         printf("%s\t%d\t.\t%c\t", h->target_name[tid], pos + 1, _ref0b);
577                         if (is_var) {
578                                 if (_ref0 == r.ref) putchar("ACGTN"[r.alt]);
579                                 else printf("%c,%c", "ACGTN"[r.ref], "ACGTN"[r.alt]);
580                         } else putchar('.');
581                         printf("\t%d\t", qref);
582                         if (!tot) printf("Q13\t");
583                         else printf(".\t");
584                         for (i = depth = 0, sqr_sum = 0; i < n; ++i) {
585                                 depth += n_plp[i];
586                                 for (j = 0; j < n_plp[i]; ++j) {
587                                         int q = plp[i][j].b->core.qual;
588                                         if (q > conf->max_mq) q = conf->max_mq;
589                                         sqr_sum += q * q;
590                                 }
591                         }
592                         rms_q = (int)(sqrt((double)sqr_sum / depth) + .499);
593                         printf("DP=%d;MQ=%d", depth, rms_q);
594                         printf(";AF=%.3lg", 1. - r.f_em);
595                         if (conf->afs)
596                                 printf(";AF0=%.3lg;AFN=%.3lg;AFE=%.3lg", 1-r.f_naive, 1-r.f_nielsen, 1-r.f_exp);
597                         printf("\tGT:GQ:DP");
598                         if (tot) {
599                                 for (i = 0; i < n; ++i) {
600                                         int x = mc_call_gt(ma, r.f_em, i);
601                                         printf("\t%c/%c:%d:%d", "10"[((x&3)==2)], "10"[((x&3)>0)], x>>2, n_plp[i]);
602                                 }
603                         } else for (i = 0; i < n; ++i) printf("\t./.:0:0");
604                         putchar('\n');
605                 } else {
606                         printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
607                         for (i = 0; i < n; ++i) {
608                                 int j;
609                                 printf("\t%d\t", n_plp[i]);
610                                 if (n_plp[i] == 0) printf("*\t*");
611                                 else {
612                                         for (j = 0; j < n_plp[i]; ++j)
613                                                 pileup_seq(plp[i] + j, pos, ref_len, ref);
614                                         putchar('\t');
615                                         for (j = 0; j < n_plp[i]; ++j) {
616                                                 const bam_pileup1_t *p = plp[i] + j;
617                                                 int c = bam1_qual(p->b)[p->qpos] + 33;
618                                                 if (c > 126) c = 126;
619                                                 putchar(c);
620                                         }
621                                 }
622                         }
623                         putchar('\n');
624                 }
625         }
626         if (hash) { // free the hash table
627                 khint_t k;
628                 for (k = kh_begin(hash); k < kh_end(hash); ++k)
629                         if (kh_exist(hash, k)) free(kh_val(hash, k));
630                 kh_destroy(64, hash);
631         }
632         mc_destroy(ma);
633         bam_mplp_destroy(iter);
634         bam_header_destroy(h);
635         for (i = 0; i < n; ++i) {
636                 bam_close(data[i]->fp);
637                 if (data[i]->iter) bam_iter_destroy(data[i]->iter);
638                 free(data[i]);
639         }
640         free(data); free(plp); free(ref); free(n_plp);
641         return 0;
642 }
643
644 int bam_mpileup(int argc, char *argv[])
645 {
646         int c;
647         mplp_conf_t mplp;
648         memset(&mplp, 0, sizeof(mplp_conf_t));
649         mplp.max_mq = 60;
650         mplp.prior_type = MC_PTYPE_FULL;
651         mplp.theta = 1e-3;
652         while ((c = getopt(argc, argv, "f:r:l:VvM:q:t:2F")) >= 0) {
653                 switch (c) {
654                 case 't': mplp.theta = atof(optarg); break;
655                 case '2': mplp.prior_type = MC_PTYPE_COND2; break;
656                 case 'f':
657                         mplp.fai = fai_load(optarg);
658                         if (mplp.fai == 0) return 1;
659                         break;
660                 case 'r': mplp.reg = strdup(optarg); break;
661                 case 'l': mplp.fn_pos = strdup(optarg); break;
662                 case 'v': mplp.var_only = 1; break;
663                 case 'V': mplp.vcf = 1; break;
664                 case 'M': mplp.max_mq = atoi(optarg); break;
665                 case 'q': mplp.min_mq = atoi(optarg); break;
666                 case 'F': mplp.afs = 1; break;
667                 }
668         }
669         if (argc == 1) {
670                 fprintf(stderr, "\n");
671                 fprintf(stderr, "Usage:   samtools mpileup [options] in1.bam [in2.bam [...]]\n\n");
672                 fprintf(stderr, "Options: -f FILE     reference sequence file [null]\n");
673                 fprintf(stderr, "         -r STR      region in which pileup is generated [null]\n");
674                 fprintf(stderr, "         -l FILE     list of positions (format: chr pos) [null]\n");
675                 fprintf(stderr, "         -M INT      cap mapping quality at INT [%d]\n", mplp.max_mq);
676                 fprintf(stderr, "         -q INT      filter out alignment with MQ smaller than INT [%d]\n", mplp.min_mq);
677                 fprintf(stderr, "         -t FLOAT    scaled mutation rate [%lg]\n", mplp.theta);
678                 fprintf(stderr, "         -V          generate VCF output (SNP calling)\n");
679                 fprintf(stderr, "         -v          show variant sites only\n");
680                 fprintf(stderr, "         -2          conditional prior\n");
681                 fprintf(stderr, "\n");
682                 fprintf(stderr, "Notes: Assuming error independency and diploid individuals.\n\n");
683                 return 1;
684         }
685         mpileup(&mplp, argc - optind, argv + optind);
686         free(mplp.reg);
687         if (mplp.fai) fai_destroy(mplp.fai);
688         return 0;
689 }