Uploaded samtools_0.1.16-1_amd64.changes.
[samtools.git] / bam_plcmd.c
1 #include <math.h>
2 #include <stdio.h>
3 #include <unistd.h>
4 #include <ctype.h>
5 #include <string.h>
6 #include <errno.h>
7 #include "sam.h"
8 #include "faidx.h"
9 #include "bam_maqcns.h"
10 #include "khash.h"
11 #include "glf.h"
12 #include "kstring.h"
13
14 typedef int *indel_list_t;
15 KHASH_MAP_INIT_INT64(64, indel_list_t)
16
17 #define BAM_PLF_SIMPLE     0x01
18 #define BAM_PLF_CNS        0x02
19 #define BAM_PLF_INDEL_ONLY 0x04
20 #define BAM_PLF_GLF        0x08
21 #define BAM_PLF_VAR_ONLY   0x10
22 #define BAM_PLF_2ND        0x20
23 #define BAM_PLF_RANBASE    0x40
24 #define BAM_PLF_1STBASE    0x80
25 #define BAM_PLF_ALLBASE    0x100
26 #define BAM_PLF_READPOS    0x200
27 #define BAM_PLF_NOBAQ      0x400
28
29 typedef struct {
30         bam_header_t *h;
31         bam_maqcns_t *c;
32         bam_maqindel_opt_t *ido;
33         faidx_t *fai;
34         khash_t(64) *hash;
35         uint32_t format;
36         int tid, len, last_pos;
37         int mask;
38         int capQ_thres, min_baseQ;
39     int max_depth;  // for indel calling, ignore reads with the depth too high. 0 for unlimited
40         char *ref;
41         glfFile fp_glf; // for glf output only
42 } pu_data_t;
43
44 char **__bam_get_lines(const char *fn, int *_n);
45 void bam_init_header_hash(bam_header_t *header);
46 int32_t bam_get_tid(const bam_header_t *header, const char *seq_name);
47
48 static khash_t(64) *load_pos(const char *fn, bam_header_t *h)
49 {
50         char **list;
51         int i, j, n, *fields, max_fields;
52         khash_t(64) *hash;
53         bam_init_header_hash(h);
54         list = __bam_get_lines(fn, &n);
55         hash = kh_init(64);
56         max_fields = 0; fields = 0;
57         for (i = 0; i < n; ++i) {
58                 char *str = list[i];
59                 int chr, n_fields, ret;
60                 khint_t k;
61                 uint64_t x;
62                 n_fields = ksplit_core(str, 0, &max_fields, &fields);
63                 if (n_fields < 2) continue;
64                 chr = bam_get_tid(h, str + fields[0]);
65                 if (chr < 0) {
66                         fprintf(stderr, "[load_pos] unknown reference sequence name: %s\n", str + fields[0]);
67                         continue;
68                 }
69                 x = (uint64_t)chr << 32 | (atoi(str + fields[1]) - 1);
70                 k = kh_put(64, hash, x, &ret);
71                 if (ret == 0) {
72                         fprintf(stderr, "[load_pos] position %s:%s has been loaded.\n", str+fields[0], str+fields[1]);
73                         continue;
74                 }
75                 kh_val(hash, k) = 0;
76                 if (n_fields > 2) {
77                         // count
78                         for (j = 2; j < n_fields; ++j) {
79                                 char *s = str + fields[j];
80                                 if ((*s != '+' && *s != '-') || !isdigit(s[1])) break;
81                         }
82                         if (j > 2) { // update kh_val()
83                                 int *q, y, z;
84                                 q = kh_val(hash, k) = (int*)calloc(j - 1, sizeof(int));
85                                 q[0] = j - 2; z = j; y = 1;
86                                 for (j = 2; j < z; ++j)
87                                         q[y++] = atoi(str + fields[j]);
88                         }
89                 }
90                 free(str);
91         }
92         free(list); free(fields);
93         return hash;
94 }
95
96 static inline int printw(int c, FILE *fp)
97 {
98         char buf[16];
99         int l, x;
100         if (c == 0) return fputc('0', fp);
101         for (l = 0, x = c < 0? -c : c; x > 0; x /= 10) buf[l++] = x%10 + '0';
102         if (c < 0) buf[l++] = '-';
103         buf[l] = 0;
104         for (x = 0; x < l/2; ++x) {
105                 int y = buf[x]; buf[x] = buf[l-1-x]; buf[l-1-x] = y;
106         }
107         fputs(buf, fp);
108         return 0;
109 }
110
111 // an analogy to pileup_func() below
112 static int glt3_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
113 {
114         pu_data_t *d = (pu_data_t*)data;
115         bam_maqindel_ret_t *r = 0;
116         int rb, *proposed_indels = 0;
117         glf1_t *g;
118         glf3_t *g3;
119
120         if (d->fai == 0) {
121                 fprintf(stderr, "[glt3_func] reference sequence is required for generating GLT. Abort!\n");
122                 exit(1);
123         }
124         if (d->hash) { // only output a list of sites
125                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
126                 if (k == kh_end(d->hash)) return 0;
127                 proposed_indels = kh_val(d->hash, k);
128         }
129         g3 = glf3_init1();
130         if (d->fai && (int)tid != d->tid) {
131                 if (d->ref) { // then write the end mark
132                         g3->rtype = GLF3_RTYPE_END;
133                         glf3_write1(d->fp_glf, g3);
134                 }
135                 glf3_ref_write(d->fp_glf, d->h->target_name[tid], d->h->target_len[tid]); // write reference
136                 free(d->ref);
137                 d->ref = fai_fetch(d->fai, d->h->target_name[tid], &d->len);
138                 d->tid = tid;
139                 d->last_pos = 0;
140         }
141         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
142         g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
143         memcpy(g3, g, sizeof(glf1_t));
144         g3->rtype = GLF3_RTYPE_SUB;
145         g3->offset = pos - d->last_pos;
146         d->last_pos = pos;
147         glf3_write1(d->fp_glf, g3);
148     if (pos < d->len) {
149         int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
150                 if (proposed_indels)
151                         r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
152                 else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
153         }
154         if (r) { // then write indel line
155                 int het = 3 * n, min;
156                 min = het;
157                 if (min > r->gl[0]) min = r->gl[0];
158                 if (min > r->gl[1]) min = r->gl[1];
159                 g3->ref_base = 0;
160                 g3->rtype = GLF3_RTYPE_INDEL;
161                 memset(g3->lk, 0, 10);
162                 g3->lk[0] = r->gl[0] - min < 255? r->gl[0] - min : 255;
163                 g3->lk[1] = r->gl[1] - min < 255? r->gl[1] - min : 255;
164                 g3->lk[2] = het - min < 255? het - min : 255;
165                 g3->offset = 0;
166                 g3->indel_len[0] = r->indel1;
167                 g3->indel_len[1] = r->indel2;
168                 g3->min_lk = min < 255? min : 255;
169                 g3->max_len = (abs(r->indel1) > abs(r->indel2)? abs(r->indel1) : abs(r->indel2)) + 1;
170                 g3->indel_seq[0] = strdup(r->s[0]+1);
171                 g3->indel_seq[1] = strdup(r->s[1]+1);
172                 glf3_write1(d->fp_glf, g3);
173                 bam_maqindel_ret_destroy(r);
174         }
175         free(g);
176         glf3_destroy1(g3);
177         return 0;
178 }
179
180 static inline void pileup_seq(const bam_pileup1_t *p, int pos, int ref_len, const char *ref)
181 {
182         int j;
183         if (p->is_head) {
184                 putchar('^');
185                 putchar(p->b->core.qual > 93? 126 : p->b->core.qual + 33);
186         }
187         if (!p->is_del) {
188                 int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)];
189                 if (ref) {
190                         int rb = pos < ref_len? ref[pos] : 'N';
191                         if (c == '=' || bam_nt16_table[c] == bam_nt16_table[rb]) c = bam1_strand(p->b)? ',' : '.';
192                         else c = bam1_strand(p->b)? tolower(c) : toupper(c);
193                 } else {
194                         if (c == '=') c = bam1_strand(p->b)? ',' : '.';
195                         else c = bam1_strand(p->b)? tolower(c) : toupper(c);
196                 }
197                 putchar(c);
198         } else putchar(p->is_refskip? (bam1_strand(p->b)? '<' : '>') : '*');
199         if (p->indel > 0) {
200                 putchar('+'); printw(p->indel, stdout);
201                 for (j = 1; j <= p->indel; ++j) {
202                         int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)];
203                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
204                 }
205         } else if (p->indel < 0) {
206                 printw(p->indel, stdout);
207                 for (j = 1; j <= -p->indel; ++j) {
208                         int c = (ref && (int)pos+j < ref_len)? ref[pos+j] : 'N';
209                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
210                 }
211         }
212         if (p->is_tail) putchar('$');
213 }
214
215 static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
216 {
217         pu_data_t *d = (pu_data_t*)data;
218         bam_maqindel_ret_t *r = 0;
219         int i, rb, rms_mapq = -1, *proposed_indels = 0;
220         uint64_t rms_aux;
221         uint32_t cns = 0;
222
223         // if GLF is required, suppress -c completely
224         if (d->format & BAM_PLF_GLF) return glt3_func(tid, pos, n, pu, data);
225         // if d->hash is initialized, only output the sites in the hash table
226         if (d->hash) {
227                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
228                 if (k == kh_end(d->hash)) return 0;
229                 proposed_indels = kh_val(d->hash, k);
230         }
231         // update d->ref if necessary
232         if (d->fai && (int)tid != d->tid) {
233                 free(d->ref);
234                 d->ref = faidx_fetch_seq(d->fai, d->h->target_name[tid], 0, 0x7fffffff, &d->len);
235                 d->tid = tid;
236         }
237         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
238         // when the indel-only mode is asked for, return if no reads mapped with indels
239         if (d->format & BAM_PLF_INDEL_ONLY) {
240                 for (i = 0; i < n; ++i)
241                         if (pu[i].indel != 0) break;
242                 if (i == n) return 0;
243         }
244         // call the consensus and indel
245         if (d->format & BAM_PLF_CNS) { // call consensus
246                 if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE)) { // use a random base or the 1st base as the consensus call
247                         const bam_pileup1_t *p = (d->format & BAM_PLF_1STBASE)? pu : pu + (int)(drand48() * n);
248                         int q = bam1_qual(p->b)[p->qpos];
249                         int mapQ = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
250                         uint32_t b = bam1_seqi(bam1_seq(p->b), p->qpos);
251                         cns = b<<28 | 0xf<<24 | mapQ<<16 | q<<8;
252                 } else if (d->format & BAM_PLF_ALLBASE) { // collapse all bases
253                         uint64_t rmsQ = 0;
254                         uint32_t b = 0;
255                         for (i = 0; i < n; ++i) {
256                                 const bam_pileup1_t *p = pu + i;
257                                 int q = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
258                                 b |= bam1_seqi(bam1_seq(p->b), p->qpos);
259                                 rmsQ += q * q;
260                         }
261                         rmsQ = (uint64_t)(sqrt((double)rmsQ / n) + .499);
262                         cns = b<<28 | 0xf<<24 | rmsQ<<16 | 60<<8;
263                 } else {
264                         glf1_t *g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
265                         cns = g->depth == 0? (0xfu<<28 | 0xf<<24) : glf2cns(g, (int)(d->c->q_r + .499));
266                         free(g);
267                 }
268         }
269     if ((d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)) && d->ref && pos < d->len) { // call indels
270         int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
271         if (proposed_indels) // the first element gives the size of the array
272             r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
273         else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
274         }
275         // when only variant sites are asked for, test if the site is a variant
276         if ((d->format & BAM_PLF_CNS) && (d->format & BAM_PLF_VAR_ONLY)) {
277                 if (!(bam_nt16_table[rb] != 15 && cns>>28 != 15 && cns>>28 != bam_nt16_table[rb])) { // not a SNP
278                         if (!(r && (r->gt == 2 || strcmp(r->s[r->gt], "*")))) { // not an indel
279                                 if (r) bam_maqindel_ret_destroy(r);
280                                 return 0;
281                         }
282                 }
283         }
284         // print the first 3 columns
285         fputs(d->h->target_name[tid], stdout); putchar('\t');
286         printw(pos+1, stdout); putchar('\t'); putchar(rb); putchar('\t');
287         // print consensus information if required
288         if (d->format & BAM_PLF_CNS) {
289                 putchar(bam_nt16_rev_table[cns>>28]); putchar('\t');
290                 printw(cns>>8&0xff, stdout); putchar('\t');
291                 printw(cns&0xff, stdout); putchar('\t');
292                 printw(cns>>16&0xff, stdout); putchar('\t');
293         }
294         // print pileup sequences
295         printw(n, stdout); putchar('\t');
296         for (i = 0; i < n; ++i)
297                 pileup_seq(pu + i, pos, d->len, d->ref);
298         // finalize rms_mapq
299         if (d->format & BAM_PLF_CNS) {
300                 for (i = rms_aux = 0; i < n; ++i) {
301                         const bam_pileup1_t *p = pu + i;
302                         int tmp = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
303                         rms_aux += tmp * tmp;
304                 }
305                 rms_aux = (uint64_t)(sqrt((double)rms_aux / n) + .499);
306                 if (rms_mapq < 0) rms_mapq = rms_aux;
307         }
308         putchar('\t');
309         // print quality
310         for (i = 0; i < n; ++i) {
311                 const bam_pileup1_t *p = pu + i;
312                 int c = bam1_qual(p->b)[p->qpos] + 33;
313                 if (c > 126) c = 126;
314                 putchar(c);
315         }
316         if (d->format & BAM_PLF_2ND) { // print 2nd calls and qualities
317                 const unsigned char *q;
318                 putchar('\t');
319                 for (i = 0; i < n; ++i) {
320                         const bam_pileup1_t *p = pu + i;
321                         q = bam_aux_get(p->b, "E2");
322                         putchar(q? q[p->qpos + 1] : 'N');
323                 }
324                 putchar('\t');
325                 for (i = 0; i < n; ++i) {
326                         const bam_pileup1_t *p = pu + i;
327                         q = bam_aux_get(p->b, "U2");
328                         putchar(q? q[p->qpos + 1] : '!');
329                 }
330         }
331         // print mapping quality if -s is flagged on the command line
332         if (d->format & BAM_PLF_SIMPLE) {
333                 putchar('\t');
334                 for (i = 0; i < n; ++i) {
335                         int c = pu[i].b->core.qual + 33;
336                         if (c > 126) c = 126;
337                         putchar(c);
338                 }
339         }
340         // print read position
341         if (d->format & BAM_PLF_READPOS) {
342                 putchar('\t');
343                 for (i = 0; i < n; ++i) {
344                         int x = pu[i].qpos;
345                         int l = pu[i].b->core.l_qseq;
346                         printw(x < l/2? x+1 : -((l-1)-x+1), stdout); putchar(',');
347                 }
348         }
349         putchar('\n');
350         // print the indel line if r has been calculated. This only happens if:
351         // a) -c or -i are flagged, AND b) the reference sequence is available
352         if (r) {
353                 printf("%s\t%d\t*\t", d->h->target_name[tid], pos + 1);
354                 if (r->gt < 2) printf("%s/%s\t", r->s[r->gt], r->s[r->gt]);
355                 else printf("%s/%s\t", r->s[0], r->s[1]);
356                 printf("%d\t%d\t", r->q_cns, r->q_ref);
357                 printf("%d\t%d\t", rms_mapq, n);
358                 printf("%s\t%s\t", r->s[0], r->s[1]);
359                 //printf("%d\t%d\t", r->gl[0], r->gl[1]);
360                 printf("%d\t%d\t%d\t", r->cnt1, r->cnt2, r->cnt_anti);
361                 printf("%d\t%d\n", r->cnt_ref, r->cnt_ambi);
362                 bam_maqindel_ret_destroy(r);
363         }
364         return 0;
365 }
366
367 int bam_pileup(int argc, char *argv[])
368 {
369         int c, is_SAM = 0;
370         char *fn_list = 0, *fn_fa = 0, *fn_pos = 0;
371         pu_data_t *d = (pu_data_t*)calloc(1, sizeof(pu_data_t));
372     d->max_depth = 1024; d->tid = -1; d->mask = BAM_DEF_MASK; d->min_baseQ = 13;
373         d->c = bam_maqcns_init();
374         d->c->errmod = BAM_ERRMOD_MAQ2; // change the default model
375         d->ido = bam_maqindel_opt_init();
376         while ((c = getopt(argc, argv, "st:f:cT:N:r:l:d:im:gI:G:vM:S2aR:PAQ:C:B")) >= 0) {
377                 switch (c) {
378                 case 'Q': d->c->min_baseQ = atoi(optarg); break;
379                 case 'C': d->capQ_thres = atoi(optarg); break;
380                 case 'B': d->format |= BAM_PLF_NOBAQ; break;
381                 case 'a': d->c->errmod = BAM_ERRMOD_SOAP; break;
382                 case 'A': d->c->errmod = BAM_ERRMOD_MAQ; break;
383                 case 's': d->format |= BAM_PLF_SIMPLE; break;
384                 case 't': fn_list = strdup(optarg); break;
385                 case 'l': fn_pos = strdup(optarg); break;
386                 case 'f': fn_fa = strdup(optarg); break;
387                 case 'T': d->c->theta = atof(optarg); break;
388                 case 'N': d->c->n_hap = atoi(optarg); break;
389                 case 'r': d->c->het_rate = atof(optarg); d->ido->r_snp = d->c->het_rate; break;
390                 case 'M': d->c->cap_mapQ = atoi(optarg); break;
391                 case 'd': d->max_depth = atoi(optarg); break;
392                 case 'c': d->format |= BAM_PLF_CNS; break;
393                 case 'i': d->format |= BAM_PLF_INDEL_ONLY; break;
394                 case 'v': d->format |= BAM_PLF_VAR_ONLY; break;
395                 case 'm': d->mask = strtol(optarg, 0, 0); break;
396                 case 'g': d->format |= BAM_PLF_GLF; break;
397                 case '2': d->format |= BAM_PLF_2ND; break;
398                 case 'P': d->format |= BAM_PLF_READPOS; break;
399                 case 'I': d->ido->q_indel = atoi(optarg); break;
400                 case 'G': d->ido->r_indel = atof(optarg); break;
401                 case 'S': is_SAM = 1; break;
402                 case 'R':
403                         if (strcmp(optarg, "random") == 0) d->format |= BAM_PLF_RANBASE;
404                         else if (strcmp(optarg, "first") == 0) d->format |= BAM_PLF_1STBASE;
405                         else if (strcmp(optarg, "all") == 0) d->format |= BAM_PLF_ALLBASE;
406                         else fprintf(stderr, "[bam_pileup] unrecognized -R\n");
407                         break;
408                 default: fprintf(stderr, "Unrecognizd option '-%c'.\n", c); return 1;
409                 }
410         }
411         if (d->c->errmod != BAM_ERRMOD_MAQ2) d->c->theta += 0.02;
412         if (d->c->theta > 1.0) d->c->theta = 1.0;
413         if (fn_list) is_SAM = 1;
414         if (optind == argc) {
415                 fprintf(stderr, "\n");
416                 fprintf(stderr, "Usage:  samtools pileup [options] <in.bam>|<in.sam>\n\n");
417                 fprintf(stderr, "Option: -s        simple (yet incomplete) pileup format\n");
418                 fprintf(stderr, "        -S        the input is in SAM\n");
419                 fprintf(stderr, "        -B        disable BAQ computation\n");
420                 fprintf(stderr, "        -A        use the original MAQ model for SNP calling (DEPRECATED)\n");
421                 fprintf(stderr, "        -2        output the 2nd best call and quality\n");
422                 fprintf(stderr, "        -i        only show lines/consensus with indels\n");
423                 fprintf(stderr, "        -Q INT    min base quality (possibly capped by BAQ) [%d]\n", d->c->min_baseQ);
424                 fprintf(stderr, "        -C INT    coefficient for adjusting mapQ of poor mappings [%d]\n", d->capQ_thres);
425                 fprintf(stderr, "        -m INT    filtering reads with bits in INT [0x%x]\n", d->mask);
426                 fprintf(stderr, "        -M INT    cap mapping quality at INT [%d]\n", d->c->cap_mapQ);
427         fprintf(stderr, "        -d INT    limit maximum depth for indels [%d]\n", d->max_depth);
428                 fprintf(stderr, "        -t FILE   list of reference sequences (force -S)\n");
429                 fprintf(stderr, "        -l FILE   list of sites at which pileup is output\n");
430                 fprintf(stderr, "        -f FILE   reference sequence in the FASTA format\n\n");
431                 fprintf(stderr, "        -c        compute the consensus sequence\n");
432                 fprintf(stderr, "        -v        print variants only (for -c)\n");
433                 fprintf(stderr, "        -g        output in the GLFv3 format (DEPRECATED)\n");
434                 fprintf(stderr, "        -T FLOAT  theta in maq consensus calling model (for -c) [%.4g]\n", d->c->theta);
435                 fprintf(stderr, "        -N INT    number of haplotypes in the sample (for -c) [%d]\n", d->c->n_hap);
436                 fprintf(stderr, "        -r FLOAT  prior of a difference between two haplotypes (for -c) [%.4g]\n", d->c->het_rate);
437                 fprintf(stderr, "        -G FLOAT  prior of an indel between two haplotypes (for -c) [%.4g]\n", d->ido->r_indel);
438                 fprintf(stderr, "        -I INT    phred prob. of an indel in sequencing/prep. (for -c) [%d]\n", d->ido->q_indel);
439                 fprintf(stderr, "\n");
440                 fprintf(stderr, "Warning: Please use the `mpileup' command instead `pileup'. `Pileup' is deprecated!\n\n");
441                 free(fn_list); free(fn_fa); free(d);
442                 return 1;
443         }
444         if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE|BAM_PLF_ALLBASE)) d->format |= BAM_PLF_CNS;
445         if (fn_fa) d->fai = fai_load(fn_fa);
446         if (d->format & (BAM_PLF_CNS|BAM_PLF_GLF)) bam_maqcns_prepare(d->c); // consensus calling
447         if (d->format & BAM_PLF_GLF) { // for glf output
448                 glf3_header_t *h;
449                 h = glf3_header_init();
450                 d->fp_glf = bgzf_fdopen(fileno(stdout), "w");
451                 glf3_header_write(d->fp_glf, h);
452                 glf3_header_destroy(h);
453         }
454         if (d->fai == 0 && (d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)))
455                 fprintf(stderr, "[bam_pileup] indels will not be called when -f is absent.\n");
456         if (fn_fa && is_SAM && fn_list == 0) fn_list = samfaipath(fn_fa);
457
458         {
459                 samfile_t *fp;
460                 fp = is_SAM? samopen(argv[optind], "r", fn_list) : samopen(argv[optind], "rb", 0);
461                 if (fp == 0 || fp->header == 0) {
462                         fprintf(stderr, "[bam_pileup] fail to read the header: non-exisiting file or wrong format.\n");
463                         return 1;
464                 }
465                 d->h = fp->header;
466                 if (fn_pos) d->hash = load_pos(fn_pos, d->h);
467                 { // run pileup
468                         extern int bam_prob_realn(bam1_t *b, const char *ref);
469                         extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
470                         bam1_t *b;
471                         int ret, tid, pos, n_plp;
472                         bam_plp_t iter;
473                         const bam_pileup1_t *plp;
474                         b = bam_init1();
475                         iter = bam_plp_init(0, 0);
476                         bam_plp_set_mask(iter, d->mask);
477                         while ((ret = samread(fp, b)) >= 0) {
478                                 int skip = 0;
479                                 if ((int)b->core.tid < 0) break;
480                                 // update d->ref if necessary
481                                 if (d->fai && (int)b->core.tid != d->tid) {
482                                         free(d->ref);
483                                         d->ref = faidx_fetch_seq(d->fai, d->h->target_name[b->core.tid], 0, 0x7fffffff, &d->len);
484                                         d->tid = b->core.tid;
485                                 }
486                                 if (d->ref && (d->format&BAM_PLF_CNS) && !(d->format&BAM_PLF_NOBAQ)) bam_prob_realn(b, d->ref);
487                                 if (d->ref && (d->format&BAM_PLF_CNS) && d->capQ_thres > 10) {
488                                         int q = bam_cap_mapQ(b, d->ref, d->capQ_thres);
489                                         if (q < 0) skip = 1;
490                                         else if (b->core.qual > q) b->core.qual = q;
491                                 } else if (b->core.flag&BAM_FUNMAP) skip = 1;
492                                 else if ((d->format&BAM_PLF_CNS) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
493                                 if (skip) continue;
494                                 bam_plp_push(iter, b);
495                                 while ((plp = bam_plp_next(iter, &tid, &pos, &n_plp)) != 0)
496                                         pileup_func(tid, pos, n_plp, plp, d);
497                         }
498                         bam_plp_push(iter, 0);
499                         while ((plp = bam_plp_next(iter, &tid, &pos, &n_plp)) != 0)
500                                 pileup_func(tid, pos, n_plp, plp, d);
501                         bam_plp_destroy(iter);
502                         bam_destroy1(b);
503                 }
504                 samclose(fp); // d->h will be destroyed here
505         }
506
507         // free
508         if (d->format & BAM_PLF_GLF) bgzf_close(d->fp_glf);
509         if (fn_pos) { // free the hash table
510                 khint_t k;
511                 for (k = kh_begin(d->hash); k < kh_end(d->hash); ++k)
512                         if (kh_exist(d->hash, k)) free(kh_val(d->hash, k));
513                 kh_destroy(64, d->hash);
514         }
515         free(fn_pos); free(fn_list); free(fn_fa);
516         if (d->fai) fai_destroy(d->fai);
517         bam_maqcns_destroy(d->c);
518         free(d->ido); free(d->ref); free(d);
519         return 0;
520 }
521
522 /***********
523  * mpileup *
524  ***********/
525
526 #include <assert.h>
527 #include "bam2bcf.h"
528 #include "sample.h"
529
530 #define MPLP_GLF   0x10
531 #define MPLP_NO_COMP 0x20
532 #define MPLP_NO_ORPHAN 0x40
533 #define MPLP_REALN   0x80
534 #define MPLP_FMT_DP 0x100
535 #define MPLP_FMT_SP 0x200
536 #define MPLP_NO_INDEL 0x400
537 #define MPLP_EXT_BAQ 0x800
538 #define MPLP_ILLUMINA13 0x1000
539
540 void *bed_read(const char *fn);
541 void bed_destroy(void *_h);
542 int bed_overlap(const void *_h, const char *chr, int beg, int end);
543
544 typedef struct {
545         int max_mq, min_mq, flag, min_baseQ, capQ_thres, max_depth, max_indel_depth;
546         int openQ, extQ, tandemQ, min_support; // for indels
547         double min_frac; // for indels
548         char *reg, *pl_list;
549         faidx_t *fai;
550         void *bed, *rghash;
551 } mplp_conf_t;
552
553 typedef struct {
554         bamFile fp;
555         bam_iter_t iter;
556         bam_header_t *h;
557         int ref_id;
558         char *ref;
559         const mplp_conf_t *conf;
560 } mplp_aux_t;
561
562 typedef struct {
563         int n;
564         int *n_plp, *m_plp;
565         bam_pileup1_t **plp;
566 } mplp_pileup_t;
567
568 static int mplp_func(void *data, bam1_t *b)
569 {
570         extern int bam_realn(bam1_t *b, const char *ref);
571         extern int bam_prob_realn_core(bam1_t *b, const char *ref, int);
572         extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
573         mplp_aux_t *ma = (mplp_aux_t*)data;
574         int ret, skip = 0;
575         do {
576                 int has_ref;
577                 ret = ma->iter? bam_iter_read(ma->fp, ma->iter, b) : bam_read1(ma->fp, b);
578                 if (ret < 0) break;
579                 if (b->core.tid < 0 || (b->core.flag&BAM_FUNMAP)) { // exclude unmapped reads
580                         skip = 1;
581                         continue;
582                 }
583                 if (ma->conf->bed) { // test overlap
584                         skip = !bed_overlap(ma->conf->bed, ma->h->target_name[b->core.tid], b->core.pos, bam_calend(&b->core, bam1_cigar(b)));
585                         if (skip) continue;
586                 }
587                 if (ma->conf->rghash) { // exclude read groups
588                         uint8_t *rg = bam_aux_get(b, "RG");
589                         skip = (rg && bcf_str2id(ma->conf->rghash, (const char*)(rg+1)) >= 0);
590                         if (skip) continue;
591                 }
592                 if (ma->conf->flag & MPLP_ILLUMINA13) {
593                         int i;
594                         uint8_t *qual = bam1_qual(b);
595                         for (i = 0; i < b->core.l_qseq; ++i)
596                                 qual[i] = qual[i] > 31? qual[i] - 31 : 0;
597                 }
598                 has_ref = (ma->ref && ma->ref_id == b->core.tid)? 1 : 0;
599                 skip = 0;
600                 if (has_ref && (ma->conf->flag&MPLP_REALN)) bam_prob_realn_core(b, ma->ref, (ma->conf->flag & MPLP_EXT_BAQ)? 3 : 1);
601                 if (has_ref && ma->conf->capQ_thres > 10) {
602                         int q = bam_cap_mapQ(b, ma->ref, ma->conf->capQ_thres);
603                         if (q < 0) skip = 1;
604                         else if (b->core.qual > q) b->core.qual = q;
605                 }
606                 else if (b->core.qual < ma->conf->min_mq) skip = 1; 
607                 else if ((ma->conf->flag&MPLP_NO_ORPHAN) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
608         } while (skip);
609         return ret;
610 }
611
612 static void group_smpl(mplp_pileup_t *m, bam_sample_t *sm, kstring_t *buf,
613                                            int n, char *const*fn, int *n_plp, const bam_pileup1_t **plp)
614 {
615         int i, j;
616         memset(m->n_plp, 0, m->n * sizeof(int));
617         for (i = 0; i < n; ++i) {
618                 for (j = 0; j < n_plp[i]; ++j) {
619                         const bam_pileup1_t *p = plp[i] + j;
620                         uint8_t *q;
621                         int id = -1;
622                         q = bam_aux_get(p->b, "RG");
623                         if (q) id = bam_smpl_rg2smid(sm, fn[i], (char*)q+1, buf);
624                         if (id < 0) id = bam_smpl_rg2smid(sm, fn[i], 0, buf);
625                         if (id < 0 || id >= m->n) {
626                                 assert(q); // otherwise a bug
627                                 fprintf(stderr, "[%s] Read group %s used in file %s but absent from the header or an alignment missing read group.\n", __func__, (char*)q+1, fn[i]);
628                                 exit(1);
629                         }
630                         if (m->n_plp[id] == m->m_plp[id]) {
631                                 m->m_plp[id] = m->m_plp[id]? m->m_plp[id]<<1 : 8;
632                                 m->plp[id] = realloc(m->plp[id], sizeof(bam_pileup1_t) * m->m_plp[id]);
633                         }
634                         m->plp[id][m->n_plp[id]++] = *p;
635                 }
636         }
637 }
638
639 static int mpileup(mplp_conf_t *conf, int n, char **fn)
640 {
641         extern void *bcf_call_add_rg(void *rghash, const char *hdtext, const char *list);
642         extern void bcf_call_del_rghash(void *rghash);
643         mplp_aux_t **data;
644         int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid, max_depth, max_indel_depth;
645         const bam_pileup1_t **plp;
646         bam_mplp_t iter;
647         bam_header_t *h = 0;
648         char *ref;
649         void *rghash = 0;
650
651         bcf_callaux_t *bca = 0;
652         bcf_callret1_t *bcr = 0;
653         bcf_call_t bc;
654         bcf_t *bp = 0;
655         bcf_hdr_t *bh = 0;
656
657         bam_sample_t *sm = 0;
658         kstring_t buf;
659         mplp_pileup_t gplp;
660
661         memset(&gplp, 0, sizeof(mplp_pileup_t));
662         memset(&buf, 0, sizeof(kstring_t));
663         memset(&bc, 0, sizeof(bcf_call_t));
664         data = calloc(n, sizeof(void*));
665         plp = calloc(n, sizeof(void*));
666         n_plp = calloc(n, sizeof(int*));
667         sm = bam_smpl_init();
668
669         // read the header and initialize data
670         for (i = 0; i < n; ++i) {
671                 bam_header_t *h_tmp;
672                 data[i] = calloc(1, sizeof(mplp_aux_t));
673                 data[i]->fp = strcmp(fn[i], "-") == 0? bam_dopen(fileno(stdin), "r") : bam_open(fn[i], "r");
674                 data[i]->conf = conf;
675                 h_tmp = bam_header_read(data[i]->fp);
676                 data[i]->h = i? h : h_tmp; // for i==0, "h" has not been set yet
677                 bam_smpl_add(sm, fn[i], h_tmp->text);
678                 rghash = bcf_call_add_rg(rghash, h_tmp->text, conf->pl_list);
679                 if (conf->reg) {
680                         int beg, end;
681                         bam_index_t *idx;
682                         idx = bam_index_load(fn[i]);
683                         if (idx == 0) {
684                                 fprintf(stderr, "[%s] fail to load index for %d-th input.\n", __func__, i+1);
685                                 exit(1);
686                         }
687                         if (bam_parse_region(h_tmp, conf->reg, &tid, &beg, &end) < 0) {
688                                 fprintf(stderr, "[%s] malformatted region or wrong seqname for %d-th input.\n", __func__, i+1);
689                                 exit(1);
690                         }
691                         if (i == 0) beg0 = beg, end0 = end;
692                         data[i]->iter = bam_iter_query(idx, tid, beg, end);
693                         bam_index_destroy(idx);
694                 }
695                 if (i == 0) h = h_tmp;
696                 else {
697                         // FIXME: to check consistency
698                         bam_header_destroy(h_tmp);
699                 }
700         }
701         gplp.n = sm->n;
702         gplp.n_plp = calloc(sm->n, sizeof(int));
703         gplp.m_plp = calloc(sm->n, sizeof(int));
704         gplp.plp = calloc(sm->n, sizeof(void*));
705
706         fprintf(stderr, "[%s] %d samples in %d input files\n", __func__, sm->n, n);
707         // write the VCF header
708         if (conf->flag & MPLP_GLF) {
709                 kstring_t s;
710                 bh = calloc(1, sizeof(bcf_hdr_t));
711                 s.l = s.m = 0; s.s = 0;
712                 bp = bcf_open("-", (conf->flag&MPLP_NO_COMP)? "wu" : "w");
713                 for (i = 0; i < h->n_targets; ++i) {
714                         kputs(h->target_name[i], &s);
715                         kputc('\0', &s);
716                 }
717                 bh->l_nm = s.l;
718                 bh->name = malloc(s.l);
719                 memcpy(bh->name, s.s, s.l);
720                 s.l = 0;
721                 for (i = 0; i < sm->n; ++i) {
722                         kputs(sm->smpl[i], &s); kputc('\0', &s);
723                 }
724                 bh->l_smpl = s.l;
725                 bh->sname = malloc(s.l);
726                 memcpy(bh->sname, s.s, s.l);
727                 bh->txt = malloc(strlen(BAM_VERSION) + 64);
728                 bh->l_txt = 1 + sprintf(bh->txt, "##samtoolsVersion=%s\n", BAM_VERSION);
729                 free(s.s);
730                 bcf_hdr_sync(bh);
731                 bcf_hdr_write(bp, bh);
732                 bca = bcf_call_init(-1., conf->min_baseQ);
733                 bcr = calloc(sm->n, sizeof(bcf_callret1_t));
734                 bca->rghash = rghash;
735                 bca->openQ = conf->openQ, bca->extQ = conf->extQ, bca->tandemQ = conf->tandemQ;
736                 bca->min_frac = conf->min_frac;
737                 bca->min_support = conf->min_support;
738         }
739         ref_tid = -1; ref = 0;
740         iter = bam_mplp_init(n, mplp_func, (void**)data);
741         max_depth = conf->max_depth;
742         if (max_depth * sm->n > 1<<20)
743                 fprintf(stderr, "(%s) Max depth is above 1M. Potential memory hog!\n", __func__);
744         if (max_depth * sm->n < 8000) {
745                 max_depth = 8000 / sm->n;
746                 fprintf(stderr, "<%s> Set max per-file depth to %d\n", __func__, max_depth);
747         }
748         max_indel_depth = conf->max_indel_depth * sm->n;
749         bam_mplp_set_maxcnt(iter, max_depth);
750         while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
751                 if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested
752                 if (conf->bed && tid >= 0 && !bed_overlap(conf->bed, h->target_name[tid], pos, pos+1)) continue;
753                 if (tid != ref_tid) {
754                         free(ref); ref = 0;
755                         if (conf->fai) ref = fai_fetch(conf->fai, h->target_name[tid], &ref_len);
756                         for (i = 0; i < n; ++i) data[i]->ref = ref, data[i]->ref_id = tid;
757                         ref_tid = tid;
758                 }
759                 if (conf->flag & MPLP_GLF) {
760                         int total_depth, _ref0, ref16;
761                         bcf1_t *b = calloc(1, sizeof(bcf1_t));
762                         for (i = total_depth = 0; i < n; ++i) total_depth += n_plp[i];
763                         group_smpl(&gplp, sm, &buf, n, fn, n_plp, plp);
764                         _ref0 = (ref && pos < ref_len)? ref[pos] : 'N';
765                         ref16 = bam_nt16_table[_ref0];
766                         for (i = 0; i < gplp.n; ++i)
767                                 bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], ref16, bca, bcr + i);
768                         bcf_call_combine(gplp.n, bcr, ref16, &bc);
769                         bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
770                                                  (conf->flag&MPLP_FMT_SP), 0, 0);
771                         bcf_write(bp, bh, b);
772                         bcf_destroy(b);
773                         // call indels
774                         if (!(conf->flag&MPLP_NO_INDEL) && total_depth < max_indel_depth && bcf_call_gap_prep(gplp.n, gplp.n_plp, gplp.plp, pos, bca, ref, rghash) >= 0) {
775                                 for (i = 0; i < gplp.n; ++i)
776                                         bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], -1, bca, bcr + i);
777                                 if (bcf_call_combine(gplp.n, bcr, -1, &bc) >= 0) {
778                                         b = calloc(1, sizeof(bcf1_t));
779                                         bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
780                                                                  (conf->flag&MPLP_FMT_SP), bca, ref);
781                                         bcf_write(bp, bh, b);
782                                         bcf_destroy(b);
783                                 }
784                         }
785                 } else {
786                         printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
787                         for (i = 0; i < n; ++i) {
788                                 int j;
789                                 printf("\t%d\t", n_plp[i]);
790                                 if (n_plp[i] == 0) printf("*\t*");
791                                 else {
792                                         for (j = 0; j < n_plp[i]; ++j)
793                                                 pileup_seq(plp[i] + j, pos, ref_len, ref);
794                                         putchar('\t');
795                                         for (j = 0; j < n_plp[i]; ++j) {
796                                                 const bam_pileup1_t *p = plp[i] + j;
797                                                 int c = bam1_qual(p->b)[p->qpos] + 33;
798                                                 if (c > 126) c = 126;
799                                                 putchar(c);
800                                         }
801                                 }
802                         }
803                         putchar('\n');
804                 }
805         }
806
807         bcf_close(bp);
808         bam_smpl_destroy(sm); free(buf.s);
809         for (i = 0; i < gplp.n; ++i) free(gplp.plp[i]);
810         free(gplp.plp); free(gplp.n_plp); free(gplp.m_plp);
811         bcf_call_del_rghash(rghash);
812         bcf_hdr_destroy(bh); bcf_call_destroy(bca); free(bc.PL); free(bcr);
813         bam_mplp_destroy(iter);
814         bam_header_destroy(h);
815         for (i = 0; i < n; ++i) {
816                 bam_close(data[i]->fp);
817                 if (data[i]->iter) bam_iter_destroy(data[i]->iter);
818                 free(data[i]);
819         }
820         free(data); free(plp); free(ref); free(n_plp);
821         return 0;
822 }
823
824 #define MAX_PATH_LEN 1024
825 static int read_file_list(const char *file_list,int *n,char **argv[])
826 {
827     char buf[MAX_PATH_LEN];
828     int len, nfiles;
829     char **files;
830
831     FILE *fh = fopen(file_list,"r");
832     if ( !fh )
833     {
834         fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
835         return 1;
836     }
837
838     // Speed is not an issue here, determine the number of files by reading the file twice
839     nfiles = 0;
840     while ( fgets(buf,MAX_PATH_LEN,fh) ) nfiles++;
841
842     if ( fseek(fh, 0L, SEEK_SET) )
843     {
844         fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
845         return 1;
846     }
847
848     files = calloc(nfiles,sizeof(char*));
849     nfiles = 0;
850     while ( fgets(buf,MAX_PATH_LEN,fh) ) 
851     {
852         len = strlen(buf);
853         while ( len>0 && isspace(buf[len-1]) ) len--;
854         if ( !len ) continue;
855
856         files[nfiles] = malloc(sizeof(char)*(len+1)); 
857         strncpy(files[nfiles],buf,len);
858         files[nfiles][len] = 0;
859         nfiles++;
860     }
861     fclose(fh);
862     if ( !nfiles )
863     {
864         fprintf(stderr,"No files read from %s\n", file_list);
865         return 1;
866     }
867     *argv = files;
868     *n    = nfiles;
869     return 0;
870 }
871 #undef MAX_PATH_LEN
872
873 int bam_mpileup(int argc, char *argv[])
874 {
875         int c;
876     const char *file_list = NULL;
877     char **fn = NULL;
878     int nfiles = 0, use_orphan = 0;
879         mplp_conf_t mplp;
880         memset(&mplp, 0, sizeof(mplp_conf_t));
881         mplp.max_mq = 60;
882         mplp.min_baseQ = 13;
883         mplp.capQ_thres = 0;
884         mplp.max_depth = 250; mplp.max_indel_depth = 250;
885         mplp.openQ = 40; mplp.extQ = 20; mplp.tandemQ = 100;
886         mplp.min_frac = 0.002; mplp.min_support = 1;
887         mplp.flag = MPLP_NO_ORPHAN | MPLP_REALN;
888         while ((c = getopt(argc, argv, "Agf:r:l:M:q:Q:uaRC:BDSd:L:b:P:o:e:h:Im:F:EG:6")) >= 0) {
889                 switch (c) {
890                 case 'f':
891                         mplp.fai = fai_load(optarg);
892                         if (mplp.fai == 0) return 1;
893                         break;
894                 case 'd': mplp.max_depth = atoi(optarg); break;
895                 case 'r': mplp.reg = strdup(optarg); break;
896                 case 'l': mplp.bed = bed_read(optarg); break;
897                 case 'P': mplp.pl_list = strdup(optarg); break;
898                 case 'g': mplp.flag |= MPLP_GLF; break;
899                 case 'u': mplp.flag |= MPLP_NO_COMP | MPLP_GLF; break;
900                 case 'a': mplp.flag |= MPLP_NO_ORPHAN | MPLP_REALN; break;
901                 case 'B': mplp.flag &= ~MPLP_REALN; break;
902                 case 'R': mplp.flag |= MPLP_REALN; break;
903                 case 'D': mplp.flag |= MPLP_FMT_DP; break;
904                 case 'S': mplp.flag |= MPLP_FMT_SP; break;
905                 case 'I': mplp.flag |= MPLP_NO_INDEL; break;
906                 case 'E': mplp.flag |= MPLP_EXT_BAQ; break;
907                 case '6': mplp.flag |= MPLP_ILLUMINA13; break;
908                 case 'C': mplp.capQ_thres = atoi(optarg); break;
909                 case 'M': mplp.max_mq = atoi(optarg); break;
910                 case 'q': mplp.min_mq = atoi(optarg); break;
911                 case 'Q': mplp.min_baseQ = atoi(optarg); break;
912         case 'b': file_list = optarg; break;
913                 case 'o': mplp.openQ = atoi(optarg); break;
914                 case 'e': mplp.extQ = atoi(optarg); break;
915                 case 'h': mplp.tandemQ = atoi(optarg); break;
916                 case 'A': use_orphan = 1; break;
917                 case 'F': mplp.min_frac = atof(optarg); break;
918                 case 'm': mplp.min_support = atoi(optarg); break;
919                 case 'L': mplp.max_indel_depth = atoi(optarg); break;
920                 case 'G': {
921                                 FILE *fp_rg;
922                                 char buf[1024];
923                                 mplp.rghash = bcf_str2id_init();
924                                 if ((fp_rg = fopen(optarg, "r")) == 0)
925                                         fprintf(stderr, "(%s) Fail to open file %s. Continue anyway.\n", __func__, optarg);
926                                 while (!feof(fp_rg) && fscanf(fp_rg, "%s", buf) > 0) // this is not a good style, but forgive me...
927                                         bcf_str2id_add(mplp.rghash, strdup(buf));
928                                 fclose(fp_rg);
929                         }
930                         break;
931                 }
932         }
933         if (use_orphan) mplp.flag &= ~MPLP_NO_ORPHAN;
934         if (argc == 1) {
935                 fprintf(stderr, "\n");
936                 fprintf(stderr, "Usage: samtools mpileup [options] in1.bam [in2.bam [...]]\n\n");
937                 fprintf(stderr, "Input options:\n\n");
938                 fprintf(stderr, "       -6           assume the quality is in the Illumina-1.3+ encoding\n");
939                 fprintf(stderr, "       -A           count anomalous read pairs\n");
940                 fprintf(stderr, "       -B           disable BAQ computation\n");
941                 fprintf(stderr, "       -b FILE      list of input BAM files [null]\n");
942                 fprintf(stderr, "       -d INT       max per-BAM depth to avoid excessive memory usage [%d]\n", mplp.max_depth);
943                 fprintf(stderr, "       -E           extended BAQ for higher sensitivity but lower specificity\n");
944                 fprintf(stderr, "       -f FILE      faidx indexed reference sequence file [null]\n");
945                 fprintf(stderr, "       -G FILE      exclude read groups listed in FILE [null]\n");
946                 fprintf(stderr, "       -l FILE      list of positions (chr pos) or regions (BED) [null]\n");
947                 fprintf(stderr, "       -M INT       cap mapping quality at INT [%d]\n", mplp.max_mq);
948                 fprintf(stderr, "       -r STR       region in which pileup is generated [null]\n");
949                 fprintf(stderr, "       -q INT       skip alignments with mapQ smaller than INT [%d]\n", mplp.min_mq);
950                 fprintf(stderr, "       -Q INT       skip bases with baseQ/BAQ smaller than INT [%d]\n", mplp.min_baseQ);
951                 fprintf(stderr, "\nOutput options:\n\n");
952                 fprintf(stderr, "       -D           output per-sample DP in BCF (require -g/-u)\n");
953                 fprintf(stderr, "       -g           generate BCF output (genotype likelihoods)\n");
954                 fprintf(stderr, "       -S           output per-sample strand bias P-value in BCF (require -g/-u)\n");
955                 fprintf(stderr, "       -u           generate uncompress BCF output\n");
956                 fprintf(stderr, "\nSNP/INDEL genotype likelihoods options (effective with `-g' or `-u'):\n\n");
957                 fprintf(stderr, "       -e INT       Phred-scaled gap extension seq error probability [%d]\n", mplp.extQ);
958                 fprintf(stderr, "       -F FLOAT     minimum fraction of gapped reads for candidates [%g]\n", mplp.min_frac);
959                 fprintf(stderr, "       -h INT       coefficient for homopolymer errors [%d]\n", mplp.tandemQ);
960                 fprintf(stderr, "       -I           do not perform indel calling\n");
961                 fprintf(stderr, "       -L INT       max per-sample depth for INDEL calling [%d]\n", mplp.max_indel_depth);
962                 fprintf(stderr, "       -m INT       minimum gapped reads for indel candidates [%d]\n", mplp.min_support);
963                 fprintf(stderr, "       -o INT       Phred-scaled gap open sequencing error probability [%d]\n", mplp.openQ);
964                 fprintf(stderr, "       -P STR       comma separated list of platforms for indels [all]\n");
965                 fprintf(stderr, "\n");
966                 fprintf(stderr, "Notes: Assuming diploid individuals.\n\n");
967                 return 1;
968         }
969     if (file_list) {
970         if ( read_file_list(file_list,&nfiles,&fn) ) return 1;
971         mpileup(&mplp,nfiles,fn);
972         for (c=0; c<nfiles; c++) free(fn[c]);
973         free(fn);
974     } else mpileup(&mplp, argc - optind, argv + optind);
975         if (mplp.rghash) bcf_str2id_thorough_destroy(mplp.rghash);
976         free(mplp.reg); free(mplp.pl_list);
977         if (mplp.fai) fai_destroy(mplp.fai);
978         if (mplp.bed) bed_destroy(mplp.bed);
979         return 0;
980 }