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