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