# Install the manual page of bcftools.
[samtools.git] / bam_pileup.c
1 #include <stdio.h>
2 #include <stdlib.h>
3 #include <ctype.h>
4 #include <assert.h>
5 #include "sam.h"
6
7 typedef struct {
8         int k, x, y, end;
9 } cstate_t;
10
11 static cstate_t g_cstate_null = { -1, 0, 0, 0 };
12
13 typedef struct __linkbuf_t {
14         bam1_t b;
15         uint32_t beg, end;
16         cstate_t s;
17         struct __linkbuf_t *next;
18 } lbnode_t;
19
20 /* --- BEGIN: Memory pool */
21
22 typedef struct {
23         int cnt, n, max;
24         lbnode_t **buf;
25 } mempool_t;
26
27 static mempool_t *mp_init()
28 {
29         mempool_t *mp;
30         mp = (mempool_t*)calloc(1, sizeof(mempool_t));
31         return mp;
32 }
33 static void mp_destroy(mempool_t *mp)
34 {
35         int k;
36         for (k = 0; k < mp->n; ++k) {
37                 free(mp->buf[k]->b.data);
38                 free(mp->buf[k]);
39         }
40         free(mp->buf);
41         free(mp);
42 }
43 static inline lbnode_t *mp_alloc(mempool_t *mp)
44 {
45         ++mp->cnt;
46         if (mp->n == 0) return (lbnode_t*)calloc(1, sizeof(lbnode_t));
47         else return mp->buf[--mp->n];
48 }
49 static inline void mp_free(mempool_t *mp, lbnode_t *p)
50 {
51         --mp->cnt; p->next = 0; // clear lbnode_t::next here
52         if (mp->n == mp->max) {
53                 mp->max = mp->max? mp->max<<1 : 256;
54                 mp->buf = (lbnode_t**)realloc(mp->buf, sizeof(lbnode_t*) * mp->max);
55         }
56         mp->buf[mp->n++] = p;
57 }
58
59 /* --- END: Memory pool */
60
61 /* --- BEGIN: Auxiliary functions */
62
63 /* s->k: the index of the CIGAR operator that has just been processed.
64    s->x: the reference coordinate of the start of s->k
65    s->y: the query coordiante of the start of s->k
66  */
67 static inline int resolve_cigar2(bam_pileup1_t *p, uint32_t pos, cstate_t *s)
68 {
69 #define _cop(c) ((c)&BAM_CIGAR_MASK)
70 #define _cln(c) ((c)>>BAM_CIGAR_SHIFT)
71
72         bam1_t *b = p->b;
73         bam1_core_t *c = &b->core;
74         uint32_t *cigar = bam1_cigar(b);
75         int k, is_head = 0;
76         // determine the current CIGAR operation
77 //      fprintf(stderr, "%s\tpos=%d\tend=%d\t(%d,%d,%d)\n", bam1_qname(b), pos, s->end, s->k, s->x, s->y);
78         if (s->k == -1) { // never processed
79                 is_head = 1;
80                 if (c->n_cigar == 1) { // just one operation, save a loop
81                         if (_cop(cigar[0]) == BAM_CMATCH) s->k = 0, s->x = c->pos, s->y = 0;
82                 } else { // find the first match or deletion
83                         for (k = 0, s->x = c->pos, s->y = 0; k < c->n_cigar; ++k) {
84                                 int op = _cop(cigar[k]);
85                                 int l = _cln(cigar[k]);
86                                 if (op == BAM_CMATCH || op == BAM_CDEL) break;
87                                 else if (op == BAM_CREF_SKIP) s->x += l;
88                                 else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) s->y += l;
89                         }
90                         assert(k < c->n_cigar);
91                         s->k = k;
92                 }
93         } else { // the read has been processed before
94                 int op, l = _cln(cigar[s->k]);
95                 if (pos - s->x >= l) { // jump to the next operation
96                         assert(s->k < c->n_cigar); // otherwise a bug: this function should not be called in this case
97                         op = _cop(cigar[s->k+1]);
98                         if (op == BAM_CMATCH || op == BAM_CDEL || op == BAM_CREF_SKIP) { // jump to the next without a loop
99                                 if (_cop(cigar[s->k]) == BAM_CMATCH) s->y += l;
100                                 s->x += l;
101                                 ++s->k;
102                         } else { // find the next M/D/N
103                                 if (_cop(cigar[s->k]) == BAM_CMATCH) s->y += l;
104                                 s->x += l;
105                                 for (k = s->k + 1; k < c->n_cigar; ++k) {
106                                         op = _cop(cigar[k]), l = _cln(cigar[k]);
107                                         if (op == BAM_CMATCH || op == BAM_CDEL || op == BAM_CREF_SKIP) break;
108                                         else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) s->y += l;
109                                 }
110                                 s->k = k;
111                         }
112                         assert(s->k < c->n_cigar); // otherwise a bug
113                 } // else, do nothing
114         }
115         { // collect pileup information
116                 int op, l;
117                 op = _cop(cigar[s->k]); l = _cln(cigar[s->k]);
118                 p->is_del = p->indel = p->is_refskip = 0;
119                 if (s->x + l - 1 == pos && s->k + 1 < c->n_cigar) { // peek the next operation
120                         int op2 = _cop(cigar[s->k+1]);
121                         int l2 = _cln(cigar[s->k+1]);
122                         if (op2 == BAM_CDEL) p->indel = -(int)l2;
123                         else if (op2 == BAM_CINS) p->indel = l2;
124                         else if (op2 == BAM_CPAD && s->k + 2 < c->n_cigar) { // no working for adjacent padding
125                                 int l3 = 0;
126                                 for (k = s->k + 2; k < c->n_cigar; ++k) {
127                                         op2 = _cop(cigar[k]); l2 = _cln(cigar[k]);
128                                         if (op2 == BAM_CINS) l3 += l2;
129                                         else if (op2 == BAM_CDEL || op2 == BAM_CMATCH || op2 == BAM_CREF_SKIP) break;
130                                 }
131                                 if (l3 > 0) p->indel = l3;
132                         }
133                 }
134                 if (op == BAM_CMATCH) {
135                         p->qpos = s->y + (pos - s->x);
136                 } else if (op == BAM_CDEL || op == BAM_CREF_SKIP) {
137                         p->is_del = 1; p->qpos = s->y; // FIXME: distinguish D and N!!!!!
138                         p->is_refskip = (op == BAM_CREF_SKIP);
139                 } // cannot be other operations; otherwise a bug
140                 p->is_head = (pos == c->pos); p->is_tail = (pos == s->end);
141         }
142         return 1;
143 }
144
145 /* --- END: Auxiliary functions */
146
147 /*******************
148  * pileup iterator *
149  *******************/
150
151 struct __bam_plp_t {
152         mempool_t *mp;
153         lbnode_t *head, *tail, *dummy;
154         int32_t tid, pos, max_tid, max_pos;
155         int is_eof, flag_mask, max_plp, error, maxcnt;
156         bam_pileup1_t *plp;
157         // for the "auto" interface only
158         bam1_t *b;
159         bam_plp_auto_f func;
160         void *data;
161 };
162
163 bam_plp_t bam_plp_init(bam_plp_auto_f func, void *data)
164 {
165         bam_plp_t iter;
166         iter = calloc(1, sizeof(struct __bam_plp_t));
167         iter->mp = mp_init();
168         iter->head = iter->tail = mp_alloc(iter->mp);
169         iter->dummy = mp_alloc(iter->mp);
170         iter->max_tid = iter->max_pos = -1;
171         iter->flag_mask = BAM_DEF_MASK;
172         iter->maxcnt = 8000;
173         if (func) {
174                 iter->func = func;
175                 iter->data = data;
176                 iter->b = bam_init1();
177         }
178         return iter;
179 }
180
181 void bam_plp_destroy(bam_plp_t iter)
182 {
183         mp_free(iter->mp, iter->dummy);
184         mp_free(iter->mp, iter->head);
185         if (iter->mp->cnt != 0)
186                 fprintf(stderr, "[bam_plp_destroy] memory leak: %d. Continue anyway.\n", iter->mp->cnt);
187         mp_destroy(iter->mp);
188         if (iter->b) bam_destroy1(iter->b);
189         free(iter->plp);
190         free(iter);
191 }
192
193 const bam_pileup1_t *bam_plp_next(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp)
194 {
195         if (iter->error) { *_n_plp = -1; return 0; }
196         *_n_plp = 0;
197         if (iter->is_eof && iter->head->next == 0) return 0;
198         while (iter->is_eof || iter->max_tid > iter->tid || (iter->max_tid == iter->tid && iter->max_pos > iter->pos)) {
199                 int n_plp = 0;
200                 lbnode_t *p, *q;
201                 // write iter->plp at iter->pos
202                 iter->dummy->next = iter->head;
203                 for (p = iter->head, q = iter->dummy; p->next; q = p, p = p->next) {
204                         if (p->b.core.tid < iter->tid || (p->b.core.tid == iter->tid && p->end <= iter->pos)) { // then remove
205                                 q->next = p->next; mp_free(iter->mp, p); p = q;
206                         } else if (p->b.core.tid == iter->tid && p->beg <= iter->pos) { // here: p->end > pos; then add to pileup
207                                 if (n_plp == iter->max_plp) { // then double the capacity
208                                         iter->max_plp = iter->max_plp? iter->max_plp<<1 : 256;
209                                         iter->plp = (bam_pileup1_t*)realloc(iter->plp, sizeof(bam_pileup1_t) * iter->max_plp);
210                                 }
211                                 iter->plp[n_plp].b = &p->b;
212                                 if (resolve_cigar2(iter->plp + n_plp, iter->pos, &p->s)) ++n_plp; // actually always true...
213                         }
214                 }
215                 iter->head = iter->dummy->next; // dummy->next may be changed
216                 *_n_plp = n_plp; *_tid = iter->tid; *_pos = iter->pos;
217                 // update iter->tid and iter->pos
218                 if (iter->head->next) {
219                         if (iter->tid > iter->head->b.core.tid) {
220                                 fprintf(stderr, "[%s] unsorted input. Pileup aborts.\n", __func__);
221                                 iter->error = 1;
222                                 *_n_plp = -1;
223                                 return 0;
224                         }
225                 }
226                 if (iter->tid < iter->head->b.core.tid) { // come to a new reference sequence
227                         iter->tid = iter->head->b.core.tid; iter->pos = iter->head->beg; // jump to the next reference
228                 } else if (iter->pos < iter->head->beg) { // here: tid == head->b.core.tid
229                         iter->pos = iter->head->beg; // jump to the next position
230                 } else ++iter->pos; // scan contiguously
231                 // return
232                 if (n_plp) return iter->plp;
233                 if (iter->is_eof && iter->head->next == 0) break;
234         }
235         return 0;
236 }
237
238 int bam_plp_push(bam_plp_t iter, const bam1_t *b)
239 {
240         if (iter->error) return -1;
241         if (b) {
242                 if (b->core.tid < 0) return 0;
243                 if (b->core.flag & iter->flag_mask) return 0;
244                 if (iter->tid == b->core.tid && iter->pos == b->core.pos && iter->mp->cnt > iter->maxcnt) return 0;
245                 bam_copy1(&iter->tail->b, b);
246                 iter->tail->beg = b->core.pos; iter->tail->end = bam_calend(&b->core, bam1_cigar(b));
247                 iter->tail->s = g_cstate_null; iter->tail->s.end = iter->tail->end - 1; // initialize cstate_t
248                 if (b->core.tid < iter->max_tid) {
249                         fprintf(stderr, "[bam_pileup_core] the input is not sorted (chromosomes out of order)\n");
250                         iter->error = 1;
251                         return -1;
252                 }
253                 if ((b->core.tid == iter->max_tid) && (iter->tail->beg < iter->max_pos)) {
254                         fprintf(stderr, "[bam_pileup_core] the input is not sorted (reads out of order)\n");
255                         iter->error = 1;
256                         return -1;
257                 }
258                 iter->max_tid = b->core.tid; iter->max_pos = iter->tail->beg;
259                 if (iter->tail->end > iter->pos || iter->tail->b.core.tid > iter->tid) {
260                         iter->tail->next = mp_alloc(iter->mp);
261                         iter->tail = iter->tail->next;
262                 }
263         } else iter->is_eof = 1;
264         return 0;
265 }
266
267 const bam_pileup1_t *bam_plp_auto(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp)
268 {
269         const bam_pileup1_t *plp;
270         if (iter->func == 0 || iter->error) { *_n_plp = -1; return 0; }
271         if ((plp = bam_plp_next(iter, _tid, _pos, _n_plp)) != 0) return plp;
272         else { // no pileup line can be obtained; read alignments
273                 *_n_plp = 0;
274                 if (iter->is_eof) return 0;
275                 while (iter->func(iter->data, iter->b) >= 0) {
276                         if (bam_plp_push(iter, iter->b) < 0) {
277                                 *_n_plp = -1;
278                                 return 0;
279                         }
280                         if ((plp = bam_plp_next(iter, _tid, _pos, _n_plp)) != 0) return plp;
281                         // otherwise no pileup line can be returned; read the next alignment.
282                 }
283                 bam_plp_push(iter, 0);
284                 if ((plp = bam_plp_next(iter, _tid, _pos, _n_plp)) != 0) return plp;
285                 return 0;
286         }
287 }
288
289 void bam_plp_reset(bam_plp_t iter)
290 {
291         lbnode_t *p, *q;
292         iter->max_tid = iter->max_pos = -1;
293         iter->tid = iter->pos = 0;
294         iter->is_eof = 0;
295         for (p = iter->head; p->next;) {
296                 q = p->next;
297                 mp_free(iter->mp, p);
298                 p = q;
299         }
300         iter->head = iter->tail;
301 }
302
303 void bam_plp_set_mask(bam_plp_t iter, int mask)
304 {
305         iter->flag_mask = mask < 0? BAM_DEF_MASK : (BAM_FUNMAP | mask);
306 }
307
308 void bam_plp_set_maxcnt(bam_plp_t iter, int maxcnt)
309 {
310         iter->maxcnt = maxcnt;
311 }
312
313 /*****************
314  * callback APIs *
315  *****************/
316
317 int bam_pileup_file(bamFile fp, int mask, bam_pileup_f func, void *func_data)
318 {
319         bam_plbuf_t *buf;
320         int ret;
321         bam1_t *b;
322         b = bam_init1();
323         buf = bam_plbuf_init(func, func_data);
324         bam_plbuf_set_mask(buf, mask);
325         while ((ret = bam_read1(fp, b)) >= 0)
326                 bam_plbuf_push(b, buf);
327         bam_plbuf_push(0, buf);
328         bam_plbuf_destroy(buf);
329         bam_destroy1(b);
330         return 0;
331 }
332
333 void bam_plbuf_set_mask(bam_plbuf_t *buf, int mask)
334 {
335         bam_plp_set_mask(buf->iter, mask);
336 }
337
338 void bam_plbuf_reset(bam_plbuf_t *buf)
339 {
340         bam_plp_reset(buf->iter);
341 }
342
343 bam_plbuf_t *bam_plbuf_init(bam_pileup_f func, void *data)
344 {
345         bam_plbuf_t *buf;
346         buf = calloc(1, sizeof(bam_plbuf_t));
347         buf->iter = bam_plp_init(0, 0);
348         buf->func = func;
349         buf->data = data;
350         return buf;
351 }
352
353 void bam_plbuf_destroy(bam_plbuf_t *buf)
354 {
355         bam_plp_destroy(buf->iter);
356         free(buf);
357 }
358
359 int bam_plbuf_push(const bam1_t *b, bam_plbuf_t *buf)
360 {
361         int ret, n_plp, tid, pos;
362         const bam_pileup1_t *plp;
363         ret = bam_plp_push(buf->iter, b);
364         if (ret < 0) return ret;
365         while ((plp = bam_plp_next(buf->iter, &tid, &pos, &n_plp)) != 0)
366                 buf->func(tid, pos, n_plp, plp, buf->data);
367         return 0;
368 }
369
370 /***********
371  * mpileup *
372  ***********/
373
374 struct __bam_mplp_t {
375         int n;
376         uint64_t min, *pos;
377         bam_plp_t *iter;
378         int *n_plp;
379         const bam_pileup1_t **plp;
380 };
381
382 bam_mplp_t bam_mplp_init(int n, bam_plp_auto_f func, void **data)
383 {
384         int i;
385         bam_mplp_t iter;
386         iter = calloc(1, sizeof(struct __bam_mplp_t));
387         iter->pos = calloc(n, 8);
388         iter->n_plp = calloc(n, sizeof(int));
389         iter->plp = calloc(n, sizeof(void*));
390         iter->iter = calloc(n, sizeof(void*));
391         iter->n = n;
392         iter->min = (uint64_t)-1;
393         for (i = 0; i < n; ++i) {
394                 iter->iter[i] = bam_plp_init(func, data[i]);
395                 iter->pos[i] = iter->min;
396         }
397         return iter;
398 }
399
400 void bam_mplp_set_maxcnt(bam_mplp_t iter, int maxcnt)
401 {
402         int i;
403         for (i = 0; i < iter->n; ++i)
404                 iter->iter[i]->maxcnt = maxcnt;
405 }
406
407 void bam_mplp_destroy(bam_mplp_t iter)
408 {
409         int i;
410         for (i = 0; i < iter->n; ++i) bam_plp_destroy(iter->iter[i]);
411         free(iter->iter); free(iter->pos); free(iter->n_plp); free(iter->plp);
412         free(iter);
413 }
414
415 int bam_mplp_auto(bam_mplp_t iter, int *_tid, int *_pos, int *n_plp, const bam_pileup1_t **plp)
416 {
417         int i, ret = 0;
418         uint64_t new_min = (uint64_t)-1;
419         for (i = 0; i < iter->n; ++i) {
420                 if (iter->pos[i] == iter->min) {
421                         int tid, pos;
422                         iter->plp[i] = bam_plp_auto(iter->iter[i], &tid, &pos, &iter->n_plp[i]);
423                         iter->pos[i] = (uint64_t)tid<<32 | pos;
424                 }
425                 if (iter->plp[i] && iter->pos[i] < new_min) new_min = iter->pos[i];
426         }
427         iter->min = new_min;
428         if (new_min == (uint64_t)-1) return 0;
429         *_tid = new_min>>32; *_pos = (uint32_t)new_min;
430         for (i = 0; i < iter->n; ++i) {
431                 if (iter->pos[i] == iter->min) { // FIXME: valgrind reports "uninitialised value(s) at this line"
432                         n_plp[i] = iter->n_plp[i], plp[i] = iter->plp[i];
433                         ++ret;
434                 } else n_plp[i] = 0, plp[i] = 0;
435         }
436         return ret;
437 }