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