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