]> git.donarmstrong.com Git - samtools.git/blob - bam_plcmd.c
bam_plcmd.c: added sanity check in read_file_list, complain if file list is not a...
[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 "sam.h"
9 #include "faidx.h"
10 #include "kstring.h"
11
12 static inline int printw(int c, FILE *fp)
13 {
14         char buf[16];
15         int l, x;
16         if (c == 0) return fputc('0', fp);
17         for (l = 0, x = c < 0? -c : c; x > 0; x /= 10) buf[l++] = x%10 + '0';
18         if (c < 0) buf[l++] = '-';
19         buf[l] = 0;
20         for (x = 0; x < l/2; ++x) {
21                 int y = buf[x]; buf[x] = buf[l-1-x]; buf[l-1-x] = y;
22         }
23         fputs(buf, fp);
24         return 0;
25 }
26
27 static inline void pileup_seq(const bam_pileup1_t *p, int pos, int ref_len, const char *ref)
28 {
29         int j;
30         if (p->is_head) {
31                 putchar('^');
32                 putchar(p->b->core.qual > 93? 126 : p->b->core.qual + 33);
33         }
34         if (!p->is_del) {
35                 int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)];
36                 if (ref) {
37                         int rb = pos < ref_len? ref[pos] : 'N';
38                         if (c == '=' || bam_nt16_table[c] == bam_nt16_table[rb]) c = bam1_strand(p->b)? ',' : '.';
39                         else c = bam1_strand(p->b)? tolower(c) : toupper(c);
40                 } else {
41                         if (c == '=') c = bam1_strand(p->b)? ',' : '.';
42                         else c = bam1_strand(p->b)? tolower(c) : toupper(c);
43                 }
44                 putchar(c);
45         } else putchar(p->is_refskip? (bam1_strand(p->b)? '<' : '>') : '*');
46         if (p->indel > 0) {
47                 putchar('+'); printw(p->indel, stdout);
48                 for (j = 1; j <= p->indel; ++j) {
49                         int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)];
50                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
51                 }
52         } else if (p->indel < 0) {
53                 printw(p->indel, stdout);
54                 for (j = 1; j <= -p->indel; ++j) {
55                         int c = (ref && (int)pos+j < ref_len)? ref[pos+j] : 'N';
56                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
57                 }
58         }
59         if (p->is_tail) putchar('$');
60 }
61
62 #include <assert.h>
63 #include "bam2bcf.h"
64 #include "sample.h"
65
66 #define MPLP_GLF   0x10
67 #define MPLP_NO_COMP 0x20
68 #define MPLP_NO_ORPHAN 0x40
69 #define MPLP_REALN   0x80
70 #define MPLP_NO_INDEL 0x400
71 #define MPLP_REDO_BAQ 0x800
72 #define MPLP_ILLUMINA13 0x1000
73 #define MPLP_IGNORE_RG 0x2000
74 #define MPLP_PRINT_POS 0x4000
75 #define MPLP_PRINT_MAPQ 0x8000
76 #define MPLP_PER_SAMPLE 0x10000
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, fmt_flag;
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_REDO_BAQ)? 7 : 3);
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         if ( !data[i]->fp )
213         {
214             fprintf(stderr, "[%s] failed to open %s: %s\n", __func__, fn[i], strerror(errno));
215             exit(1);
216         }
217                 data[i]->conf = conf;
218                 h_tmp = bam_header_read(data[i]->fp);
219                 data[i]->h = i? h : h_tmp; // for i==0, "h" has not been set yet
220                 bam_smpl_add(sm, fn[i], (conf->flag&MPLP_IGNORE_RG)? 0 : h_tmp->text);
221                 rghash = bcf_call_add_rg(rghash, h_tmp->text, conf->pl_list);
222                 if (conf->reg) {
223                         int beg, end;
224                         bam_index_t *idx;
225                         idx = bam_index_load(fn[i]);
226                         if (idx == 0) {
227                                 fprintf(stderr, "[%s] fail to load index for %s\n", __func__, fn[i]);
228                                 exit(1);
229                         }
230                         if (bam_parse_region(h_tmp, conf->reg, &tid, &beg, &end) < 0) {
231                                 fprintf(stderr, "[%s] malformatted region or wrong seqname for %s\n", __func__, fn[i]);
232                                 exit(1);
233                         }
234                         if (i == 0) tid0 = tid, beg0 = beg, end0 = end;
235                         data[i]->iter = bam_iter_query(idx, tid, beg, end);
236                         bam_index_destroy(idx);
237                 }
238                 if (i == 0) h = h_tmp;
239                 else {
240                         // FIXME: to check consistency
241                         bam_header_destroy(h_tmp);
242                 }
243         }
244         gplp.n = sm->n;
245         gplp.n_plp = calloc(sm->n, sizeof(int));
246         gplp.m_plp = calloc(sm->n, sizeof(int));
247         gplp.plp = calloc(sm->n, sizeof(void*));
248
249         fprintf(stderr, "[%s] %d samples in %d input files\n", __func__, sm->n, n);
250         // write the VCF header
251         if (conf->flag & MPLP_GLF) {
252                 kstring_t s;
253                 bh = calloc(1, sizeof(bcf_hdr_t));
254                 s.l = s.m = 0; s.s = 0;
255                 bp = bcf_open("-", (conf->flag&MPLP_NO_COMP)? "wu" : "w");
256                 for (i = 0; i < h->n_targets; ++i) {
257                         kputs(h->target_name[i], &s);
258                         kputc('\0', &s);
259                 }
260                 bh->l_nm = s.l;
261                 bh->name = malloc(s.l);
262                 memcpy(bh->name, s.s, s.l);
263                 s.l = 0;
264                 for (i = 0; i < sm->n; ++i) {
265                         kputs(sm->smpl[i], &s); kputc('\0', &s);
266                 }
267                 bh->l_smpl = s.l;
268                 bh->sname = malloc(s.l);
269                 memcpy(bh->sname, s.s, s.l);
270                 bh->txt = malloc(strlen(BAM_VERSION) + 64);
271                 bh->l_txt = 1 + sprintf(bh->txt, "##samtoolsVersion=%s\n", BAM_VERSION);
272                 free(s.s);
273                 bcf_hdr_sync(bh);
274                 bcf_hdr_write(bp, bh);
275                 bca = bcf_call_init(-1., conf->min_baseQ);
276                 bcr = calloc(sm->n, sizeof(bcf_callret1_t));
277                 bca->rghash = rghash;
278                 bca->openQ = conf->openQ, bca->extQ = conf->extQ, bca->tandemQ = conf->tandemQ;
279                 bca->min_frac = conf->min_frac;
280                 bca->min_support = conf->min_support;
281         bca->per_sample_flt = conf->flag & MPLP_PER_SAMPLE;
282         }
283         if (tid0 >= 0 && conf->fai) { // region is set
284                 ref = faidx_fetch_seq(conf->fai, h->target_name[tid0], 0, 0x7fffffff, &ref_len);
285                 ref_tid = tid0;
286                 for (i = 0; i < n; ++i) data[i]->ref = ref, data[i]->ref_id = tid0;
287         } else ref_tid = -1, ref = 0;
288         iter = bam_mplp_init(n, mplp_func, (void**)data);
289         max_depth = conf->max_depth;
290         if (max_depth * sm->n > 1<<20)
291                 fprintf(stderr, "(%s) Max depth is above 1M. Potential memory hog!\n", __func__);
292         if (max_depth * sm->n < 8000) {
293                 max_depth = 8000 / sm->n;
294                 fprintf(stderr, "<%s> Set max per-file depth to %d\n", __func__, max_depth);
295         }
296         max_indel_depth = conf->max_indel_depth * sm->n;
297         bam_mplp_set_maxcnt(iter, max_depth);
298         while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
299                 if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested
300                 if (conf->bed && tid >= 0 && !bed_overlap(conf->bed, h->target_name[tid], pos, pos+1)) continue;
301                 if (tid != ref_tid) {
302                         free(ref); ref = 0;
303                         if (conf->fai) ref = faidx_fetch_seq(conf->fai, h->target_name[tid], 0, 0x7fffffff, &ref_len);
304                         for (i = 0; i < n; ++i) data[i]->ref = ref, data[i]->ref_id = tid;
305                         ref_tid = tid;
306                 }
307                 if (conf->flag & MPLP_GLF) {
308                         int total_depth, _ref0, ref16;
309                         bcf1_t *b = calloc(1, sizeof(bcf1_t));
310                         for (i = total_depth = 0; i < n; ++i) total_depth += n_plp[i];
311                         group_smpl(&gplp, sm, &buf, n, fn, n_plp, plp, conf->flag & MPLP_IGNORE_RG);
312                         _ref0 = (ref && pos < ref_len)? ref[pos] : 'N';
313                         ref16 = bam_nt16_table[_ref0];
314                         for (i = 0; i < gplp.n; ++i)
315                                 bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], ref16, bca, bcr + i);
316                         bcf_call_combine(gplp.n, bcr, ref16, &bc);
317                         bcf_call2bcf(tid, pos, &bc, b, bcr, conf->fmt_flag, 0, 0);
318                         bcf_write(bp, bh, b);
319                         bcf_destroy(b);
320                         // call indels
321                         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) {
322                                 for (i = 0; i < gplp.n; ++i)
323                                         bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], -1, bca, bcr + i);
324                                 if (bcf_call_combine(gplp.n, bcr, -1, &bc) >= 0) {
325                                         b = calloc(1, sizeof(bcf1_t));
326                                         bcf_call2bcf(tid, pos, &bc, b, bcr, conf->fmt_flag, bca, ref);
327                                         bcf_write(bp, bh, b);
328                                         bcf_destroy(b);
329                                 }
330                         }
331                 } else {
332                         printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
333                         for (i = 0; i < n; ++i) {
334                                 int j, cnt;
335                                 for (j = cnt = 0; j < n_plp[i]; ++j) {
336                                         const bam_pileup1_t *p = plp[i] + j;
337                                         if (bam1_qual(p->b)[p->qpos] >= conf->min_baseQ) ++cnt;
338                                 }
339                                 printf("\t%d\t", cnt);
340                                 if (n_plp[i] == 0) {
341                                         printf("*\t*"); // FIXME: printf() is very slow...
342                                         if (conf->flag & MPLP_PRINT_POS) printf("\t*");
343                                 } else {
344                                         for (j = 0; j < n_plp[i]; ++j) {
345                                                 const bam_pileup1_t *p = plp[i] + j;
346                                                 if (bam1_qual(p->b)[p->qpos] >= conf->min_baseQ)
347                                                         pileup_seq(plp[i] + j, pos, ref_len, ref);
348                                         }
349                                         putchar('\t');
350                                         for (j = 0; j < n_plp[i]; ++j) {
351                                                 const bam_pileup1_t *p = plp[i] + j;
352                                                 int c = bam1_qual(p->b)[p->qpos];
353                                                 if (c >= conf->min_baseQ) {
354                                                         c = c + 33 < 126? c + 33 : 126;
355                                                         putchar(c);
356                                                 }
357                                         }
358                                         if (conf->flag & MPLP_PRINT_MAPQ) {
359                                                 putchar('\t');
360                                                 for (j = 0; j < n_plp[i]; ++j) {
361                                                         int c = plp[i][j].b->core.qual + 33;
362                                                         if (c > 126) c = 126;
363                                                         putchar(c);
364                                                 }
365                                         }
366                                         if (conf->flag & MPLP_PRINT_POS) {
367                                                 putchar('\t');
368                                                 for (j = 0; j < n_plp[i]; ++j) {
369                                                         if (j > 0) putchar(',');
370                                                         printf("%d", plp[i][j].qpos + 1); // FIXME: printf() is very slow...
371                                                 }
372                                         }
373                                 }
374                         }
375                         putchar('\n');
376                 }
377         }
378
379         bcf_close(bp);
380         bam_smpl_destroy(sm); free(buf.s);
381         for (i = 0; i < gplp.n; ++i) free(gplp.plp[i]);
382         free(gplp.plp); free(gplp.n_plp); free(gplp.m_plp);
383         bcf_call_del_rghash(rghash);
384         bcf_hdr_destroy(bh); bcf_call_destroy(bca); free(bc.PL); free(bcr);
385         bam_mplp_destroy(iter);
386         bam_header_destroy(h);
387         for (i = 0; i < n; ++i) {
388                 bam_close(data[i]->fp);
389                 if (data[i]->iter) bam_iter_destroy(data[i]->iter);
390                 free(data[i]);
391         }
392         free(data); free(plp); free(ref); free(n_plp);
393         return 0;
394 }
395
396 #define MAX_PATH_LEN 1024
397 static int read_file_list(const char *file_list,int *n,char **argv[])
398 {
399     char buf[MAX_PATH_LEN];
400     int len, nfiles = 0;
401     char **files = NULL;
402     struct stat sb;
403
404     *n = 0;
405     *argv = NULL;
406
407     FILE *fh = fopen(file_list,"r");
408     if ( !fh )
409     {
410         fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
411         return 1;
412     }
413
414     files = calloc(nfiles,sizeof(char*));
415     nfiles = 0;
416     while ( fgets(buf,MAX_PATH_LEN,fh) ) 
417     {
418         // allow empty lines and trailing spaces
419         len = strlen(buf);
420         while ( len>0 && isspace(buf[len-1]) ) len--;
421         if ( !len ) continue;
422
423         // check sanity of the file list
424         buf[len] = 0;
425         if (stat(buf, &sb) != 0)
426         {
427             // no such file, check if it is safe to print its name
428             int i, safe_to_print = 1;
429             for (i=0; i<len; i++)
430                 if (!isprint(buf[i])) { safe_to_print = 0; break; } 
431             if ( safe_to_print )
432                 fprintf(stderr,"The file list \"%s\" appears broken, could not locate: %s\n", file_list,buf);
433             else
434                 fprintf(stderr,"Does the file \"%s\" really contain a list of files and do all exist?\n", file_list);
435             return 1;
436         }
437
438         nfiles++;
439         files = realloc(files,nfiles*sizeof(char*));
440         files[nfiles-1] = strdup(buf);
441     }
442     fclose(fh);
443     if ( !nfiles )
444     {
445         fprintf(stderr,"No files read from %s\n", file_list);
446         return 1;
447     }
448     *argv = files;
449     *n    = nfiles;
450     return 0;
451 }
452 #undef MAX_PATH_LEN
453
454 int bam_mpileup(int argc, char *argv[])
455 {
456         int c;
457     const char *file_list = NULL;
458     char **fn = NULL;
459     int nfiles = 0, use_orphan = 0;
460         mplp_conf_t mplp;
461         memset(&mplp, 0, sizeof(mplp_conf_t));
462         mplp.max_mq = 60;
463         mplp.min_baseQ = 13;
464         mplp.capQ_thres = 0;
465         mplp.max_depth = 250; mplp.max_indel_depth = 250;
466         mplp.openQ = 40; mplp.extQ = 20; mplp.tandemQ = 100;
467         mplp.min_frac = 0.002; mplp.min_support = 1;
468         mplp.flag = MPLP_NO_ORPHAN | MPLP_REALN;
469         while ((c = getopt(argc, argv, "Agf:r:l:M:q:Q:uaRC:BDSd:L:b:P:po:e:h:Im:F:EG:6OsV")) >= 0) {
470                 switch (c) {
471                 case 'f':
472                         mplp.fai = fai_load(optarg);
473                         if (mplp.fai == 0) return 1;
474                         break;
475                 case 'd': mplp.max_depth = atoi(optarg); break;
476                 case 'r': mplp.reg = strdup(optarg); break;
477                 case 'l': mplp.bed = bed_read(optarg); break;
478                 case 'P': mplp.pl_list = strdup(optarg); break;
479                 case 'p': mplp.flag |= MPLP_PER_SAMPLE; break;
480                 case 'g': mplp.flag |= MPLP_GLF; break;
481                 case 'u': mplp.flag |= MPLP_NO_COMP | MPLP_GLF; break;
482                 case 'a': mplp.flag |= MPLP_NO_ORPHAN | MPLP_REALN; break;
483                 case 'B': mplp.flag &= ~MPLP_REALN; break;
484                 case 'D': mplp.fmt_flag |= B2B_FMT_DP; break;
485                 case 'S': mplp.fmt_flag |= B2B_FMT_SP; break;
486                 case 'V': mplp.fmt_flag |= B2B_FMT_DV; break;
487                 case 'I': mplp.flag |= MPLP_NO_INDEL; break;
488                 case 'E': mplp.flag |= MPLP_REDO_BAQ; break;
489                 case '6': mplp.flag |= MPLP_ILLUMINA13; break;
490                 case 'R': mplp.flag |= MPLP_IGNORE_RG; break;
491                 case 's': mplp.flag |= MPLP_PRINT_MAPQ; break;
492                 case 'O': mplp.flag |= MPLP_PRINT_POS; break;
493                 case 'C': mplp.capQ_thres = atoi(optarg); break;
494                 case 'M': mplp.max_mq = atoi(optarg); break;
495                 case 'q': mplp.min_mq = atoi(optarg); break;
496                 case 'Q': mplp.min_baseQ = atoi(optarg); break;
497         case 'b': file_list = optarg; break;
498                 case 'o': mplp.openQ = atoi(optarg); break;
499                 case 'e': mplp.extQ = atoi(optarg); break;
500                 case 'h': mplp.tandemQ = atoi(optarg); break;
501                 case 'A': use_orphan = 1; break;
502                 case 'F': mplp.min_frac = atof(optarg); break;
503                 case 'm': mplp.min_support = atoi(optarg); break;
504                 case 'L': mplp.max_indel_depth = atoi(optarg); break;
505                 case 'G': {
506                                 FILE *fp_rg;
507                                 char buf[1024];
508                                 mplp.rghash = bcf_str2id_init();
509                                 if ((fp_rg = fopen(optarg, "r")) == 0)
510                                         fprintf(stderr, "(%s) Fail to open file %s. Continue anyway.\n", __func__, optarg);
511                                 while (!feof(fp_rg) && fscanf(fp_rg, "%s", buf) > 0) // this is not a good style, but forgive me...
512                                         bcf_str2id_add(mplp.rghash, strdup(buf));
513                                 fclose(fp_rg);
514                         }
515                         break;
516                 }
517         }
518         if (use_orphan) mplp.flag &= ~MPLP_NO_ORPHAN;
519         if (argc == 1) {
520                 fprintf(stderr, "\n");
521                 fprintf(stderr, "Usage: samtools mpileup [options] in1.bam [in2.bam [...]]\n\n");
522                 fprintf(stderr, "Input options:\n\n");
523                 fprintf(stderr, "       -6           assume the quality is in the Illumina-1.3+ encoding\n");
524                 fprintf(stderr, "       -A           count anomalous read pairs\n");
525                 fprintf(stderr, "       -B           disable BAQ computation\n");
526                 fprintf(stderr, "       -b FILE      list of input BAM files [null]\n");
527                 fprintf(stderr, "       -C INT       parameter for adjusting mapQ; 0 to disable [0]\n");
528                 fprintf(stderr, "       -d INT       max per-BAM depth to avoid excessive memory usage [%d]\n", mplp.max_depth);
529                 fprintf(stderr, "       -E           recalculate extended BAQ on the fly thus ignoring existing BQs\n");
530                 fprintf(stderr, "       -f FILE      faidx indexed reference sequence file [null]\n");
531                 fprintf(stderr, "       -G FILE      exclude read groups listed in FILE [null]\n");
532                 fprintf(stderr, "       -l FILE      list of positions (chr pos) or regions (BED) [null]\n");
533                 fprintf(stderr, "       -M INT       cap mapping quality at INT [%d]\n", mplp.max_mq);
534                 fprintf(stderr, "       -r STR       region in which pileup is generated [null]\n");
535                 fprintf(stderr, "       -R           ignore RG tags\n");
536                 fprintf(stderr, "       -q INT       skip alignments with mapQ smaller than INT [%d]\n", mplp.min_mq);
537                 fprintf(stderr, "       -Q INT       skip bases with baseQ/BAQ smaller than INT [%d]\n", mplp.min_baseQ);
538                 fprintf(stderr, "\nOutput options:\n\n");
539                 fprintf(stderr, "       -D           output per-sample DP in BCF (require -g/-u)\n");
540                 fprintf(stderr, "       -g           generate BCF output (genotype likelihoods)\n");
541                 fprintf(stderr, "       -O           output base positions on reads (disabled by -g/-u)\n");
542                 fprintf(stderr, "       -s           output mapping quality (disabled by -g/-u)\n");
543                 fprintf(stderr, "       -S           output per-sample strand bias P-value in BCF (require -g/-u)\n");
544                 fprintf(stderr, "       -u           generate uncompress BCF output\n");
545                 fprintf(stderr, "\nSNP/INDEL genotype likelihoods options (effective with `-g' or `-u'):\n\n");
546                 fprintf(stderr, "       -e INT       Phred-scaled gap extension seq error probability [%d]\n", mplp.extQ);
547                 fprintf(stderr, "       -F FLOAT     minimum fraction of gapped reads for candidates [%g]\n", mplp.min_frac);
548                 fprintf(stderr, "       -h INT       coefficient for homopolymer errors [%d]\n", mplp.tandemQ);
549                 fprintf(stderr, "       -I           do not perform indel calling\n");
550                 fprintf(stderr, "       -L INT       max per-sample depth for INDEL calling [%d]\n", mplp.max_indel_depth);
551                 fprintf(stderr, "       -m INT       minimum gapped reads for indel candidates [%d]\n", mplp.min_support);
552                 fprintf(stderr, "       -o INT       Phred-scaled gap open sequencing error probability [%d]\n", mplp.openQ);
553                 fprintf(stderr, "       -p           apply -m and -F per-sample to increase sensitivity\n");
554                 fprintf(stderr, "       -P STR       comma separated list of platforms for indels [all]\n");
555                 fprintf(stderr, "\n");
556                 fprintf(stderr, "Notes: Assuming diploid individuals.\n\n");
557                 return 1;
558         }
559         bam_no_B = 1;
560     if (file_list) {
561         if ( read_file_list(file_list,&nfiles,&fn) ) return 1;
562         mpileup(&mplp,nfiles,fn);
563         for (c=0; c<nfiles; c++) free(fn[c]);
564         free(fn);
565     } else mpileup(&mplp, argc - optind, argv + optind);
566         if (mplp.rghash) bcf_str2id_thorough_destroy(mplp.rghash);
567         free(mplp.reg); free(mplp.pl_list);
568         if (mplp.fai) fai_destroy(mplp.fai);
569         if (mplp.bed) bed_destroy(mplp.bed);
570         return 0;
571 }