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