]> git.donarmstrong.com Git - samtools.git/blob - bam_plcmd.c
* samtools-0.1.8-20 (r778)
[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 "khash.h"
9 #include "glf.h"
10 #include "kstring.h"
11
12 typedef int *indel_list_t;
13 KHASH_MAP_INIT_INT64(64, indel_list_t)
14
15 #define BAM_PLF_SIMPLE     0x01
16 #define BAM_PLF_CNS        0x02
17 #define BAM_PLF_INDEL_ONLY 0x04
18 #define BAM_PLF_GLF        0x08
19 #define BAM_PLF_VAR_ONLY   0x10
20 #define BAM_PLF_2ND        0x20
21 #define BAM_PLF_RANBASE    0x40
22 #define BAM_PLF_1STBASE    0x80
23 #define BAM_PLF_ALLBASE    0x100
24 #define BAM_PLF_READPOS    0x200
25 #define BAM_PLF_NOBAQ      0x400
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 capQ_thres, min_baseQ;
37     int max_depth;  // for indel calling, ignore reads with the depth too high. 0 for unlimited
38         char *ref;
39         glfFile fp_glf; // for glf output only
40 } pu_data_t;
41
42 char **__bam_get_lines(const char *fn, int *_n);
43 void bam_init_header_hash(bam_header_t *header);
44 int32_t bam_get_tid(const bam_header_t *header, const char *seq_name);
45
46 static khash_t(64) *load_pos(const char *fn, bam_header_t *h)
47 {
48         char **list;
49         int i, j, n, *fields, max_fields;
50         khash_t(64) *hash;
51         bam_init_header_hash(h);
52         list = __bam_get_lines(fn, &n);
53         hash = kh_init(64);
54         max_fields = 0; fields = 0;
55         for (i = 0; i < n; ++i) {
56                 char *str = list[i];
57                 int chr, n_fields, ret;
58                 khint_t k;
59                 uint64_t x;
60                 n_fields = ksplit_core(str, 0, &max_fields, &fields);
61                 if (n_fields < 2) continue;
62                 chr = bam_get_tid(h, str + fields[0]);
63                 if (chr < 0) {
64                         fprintf(stderr, "[load_pos] unknown reference sequence name: %s\n", str + fields[0]);
65                         continue;
66                 }
67                 x = (uint64_t)chr << 32 | (atoi(str + fields[1]) - 1);
68                 k = kh_put(64, hash, x, &ret);
69                 if (ret == 0) {
70                         fprintf(stderr, "[load_pos] position %s:%s has been loaded.\n", str+fields[0], str+fields[1]);
71                         continue;
72                 }
73                 kh_val(hash, k) = 0;
74                 if (n_fields > 2) {
75                         // count
76                         for (j = 2; j < n_fields; ++j) {
77                                 char *s = str + fields[j];
78                                 if ((*s != '+' && *s != '-') || !isdigit(s[1])) break;
79                         }
80                         if (j > 2) { // update kh_val()
81                                 int *q, y, z;
82                                 q = kh_val(hash, k) = (int*)calloc(j - 1, sizeof(int));
83                                 q[0] = j - 2; z = j; y = 1;
84                                 for (j = 2; j < z; ++j)
85                                         q[y++] = atoi(str + fields[j]);
86                         }
87                 }
88                 free(str);
89         }
90         free(list); free(fields);
91         return hash;
92 }
93
94 static inline int printw(int c, FILE *fp)
95 {
96         char buf[16];
97         int l, x;
98         if (c == 0) return fputc('0', fp);
99         for (l = 0, x = c < 0? -c : c; x > 0; x /= 10) buf[l++] = x%10 + '0';
100         if (c < 0) buf[l++] = '-';
101         buf[l] = 0;
102         for (x = 0; x < l/2; ++x) {
103                 int y = buf[x]; buf[x] = buf[l-1-x]; buf[l-1-x] = y;
104         }
105         fputs(buf, fp);
106         return 0;
107 }
108
109 // an analogy to pileup_func() below
110 static int glt3_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
111 {
112         pu_data_t *d = (pu_data_t*)data;
113         bam_maqindel_ret_t *r = 0;
114         int rb, *proposed_indels = 0;
115         glf1_t *g;
116         glf3_t *g3;
117
118         if (d->fai == 0) {
119                 fprintf(stderr, "[glt3_func] reference sequence is required for generating GLT. Abort!\n");
120                 exit(1);
121         }
122         if (d->hash) { // only output a list of sites
123                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
124                 if (k == kh_end(d->hash)) return 0;
125                 proposed_indels = kh_val(d->hash, k);
126         }
127         g3 = glf3_init1();
128         if (d->fai && (int)tid != d->tid) {
129                 if (d->ref) { // then write the end mark
130                         g3->rtype = GLF3_RTYPE_END;
131                         glf3_write1(d->fp_glf, g3);
132                 }
133                 glf3_ref_write(d->fp_glf, d->h->target_name[tid], d->h->target_len[tid]); // write reference
134                 free(d->ref);
135                 d->ref = fai_fetch(d->fai, d->h->target_name[tid], &d->len);
136                 d->tid = tid;
137                 d->last_pos = 0;
138         }
139         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
140         g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
141         memcpy(g3, g, sizeof(glf1_t));
142         g3->rtype = GLF3_RTYPE_SUB;
143         g3->offset = pos - d->last_pos;
144         d->last_pos = pos;
145         glf3_write1(d->fp_glf, g3);
146     if (pos < d->len) {
147         int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
148                 if (proposed_indels)
149                         r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
150                 else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
151         }
152         if (r) { // then write indel line
153                 int het = 3 * n, min;
154                 min = het;
155                 if (min > r->gl[0]) min = r->gl[0];
156                 if (min > r->gl[1]) min = r->gl[1];
157                 g3->ref_base = 0;
158                 g3->rtype = GLF3_RTYPE_INDEL;
159                 memset(g3->lk, 0, 10);
160                 g3->lk[0] = r->gl[0] - min < 255? r->gl[0] - min : 255;
161                 g3->lk[1] = r->gl[1] - min < 255? r->gl[1] - min : 255;
162                 g3->lk[2] = het - min < 255? het - min : 255;
163                 g3->offset = 0;
164                 g3->indel_len[0] = r->indel1;
165                 g3->indel_len[1] = r->indel2;
166                 g3->min_lk = min < 255? min : 255;
167                 g3->max_len = (abs(r->indel1) > abs(r->indel2)? abs(r->indel1) : abs(r->indel2)) + 1;
168                 g3->indel_seq[0] = strdup(r->s[0]+1);
169                 g3->indel_seq[1] = strdup(r->s[1]+1);
170                 glf3_write1(d->fp_glf, g3);
171                 bam_maqindel_ret_destroy(r);
172         }
173         free(g);
174         glf3_destroy1(g3);
175         return 0;
176 }
177
178 static void pileup_seq(const bam_pileup1_t *p, int pos, int ref_len, const char *ref)
179 {
180         if (p->is_head) {
181                 putchar('^');
182                 putchar(p->b->core.qual > 93? 126 : p->b->core.qual + 33);
183         }
184         if (!p->is_del) {
185                 int j, rb, c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)];
186                 rb = (ref && pos < ref_len)? ref[pos] : 'N';
187                 if (c == '=' || toupper(c) == toupper(rb)) c = bam1_strand(p->b)? ',' : '.';
188                 else c = bam1_strand(p->b)? tolower(c) : toupper(c);
189                 putchar(c);
190                 if (p->indel > 0) {
191                         putchar('+'); printw(p->indel, stdout);
192                         for (j = 1; j <= p->indel; ++j) {
193                                 c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)];
194                                 putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
195                         }
196                 } else if (p->indel < 0) {
197                         printw(p->indel, stdout);
198                         for (j = 1; j <= -p->indel; ++j) {
199                                 c = (ref && (int)pos+j < ref_len)? ref[pos+j] : 'N';
200                                 putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
201                         }
202                 }
203         } else putchar(p->is_refskip? (bam1_strand(p->b)? '<' : '>') : '*');
204         if (p->is_tail) putchar('$');
205 }
206
207 static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
208 {
209         pu_data_t *d = (pu_data_t*)data;
210         bam_maqindel_ret_t *r = 0;
211         int i, rb, rms_mapq = -1, *proposed_indels = 0;
212         uint64_t rms_aux;
213         uint32_t cns = 0;
214
215         // if GLF is required, suppress -c completely
216         if (d->format & BAM_PLF_GLF) return glt3_func(tid, pos, n, pu, data);
217         // if d->hash is initialized, only output the sites in the hash table
218         if (d->hash) {
219                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
220                 if (k == kh_end(d->hash)) return 0;
221                 proposed_indels = kh_val(d->hash, k);
222         }
223         // update d->ref if necessary
224         if (d->fai && (int)tid != d->tid) {
225                 free(d->ref);
226                 d->ref = faidx_fetch_seq(d->fai, d->h->target_name[tid], 0, 0x7fffffff, &d->len);
227                 d->tid = tid;
228         }
229         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
230         // when the indel-only mode is asked for, return if no reads mapped with indels
231         if (d->format & BAM_PLF_INDEL_ONLY) {
232                 for (i = 0; i < n; ++i)
233                         if (pu[i].indel != 0) break;
234                 if (i == n) return 0;
235         }
236         // call the consensus and indel
237         if (d->format & BAM_PLF_CNS) { // call consensus
238                 if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE)) { // use a random base or the 1st base as the consensus call
239                         const bam_pileup1_t *p = (d->format & BAM_PLF_1STBASE)? pu : pu + (int)(drand48() * n);
240                         int q = bam1_qual(p->b)[p->qpos];
241                         int mapQ = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
242                         uint32_t b = bam1_seqi(bam1_seq(p->b), p->qpos);
243                         cns = b<<28 | 0xf<<24 | mapQ<<16 | q<<8;
244                 } else if (d->format & BAM_PLF_ALLBASE) { // collapse all bases
245                         uint64_t rmsQ = 0;
246                         uint32_t b = 0;
247                         for (i = 0; i < n; ++i) {
248                                 const bam_pileup1_t *p = pu + i;
249                                 int q = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
250                                 b |= bam1_seqi(bam1_seq(p->b), p->qpos);
251                                 rmsQ += q * q;
252                         }
253                         rmsQ = (uint64_t)(sqrt((double)rmsQ / n) + .499);
254                         cns = b<<28 | 0xf<<24 | rmsQ<<16 | 60<<8;
255                 } else {
256                         glf1_t *g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
257                         cns = g->depth == 0? (0xfu<<28 | 0xf<<24) : glf2cns(g, (int)(d->c->q_r + .499));
258                         free(g);
259                 }
260         }
261     if ((d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)) && d->ref && pos < d->len) { // call indels
262         int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
263         if (proposed_indels) // the first element gives the size of the array
264             r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
265         else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
266         }
267         // when only variant sites are asked for, test if the site is a variant
268         if ((d->format & BAM_PLF_CNS) && (d->format & BAM_PLF_VAR_ONLY)) {
269                 if (!(bam_nt16_table[rb] != 15 && cns>>28 != 15 && cns>>28 != bam_nt16_table[rb])) { // not a SNP
270                         if (!(r && (r->gt == 2 || strcmp(r->s[r->gt], "*")))) { // not an indel
271                                 if (r) bam_maqindel_ret_destroy(r);
272                                 return 0;
273                         }
274                 }
275         }
276         // print the first 3 columns
277         fputs(d->h->target_name[tid], stdout); putchar('\t');
278         printw(pos+1, stdout); putchar('\t'); putchar(rb); putchar('\t');
279         // print consensus information if required
280         if (d->format & BAM_PLF_CNS) {
281                 putchar(bam_nt16_rev_table[cns>>28]); putchar('\t');
282                 printw(cns>>8&0xff, stdout); putchar('\t');
283                 printw(cns&0xff, stdout); putchar('\t');
284                 printw(cns>>16&0xff, stdout); putchar('\t');
285         }
286         // print pileup sequences
287         printw(n, stdout); putchar('\t');
288         rms_aux = 0; // we need to recalculate rms_mapq when -c is not flagged on the command line
289         for (i = 0; i < n; ++i) {
290                 const bam_pileup1_t *p = pu + i;
291                 int tmp = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
292                 rms_aux += tmp * tmp;
293                 pileup_seq(p, pos, d->len, d->ref);
294         }
295         // finalize rms_mapq
296         rms_aux = (uint64_t)(sqrt((double)rms_aux / n) + .499);
297         if (rms_mapq < 0) rms_mapq = rms_aux;
298         putchar('\t');
299         // print quality
300         for (i = 0; i < n; ++i) {
301                 const bam_pileup1_t *p = pu + i;
302                 int c = bam1_qual(p->b)[p->qpos] + 33;
303                 if (c > 126) c = 126;
304                 putchar(c);
305         }
306         if (d->format & BAM_PLF_2ND) { // print 2nd calls and qualities
307                 const unsigned char *q;
308                 putchar('\t');
309                 for (i = 0; i < n; ++i) {
310                         const bam_pileup1_t *p = pu + i;
311                         q = bam_aux_get(p->b, "E2");
312                         putchar(q? q[p->qpos + 1] : 'N');
313                 }
314                 putchar('\t');
315                 for (i = 0; i < n; ++i) {
316                         const bam_pileup1_t *p = pu + i;
317                         q = bam_aux_get(p->b, "U2");
318                         putchar(q? q[p->qpos + 1] : '!');
319                 }
320         }
321         // print mapping quality if -s is flagged on the command line
322         if (d->format & BAM_PLF_SIMPLE) {
323                 putchar('\t');
324                 for (i = 0; i < n; ++i) {
325                         int c = pu[i].b->core.qual + 33;
326                         if (c > 126) c = 126;
327                         putchar(c);
328                 }
329         }
330         // print read position
331         if (d->format & BAM_PLF_READPOS) {
332                 putchar('\t');
333                 for (i = 0; i < n; ++i) {
334                         int x = pu[i].qpos;
335                         int l = pu[i].b->core.l_qseq;
336                         printw(x < l/2? x+1 : -((l-1)-x+1), stdout); putchar(',');
337                 }
338         }
339         putchar('\n');
340         // print the indel line if r has been calculated. This only happens if:
341         // a) -c or -i are flagged, AND b) the reference sequence is available
342         if (r) {
343                 printf("%s\t%d\t*\t", d->h->target_name[tid], pos + 1);
344                 if (r->gt < 2) printf("%s/%s\t", r->s[r->gt], r->s[r->gt]);
345                 else printf("%s/%s\t", r->s[0], r->s[1]);
346                 printf("%d\t%d\t", r->q_cns, r->q_ref);
347                 printf("%d\t%d\t", rms_mapq, n);
348                 printf("%s\t%s\t", r->s[0], r->s[1]);
349                 //printf("%d\t%d\t", r->gl[0], r->gl[1]);
350                 printf("%d\t%d\t%d\t", r->cnt1, r->cnt2, r->cnt_anti);
351                 printf("%d\t%d\n", r->cnt_ref, r->cnt_ambi);
352                 bam_maqindel_ret_destroy(r);
353         }
354         return 0;
355 }
356
357 int bam_pileup(int argc, char *argv[])
358 {
359         int c, is_SAM = 0;
360         char *fn_list = 0, *fn_fa = 0, *fn_pos = 0;
361         pu_data_t *d = (pu_data_t*)calloc(1, sizeof(pu_data_t));
362     d->max_depth = 1024; d->tid = -1; d->mask = BAM_DEF_MASK; d->min_baseQ = 13;
363         d->c = bam_maqcns_init();
364         d->c->errmod = BAM_ERRMOD_MAQ2; // change the default model
365         d->ido = bam_maqindel_opt_init();
366         while ((c = getopt(argc, argv, "st:f:cT:N:r:l:d:im:gI:G:vM:S2aR:PAQ:C:B")) >= 0) {
367                 switch (c) {
368                 case 'Q': d->c->min_baseQ = atoi(optarg); break;
369                 case 'C': d->capQ_thres = atoi(optarg); break;
370                 case 'B': d->format |= BAM_PLF_NOBAQ; break;
371                 case 'a': d->c->errmod = BAM_ERRMOD_SOAP; break;
372                 case 'A': d->c->errmod = BAM_ERRMOD_MAQ; break;
373                 case 's': d->format |= BAM_PLF_SIMPLE; break;
374                 case 't': fn_list = strdup(optarg); break;
375                 case 'l': fn_pos = strdup(optarg); break;
376                 case 'f': fn_fa = strdup(optarg); break;
377                 case 'T': d->c->theta = atof(optarg); break;
378                 case 'N': d->c->n_hap = atoi(optarg); break;
379                 case 'r': d->c->het_rate = atof(optarg); d->ido->r_snp = d->c->het_rate; break;
380                 case 'M': d->c->cap_mapQ = atoi(optarg); break;
381                 case 'd': d->max_depth = atoi(optarg); break;
382                 case 'c': d->format |= BAM_PLF_CNS; break;
383                 case 'i': d->format |= BAM_PLF_INDEL_ONLY; break;
384                 case 'v': d->format |= BAM_PLF_VAR_ONLY; break;
385                 case 'm': d->mask = strtol(optarg, 0, 0); break;
386                 case 'g': d->format |= BAM_PLF_GLF; break;
387                 case '2': d->format |= BAM_PLF_2ND; break;
388                 case 'P': d->format |= BAM_PLF_READPOS; break;
389                 case 'I': d->ido->q_indel = atoi(optarg); break;
390                 case 'G': d->ido->r_indel = atof(optarg); break;
391                 case 'S': is_SAM = 1; break;
392                 case 'R':
393                         if (strcmp(optarg, "random") == 0) d->format |= BAM_PLF_RANBASE;
394                         else if (strcmp(optarg, "first") == 0) d->format |= BAM_PLF_1STBASE;
395                         else if (strcmp(optarg, "all") == 0) d->format |= BAM_PLF_ALLBASE;
396                         else fprintf(stderr, "[bam_pileup] unrecognized -R\n");
397                         break;
398                 default: fprintf(stderr, "Unrecognizd option '-%c'.\n", c); return 1;
399                 }
400         }
401         if (d->c->errmod != BAM_ERRMOD_MAQ2) d->c->theta += 0.02;
402         if (d->c->theta > 1.0) d->c->theta = 1.0;
403         if (fn_list) is_SAM = 1;
404         if (optind == argc) {
405                 fprintf(stderr, "\n");
406                 fprintf(stderr, "Usage:  samtools pileup [options] <in.bam>|<in.sam>\n\n");
407                 fprintf(stderr, "Option: -s        simple (yet incomplete) pileup format\n");
408                 fprintf(stderr, "        -S        the input is in SAM\n");
409                 fprintf(stderr, "        -B        disable BAQ computation\n");
410                 fprintf(stderr, "        -A        use the original MAQ model for SNP calling (DEPRECATED)\n");
411                 fprintf(stderr, "        -2        output the 2nd best call and quality\n");
412                 fprintf(stderr, "        -i        only show lines/consensus with indels\n");
413                 fprintf(stderr, "        -Q INT    min base quality (possibly capped by BAQ) [%d]\n", d->c->min_baseQ);
414                 fprintf(stderr, "        -C INT    coefficient for adjusting mapQ of poor mappings [%d]\n", d->capQ_thres);
415                 fprintf(stderr, "        -m INT    filtering reads with bits in INT [0x%x]\n", d->mask);
416                 fprintf(stderr, "        -M INT    cap mapping quality at INT [%d]\n", d->c->cap_mapQ);
417         fprintf(stderr, "        -d INT    limit maximum depth for indels [%d]\n", d->max_depth);
418                 fprintf(stderr, "        -t FILE   list of reference sequences (force -S)\n");
419                 fprintf(stderr, "        -l FILE   list of sites at which pileup is output\n");
420                 fprintf(stderr, "        -f FILE   reference sequence in the FASTA format\n\n");
421                 fprintf(stderr, "        -c        compute the consensus sequence\n");
422                 fprintf(stderr, "        -v        print variants only (for -c)\n");
423                 fprintf(stderr, "        -g        output in the GLFv3 format (DEPRECATED)\n");
424                 fprintf(stderr, "        -T FLOAT  theta in maq consensus calling model (for -c) [%.4g]\n", d->c->theta);
425                 fprintf(stderr, "        -N INT    number of haplotypes in the sample (for -c) [%d]\n", d->c->n_hap);
426                 fprintf(stderr, "        -r FLOAT  prior of a difference between two haplotypes (for -c) [%.4g]\n", d->c->het_rate);
427                 fprintf(stderr, "        -G FLOAT  prior of an indel between two haplotypes (for -c) [%.4g]\n", d->ido->r_indel);
428                 fprintf(stderr, "        -I INT    phred prob. of an indel in sequencing/prep. (for -c) [%d]\n", d->ido->q_indel);
429                 fprintf(stderr, "\n");
430                 free(fn_list); free(fn_fa); free(d);
431                 return 1;
432         }
433         if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE|BAM_PLF_ALLBASE)) d->format |= BAM_PLF_CNS;
434         if (fn_fa) d->fai = fai_load(fn_fa);
435         if (d->format & (BAM_PLF_CNS|BAM_PLF_GLF)) bam_maqcns_prepare(d->c); // consensus calling
436         if (d->format & BAM_PLF_GLF) { // for glf output
437                 glf3_header_t *h;
438                 h = glf3_header_init();
439                 d->fp_glf = bgzf_fdopen(fileno(stdout), "w");
440                 glf3_header_write(d->fp_glf, h);
441                 glf3_header_destroy(h);
442         }
443         if (d->fai == 0 && (d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)))
444                 fprintf(stderr, "[bam_pileup] indels will not be called when -f is absent.\n");
445         if (fn_fa && is_SAM && fn_list == 0) fn_list = samfaipath(fn_fa);
446
447         {
448                 samfile_t *fp;
449                 fp = is_SAM? samopen(argv[optind], "r", fn_list) : samopen(argv[optind], "rb", 0);
450                 if (fp == 0 || fp->header == 0) {
451                         fprintf(stderr, "[bam_pileup] fail to read the header: non-exisiting file or wrong format.\n");
452                         return 1;
453                 }
454                 d->h = fp->header;
455                 if (fn_pos) d->hash = load_pos(fn_pos, d->h);
456                 { // run pileup
457                         extern int bam_prob_realn(bam1_t *b, const char *ref);
458                         extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
459                         bam1_t *b;
460                         int ret, tid, pos, n_plp;
461                         bam_plp_t iter;
462                         const bam_pileup1_t *plp;
463                         b = bam_init1();
464                         iter = bam_plp_init(0, 0);
465                         bam_plp_set_mask(iter, d->mask);
466                         while ((ret = samread(fp, b)) >= 0) {
467                                 int skip = 0;
468                                 // update d->ref if necessary
469                                 if (d->fai && (int)b->core.tid != d->tid) {
470                                         free(d->ref);
471                                         d->ref = faidx_fetch_seq(d->fai, d->h->target_name[b->core.tid], 0, 0x7fffffff, &d->len);
472                                         d->tid = b->core.tid;
473                                 }
474                                 if (d->ref && (d->format&BAM_PLF_CNS) && !(d->format&BAM_PLF_NOBAQ)) bam_prob_realn(b, d->ref);
475                                 if (d->ref && (d->format&BAM_PLF_CNS) && d->capQ_thres > 10) {
476                                         int q = bam_cap_mapQ(b, d->ref, d->capQ_thres);
477                                         if (q < 0) skip = 1;
478                                         else if (b->core.qual > q) b->core.qual = q;
479                                 } else if (b->core.flag&BAM_FUNMAP) skip = 1;
480                                 else if ((d->format&BAM_PLF_CNS) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
481                                 if (skip) continue;
482                                 bam_plp_push(iter, b);
483                                 while ((plp = bam_plp_next(iter, &tid, &pos, &n_plp)) != 0)
484                                         pileup_func(tid, pos, n_plp, plp, d);
485                         }
486                         bam_plp_push(iter, 0);
487                         while ((plp = bam_plp_next(iter, &tid, &pos, &n_plp)) != 0)
488                                 pileup_func(tid, pos, n_plp, plp, d);
489                         bam_plp_destroy(iter);
490                         bam_destroy1(b);
491                 }
492                 samclose(fp); // d->h will be destroyed here
493         }
494
495         // free
496         if (d->format & BAM_PLF_GLF) bgzf_close(d->fp_glf);
497         if (fn_pos) { // free the hash table
498                 khint_t k;
499                 for (k = kh_begin(d->hash); k < kh_end(d->hash); ++k)
500                         if (kh_exist(d->hash, k)) free(kh_val(d->hash, k));
501                 kh_destroy(64, d->hash);
502         }
503         free(fn_pos); free(fn_list); free(fn_fa);
504         if (d->fai) fai_destroy(d->fai);
505         bam_maqcns_destroy(d->c);
506         free(d->ido); free(d->ref); free(d);
507         return 0;
508 }
509
510 /***********
511  * mpileup *
512  ***********/
513
514 #include <assert.h>
515 #include "bam2bcf.h"
516 #include "sample.h"
517
518 #define MPLP_GLF   0x10
519 #define MPLP_NO_COMP 0x20
520 #define MPLP_NO_ORPHAN 0x40
521 #define MPLP_REALN   0x80
522
523 typedef struct {
524         int max_mq, min_mq, flag, min_baseQ, capQ_thres;
525         char *reg, *fn_pos;
526         faidx_t *fai;
527         kh_64_t *hash;
528 } mplp_conf_t;
529
530 typedef struct {
531         bamFile fp;
532         bam_iter_t iter;
533         int min_mq, flag, ref_id, capQ_thres;
534         char *ref;
535 } mplp_aux_t;
536
537 typedef struct {
538         int n;
539         int *n_plp, *m_plp;
540         bam_pileup1_t **plp;
541 } mplp_pileup_t;
542
543 static int mplp_func(void *data, bam1_t *b)
544 {
545         extern int bam_realn(bam1_t *b, const char *ref);
546         extern int bam_prob_realn(bam1_t *b, const char *ref);
547         extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
548         mplp_aux_t *ma = (mplp_aux_t*)data;
549         int ret, skip = 0;
550         do {
551                 int has_ref = (ma->ref && ma->ref_id == b->core.tid)? 1 : 0;
552                 ret = ma->iter? bam_iter_read(ma->fp, ma->iter, b) : bam_read1(ma->fp, b);
553                 if (ret < 0) break;
554                 skip = 0;
555                 if (has_ref && (ma->flag&MPLP_REALN)) bam_prob_realn(b, ma->ref);
556                 if (has_ref && ma->capQ_thres > 10) {
557                         int q = bam_cap_mapQ(b, ma->ref, ma->capQ_thres);
558                         if (q < 0) skip = 1;
559                         else if (b->core.qual > q) b->core.qual = q;
560                 } else if (b->core.flag&BAM_FUNMAP) skip = 1;
561                 else if (b->core.qual < ma->min_mq) skip = 1; 
562                 else if ((ma->flag&MPLP_NO_ORPHAN) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
563         } while (skip);
564         return ret;
565 }
566
567 static void group_smpl(mplp_pileup_t *m, bam_sample_t *sm, kstring_t *buf,
568                                            int n, char *const*fn, int *n_plp, const bam_pileup1_t **plp)
569 {
570         int i, j;
571         memset(m->n_plp, 0, m->n * sizeof(int));
572         for (i = 0; i < n; ++i) {
573                 for (j = 0; j < n_plp[i]; ++j) {
574                         const bam_pileup1_t *p = plp[i] + j;
575                         uint8_t *q;
576                         int id = -1;
577                         q = bam_aux_get(p->b, "RG");
578                         if (q) id = bam_smpl_rg2smid(sm, fn[i], (char*)q+1, buf);
579                         if (id < 0) id = bam_smpl_rg2smid(sm, fn[i], 0, buf);
580                         assert(id >= 0 && id < m->n);
581                         if (m->n_plp[id] == m->m_plp[id]) {
582                                 m->m_plp[id] = m->m_plp[id]? m->m_plp[id]<<1 : 8;
583                                 m->plp[id] = realloc(m->plp[id], sizeof(bam_pileup1_t) * m->m_plp[id]);
584                         }
585                         m->plp[id][m->n_plp[id]++] = *p;
586                 }
587         }
588 }
589
590 static int mpileup(mplp_conf_t *conf, int n, char **fn)
591 {
592         mplp_aux_t **data;
593         int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid;
594         const bam_pileup1_t **plp;
595         bam_mplp_t iter;
596         bam_header_t *h = 0;
597         char *ref;
598         khash_t(64) *hash = 0;
599
600         bcf_callaux_t *bca = 0;
601         bcf_callret1_t *bcr = 0;
602         bcf_call_t bc;
603         bcf_t *bp = 0;
604         bcf_hdr_t *bh = 0;
605
606         bam_sample_t *sm = 0;
607         kstring_t buf;
608         mplp_pileup_t gplp;
609
610         memset(&gplp, 0, sizeof(mplp_pileup_t));
611         memset(&buf, 0, sizeof(kstring_t));
612         memset(&bc, 0, sizeof(bcf_call_t));
613         data = calloc(n, sizeof(void*));
614         plp = calloc(n, sizeof(void*));
615         n_plp = calloc(n, sizeof(int*));
616         sm = bam_smpl_init();
617
618         // read the header and initialize data
619         for (i = 0; i < n; ++i) {
620                 bam_header_t *h_tmp;
621                 data[i] = calloc(1, sizeof(mplp_aux_t));
622                 data[i]->min_mq = conf->min_mq;
623                 data[i]->flag = conf->flag;
624                 data[i]->capQ_thres = conf->capQ_thres;
625                 data[i]->fp = strcmp(fn[i], "-") == 0? bam_dopen(fileno(stdin), "r") : bam_open(fn[i], "r");
626                 h_tmp = bam_header_read(data[i]->fp);
627                 bam_smpl_add(sm, fn[i], h_tmp->text);
628                 if (conf->reg) {
629                         int beg, end;
630                         bam_index_t *idx;
631                         idx = bam_index_load(fn[i]);
632                         if (idx == 0) {
633                                 fprintf(stderr, "[%s] fail to load index for %d-th input.\n", __func__, i+1);
634                                 exit(1);
635                         }
636                         if (bam_parse_region(h_tmp, conf->reg, &tid, &beg, &end) < 0) {
637                                 fprintf(stderr, "[%s] malformatted region or wrong seqname for %d-th input.\n", __func__, i+1);
638                                 exit(1);
639                         }
640                         if (i == 0) beg0 = beg, end0 = end;
641                         data[i]->iter = bam_iter_query(idx, tid, beg, end);
642                         bam_index_destroy(idx);
643                 }
644                 if (i == 0) h = h_tmp;
645                 else {
646                         // FIXME: to check consistency
647                         bam_header_destroy(h_tmp);
648                 }
649         }
650         gplp.n = sm->n;
651         gplp.n_plp = calloc(sm->n, sizeof(int));
652         gplp.m_plp = calloc(sm->n, sizeof(int));
653         gplp.plp = calloc(sm->n, sizeof(void*));
654
655         fprintf(stderr, "[%s] %d samples in %d input files\n", __func__, sm->n, n);
656         if (conf->fn_pos) hash = load_pos(conf->fn_pos, h);
657         // write the VCF header
658         if (conf->flag & MPLP_GLF) {
659                 kstring_t s;
660                 bh = calloc(1, sizeof(bcf_hdr_t));
661                 s.l = s.m = 0; s.s = 0;
662                 bp = bcf_open("-", (conf->flag&MPLP_NO_COMP)? "wu" : "w");
663                 for (i = 0; i < h->n_targets; ++i) {
664                         kputs(h->target_name[i], &s);
665                         kputc('\0', &s);
666                 }
667                 bh->l_nm = s.l;
668                 bh->name = malloc(s.l);
669                 memcpy(bh->name, s.s, s.l);
670                 s.l = 0;
671                 for (i = 0; i < sm->n; ++i) {
672                         kputs(sm->smpl[i], &s); kputc('\0', &s);
673                 }
674                 bh->l_smpl = s.l;
675                 bh->sname = malloc(s.l);
676                 memcpy(bh->sname, s.s, s.l);
677                 bh->l_txt = 0;
678                 free(s.s);
679                 bcf_hdr_sync(bh);
680                 bcf_hdr_write(bp, bh);
681                 bca = bcf_call_init(-1., conf->min_baseQ);
682                 bcr = calloc(sm->n, sizeof(bcf_callret1_t));
683         }
684         ref_tid = -1; ref = 0;
685         iter = bam_mplp_init(n, mplp_func, (void**)data);
686         while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
687                 if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested
688                 if (hash) {
689                         khint_t k;
690                         k = kh_get(64, hash, (uint64_t)tid<<32 | pos);
691                         if (k == kh_end(hash)) continue;
692                 }
693                 if (tid != ref_tid) {
694                         free(ref); ref = 0;
695                         if (conf->fai) ref = fai_fetch(conf->fai, h->target_name[tid], &ref_len);
696                         for (i = 0; i < n; ++i) data[i]->ref = ref, data[i]->ref_id = tid;
697                         ref_tid = tid;
698                 }
699                 if (conf->flag & MPLP_GLF) {
700                         int _ref0, ref16;
701                         bcf1_t *b = calloc(1, sizeof(bcf1_t));
702                         group_smpl(&gplp, sm, &buf, n, fn, n_plp, plp);
703                         _ref0 = (ref && pos < ref_len)? ref[pos] : 'N';
704                         ref16 = bam_nt16_table[_ref0];
705                         for (i = 0; i < gplp.n; ++i)
706                                 bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], ref16, bca, bcr + i);
707                         bcf_call_combine(gplp.n, bcr, ref16, &bc);
708                         bcf_call2bcf(tid, pos, &bc, b);
709                         bcf_write(bp, bh, b);
710                         bcf_destroy(b);
711                 } else {
712                         printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
713                         for (i = 0; i < n; ++i) {
714                                 int j;
715                                 printf("\t%d\t", n_plp[i]);
716                                 if (n_plp[i] == 0) printf("*\t*");
717                                 else {
718                                         for (j = 0; j < n_plp[i]; ++j)
719                                                 pileup_seq(plp[i] + j, pos, ref_len, ref);
720                                         putchar('\t');
721                                         for (j = 0; j < n_plp[i]; ++j) {
722                                                 const bam_pileup1_t *p = plp[i] + j;
723                                                 int c = bam1_qual(p->b)[p->qpos] + 33;
724                                                 if (c > 126) c = 126;
725                                                 putchar(c);
726                                         }
727                                 }
728                         }
729                         putchar('\n');
730                 }
731         }
732
733         bcf_close(bp);
734         bam_smpl_destroy(sm); free(buf.s);
735         for (i = 0; i < gplp.n; ++i) free(gplp.plp[i]);
736         free(gplp.plp); free(gplp.n_plp); free(gplp.m_plp);
737         if (hash) { // free the hash table
738                 khint_t k;
739                 for (k = kh_begin(hash); k < kh_end(hash); ++k)
740                         if (kh_exist(hash, k)) free(kh_val(hash, k));
741                 kh_destroy(64, hash);
742         }
743         bcf_hdr_destroy(bh); bcf_call_destroy(bca); free(bc.PL); free(bcr);
744         bam_mplp_destroy(iter);
745         bam_header_destroy(h);
746         for (i = 0; i < n; ++i) {
747                 bam_close(data[i]->fp);
748                 if (data[i]->iter) bam_iter_destroy(data[i]->iter);
749                 free(data[i]);
750         }
751         free(data); free(plp); free(ref); free(n_plp);
752         return 0;
753 }
754
755 int bam_mpileup(int argc, char *argv[])
756 {
757         int c;
758         mplp_conf_t mplp;
759         memset(&mplp, 0, sizeof(mplp_conf_t));
760         mplp.max_mq = 60;
761         mplp.min_baseQ = 13;
762         mplp.capQ_thres = 0;
763         while ((c = getopt(argc, argv, "gf:r:l:M:q:Q:uaORC:")) >= 0) {
764                 switch (c) {
765                 case 'f':
766                         mplp.fai = fai_load(optarg);
767                         if (mplp.fai == 0) return 1;
768                         break;
769                 case 'r': mplp.reg = strdup(optarg); break;
770                 case 'l': mplp.fn_pos = strdup(optarg); break;
771                 case 'g': mplp.flag |= MPLP_GLF; break;
772                 case 'u': mplp.flag |= MPLP_NO_COMP | MPLP_GLF; break;
773                 case 'a': mplp.flag |= MPLP_NO_ORPHAN | MPLP_REALN; break;
774                 case 'O': mplp.flag |= MPLP_NO_ORPHAN; break;
775                 case 'R': mplp.flag |= MPLP_REALN; break;
776                 case 'C': mplp.capQ_thres = atoi(optarg); break;
777                 case 'M': mplp.max_mq = atoi(optarg); break;
778                 case 'q': mplp.min_mq = atoi(optarg); break;
779                 case 'Q': mplp.min_baseQ = atoi(optarg); break;
780                 }
781         }
782         if (argc == 1) {
783                 fprintf(stderr, "\n");
784                 fprintf(stderr, "Usage:   samtools mpileup [options] in1.bam [in2.bam [...]]\n\n");
785                 fprintf(stderr, "Options: -f FILE     reference sequence file [null]\n");
786                 fprintf(stderr, "         -r STR      region in which pileup is generated [null]\n");
787                 fprintf(stderr, "         -l FILE     list of positions (format: chr pos) [null]\n");
788                 fprintf(stderr, "         -M INT      cap mapping quality at INT [%d]\n", mplp.max_mq);
789                 fprintf(stderr, "         -Q INT      min base quality [%d]\n", mplp.min_baseQ);
790                 fprintf(stderr, "         -q INT      filter out alignment with MQ smaller than INT [%d]\n", mplp.min_mq);
791                 fprintf(stderr, "         -g          generate BCF output\n");
792                 fprintf(stderr, "         -u          do not compress BCF output\n");
793                 fprintf(stderr, "\n");
794                 fprintf(stderr, "Notes: Assuming diploid individuals.\n\n");
795                 return 1;
796         }
797         mpileup(&mplp, argc - optind, argv + optind);
798         free(mplp.reg);
799         if (mplp.fai) fai_destroy(mplp.fai);
800         return 0;
801 }