New upstream release.
[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 "kstring.h"
10
11 static inline int printw(int c, FILE *fp)
12 {
13         char buf[16];
14         int l, x;
15         if (c == 0) return fputc('0', fp);
16         for (l = 0, x = c < 0? -c : c; x > 0; x /= 10) buf[l++] = x%10 + '0';
17         if (c < 0) buf[l++] = '-';
18         buf[l] = 0;
19         for (x = 0; x < l/2; ++x) {
20                 int y = buf[x]; buf[x] = buf[l-1-x]; buf[l-1-x] = y;
21         }
22         fputs(buf, fp);
23         return 0;
24 }
25
26 static inline void pileup_seq(const bam_pileup1_t *p, int pos, int ref_len, const char *ref)
27 {
28         int j;
29         if (p->is_head) {
30                 putchar('^');
31                 putchar(p->b->core.qual > 93? 126 : p->b->core.qual + 33);
32         }
33         if (!p->is_del) {
34                 int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)];
35                 if (ref) {
36                         int rb = pos < ref_len? ref[pos] : 'N';
37                         if (c == '=' || bam_nt16_table[c] == bam_nt16_table[rb]) c = bam1_strand(p->b)? ',' : '.';
38                         else c = bam1_strand(p->b)? tolower(c) : toupper(c);
39                 } else {
40                         if (c == '=') c = bam1_strand(p->b)? ',' : '.';
41                         else c = bam1_strand(p->b)? tolower(c) : toupper(c);
42                 }
43                 putchar(c);
44         } else putchar(p->is_refskip? (bam1_strand(p->b)? '<' : '>') : '*');
45         if (p->indel > 0) {
46                 putchar('+'); printw(p->indel, stdout);
47                 for (j = 1; j <= p->indel; ++j) {
48                         int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)];
49                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
50                 }
51         } else if (p->indel < 0) {
52                 printw(p->indel, stdout);
53                 for (j = 1; j <= -p->indel; ++j) {
54                         int c = (ref && (int)pos+j < ref_len)? ref[pos+j] : 'N';
55                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
56                 }
57         }
58         if (p->is_tail) putchar('$');
59 }
60
61 #include <assert.h>
62 #include "bam2bcf.h"
63 #include "sample.h"
64
65 #define MPLP_GLF   0x10
66 #define MPLP_NO_COMP 0x20
67 #define MPLP_NO_ORPHAN 0x40
68 #define MPLP_REALN   0x80
69 #define MPLP_FMT_DP 0x100
70 #define MPLP_FMT_SP 0x200
71 #define MPLP_NO_INDEL 0x400
72 #define MPLP_EXT_BAQ 0x800
73 #define MPLP_ILLUMINA13 0x1000
74 #define MPLP_IGNORE_RG 0x2000
75 #define MPLP_PRINT_POS 0x4000
76 #define MPLP_PRINT_MAPQ 0x8000
77
78 void *bed_read(const char *fn);
79 void bed_destroy(void *_h);
80 int bed_overlap(const void *_h, const char *chr, int beg, int end);
81
82 typedef struct {
83         int max_mq, min_mq, flag, min_baseQ, capQ_thres, max_depth, max_indel_depth;
84         int openQ, extQ, tandemQ, min_support; // for indels
85         double min_frac; // for indels
86         char *reg, *pl_list;
87         faidx_t *fai;
88         void *bed, *rghash;
89 } mplp_conf_t;
90
91 typedef struct {
92         bamFile fp;
93         bam_iter_t iter;
94         bam_header_t *h;
95         int ref_id;
96         char *ref;
97         const mplp_conf_t *conf;
98 } mplp_aux_t;
99
100 typedef struct {
101         int n;
102         int *n_plp, *m_plp;
103         bam_pileup1_t **plp;
104 } mplp_pileup_t;
105
106 static int mplp_func(void *data, bam1_t *b)
107 {
108         extern int bam_realn(bam1_t *b, const char *ref);
109         extern int bam_prob_realn_core(bam1_t *b, const char *ref, int);
110         extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
111         mplp_aux_t *ma = (mplp_aux_t*)data;
112         int ret, skip = 0;
113         do {
114                 int has_ref;
115                 ret = ma->iter? bam_iter_read(ma->fp, ma->iter, b) : bam_read1(ma->fp, b);
116                 if (ret < 0) break;
117                 if (b->core.tid < 0 || (b->core.flag&BAM_FUNMAP)) { // exclude unmapped reads
118                         skip = 1;
119                         continue;
120                 }
121                 if (ma->conf->bed) { // test overlap
122                         skip = !bed_overlap(ma->conf->bed, ma->h->target_name[b->core.tid], b->core.pos, bam_calend(&b->core, bam1_cigar(b)));
123                         if (skip) continue;
124                 }
125                 if (ma->conf->rghash) { // exclude read groups
126                         uint8_t *rg = bam_aux_get(b, "RG");
127                         skip = (rg && bcf_str2id(ma->conf->rghash, (const char*)(rg+1)) >= 0);
128                         if (skip) continue;
129                 }
130                 if (ma->conf->flag & MPLP_ILLUMINA13) {
131                         int i;
132                         uint8_t *qual = bam1_qual(b);
133                         for (i = 0; i < b->core.l_qseq; ++i)
134                                 qual[i] = qual[i] > 31? qual[i] - 31 : 0;
135                 }
136                 has_ref = (ma->ref && ma->ref_id == b->core.tid)? 1 : 0;
137                 skip = 0;
138                 if (has_ref && (ma->conf->flag&MPLP_REALN)) bam_prob_realn_core(b, ma->ref, (ma->conf->flag & MPLP_EXT_BAQ)? 3 : 1);
139                 if (has_ref && ma->conf->capQ_thres > 10) {
140                         int q = bam_cap_mapQ(b, ma->ref, ma->conf->capQ_thres);
141                         if (q < 0) skip = 1;
142                         else if (b->core.qual > q) b->core.qual = q;
143                 }
144                 else if (b->core.qual < ma->conf->min_mq) skip = 1; 
145                 else if ((ma->conf->flag&MPLP_NO_ORPHAN) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
146         } while (skip);
147         return ret;
148 }
149
150 static void group_smpl(mplp_pileup_t *m, bam_sample_t *sm, kstring_t *buf,
151                                            int n, char *const*fn, int *n_plp, const bam_pileup1_t **plp, int ignore_rg)
152 {
153         int i, j;
154         memset(m->n_plp, 0, m->n * sizeof(int));
155         for (i = 0; i < n; ++i) {
156                 for (j = 0; j < n_plp[i]; ++j) {
157                         const bam_pileup1_t *p = plp[i] + j;
158                         uint8_t *q;
159                         int id = -1;
160                         q = ignore_rg? 0 : bam_aux_get(p->b, "RG");
161                         if (q) id = bam_smpl_rg2smid(sm, fn[i], (char*)q+1, buf);
162                         if (id < 0) id = bam_smpl_rg2smid(sm, fn[i], 0, buf);
163                         if (id < 0 || id >= m->n) {
164                                 assert(q); // otherwise a bug
165                                 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]);
166                                 exit(1);
167                         }
168                         if (m->n_plp[id] == m->m_plp[id]) {
169                                 m->m_plp[id] = m->m_plp[id]? m->m_plp[id]<<1 : 8;
170                                 m->plp[id] = realloc(m->plp[id], sizeof(bam_pileup1_t) * m->m_plp[id]);
171                         }
172                         m->plp[id][m->n_plp[id]++] = *p;
173                 }
174         }
175 }
176
177 static int mpileup(mplp_conf_t *conf, int n, char **fn)
178 {
179         extern void *bcf_call_add_rg(void *rghash, const char *hdtext, const char *list);
180         extern void bcf_call_del_rghash(void *rghash);
181         mplp_aux_t **data;
182         int i, tid, pos, *n_plp, tid0 = -1, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid = -1, max_depth, max_indel_depth;
183         const bam_pileup1_t **plp;
184         bam_mplp_t iter;
185         bam_header_t *h = 0;
186         char *ref;
187         void *rghash = 0;
188
189         bcf_callaux_t *bca = 0;
190         bcf_callret1_t *bcr = 0;
191         bcf_call_t bc;
192         bcf_t *bp = 0;
193         bcf_hdr_t *bh = 0;
194
195         bam_sample_t *sm = 0;
196         kstring_t buf;
197         mplp_pileup_t gplp;
198
199         memset(&gplp, 0, sizeof(mplp_pileup_t));
200         memset(&buf, 0, sizeof(kstring_t));
201         memset(&bc, 0, sizeof(bcf_call_t));
202         data = calloc(n, sizeof(void*));
203         plp = calloc(n, sizeof(void*));
204         n_plp = calloc(n, sizeof(int*));
205         sm = bam_smpl_init();
206
207         // read the header and initialize data
208         for (i = 0; i < n; ++i) {
209                 bam_header_t *h_tmp;
210                 data[i] = calloc(1, sizeof(mplp_aux_t));
211                 data[i]->fp = strcmp(fn[i], "-") == 0? bam_dopen(fileno(stdin), "r") : bam_open(fn[i], "r");
212                 data[i]->conf = conf;
213                 h_tmp = bam_header_read(data[i]->fp);
214                 data[i]->h = i? h : h_tmp; // for i==0, "h" has not been set yet
215                 bam_smpl_add(sm, fn[i], (conf->flag&MPLP_IGNORE_RG)? 0 : h_tmp->text);
216                 rghash = bcf_call_add_rg(rghash, h_tmp->text, conf->pl_list);
217                 if (conf->reg) {
218                         int beg, end;
219                         bam_index_t *idx;
220                         idx = bam_index_load(fn[i]);
221                         if (idx == 0) {
222                                 fprintf(stderr, "[%s] fail to load index for %d-th input.\n", __func__, i+1);
223                                 exit(1);
224                         }
225                         if (bam_parse_region(h_tmp, conf->reg, &tid, &beg, &end) < 0) {
226                                 fprintf(stderr, "[%s] malformatted region or wrong seqname for %d-th input.\n", __func__, i+1);
227                                 exit(1);
228                         }
229                         if (i == 0) tid0 = tid, beg0 = beg, end0 = end;
230                         data[i]->iter = bam_iter_query(idx, tid, beg, end);
231                         bam_index_destroy(idx);
232                 }
233                 if (i == 0) h = h_tmp;
234                 else {
235                         // FIXME: to check consistency
236                         bam_header_destroy(h_tmp);
237                 }
238         }
239         gplp.n = sm->n;
240         gplp.n_plp = calloc(sm->n, sizeof(int));
241         gplp.m_plp = calloc(sm->n, sizeof(int));
242         gplp.plp = calloc(sm->n, sizeof(void*));
243
244         fprintf(stderr, "[%s] %d samples in %d input files\n", __func__, sm->n, n);
245         // write the VCF header
246         if (conf->flag & MPLP_GLF) {
247                 kstring_t s;
248                 bh = calloc(1, sizeof(bcf_hdr_t));
249                 s.l = s.m = 0; s.s = 0;
250                 bp = bcf_open("-", (conf->flag&MPLP_NO_COMP)? "wu" : "w");
251                 for (i = 0; i < h->n_targets; ++i) {
252                         kputs(h->target_name[i], &s);
253                         kputc('\0', &s);
254                 }
255                 bh->l_nm = s.l;
256                 bh->name = malloc(s.l);
257                 memcpy(bh->name, s.s, s.l);
258                 s.l = 0;
259                 for (i = 0; i < sm->n; ++i) {
260                         kputs(sm->smpl[i], &s); kputc('\0', &s);
261                 }
262                 bh->l_smpl = s.l;
263                 bh->sname = malloc(s.l);
264                 memcpy(bh->sname, s.s, s.l);
265                 bh->txt = malloc(strlen(BAM_VERSION) + 64);
266                 bh->l_txt = 1 + sprintf(bh->txt, "##samtoolsVersion=%s\n", BAM_VERSION);
267                 free(s.s);
268                 bcf_hdr_sync(bh);
269                 bcf_hdr_write(bp, bh);
270                 bca = bcf_call_init(-1., conf->min_baseQ);
271                 bcr = calloc(sm->n, sizeof(bcf_callret1_t));
272                 bca->rghash = rghash;
273                 bca->openQ = conf->openQ, bca->extQ = conf->extQ, bca->tandemQ = conf->tandemQ;
274                 bca->min_frac = conf->min_frac;
275                 bca->min_support = conf->min_support;
276         }
277         if (tid0 >= 0 && conf->fai) { // region is set
278                 ref = faidx_fetch_seq(conf->fai, h->target_name[tid0], 0, 0x7fffffff, &ref_len);
279                 ref_tid = tid0;
280                 for (i = 0; i < n; ++i) data[i]->ref = ref, data[i]->ref_id = tid0;
281         } else ref_tid = -1, ref = 0;
282         iter = bam_mplp_init(n, mplp_func, (void**)data);
283         max_depth = conf->max_depth;
284         if (max_depth * sm->n > 1<<20)
285                 fprintf(stderr, "(%s) Max depth is above 1M. Potential memory hog!\n", __func__);
286         if (max_depth * sm->n < 8000) {
287                 max_depth = 8000 / sm->n;
288                 fprintf(stderr, "<%s> Set max per-file depth to %d\n", __func__, max_depth);
289         }
290         max_indel_depth = conf->max_indel_depth * sm->n;
291         bam_mplp_set_maxcnt(iter, max_depth);
292         while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
293                 if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested
294                 if (conf->bed && tid >= 0 && !bed_overlap(conf->bed, h->target_name[tid], pos, pos+1)) continue;
295                 if (tid != ref_tid) {
296                         free(ref); ref = 0;
297                         if (conf->fai) ref = faidx_fetch_seq(conf->fai, h->target_name[tid], 0, 0x7fffffff, &ref_len);
298                         for (i = 0; i < n; ++i) data[i]->ref = ref, data[i]->ref_id = tid;
299                         ref_tid = tid;
300                 }
301                 if (conf->flag & MPLP_GLF) {
302                         int total_depth, _ref0, ref16;
303                         bcf1_t *b = calloc(1, sizeof(bcf1_t));
304                         for (i = total_depth = 0; i < n; ++i) total_depth += n_plp[i];
305                         group_smpl(&gplp, sm, &buf, n, fn, n_plp, plp, conf->flag & MPLP_IGNORE_RG);
306                         _ref0 = (ref && pos < ref_len)? ref[pos] : 'N';
307                         ref16 = bam_nt16_table[_ref0];
308                         for (i = 0; i < gplp.n; ++i)
309                                 bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], ref16, bca, bcr + i);
310                         bcf_call_combine(gplp.n, bcr, ref16, &bc);
311                         bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
312                                                  (conf->flag&MPLP_FMT_SP), 0, 0);
313                         bcf_write(bp, bh, b);
314                         bcf_destroy(b);
315                         // call indels
316                         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) {
317                                 for (i = 0; i < gplp.n; ++i)
318                                         bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], -1, bca, bcr + i);
319                                 if (bcf_call_combine(gplp.n, bcr, -1, &bc) >= 0) {
320                                         b = calloc(1, sizeof(bcf1_t));
321                                         bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
322                                                                  (conf->flag&MPLP_FMT_SP), bca, ref);
323                                         bcf_write(bp, bh, b);
324                                         bcf_destroy(b);
325                                 }
326                         }
327                 } else {
328                         printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
329                         for (i = 0; i < n; ++i) {
330                                 int j;
331                                 printf("\t%d\t", n_plp[i]);
332                                 if (n_plp[i] == 0) {
333                                         printf("*\t*"); // FIXME: printf() is very slow...
334                                         if (conf->flag & MPLP_PRINT_POS) printf("\t*");
335                                 } else {
336                                         for (j = 0; j < n_plp[i]; ++j)
337                                                 pileup_seq(plp[i] + j, pos, ref_len, ref);
338                                         putchar('\t');
339                                         for (j = 0; j < n_plp[i]; ++j) {
340                                                 const bam_pileup1_t *p = plp[i] + j;
341                                                 int c = bam1_qual(p->b)[p->qpos] + 33;
342                                                 if (c > 126) c = 126;
343                                                 putchar(c);
344                                         }
345                                         if (conf->flag & MPLP_PRINT_MAPQ) {
346                                                 putchar('\t');
347                                                 for (j = 0; j < n_plp[i]; ++j) {
348                                                         int c = plp[i][j].b->core.qual + 33;
349                                                         if (c > 126) c = 126;
350                                                         putchar(c);
351                                                 }
352                                         }
353                                         if (conf->flag & MPLP_PRINT_POS) {
354                                                 putchar('\t');
355                                                 for (j = 0; j < n_plp[i]; ++j) {
356                                                         if (j > 0) putchar(',');
357                                                         printf("%d", plp[i][j].qpos + 1); // FIXME: printf() is very slow...
358                                                 }
359                                         }
360                                 }
361                         }
362                         putchar('\n');
363                 }
364         }
365
366         bcf_close(bp);
367         bam_smpl_destroy(sm); free(buf.s);
368         for (i = 0; i < gplp.n; ++i) free(gplp.plp[i]);
369         free(gplp.plp); free(gplp.n_plp); free(gplp.m_plp);
370         bcf_call_del_rghash(rghash);
371         bcf_hdr_destroy(bh); bcf_call_destroy(bca); free(bc.PL); free(bcr);
372         bam_mplp_destroy(iter);
373         bam_header_destroy(h);
374         for (i = 0; i < n; ++i) {
375                 bam_close(data[i]->fp);
376                 if (data[i]->iter) bam_iter_destroy(data[i]->iter);
377                 free(data[i]);
378         }
379         free(data); free(plp); free(ref); free(n_plp);
380         return 0;
381 }
382
383 #define MAX_PATH_LEN 1024
384 static int read_file_list(const char *file_list,int *n,char **argv[])
385 {
386     char buf[MAX_PATH_LEN];
387     int len, nfiles;
388     char **files;
389
390     FILE *fh = fopen(file_list,"r");
391     if ( !fh )
392     {
393         fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
394         return 1;
395     }
396
397     // Speed is not an issue here, determine the number of files by reading the file twice
398     nfiles = 0;
399     while ( fgets(buf,MAX_PATH_LEN,fh) ) nfiles++;
400
401     if ( fseek(fh, 0L, SEEK_SET) )
402     {
403         fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
404         return 1;
405     }
406
407     files = calloc(nfiles,sizeof(char*));
408     nfiles = 0;
409     while ( fgets(buf,MAX_PATH_LEN,fh) ) 
410     {
411         len = strlen(buf);
412         while ( len>0 && isspace(buf[len-1]) ) len--;
413         if ( !len ) continue;
414
415         files[nfiles] = malloc(sizeof(char)*(len+1)); 
416         strncpy(files[nfiles],buf,len);
417         files[nfiles][len] = 0;
418         nfiles++;
419     }
420     fclose(fh);
421     if ( !nfiles )
422     {
423         fprintf(stderr,"No files read from %s\n", file_list);
424         return 1;
425     }
426     *argv = files;
427     *n    = nfiles;
428     return 0;
429 }
430 #undef MAX_PATH_LEN
431
432 int bam_mpileup(int argc, char *argv[])
433 {
434         int c;
435     const char *file_list = NULL;
436     char **fn = NULL;
437     int nfiles = 0, use_orphan = 0;
438         mplp_conf_t mplp;
439         memset(&mplp, 0, sizeof(mplp_conf_t));
440         #define MPLP_PRINT_POS 0x4000
441         mplp.max_mq = 60;
442         mplp.min_baseQ = 13;
443         mplp.capQ_thres = 0;
444         mplp.max_depth = 250; mplp.max_indel_depth = 250;
445         mplp.openQ = 40; mplp.extQ = 20; mplp.tandemQ = 100;
446         mplp.min_frac = 0.002; mplp.min_support = 1;
447         mplp.flag = MPLP_NO_ORPHAN | MPLP_REALN;
448         while ((c = getopt(argc, argv, "Agf:r:l:M:q:Q:uaRC:BDSd:L:b:P:o:e:h:Im:F:EG:6Os")) >= 0) {
449                 switch (c) {
450                 case 'f':
451                         mplp.fai = fai_load(optarg);
452                         if (mplp.fai == 0) return 1;
453                         break;
454                 case 'd': mplp.max_depth = atoi(optarg); break;
455                 case 'r': mplp.reg = strdup(optarg); break;
456                 case 'l': mplp.bed = bed_read(optarg); break;
457                 case 'P': mplp.pl_list = strdup(optarg); break;
458                 case 'g': mplp.flag |= MPLP_GLF; break;
459                 case 'u': mplp.flag |= MPLP_NO_COMP | MPLP_GLF; break;
460                 case 'a': mplp.flag |= MPLP_NO_ORPHAN | MPLP_REALN; break;
461                 case 'B': mplp.flag &= ~MPLP_REALN; break;
462                 case 'D': mplp.flag |= MPLP_FMT_DP; break;
463                 case 'S': mplp.flag |= MPLP_FMT_SP; break;
464                 case 'I': mplp.flag |= MPLP_NO_INDEL; break;
465                 case 'E': mplp.flag |= MPLP_EXT_BAQ; break;
466                 case '6': mplp.flag |= MPLP_ILLUMINA13; break;
467                 case 'R': mplp.flag |= MPLP_IGNORE_RG; break;
468                 case 's': mplp.flag |= MPLP_PRINT_MAPQ; break;
469                 case 'O': mplp.flag |= MPLP_PRINT_POS; break;
470                 case 'C': mplp.capQ_thres = atoi(optarg); break;
471                 case 'M': mplp.max_mq = atoi(optarg); break;
472                 case 'q': mplp.min_mq = atoi(optarg); break;
473                 case 'Q': mplp.min_baseQ = atoi(optarg); break;
474         case 'b': file_list = optarg; break;
475                 case 'o': mplp.openQ = atoi(optarg); break;
476                 case 'e': mplp.extQ = atoi(optarg); break;
477                 case 'h': mplp.tandemQ = atoi(optarg); break;
478                 case 'A': use_orphan = 1; break;
479                 case 'F': mplp.min_frac = atof(optarg); break;
480                 case 'm': mplp.min_support = atoi(optarg); break;
481                 case 'L': mplp.max_indel_depth = atoi(optarg); break;
482                 case 'G': {
483                                 FILE *fp_rg;
484                                 char buf[1024];
485                                 mplp.rghash = bcf_str2id_init();
486                                 if ((fp_rg = fopen(optarg, "r")) == 0)
487                                         fprintf(stderr, "(%s) Fail to open file %s. Continue anyway.\n", __func__, optarg);
488                                 while (!feof(fp_rg) && fscanf(fp_rg, "%s", buf) > 0) // this is not a good style, but forgive me...
489                                         bcf_str2id_add(mplp.rghash, strdup(buf));
490                                 fclose(fp_rg);
491                         }
492                         break;
493                 }
494         }
495         if (use_orphan) mplp.flag &= ~MPLP_NO_ORPHAN;
496         if (argc == 1) {
497                 fprintf(stderr, "\n");
498                 fprintf(stderr, "Usage: samtools mpileup [options] in1.bam [in2.bam [...]]\n\n");
499                 fprintf(stderr, "Input options:\n\n");
500                 fprintf(stderr, "       -6           assume the quality is in the Illumina-1.3+ encoding\n");
501                 fprintf(stderr, "       -A           count anomalous read pairs\n");
502                 fprintf(stderr, "       -B           disable BAQ computation\n");
503                 fprintf(stderr, "       -b FILE      list of input BAM files [null]\n");
504                 fprintf(stderr, "       -C INT       parameter for adjusting mapQ; 0 to disable [0]\n");
505                 fprintf(stderr, "       -d INT       max per-BAM depth to avoid excessive memory usage [%d]\n", mplp.max_depth);
506                 fprintf(stderr, "       -E           extended BAQ for higher sensitivity but lower specificity\n");
507                 fprintf(stderr, "       -f FILE      faidx indexed reference sequence file [null]\n");
508                 fprintf(stderr, "       -G FILE      exclude read groups listed in FILE [null]\n");
509                 fprintf(stderr, "       -l FILE      list of positions (chr pos) or regions (BED) [null]\n");
510                 fprintf(stderr, "       -M INT       cap mapping quality at INT [%d]\n", mplp.max_mq);
511                 fprintf(stderr, "       -r STR       region in which pileup is generated [null]\n");
512                 fprintf(stderr, "       -R           ignore RG tags\n");
513                 fprintf(stderr, "       -q INT       skip alignments with mapQ smaller than INT [%d]\n", mplp.min_mq);
514                 fprintf(stderr, "       -Q INT       skip bases with baseQ/BAQ smaller than INT [%d]\n", mplp.min_baseQ);
515                 fprintf(stderr, "\nOutput options:\n\n");
516                 fprintf(stderr, "       -D           output per-sample DP in BCF (require -g/-u)\n");
517                 fprintf(stderr, "       -g           generate BCF output (genotype likelihoods)\n");
518                 fprintf(stderr, "       -O           output base positions on reads (disabled by -g/-u)\n");
519                 fprintf(stderr, "       -s           output mapping quality (disabled by -g/-u)\n");
520                 fprintf(stderr, "       -S           output per-sample strand bias P-value in BCF (require -g/-u)\n");
521                 fprintf(stderr, "       -u           generate uncompress BCF output\n");
522                 fprintf(stderr, "\nSNP/INDEL genotype likelihoods options (effective with `-g' or `-u'):\n\n");
523                 fprintf(stderr, "       -e INT       Phred-scaled gap extension seq error probability [%d]\n", mplp.extQ);
524                 fprintf(stderr, "       -F FLOAT     minimum fraction of gapped reads for candidates [%g]\n", mplp.min_frac);
525                 fprintf(stderr, "       -h INT       coefficient for homopolymer errors [%d]\n", mplp.tandemQ);
526                 fprintf(stderr, "       -I           do not perform indel calling\n");
527                 fprintf(stderr, "       -L INT       max per-sample depth for INDEL calling [%d]\n", mplp.max_indel_depth);
528                 fprintf(stderr, "       -m INT       minimum gapped reads for indel candidates [%d]\n", mplp.min_support);
529                 fprintf(stderr, "       -o INT       Phred-scaled gap open sequencing error probability [%d]\n", mplp.openQ);
530                 fprintf(stderr, "       -P STR       comma separated list of platforms for indels [all]\n");
531                 fprintf(stderr, "\n");
532                 fprintf(stderr, "Notes: Assuming diploid individuals.\n\n");
533                 return 1;
534         }
535     if (file_list) {
536         if ( read_file_list(file_list,&nfiles,&fn) ) return 1;
537         mpileup(&mplp,nfiles,fn);
538         for (c=0; c<nfiles; c++) free(fn[c]);
539         free(fn);
540     } else mpileup(&mplp, argc - optind, argv + optind);
541         if (mplp.rghash) bcf_str2id_thorough_destroy(mplp.rghash);
542         free(mplp.reg); free(mplp.pl_list);
543         if (mplp.fai) fai_destroy(mplp.fai);
544         if (mplp.bed) bed_destroy(mplp.bed);
545         return 0;
546 }