Copied from SVN revision 4710 and adapted to Git.
[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 __linkbuf_t {
8         bam1_t b;
9         uint32_t beg, end;
10         struct __linkbuf_t *next;
11 } lbnode_t;
12
13 /* --- BEGIN: Memory pool */
14
15 typedef struct {
16         int cnt, n, max;
17         lbnode_t **buf;
18 } mempool_t;
19
20 static mempool_t *mp_init()
21 {
22         mempool_t *mp;
23         mp = (mempool_t*)calloc(1, sizeof(mempool_t));
24         return mp;
25 }
26 static void mp_destroy(mempool_t *mp)
27 {
28         int k;
29         for (k = 0; k < mp->n; ++k) {
30                 free(mp->buf[k]->b.data);
31                 free(mp->buf[k]);
32         }
33         free(mp->buf);
34         free(mp);
35 }
36 static inline lbnode_t *mp_alloc(mempool_t *mp)
37 {
38         ++mp->cnt;
39         if (mp->n == 0) return (lbnode_t*)calloc(1, sizeof(lbnode_t));
40         else return mp->buf[--mp->n];
41 }
42 static inline void mp_free(mempool_t *mp, lbnode_t *p)
43 {
44         --mp->cnt; p->next = 0; // clear lbnode_t::next here
45         if (mp->n == mp->max) {
46                 mp->max = mp->max? mp->max<<1 : 256;
47                 mp->buf = (lbnode_t**)realloc(mp->buf, sizeof(lbnode_t*) * mp->max);
48         }
49         mp->buf[mp->n++] = p;
50 }
51
52 /* --- END: Memory pool */
53
54 /* --- BEGIN: Auxiliary functions */
55
56 static inline int resolve_cigar(bam_pileup1_t *p, uint32_t pos)
57 {
58         unsigned k;
59         bam1_t *b = p->b;
60         bam1_core_t *c = &b->core;
61         uint32_t x = c->pos, y = 0;
62         int ret = 1, is_restart = 1;
63
64         if (c->flag&BAM_FUNMAP) return 0; // unmapped read
65         assert(x <= pos); // otherwise a bug
66         p->qpos = -1; p->indel = 0; p->is_del = p->is_head = p->is_tail = 0;
67         for (k = 0; k < c->n_cigar; ++k) {
68                 int op = bam1_cigar(b)[k] & BAM_CIGAR_MASK; // operation
69                 int l = bam1_cigar(b)[k] >> BAM_CIGAR_SHIFT; // length
70                 if (op == BAM_CMATCH) { // NOTE: this assumes the first and the last operation MUST BE a match or a clip
71                         if (x + l > pos) { // overlap with pos
72                                 p->indel = p->is_del = 0;
73                                 p->qpos = y + (pos - x);
74                                 if (x == pos && is_restart) p->is_head = 1;
75                                 if (x + l - 1 == pos) { // come to the end of a match
76                                         int has_next_match = 0;
77                                         unsigned i;
78                                         for (i = k + 1; i < c->n_cigar; ++i) {
79                                                 uint32_t cigar = bam1_cigar(b)[i];
80                                                 int opi = cigar&BAM_CIGAR_MASK;
81                                                 if (opi == BAM_CMATCH) {
82                                                         has_next_match = 1;
83                                                         break;
84                                                 } else if (opi == BAM_CSOFT_CLIP || opi == BAM_CREF_SKIP || opi == BAM_CHARD_CLIP) break;
85                                         }
86                                         if (!has_next_match) p->is_tail = 1;
87                                         if (k < c->n_cigar - 1 && has_next_match) { // there are additional operation(s)
88                                                 uint32_t cigar = bam1_cigar(b)[k+1]; // next CIGAR
89                                                 int op_next = cigar&BAM_CIGAR_MASK; // next CIGAR operation
90                                                 if (op_next == BAM_CDEL) p->indel = -(int32_t)(cigar>>BAM_CIGAR_SHIFT); // del
91                                                 else if (op_next == BAM_CINS) p->indel = cigar>>BAM_CIGAR_SHIFT; // ins
92                                                 else if (op_next == BAM_CPAD && k + 2 < c->n_cigar) { // no working for adjacent padding
93                                                         cigar = bam1_cigar(b)[k+2]; op_next = cigar&BAM_CIGAR_MASK;
94                                                         if (op_next == BAM_CDEL) p->indel = -(int32_t)(cigar>>BAM_CIGAR_SHIFT); // del
95                                                         else if (op_next == BAM_CINS) p->indel = cigar>>BAM_CIGAR_SHIFT; // ins
96                                                 }
97                                         }
98                                 }
99                         }
100                         x += l; y += l;
101                 } else if (op == BAM_CDEL) { // then set ->is_del
102                         if (x + l > pos) {
103                                 p->indel = 0; p->is_del = 1;
104                                 p->qpos = y + (pos - x);
105                         }
106                         x += l;
107                 } else if (op == BAM_CREF_SKIP) x += l;
108                 else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += l;
109                 if (is_restart) is_restart ^= (op == BAM_CMATCH);
110                 else is_restart ^= (op == BAM_CREF_SKIP || op == BAM_CSOFT_CLIP || op == BAM_CHARD_CLIP);
111                 if (x > pos) {
112                         if (op == BAM_CREF_SKIP) ret = 0; // then do not put it into pileup at all
113                         break;
114                 }
115         }
116         assert(x > pos); // otherwise a bug
117         return ret;
118 }
119
120 /* --- END: Auxiliary functions */
121
122 /*******************
123  * pileup iterator *
124  *******************/
125
126 struct __bam_plp_t {
127         mempool_t *mp;
128         lbnode_t *head, *tail, *dummy;
129         int32_t tid, pos, max_tid, max_pos;
130         int is_eof, flag_mask, max_plp, error;
131         bam_pileup1_t *plp;
132         // for the "auto" interface only
133         bam1_t *b;
134         bam_plp_auto_f func;
135         void *data;
136 };
137
138 bam_plp_t bam_plp_init(bam_plp_auto_f func, void *data)
139 {
140         bam_plp_t iter;
141         iter = calloc(1, sizeof(struct __bam_plp_t));
142         iter->mp = mp_init();
143         iter->head = iter->tail = mp_alloc(iter->mp);
144         iter->dummy = mp_alloc(iter->mp);
145         iter->max_tid = iter->max_pos = -1;
146         iter->flag_mask = BAM_DEF_MASK;
147         if (func) {
148                 iter->func = func;
149                 iter->data = data;
150                 iter->b = bam_init1();
151         }
152         return iter;
153 }
154
155 void bam_plp_destroy(bam_plp_t iter)
156 {
157         mp_free(iter->mp, iter->dummy);
158         mp_free(iter->mp, iter->head);
159         if (iter->mp->cnt != 0)
160                 fprintf(stderr, "[bam_plp_destroy] memory leak: %d. Continue anyway.\n", iter->mp->cnt);
161         mp_destroy(iter->mp);
162         if (iter->b) bam_destroy1(iter->b);
163         free(iter->plp);
164         free(iter);
165 }
166
167 const bam_pileup1_t *bam_plp_next(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp)
168 {
169         if (iter->error) { *_n_plp = -1; return 0; }
170         *_n_plp = 0;
171         if (iter->is_eof && iter->head->next == 0) return 0;
172         while (iter->is_eof || iter->max_tid > iter->tid || (iter->max_tid == iter->tid && iter->max_pos > iter->pos)) {
173                 int n_plp = 0;
174                 lbnode_t *p, *q;
175                 // write iter->plp at iter->pos
176                 iter->dummy->next = iter->head;
177                 for (p = iter->head, q = iter->dummy; p->next; q = p, p = p->next) {
178                         if (p->b.core.tid < iter->tid || (p->b.core.tid == iter->tid && p->end <= iter->pos)) { // then remove
179                                 q->next = p->next; mp_free(iter->mp, p); p = q;
180                         } else if (p->b.core.tid == iter->tid && p->beg <= iter->pos) { // here: p->end > pos; then add to pileup
181                                 if (n_plp == iter->max_plp) { // then double the capacity
182                                         iter->max_plp = iter->max_plp? iter->max_plp<<1 : 256;
183                                         iter->plp = (bam_pileup1_t*)realloc(iter->plp, sizeof(bam_pileup1_t) * iter->max_plp);
184                                 }
185                                 iter->plp[n_plp].b = &p->b;
186                                 if (resolve_cigar(iter->plp + n_plp, iter->pos)) ++n_plp; // skip the read if we are looking at ref-skip
187                         }
188                 }
189                 iter->head = iter->dummy->next; // dummy->next may be changed
190                 *_n_plp = n_plp; *_tid = iter->tid; *_pos = iter->pos;
191                 // update iter->tid and iter->pos
192                 if (iter->head->next) {
193                         if (iter->tid > iter->head->b.core.tid) {
194                                 fprintf(stderr, "[%s] unsorted input. Pileup aborts.\n", __func__);
195                                 iter->error = 1;
196                                 *_n_plp = -1;
197                                 return 0;
198                         }
199                 }
200                 if (iter->tid < iter->head->b.core.tid) { // come to a new reference sequence
201                         iter->tid = iter->head->b.core.tid; iter->pos = iter->head->beg; // jump to the next reference
202                 } else if (iter->pos < iter->head->beg) { // here: tid == head->b.core.tid
203                         iter->pos = iter->head->beg; // jump to the next position
204                 } else ++iter->pos; // scan contiguously
205                 // return
206                 if (n_plp) return iter->plp;
207                 if (iter->is_eof && iter->head->next == 0) break;
208         }
209         return 0;
210 }
211
212 int bam_plp_push(bam_plp_t iter, const bam1_t *b)
213 {
214         if (iter->error) return -1;
215         if (b) {
216                 if (b->core.tid < 0) return 0;
217                 if (b->core.flag & iter->flag_mask) return 0;
218                 bam_copy1(&iter->tail->b, b);
219                 iter->tail->beg = b->core.pos; iter->tail->end = bam_calend(&b->core, bam1_cigar(b));
220                 if (b->core.tid < iter->max_tid) {
221                         fprintf(stderr, "[bam_pileup_core] the input is not sorted (chromosomes out of order)\n");
222                         iter->error = 1;
223                         return -1;
224                 }
225                 if ((b->core.tid == iter->max_tid) && (iter->tail->beg < iter->max_pos)) {
226                         fprintf(stderr, "[bam_pileup_core] the input is not sorted (reads out of order)\n");
227                         iter->error = 1;
228                         return -1;
229                 }
230                 iter->max_tid = b->core.tid; iter->max_pos = iter->tail->beg;
231                 if (iter->tail->end > iter->pos || iter->tail->b.core.tid > iter->tid) {
232                         iter->tail->next = mp_alloc(iter->mp);
233                         iter->tail = iter->tail->next;
234                 }
235         } else iter->is_eof = 1;
236         return 0;
237 }
238
239 const bam_pileup1_t *bam_plp_auto(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp)
240 {
241         const bam_pileup1_t *plp;
242         if (iter->func == 0 || iter->error) { *_n_plp = -1; return 0; }
243         if ((plp = bam_plp_next(iter, _tid, _pos, _n_plp)) != 0) return plp;
244         else {
245                 *_n_plp = 0;
246                 if (iter->is_eof) return 0;
247                 while (iter->func(iter->data, iter->b) >= 0) {
248                         if (bam_plp_push(iter, iter->b) < 0) {
249                                 *_n_plp = -1;
250                                 return 0;
251                         }
252                         if ((plp = bam_plp_next(iter, _tid, _pos, _n_plp)) != 0) return plp;
253                 }
254                 bam_plp_push(iter, 0);
255                 if ((plp = bam_plp_next(iter, _tid, _pos, _n_plp)) != 0) return plp;
256                 return 0;
257         }
258 }
259
260 void bam_plp_reset(bam_plp_t iter)
261 {
262         lbnode_t *p, *q;
263         iter->max_tid = iter->max_pos = -1;
264         iter->tid = iter->pos = 0;
265         iter->is_eof = 0;
266         for (p = iter->head; p->next;) {
267                 q = p->next;
268                 mp_free(iter->mp, p);
269                 p = q;
270         }
271         iter->head = iter->tail;
272 }
273
274 void bam_plp_set_mask(bam_plp_t iter, int mask)
275 {
276         iter->flag_mask = mask < 0? BAM_DEF_MASK : (BAM_FUNMAP | mask);
277 }
278
279 /*****************
280  * callback APIs *
281  *****************/
282
283 int bam_pileup_file(bamFile fp, int mask, bam_pileup_f func, void *func_data)
284 {
285         bam_plbuf_t *buf;
286         int ret;
287         bam1_t *b;
288         b = bam_init1();
289         buf = bam_plbuf_init(func, func_data);
290         bam_plbuf_set_mask(buf, mask);
291         while ((ret = bam_read1(fp, b)) >= 0)
292                 bam_plbuf_push(b, buf);
293         bam_plbuf_push(0, buf);
294         bam_plbuf_destroy(buf);
295         bam_destroy1(b);
296         return 0;
297 }
298
299 void bam_plbuf_set_mask(bam_plbuf_t *buf, int mask)
300 {
301         bam_plp_set_mask(buf->iter, mask);
302 }
303
304 void bam_plbuf_reset(bam_plbuf_t *buf)
305 {
306         bam_plp_reset(buf->iter);
307 }
308
309 bam_plbuf_t *bam_plbuf_init(bam_pileup_f func, void *data)
310 {
311         bam_plbuf_t *buf;
312         buf = calloc(1, sizeof(bam_plbuf_t));
313         buf->iter = bam_plp_init(0, 0);
314         buf->func = func;
315         buf->data = data;
316         return buf;
317 }
318
319 void bam_plbuf_destroy(bam_plbuf_t *buf)
320 {
321         bam_plp_destroy(buf->iter);
322         free(buf);
323 }
324
325 int bam_plbuf_push(const bam1_t *b, bam_plbuf_t *buf)
326 {
327         int ret, n_plp, tid, pos;
328         const bam_pileup1_t *plp;
329         ret = bam_plp_push(buf->iter, b);
330         if (ret < 0) return ret;
331         while ((plp = bam_plp_next(buf->iter, &tid, &pos, &n_plp)) != 0)
332                 buf->func(tid, pos, n_plp, plp, buf->data);
333         return 0;
334 }
335
336 /***********
337  * mpileup *
338  ***********/
339
340 struct __bam_mplp_t {
341         int n;
342         uint64_t min, *pos;
343         bam_plp_t *iter;
344         int *n_plp;
345         const bam_pileup1_t **plp;
346 };
347
348 bam_mplp_t bam_mplp_init(int n, bam_plp_auto_f func, void **data)
349 {
350         int i;
351         bam_mplp_t iter;
352         iter = calloc(1, sizeof(struct __bam_mplp_t));
353         iter->pos = calloc(n, 8);
354         iter->n_plp = calloc(n, sizeof(int));
355         iter->plp = calloc(n, sizeof(void*));
356         iter->iter = calloc(n, sizeof(void*));
357         iter->n = n;
358         iter->min = (uint64_t)-1;
359         for (i = 0; i < n; ++i) {
360                 iter->iter[i] = bam_plp_init(func, data[i]);
361                 iter->pos[i] = iter->min;
362         }
363         return iter;
364 }
365
366 void bam_mplp_destroy(bam_mplp_t iter)
367 {
368         int i;
369         for (i = 0; i < iter->n; ++i) bam_plp_destroy(iter->iter[i]);
370         free(iter->iter); free(iter->pos); free(iter->n_plp); free(iter->plp);
371         free(iter);
372 }
373
374 int bam_mplp_auto(bam_mplp_t iter, int *_tid, int *_pos, int *n_plp, const bam_pileup1_t **plp)
375 {
376         int i, ret = 0;
377         uint64_t new_min = (uint64_t)-1;
378         for (i = 0; i < iter->n; ++i) {
379                 if (iter->pos[i] == iter->min) {
380                         int tid, pos;
381                         iter->plp[i] = bam_plp_auto(iter->iter[i], &tid, &pos, &iter->n_plp[i]);
382                         iter->pos[i] = (uint64_t)tid<<32 | pos;
383                 }
384                 if (iter->plp[i] && iter->pos[i] < new_min) new_min = iter->pos[i];
385         }
386         iter->min = new_min;
387         if (new_min == (uint64_t)-1) return 0;
388         *_tid = new_min>>32; *_pos = (uint32_t)new_min;
389         for (i = 0; i < iter->n; ++i) {
390                 if (iter->pos[i] == iter->min) {
391                         n_plp[i] = iter->n_plp[i], plp[i] = iter->plp[i];
392                         ++ret;
393                 } else n_plp[i] = 0, plp[i] = 0;
394         }
395         return ret;
396 }