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 "bam_maqcns.h"
10 #include "khash.h"
11 #include "glf.h"
12 #include "kstring.h"
13
14 typedef int *indel_list_t;
15 KHASH_MAP_INIT_INT64(64, indel_list_t)
16
17 #define BAM_PLF_SIMPLE 0x01
18 #define BAM_PLF_CNS 0x02
19 #define BAM_PLF_INDEL_ONLY 0x04
20 #define BAM_PLF_GLF 0x08
21 #define BAM_PLF_VAR_ONLY 0x10
22 #define BAM_PLF_2ND 0x20
23 #define BAM_PLF_RANBASE 0x40
24 #define BAM_PLF_1STBASE 0x80
25 #define BAM_PLF_ALLBASE 0x100
26 #define BAM_PLF_READPOS 0x200
27 #define BAM_PLF_NOBAQ 0x400
28
29 typedef struct {
30 bam_header_t *h;
31 bam_maqcns_t *c;
32 bam_maqindel_opt_t *ido;
33 faidx_t *fai;
34 khash_t(64) *hash;
35 uint32_t format;
36 int tid, len, last_pos;
37 int mask;
38 int capQ_thres, min_baseQ;
39 int max_depth; // for indel calling, ignore reads with the depth too high. 0 for unlimited
40 char *ref;
41 glfFile fp_glf; // for glf output only
42 } pu_data_t;
43
44 char **__bam_get_lines(const char *fn, int *_n);
45 void bam_init_header_hash(bam_header_t *header);
46 int32_t bam_get_tid(const bam_header_t *header, const char *seq_name);
47
load_pos(const char * fn,bam_header_t * h)48 static khash_t(64) *load_pos(const char *fn, bam_header_t *h)
49 {
50 char **list;
51 int i, j, n, *fields, max_fields;
52 khash_t(64) *hash;
53 bam_init_header_hash(h);
54 list = __bam_get_lines(fn, &n);
55 hash = kh_init(64);
56 max_fields = 0; fields = 0;
57 for (i = 0; i < n; ++i) {
58 char *str = list[i];
59 int chr, n_fields, ret;
60 khint_t k;
61 uint64_t x;
62 n_fields = ksplit_core(str, 0, &max_fields, &fields);
63 if (n_fields < 2) continue;
64 chr = bam_get_tid(h, str + fields[0]);
65 if (chr < 0) {
66 fprintf(stderr, "[load_pos] unknown reference sequence name: %s\n", str + fields[0]);
67 continue;
68 }
69 x = (uint64_t)chr << 32 | (atoi(str + fields[1]) - 1);
70 k = kh_put(64, hash, x, &ret);
71 if (ret == 0) {
72 fprintf(stderr, "[load_pos] position %s:%s has been loaded.\n", str+fields[0], str+fields[1]);
73 continue;
74 }
75 kh_val(hash, k) = 0;
76 if (n_fields > 2) {
77 // count
78 for (j = 2; j < n_fields; ++j) {
79 char *s = str + fields[j];
80 if ((*s != '+' && *s != '-') || !isdigit(s[1])) break;
81 }
82 if (j > 2) { // update kh_val()
83 int *q, y, z;
84 q = kh_val(hash, k) = (int*)calloc(j - 1, sizeof(int));
85 q[0] = j - 2; z = j; y = 1;
86 for (j = 2; j < z; ++j)
87 q[y++] = atoi(str + fields[j]);
88 }
89 }
90 free(str);
91 }
92 free(list); free(fields);
93 return hash;
94 }
95
printw(int c,FILE * fp)96 static inline int printw(int c, FILE *fp)
97 {
98 char buf[16];
99 int l, x;
100 if (c == 0) return fputc('0', fp);
101 for (l = 0, x = c < 0? -c : c; x > 0; x /= 10) buf[l++] = x%10 + '0';
102 if (c < 0) buf[l++] = '-';
103 buf[l] = 0;
104 for (x = 0; x < l/2; ++x) {
105 int y = buf[x]; buf[x] = buf[l-1-x]; buf[l-1-x] = y;
106 }
107 fputs(buf, fp);
108 return 0;
109 }
110
111 // an analogy to pileup_func() below
glt3_func(uint32_t tid,uint32_t pos,int n,const bam_pileup1_t * pu,void * data)112 static int glt3_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
113 {
114 pu_data_t *d = (pu_data_t*)data;
115 bam_maqindel_ret_t *r = 0;
116 int rb, *proposed_indels = 0;
117 glf1_t *g;
118 glf3_t *g3;
119
120 if (d->fai == 0) {
121 fprintf(stderr, "[glt3_func] reference sequence is required for generating GLT. Abort!\n");
122 exit(1);
123 }
124 if (d->hash) { // only output a list of sites
125 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
126 if (k == kh_end(d->hash)) return 0;
127 proposed_indels = kh_val(d->hash, k);
128 }
129 g3 = glf3_init1();
130 if (d->fai && (int)tid != d->tid) {
131 if (d->ref) { // then write the end mark
132 g3->rtype = GLF3_RTYPE_END;
133 glf3_write1(d->fp_glf, g3);
134 }
135 glf3_ref_write(d->fp_glf, d->h->target_name[tid], d->h->target_len[tid]); // write reference
136 free(d->ref);
137 d->ref = fai_fetch(d->fai, d->h->target_name[tid], &d->len);
138 d->tid = tid;
139 d->last_pos = 0;
140 }
141 rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
142 g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
143 memcpy(g3, g, sizeof(glf1_t));
144 g3->rtype = GLF3_RTYPE_SUB;
145 g3->offset = pos - d->last_pos;
146 d->last_pos = pos;
147 glf3_write1(d->fp_glf, g3);
148 if (pos < d->len) {
149 int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
150 if (proposed_indels)
151 r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
152 else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
153 }
154 if (r) { // then write indel line
155 int het = 3 * n, min;
156 min = het;
157 if (min > r->gl[0]) min = r->gl[0];
158 if (min > r->gl[1]) min = r->gl[1];
159 g3->ref_base = 0;
160 g3->rtype = GLF3_RTYPE_INDEL;
161 memset(g3->lk, 0, 10);
162 g3->lk[0] = r->gl[0] - min < 255? r->gl[0] - min : 255;
163 g3->lk[1] = r->gl[1] - min < 255? r->gl[1] - min : 255;
164 g3->lk[2] = het - min < 255? het - min : 255;
165 g3->offset = 0;
166 g3->indel_len[0] = r->indel1;
167 g3->indel_len[1] = r->indel2;
168 g3->min_lk = min < 255? min : 255;
169 g3->max_len = (abs(r->indel1) > abs(r->indel2)? abs(r->indel1) : abs(r->indel2)) + 1;
170 g3->indel_seq[0] = strdup(r->s[0]+1);
171 g3->indel_seq[1] = strdup(r->s[1]+1);
172 glf3_write1(d->fp_glf, g3);
173 bam_maqindel_ret_destroy(r);
174 }
175 free(g);
176 glf3_destroy1(g3);
177 return 0;
178 }
179
pileup_seq(const bam_pileup1_t * p,int pos,int ref_len,const char * ref)180 static inline void pileup_seq(const bam_pileup1_t *p, int pos, int ref_len, const char *ref)
181 {
182 int j;
183 if (p->is_head) {
184 putchar('^');
185 putchar(p->b->core.qual > 93? 126 : p->b->core.qual + 33);
186 }
187 if (!p->is_del) {
188 int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)];
189 if (ref) {
190 int rb = pos < ref_len? ref[pos] : 'N';
191 if (c == '=' || bam_nt16_table[c] == bam_nt16_table[rb]) c = bam1_strand(p->b)? ',' : '.';
192 else c = bam1_strand(p->b)? tolower(c) : toupper(c);
193 } else {
194 if (c == '=') c = bam1_strand(p->b)? ',' : '.';
195 else c = bam1_strand(p->b)? tolower(c) : toupper(c);
196 }
197 putchar(c);
198 } else putchar(p->is_refskip? (bam1_strand(p->b)? '<' : '>') : '*');
199 if (p->indel > 0) {
200 putchar('+'); printw(p->indel, stdout);
201 for (j = 1; j <= p->indel; ++j) {
202 int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)];
203 putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
204 }
205 } else if (p->indel < 0) {
206 printw(p->indel, stdout);
207 for (j = 1; j <= -p->indel; ++j) {
208 int c = (ref && (int)pos+j < ref_len)? ref[pos+j] : 'N';
209 putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
210 }
211 }
212 if (p->is_tail) putchar('$');
213 }
214
pileup_func(uint32_t tid,uint32_t pos,int n,const bam_pileup1_t * pu,void * data)215 static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
216 {
217 pu_data_t *d = (pu_data_t*)data;
218 bam_maqindel_ret_t *r = 0;
219 int i, rb, rms_mapq = -1, *proposed_indels = 0;
220 uint64_t rms_aux;
221 uint32_t cns = 0;
222
223 // if GLF is required, suppress -c completely
224 if (d->format & BAM_PLF_GLF) return glt3_func(tid, pos, n, pu, data);
225 // if d->hash is initialized, only output the sites in the hash table
226 if (d->hash) {
227 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
228 if (k == kh_end(d->hash)) return 0;
229 proposed_indels = kh_val(d->hash, k);
230 }
231 // update d->ref if necessary
232 if (d->fai && (int)tid != d->tid) {
233 free(d->ref);
234 d->ref = faidx_fetch_seq(d->fai, d->h->target_name[tid], 0, 0x7fffffff, &d->len);
235 d->tid = tid;
236 }
237 rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
238 // when the indel-only mode is asked for, return if no reads mapped with indels
239 if (d->format & BAM_PLF_INDEL_ONLY) {
240 for (i = 0; i < n; ++i)
241 if (pu[i].indel != 0) break;
242 if (i == n) return 0;
243 }
244 // call the consensus and indel
245 if (d->format & BAM_PLF_CNS) { // call consensus
246 if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE)) { // use a random base or the 1st base as the consensus call
247 const bam_pileup1_t *p = (d->format & BAM_PLF_1STBASE)? pu : pu + (int)(drand48() * n);
248 int q = bam1_qual(p->b)[p->qpos];
249 int mapQ = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
250 uint32_t b = bam1_seqi(bam1_seq(p->b), p->qpos);
251 cns = b<<28 | 0xf<<24 | mapQ<<16 | q<<8;
252 } else if (d->format & BAM_PLF_ALLBASE) { // collapse all bases
253 uint64_t rmsQ = 0;
254 uint32_t b = 0;
255 for (i = 0; i < n; ++i) {
256 const bam_pileup1_t *p = pu + i;
257 int q = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
258 b |= bam1_seqi(bam1_seq(p->b), p->qpos);
259 rmsQ += q * q;
260 }
261 rmsQ = (uint64_t)(sqrt((double)rmsQ / n) + .499);
262 cns = b<<28 | 0xf<<24 | rmsQ<<16 | 60<<8;
263 } else {
264 glf1_t *g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
265 cns = g->depth == 0? (0xfu<<28 | 0xf<<24) : glf2cns(g, (int)(d->c->q_r + .499));
266 free(g);
267 }
268 }
269 if ((d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)) && d->ref && pos < d->len) { // call indels
270 int m = (!d->max_depth || d->max_depth>n) ? n : d->max_depth;
271 if (proposed_indels) // the first element gives the size of the array
272 r = bam_maqindel(m, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
273 else r = bam_maqindel(m, pos, d->ido, pu, d->ref, 0, 0);
274 }
275 // when only variant sites are asked for, test if the site is a variant
276 if ((d->format & BAM_PLF_CNS) && (d->format & BAM_PLF_VAR_ONLY)) {
277 if (!(bam_nt16_table[rb] != 15 && cns>>28 != 15 && cns>>28 != bam_nt16_table[rb])) { // not a SNP
278 if (!(r && (r->gt == 2 || strcmp(r->s[r->gt], "*")))) { // not an indel
279 if (r) bam_maqindel_ret_destroy(r);
280 return 0;
281 }
282 }
283 }
284 // print the first 3 columns
285 fputs(d->h->target_name[tid], stdout); putchar('\t');
286 printw(pos+1, stdout); putchar('\t'); putchar(rb); putchar('\t');
287 // print consensus information if required
288 if (d->format & BAM_PLF_CNS) {
289 putchar(bam_nt16_rev_table[cns>>28]); putchar('\t');
290 printw(cns>>8&0xff, stdout); putchar('\t');
291 printw(cns&0xff, stdout); putchar('\t');
292 printw(cns>>16&0xff, stdout); putchar('\t');
293 }
294 // print pileup sequences
295 printw(n, stdout); putchar('\t');
296 for (i = 0; i < n; ++i)
297 pileup_seq(pu + i, pos, d->len, d->ref);
298 // finalize rms_mapq
299 if (d->format & BAM_PLF_CNS) {
300 for (i = rms_aux = 0; i < n; ++i) {
301 const bam_pileup1_t *p = pu + i;
302 int tmp = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ;
303 rms_aux += tmp * tmp;
304 }
305 rms_aux = (uint64_t)(sqrt((double)rms_aux / n) + .499);
306 if (rms_mapq < 0) rms_mapq = rms_aux;
307 }
308 putchar('\t');
309 // print quality
310 for (i = 0; i < n; ++i) {
311 const bam_pileup1_t *p = pu + i;
312 int c = bam1_qual(p->b)[p->qpos] + 33;
313 if (c > 126) c = 126;
314 putchar(c);
315 }
316 if (d->format & BAM_PLF_2ND) { // print 2nd calls and qualities
317 const unsigned char *q;
318 putchar('\t');
319 for (i = 0; i < n; ++i) {
320 const bam_pileup1_t *p = pu + i;
321 q = bam_aux_get(p->b, "E2");
322 putchar(q? q[p->qpos + 1] : 'N');
323 }
324 putchar('\t');
325 for (i = 0; i < n; ++i) {
326 const bam_pileup1_t *p = pu + i;
327 q = bam_aux_get(p->b, "U2");
328 putchar(q? q[p->qpos + 1] : '!');
329 }
330 }
331 // print mapping quality if -s is flagged on the command line
332 if (d->format & BAM_PLF_SIMPLE) {
333 putchar('\t');
334 for (i = 0; i < n; ++i) {
335 int c = pu[i].b->core.qual + 33;
336 if (c > 126) c = 126;
337 putchar(c);
338 }
339 }
340 // print read position
341 if (d->format & BAM_PLF_READPOS) {
342 putchar('\t');
343 for (i = 0; i < n; ++i) {
344 int x = pu[i].qpos;
345 int l = pu[i].b->core.l_qseq;
346 printw(x < l/2? x+1 : -((l-1)-x+1), stdout); putchar(',');
347 }
348 }
349 putchar('\n');
350 // print the indel line if r has been calculated. This only happens if:
351 // a) -c or -i are flagged, AND b) the reference sequence is available
352 if (r) {
353 printf("%s\t%d\t*\t", d->h->target_name[tid], pos + 1);
354 if (r->gt < 2) printf("%s/%s\t", r->s[r->gt], r->s[r->gt]);
355 else printf("%s/%s\t", r->s[0], r->s[1]);
356 printf("%d\t%d\t", r->q_cns, r->q_ref);
357 printf("%d\t%d\t", rms_mapq, n);
358 printf("%s\t%s\t", r->s[0], r->s[1]);
359 //printf("%d\t%d\t", r->gl[0], r->gl[1]);
360 printf("%d\t%d\t%d\t", r->cnt1, r->cnt2, r->cnt_anti);
361 printf("%d\t%d\n", r->cnt_ref, r->cnt_ambi);
362 bam_maqindel_ret_destroy(r);
363 }
364 return 0;
365 }
366
bam_pileup(int argc,char * argv[])367 int bam_pileup(int argc, char *argv[])
368 {
369 int c, is_SAM = 0;
370 char *fn_list = 0, *fn_fa = 0, *fn_pos = 0;
371 pu_data_t *d = (pu_data_t*)calloc(1, sizeof(pu_data_t));
372 d->max_depth = 1024; d->tid = -1; d->mask = BAM_DEF_MASK; d->min_baseQ = 13;
373 d->c = bam_maqcns_init();
374 d->c->errmod = BAM_ERRMOD_MAQ2; // change the default model
375 d->ido = bam_maqindel_opt_init();
376 while ((c = getopt(argc, argv, "st:f:cT:N:r:l:d:im:gI:G:vM:S2aR:PAQ:C:B")) >= 0) {
377 switch (c) {
378 case 'Q': d->c->min_baseQ = atoi(optarg); break;
379 case 'C': d->capQ_thres = atoi(optarg); break;
380 case 'B': d->format |= BAM_PLF_NOBAQ; break;
381 case 'a': d->c->errmod = BAM_ERRMOD_SOAP; break;
382 case 'A': d->c->errmod = BAM_ERRMOD_MAQ; break;
383 case 's': d->format |= BAM_PLF_SIMPLE; break;
384 case 't': fn_list = strdup(optarg); break;
385 case 'l': fn_pos = strdup(optarg); break;
386 case 'f': fn_fa = strdup(optarg); break;
387 case 'T': d->c->theta = atof(optarg); break;
388 case 'N': d->c->n_hap = atoi(optarg); break;
389 case 'r': d->c->het_rate = atof(optarg); d->ido->r_snp = d->c->het_rate; break;
390 case 'M': d->c->cap_mapQ = atoi(optarg); break;
391 case 'd': d->max_depth = atoi(optarg); break;
392 case 'c': d->format |= BAM_PLF_CNS; break;
393 case 'i': d->format |= BAM_PLF_INDEL_ONLY; break;
394 case 'v': d->format |= BAM_PLF_VAR_ONLY; break;
395 case 'm': d->mask = strtol(optarg, 0, 0); break;
396 case 'g': d->format |= BAM_PLF_GLF; break;
397 case '2': d->format |= BAM_PLF_2ND; break;
398 case 'P': d->format |= BAM_PLF_READPOS; break;
399 case 'I': d->ido->q_indel = atoi(optarg); break;
400 case 'G': d->ido->r_indel = atof(optarg); break;
401 case 'S': is_SAM = 1; break;
402 case 'R':
403 if (strcmp(optarg, "random") == 0) d->format |= BAM_PLF_RANBASE;
404 else if (strcmp(optarg, "first") == 0) d->format |= BAM_PLF_1STBASE;
405 else if (strcmp(optarg, "all") == 0) d->format |= BAM_PLF_ALLBASE;
406 else fprintf(stderr, "[bam_pileup] unrecognized -R\n");
407 break;
408 default: fprintf(stderr, "Unrecognizd option '-%c'.\n", c); return 1;
409 }
410 }
411 if (d->c->errmod != BAM_ERRMOD_MAQ2) d->c->theta += 0.02;
412 if (d->c->theta > 1.0) d->c->theta = 1.0;
413 if (fn_list) is_SAM = 1;
414 if (optind == argc) {
415 fprintf(stderr, "\n");
416 fprintf(stderr, "Usage: samtools pileup [options] <in.bam>|<in.sam>\n\n");
417 fprintf(stderr, "Option: -s simple (yet incomplete) pileup format\n");
418 fprintf(stderr, " -S the input is in SAM\n");
419 fprintf(stderr, " -B disable BAQ computation\n");
420 fprintf(stderr, " -A use the original MAQ model for SNP calling (DEPRECATED)\n");
421 fprintf(stderr, " -2 output the 2nd best call and quality\n");
422 fprintf(stderr, " -i only show lines/consensus with indels\n");
423 fprintf(stderr, " -Q INT min base quality (possibly capped by BAQ) [%d]\n", d->c->min_baseQ);
424 fprintf(stderr, " -C INT coefficient for adjusting mapQ of poor mappings [%d]\n", d->capQ_thres);
425 fprintf(stderr, " -m INT filtering reads with bits in INT [0x%x]\n", d->mask);
426 fprintf(stderr, " -M INT cap mapping quality at INT [%d]\n", d->c->cap_mapQ);
427 fprintf(stderr, " -d INT limit maximum depth for indels [%d]\n", d->max_depth);
428 fprintf(stderr, " -t FILE list of reference sequences (force -S)\n");
429 fprintf(stderr, " -l FILE list of sites at which pileup is output\n");
430 fprintf(stderr, " -f FILE reference sequence in the FASTA format\n\n");
431 fprintf(stderr, " -c compute the consensus sequence\n");
432 fprintf(stderr, " -v print variants only (for -c)\n");
433 fprintf(stderr, " -g output in the GLFv3 format (DEPRECATED)\n");
434 fprintf(stderr, " -T FLOAT theta in maq consensus calling model (for -c) [%.4g]\n", d->c->theta);
435 fprintf(stderr, " -N INT number of haplotypes in the sample (for -c) [%d]\n", d->c->n_hap);
436 fprintf(stderr, " -r FLOAT prior of a difference between two haplotypes (for -c) [%.4g]\n", d->c->het_rate);
437 fprintf(stderr, " -G FLOAT prior of an indel between two haplotypes (for -c) [%.4g]\n", d->ido->r_indel);
438 fprintf(stderr, " -I INT phred prob. of an indel in sequencing/prep. (for -c) [%d]\n", d->ido->q_indel);
439 fprintf(stderr, "\n");
440 free(fn_list); free(fn_fa); free(d);
441 return 1;
442 }
443 if (d->format & (BAM_PLF_RANBASE|BAM_PLF_1STBASE|BAM_PLF_ALLBASE)) d->format |= BAM_PLF_CNS;
444 if (fn_fa) d->fai = fai_load(fn_fa);
445 if (d->format & (BAM_PLF_CNS|BAM_PLF_GLF)) bam_maqcns_prepare(d->c); // consensus calling
446 if (d->format & BAM_PLF_GLF) { // for glf output
447 glf3_header_t *h;
448 h = glf3_header_init();
449 d->fp_glf = bgzf_fdopen(fileno(stdout), "w");
450 glf3_header_write(d->fp_glf, h);
451 glf3_header_destroy(h);
452 }
453 if (d->fai == 0 && (d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)))
454 fprintf(stderr, "[bam_pileup] indels will not be called when -f is absent.\n");
455 if (fn_fa && is_SAM && fn_list == 0) fn_list = samfaipath(fn_fa);
456
457 {
458 samfile_t *fp;
459 fp = is_SAM? samopen(argv[optind], "r", fn_list) : samopen(argv[optind], "rb", 0);
460 if (fp == 0 || fp->header == 0) {
461 fprintf(stderr, "[bam_pileup] fail to read the header: non-exisiting file or wrong format.\n");
462 return 1;
463 }
464 d->h = fp->header;
465 if (fn_pos) d->hash = load_pos(fn_pos, d->h);
466 { // run pileup
467 extern int bam_prob_realn(bam1_t *b, const char *ref);
468 extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
469 bam1_t *b;
470 int ret, tid, pos, n_plp;
471 bam_plp_t iter;
472 const bam_pileup1_t *plp;
473 b = bam_init1();
474 iter = bam_plp_init(0, 0);
475 bam_plp_set_mask(iter, d->mask);
476 while ((ret = samread(fp, b)) >= 0) {
477 int skip = 0;
478 if ((int)b->core.tid < 0) break;
479 // update d->ref if necessary
480 if (d->fai && (int)b->core.tid != d->tid) {
481 free(d->ref);
482 d->ref = faidx_fetch_seq(d->fai, d->h->target_name[b->core.tid], 0, 0x7fffffff, &d->len);
483 d->tid = b->core.tid;
484 }
485 if (d->ref && (d->format&BAM_PLF_CNS) && !(d->format&BAM_PLF_NOBAQ)) bam_prob_realn(b, d->ref);
486 if (d->ref && (d->format&BAM_PLF_CNS) && d->capQ_thres > 10) {
487 int q = bam_cap_mapQ(b, d->ref, d->capQ_thres);
488 if (q < 0) skip = 1;
489 else if (b->core.qual > q) b->core.qual = q;
490 } else if (b->core.flag&BAM_FUNMAP) skip = 1;
491 else if ((d->format&BAM_PLF_CNS) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
492 if (skip) continue;
493 bam_plp_push(iter, b);
494 while ((plp = bam_plp_next(iter, &tid, &pos, &n_plp)) != 0)
495 pileup_func(tid, pos, n_plp, plp, d);
496 }
497 bam_plp_push(iter, 0);
498 while ((plp = bam_plp_next(iter, &tid, &pos, &n_plp)) != 0)
499 pileup_func(tid, pos, n_plp, plp, d);
500 bam_plp_destroy(iter);
501 bam_destroy1(b);
502 }
503 samclose(fp); // d->h will be destroyed here
504 }
505
506 // free
507 if (d->format & BAM_PLF_GLF) bgzf_close(d->fp_glf);
508 if (fn_pos) { // free the hash table
509 khint_t k;
510 for (k = kh_begin(d->hash); k < kh_end(d->hash); ++k)
511 if (kh_exist(d->hash, k)) free(kh_val(d->hash, k));
512 kh_destroy(64, d->hash);
513 }
514 free(fn_pos); free(fn_list); free(fn_fa);
515 if (d->fai) fai_destroy(d->fai);
516 bam_maqcns_destroy(d->c);
517 free(d->ido); free(d->ref); free(d);
518 return 0;
519 }
520
521 /***********
522 * mpileup *
523 ***********/
524
525 #include <assert.h>
526 #include "bam2bcf.h"
527 #include "sample.h"
528
529 #define MPLP_GLF 0x10
530 #define MPLP_NO_COMP 0x20
531 #define MPLP_NO_ORPHAN 0x40
532 #define MPLP_REALN 0x80
533 #define MPLP_FMT_DP 0x100
534 #define MPLP_FMT_SP 0x200
535 #define MPLP_NO_INDEL 0x400
536
537 typedef struct {
538 int max_mq, min_mq, flag, min_baseQ, capQ_thres, max_depth;
539 int openQ, extQ, tandemQ;
540 char *reg, *fn_pos, *pl_list;
541 faidx_t *fai;
542 kh_64_t *hash;
543 } mplp_conf_t;
544
545 typedef struct {
546 bamFile fp;
547 bam_iter_t iter;
548 int min_mq, flag, ref_id, capQ_thres;
549 char *ref;
550 } mplp_aux_t;
551
552 typedef struct {
553 int n;
554 int *n_plp, *m_plp;
555 bam_pileup1_t **plp;
556 } mplp_pileup_t;
557
mplp_func(void * data,bam1_t * b)558 static int mplp_func(void *data, bam1_t *b)
559 {
560 extern int bam_realn(bam1_t *b, const char *ref);
561 extern int bam_prob_realn_core(bam1_t *b, const char *ref, int);
562 extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
563 mplp_aux_t *ma = (mplp_aux_t*)data;
564 int ret, skip = 0;
565 do {
566 int has_ref;
567 ret = ma->iter? bam_iter_read(ma->fp, ma->iter, b) : bam_read1(ma->fp, b);
568 if (ret < 0) break;
569 has_ref = (ma->ref && ma->ref_id == b->core.tid)? 1 : 0;
570 skip = 0;
571 if (has_ref && (ma->flag&MPLP_REALN)) bam_prob_realn_core(b, ma->ref, 1);
572 if (has_ref && ma->capQ_thres > 10) {
573 int q = bam_cap_mapQ(b, ma->ref, ma->capQ_thres);
574 if (q < 0) skip = 1;
575 else if (b->core.qual > q) b->core.qual = q;
576 } else if (b->core.flag&BAM_FUNMAP) skip = 1;
577 else if (b->core.qual < ma->min_mq) skip = 1;
578 else if ((ma->flag&MPLP_NO_ORPHAN) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
579 } while (skip);
580 return ret;
581 }
582
group_smpl(mplp_pileup_t * m,bam_sample_t * sm,kstring_t * buf,int n,char * const * fn,int * n_plp,const bam_pileup1_t ** plp)583 static void group_smpl(mplp_pileup_t *m, bam_sample_t *sm, kstring_t *buf,
584 int n, char *const*fn, int *n_plp, const bam_pileup1_t **plp)
585 {
586 int i, j;
587 memset(m->n_plp, 0, m->n * sizeof(int));
588 for (i = 0; i < n; ++i) {
589 for (j = 0; j < n_plp[i]; ++j) {
590 const bam_pileup1_t *p = plp[i] + j;
591 uint8_t *q;
592 int id = -1;
593 q = bam_aux_get(p->b, "RG");
594 if (q) id = bam_smpl_rg2smid(sm, fn[i], (char*)q+1, buf);
595 if (id < 0) id = bam_smpl_rg2smid(sm, fn[i], 0, buf);
596 assert(id >= 0 && id < m->n);
597 if (m->n_plp[id] == m->m_plp[id]) {
598 m->m_plp[id] = m->m_plp[id]? m->m_plp[id]<<1 : 8;
599 m->plp[id] = realloc(m->plp[id], sizeof(bam_pileup1_t) * m->m_plp[id]);
600 }
601 m->plp[id][m->n_plp[id]++] = *p;
602 }
603 }
604 }
605
mpileup(mplp_conf_t * conf,int n,char ** fn)606 static int mpileup(mplp_conf_t *conf, int n, char **fn)
607 {
608 extern void *bcf_call_add_rg(void *rghash, const char *hdtext, const char *list);
609 extern void bcf_call_del_rghash(void *rghash);
610 mplp_aux_t **data;
611 int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid, max_depth;
612 const bam_pileup1_t **plp;
613 bam_mplp_t iter;
614 bam_header_t *h = 0;
615 char *ref;
616 khash_t(64) *hash = 0;
617 void *rghash = 0;
618
619 bcf_callaux_t *bca = 0;
620 bcf_callret1_t *bcr = 0;
621 bcf_call_t bc;
622 bcf_t *bp = 0;
623 bcf_hdr_t *bh = 0;
624
625 bam_sample_t *sm = 0;
626 kstring_t buf;
627 mplp_pileup_t gplp;
628
629 memset(&gplp, 0, sizeof(mplp_pileup_t));
630 memset(&buf, 0, sizeof(kstring_t));
631 memset(&bc, 0, sizeof(bcf_call_t));
632 data = calloc(n, sizeof(void*));
633 plp = calloc(n, sizeof(void*));
634 n_plp = calloc(n, sizeof(int*));
635 sm = bam_smpl_init();
636
637 // read the header and initialize data
638 for (i = 0; i < n; ++i) {
639 bam_header_t *h_tmp;
640 data[i] = calloc(1, sizeof(mplp_aux_t));
641 data[i]->min_mq = conf->min_mq;
642 data[i]->flag = conf->flag;
643 data[i]->capQ_thres = conf->capQ_thres;
644 data[i]->fp = strcmp(fn[i], "-") == 0? bam_dopen(fileno(stdin), "r") : bam_open(fn[i], "r");
645 h_tmp = bam_header_read(data[i]->fp);
646 bam_smpl_add(sm, fn[i], h_tmp->text);
647 rghash = bcf_call_add_rg(rghash, h_tmp->text, conf->pl_list);
648 if (conf->reg) {
649 int beg, end;
650 bam_index_t *idx;
651 idx = bam_index_load(fn[i]);
652 if (idx == 0) {
653 fprintf(stderr, "[%s] fail to load index for %d-th input.\n", __func__, i+1);
654 exit(1);
655 }
656 if (bam_parse_region(h_tmp, conf->reg, &tid, &beg, &end) < 0) {
657 fprintf(stderr, "[%s] malformatted region or wrong seqname for %d-th input.\n", __func__, i+1);
658 exit(1);
659 }
660 if (i == 0) beg0 = beg, end0 = end;
661 data[i]->iter = bam_iter_query(idx, tid, beg, end);
662 bam_index_destroy(idx);
663 }
664 if (i == 0) h = h_tmp;
665 else {
666 // FIXME: to check consistency
667 bam_header_destroy(h_tmp);
668 }
669 }
670 gplp.n = sm->n;
671 gplp.n_plp = calloc(sm->n, sizeof(int));
672 gplp.m_plp = calloc(sm->n, sizeof(int));
673 gplp.plp = calloc(sm->n, sizeof(void*));
674
675 fprintf(stderr, "[%s] %d samples in %d input files\n", __func__, sm->n, n);
676 if (conf->fn_pos) hash = load_pos(conf->fn_pos, h);
677 // write the VCF header
678 if (conf->flag & MPLP_GLF) {
679 kstring_t s;
680 bh = calloc(1, sizeof(bcf_hdr_t));
681 s.l = s.m = 0; s.s = 0;
682 bp = bcf_open("-", (conf->flag&MPLP_NO_COMP)? "wu" : "w");
683 for (i = 0; i < h->n_targets; ++i) {
684 kputs(h->target_name[i], &s);
685 kputc('\0', &s);
686 }
687 bh->l_nm = s.l;
688 bh->name = malloc(s.l);
689 memcpy(bh->name, s.s, s.l);
690 s.l = 0;
691 for (i = 0; i < sm->n; ++i) {
692 kputs(sm->smpl[i], &s); kputc('\0', &s);
693 }
694 bh->l_smpl = s.l;
695 bh->sname = malloc(s.l);
696 memcpy(bh->sname, s.s, s.l);
697 bh->l_txt = 0;
698 free(s.s);
699 bcf_hdr_sync(bh);
700 bcf_hdr_write(bp, bh);
701 bca = bcf_call_init(-1., conf->min_baseQ);
702 bcr = calloc(sm->n, sizeof(bcf_callret1_t));
703 bca->rghash = rghash;
704 bca->openQ = conf->openQ, bca->extQ = conf->extQ, bca->tandemQ = conf->tandemQ;
705 }
706 ref_tid = -1; ref = 0;
707 iter = bam_mplp_init(n, mplp_func, (void**)data);
708 max_depth = conf->max_depth;
709 if (max_depth * sm->n > 1<<20)
710 fprintf(stderr, "(%s) Max depth is above 1M. Potential memory hog!\n", __func__);
711 if (max_depth * sm->n < 8000) {
712 max_depth = 8000 / sm->n;
713 fprintf(stderr, "<%s> Set max per-sample depth to %d\n", __func__, max_depth);
714 }
715 bam_mplp_set_maxcnt(iter, max_depth);
716 while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
717 if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested
718 if (hash) {
719 khint_t k;
720 k = kh_get(64, hash, (uint64_t)tid<<32 | pos);
721 if (k == kh_end(hash)) continue;
722 }
723 if (tid != ref_tid) {
724 free(ref); ref = 0;
725 if (conf->fai) ref = fai_fetch(conf->fai, h->target_name[tid], &ref_len);
726 for (i = 0; i < n; ++i) data[i]->ref = ref, data[i]->ref_id = tid;
727 ref_tid = tid;
728 }
729 if (conf->flag & MPLP_GLF) {
730 int _ref0, ref16;
731 bcf1_t *b = calloc(1, sizeof(bcf1_t));
732 group_smpl(&gplp, sm, &buf, n, fn, n_plp, plp);
733 _ref0 = (ref && pos < ref_len)? ref[pos] : 'N';
734 ref16 = bam_nt16_table[_ref0];
735 for (i = 0; i < gplp.n; ++i)
736 bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], ref16, bca, bcr + i);
737 bcf_call_combine(gplp.n, bcr, ref16, &bc);
738 bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
739 (conf->flag&MPLP_FMT_SP), 0, 0);
740 bcf_write(bp, bh, b);
741 bcf_destroy(b);
742 // call indels
743 if (!(conf->flag&MPLP_NO_INDEL) && bcf_call_gap_prep(gplp.n, gplp.n_plp, gplp.plp, pos, bca, ref, rghash) >= 0) {
744 for (i = 0; i < gplp.n; ++i)
745 bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], -1, bca, bcr + i);
746 if (bcf_call_combine(gplp.n, bcr, -1, &bc) >= 0) {
747 b = calloc(1, sizeof(bcf1_t));
748 bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
749 (conf->flag&MPLP_FMT_SP), bca, ref);
750 bcf_write(bp, bh, b);
751 bcf_destroy(b);
752 }
753 }
754 } else {
755 printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
756 for (i = 0; i < n; ++i) {
757 int j;
758 printf("\t%d\t", n_plp[i]);
759 if (n_plp[i] == 0) printf("*\t*");
760 else {
761 for (j = 0; j < n_plp[i]; ++j)
762 pileup_seq(plp[i] + j, pos, ref_len, ref);
763 putchar('\t');
764 for (j = 0; j < n_plp[i]; ++j) {
765 const bam_pileup1_t *p = plp[i] + j;
766 int c = bam1_qual(p->b)[p->qpos] + 33;
767 if (c > 126) c = 126;
768 putchar(c);
769 }
770 }
771 }
772 putchar('\n');
773 }
774 }
775
776 bcf_close(bp);
777 bam_smpl_destroy(sm); free(buf.s);
778 for (i = 0; i < gplp.n; ++i) free(gplp.plp[i]);
779 free(gplp.plp); free(gplp.n_plp); free(gplp.m_plp);
780 bcf_call_del_rghash(rghash);
781 if (hash) { // free the hash table
782 khint_t k;
783 for (k = kh_begin(hash); k < kh_end(hash); ++k)
784 if (kh_exist(hash, k)) free(kh_val(hash, k));
785 kh_destroy(64, hash);
786 }
787 bcf_hdr_destroy(bh); bcf_call_destroy(bca); free(bc.PL); free(bcr);
788 bam_mplp_destroy(iter);
789 bam_header_destroy(h);
790 for (i = 0; i < n; ++i) {
791 bam_close(data[i]->fp);
792 if (data[i]->iter) bam_iter_destroy(data[i]->iter);
793 free(data[i]);
794 }
795 free(data); free(plp); free(ref); free(n_plp);
796 return 0;
797 }
798
799 #define MAX_PATH_LEN 1024
read_file_list(const char * file_list,int * n,char ** argv[])800 int read_file_list(const char *file_list,int *n,char **argv[])
801 {
802 char buf[MAX_PATH_LEN];
803 int len, nfiles;
804 char **files;
805
806 FILE *fh = fopen(file_list,"r");
807 if ( !fh )
808 {
809 fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
810 return 1;
811 }
812
813 // Speed is not an issue here, determine the number of files by reading the file twice
814 nfiles = 0;
815 while ( fgets(buf,MAX_PATH_LEN,fh) ) nfiles++;
816
817 if ( fseek(fh, 0L, SEEK_SET) )
818 {
819 fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
820 return 1;
821 }
822
823 files = calloc(nfiles,sizeof(char*));
824 nfiles = 0;
825 while ( fgets(buf,MAX_PATH_LEN,fh) )
826 {
827 len = strlen(buf);
828 while ( len>0 && isspace(buf[len-1]) ) len--;
829 if ( !len ) continue;
830
831 files[nfiles] = malloc(sizeof(char)*(len+1));
832 strncpy(files[nfiles],buf,len);
833 files[nfiles][len] = 0;
834 nfiles++;
835 }
836 fclose(fh);
837 if ( !nfiles )
838 {
839 fprintf(stderr,"No files read from %s\n", file_list);
840 return 1;
841 }
842 *argv = files;
843 *n = nfiles;
844 return 0;
845 }
846 #undef MAX_PATH_LEN
847
bam_mpileup(int argc,char * argv[])848 int bam_mpileup(int argc, char *argv[])
849 {
850 int c;
851 const char *file_list = NULL;
852 char **fn = NULL;
853 int nfiles = 0;
854 mplp_conf_t mplp;
855 memset(&mplp, 0, sizeof(mplp_conf_t));
856 mplp.max_mq = 60;
857 mplp.min_baseQ = 13;
858 mplp.capQ_thres = 0;
859 mplp.max_depth = 250;
860 mplp.openQ = 40; mplp.extQ = 20; mplp.tandemQ = 100;
861 mplp.flag = MPLP_NO_ORPHAN | MPLP_REALN;
862 while ((c = getopt(argc, argv, "gf:r:l:M:q:Q:uaORC:BDSd:b:P:o:e:h:I")) >= 0) {
863 switch (c) {
864 case 'f':
865 mplp.fai = fai_load(optarg);
866 if (mplp.fai == 0) return 1;
867 break;
868 case 'd': mplp.max_depth = atoi(optarg); break;
869 case 'r': mplp.reg = strdup(optarg); break;
870 case 'l': mplp.fn_pos = strdup(optarg); break;
871 case 'P': mplp.pl_list = strdup(optarg); break;
872 case 'g': mplp.flag |= MPLP_GLF; break;
873 case 'u': mplp.flag |= MPLP_NO_COMP | MPLP_GLF; break;
874 case 'a': mplp.flag |= MPLP_NO_ORPHAN | MPLP_REALN; break;
875 case 'B': mplp.flag &= ~MPLP_REALN & ~MPLP_NO_ORPHAN; break;
876 case 'O': mplp.flag |= MPLP_NO_ORPHAN; break;
877 case 'R': mplp.flag |= MPLP_REALN; break;
878 case 'D': mplp.flag |= MPLP_FMT_DP; break;
879 case 'S': mplp.flag |= MPLP_FMT_SP; break;
880 case 'I': mplp.flag |= MPLP_NO_INDEL; break;
881 case 'C': mplp.capQ_thres = atoi(optarg); break;
882 case 'M': mplp.max_mq = atoi(optarg); break;
883 case 'q': mplp.min_mq = atoi(optarg); break;
884 case 'Q': mplp.min_baseQ = atoi(optarg); break;
885 case 'b': file_list = optarg; break;
886 case 'o': mplp.openQ = atoi(optarg); break;
887 case 'e': mplp.extQ = atoi(optarg); break;
888 case 'h': mplp.tandemQ = atoi(optarg); break;
889 }
890 }
891 if (argc == 1) {
892 fprintf(stderr, "\n");
893 fprintf(stderr, "Usage: samtools mpileup [options] in1.bam [in2.bam [...]]\n\n");
894 fprintf(stderr, "Options: -f FILE reference sequence file [null]\n");
895 fprintf(stderr, " -r STR region in which pileup is generated [null]\n");
896 fprintf(stderr, " -l FILE list of positions (format: chr pos) [null]\n");
897 fprintf(stderr, " -b FILE list of input BAM files [null]\n");
898 fprintf(stderr, " -M INT cap mapping quality at INT [%d]\n", mplp.max_mq);
899 fprintf(stderr, " -Q INT min base quality [%d]\n", mplp.min_baseQ);
900 fprintf(stderr, " -q INT filter out alignment with MQ smaller than INT [%d]\n", mplp.min_mq);
901 fprintf(stderr, " -d INT max per-sample depth [%d]\n", mplp.max_depth);
902 fprintf(stderr, " -P STR comma separated list of platforms for indels [all]\n");
903 fprintf(stderr, " -o INT Phred-scaled gap open sequencing error probability [%d]\n", mplp.openQ);
904 fprintf(stderr, " -e INT Phred-scaled gap extension seq error probability [%d]\n", mplp.extQ);
905 fprintf(stderr, " -h INT coefficient for homopolyer errors [%d]\n", mplp.tandemQ);
906 fprintf(stderr, " -g generate BCF output\n");
907 fprintf(stderr, " -u do not compress BCF output\n");
908 fprintf(stderr, " -B disable BAQ computation\n");
909 fprintf(stderr, " -D output per-sample DP\n");
910 fprintf(stderr, " -S output per-sample SP (strand bias P-value, slow)\n");
911 fprintf(stderr, " -I do not perform indel calling\n");
912 fprintf(stderr, "\n");
913 fprintf(stderr, "Notes: Assuming diploid individuals.\n\n");
914 return 1;
915 }
916 if ( file_list )
917 {
918 if ( read_file_list(file_list,&nfiles,&fn) ) return 1;
919 mpileup(&mplp,nfiles,fn);
920 for (c=0; c<nfiles; c++) free(fn[c]);
921 free(fn);
922 }
923 else
924 mpileup(&mplp, argc - optind, argv + optind);
925 free(mplp.reg); free(mplp.pl_list);
926 if (mplp.fai) fai_destroy(mplp.fai);
927 return 0;
928 }
929