Imported Upstream version 0.2
[pysam.git] / pysam / pysam_util.c
1 #include <ctype.h>
2 #include <assert.h>
3 #include "bam.h"
4 #include "khash.h"
5 #include "ksort.h"
6 #include "bam_endian.h"
7 #include "knetfile.h"
8 #include "pysam_util.h"
9
10 // #######################################################
11 // utility routines to avoid using callbacks in bam_fetch
12 // taken from bam_index.c
13 // The order of the following declarations is important.
14 // #######################################################
15
16 typedef struct
17 {
18   uint64_t u, v;
19 } pair64_t;
20
21 #define pair64_lt(a,b) ((a).u < (b).u)
22
23 typedef struct {
24         uint32_t m, n;
25         pair64_t *list;
26 } bam_binlist_t;
27
28 typedef struct {
29         int32_t n, m;
30         uint64_t *offset;
31 } bam_lidx_t;
32
33 KSORT_INIT(my_off, pair64_t, pair64_lt);
34 KHASH_MAP_INIT_INT(my_i, bam_binlist_t);
35
36 struct __bam_index_t
37 {
38   int32_t n;
39   khash_t(my_i) **index;
40   bam_lidx_t *index2;
41 };
42
43 typedef struct __linkbuf_t {
44         bam1_t b;
45         uint32_t beg, end;
46         struct __linkbuf_t *next;
47 } lbnode_t;
48
49 typedef struct {
50         int cnt, n, max;
51         lbnode_t **buf;
52 } mempool_t;
53
54 struct __bam_plbuf_t {
55         mempool_t *mp;
56         lbnode_t *head, *tail, *dummy;
57         bam_pileup_f func;
58         void *func_data;
59         int32_t tid, pos, max_tid, max_pos;
60         int max_pu, is_eof;
61         bam_pileup1_t *pu;
62         int flag_mask;
63 };
64
65 static mempool_t *mp_init()
66 {
67         mempool_t *mp;
68         mp = (mempool_t*)calloc(1, sizeof(mempool_t));
69         return mp;
70 }
71 static void mp_destroy(mempool_t *mp)
72 {
73         int k;
74         for (k = 0; k < mp->n; ++k) {
75                 free(mp->buf[k]->b.data);
76                 free(mp->buf[k]);
77         }
78         free(mp->buf);
79         free(mp);
80 }
81 static inline lbnode_t *mp_alloc(mempool_t *mp)
82 {
83         ++mp->cnt;
84         if (mp->n == 0) return (lbnode_t*)calloc(1, sizeof(lbnode_t));
85         else return mp->buf[--mp->n];
86 }
87 static inline void mp_free(mempool_t *mp, lbnode_t *p)
88 {
89         --mp->cnt; p->next = 0; // clear lbnode_t::next here
90         if (mp->n == mp->max) {
91                 mp->max = mp->max? mp->max<<1 : 256;
92                 mp->buf = (lbnode_t**)realloc(mp->buf, sizeof(lbnode_t*) * mp->max);
93         }
94         mp->buf[mp->n++] = p;
95 }
96
97 static inline int resolve_cigar(bam_pileup1_t *p, uint32_t pos)
98 {
99         unsigned k;
100         bam1_t *b = p->b;
101         bam1_core_t *c = &b->core;
102         uint32_t x = c->pos, y = 0;
103         int ret = 1, is_restart = 1;
104
105         if (c->flag&BAM_FUNMAP) return 0; // unmapped read
106         assert(x <= pos); // otherwise a bug
107         p->qpos = -1; p->indel = 0; p->is_del = p->is_head = p->is_tail = 0;
108         for (k = 0; k < c->n_cigar; ++k) {
109                 int op = bam1_cigar(b)[k] & BAM_CIGAR_MASK; // operation
110                 int l = bam1_cigar(b)[k] >> BAM_CIGAR_SHIFT; // length
111                 if (op == BAM_CMATCH) { // NOTE: this assumes the first and the last operation MUST BE a match or a clip
112                         if (x + l > pos) { // overlap with pos
113                                 p->indel = p->is_del = 0;
114                                 p->qpos = y + (pos - x);
115                                 if (x == pos && is_restart) p->is_head = 1;
116                                 if (x + l - 1 == pos) { // come to the end of a match
117                                         if (k < c->n_cigar - 1) { // there are additional operation(s)
118                                                 uint32_t cigar = bam1_cigar(b)[k+1]; // next CIGAR
119                                                 int op_next = cigar&BAM_CIGAR_MASK; // next CIGAR operation
120                                                 if (op_next == BAM_CDEL) p->indel = -(int32_t)(cigar>>BAM_CIGAR_SHIFT); // del
121                                                 else if (op_next == BAM_CINS) p->indel = cigar>>BAM_CIGAR_SHIFT; // ins
122                                                 if (op_next == BAM_CSOFT_CLIP || op_next == BAM_CREF_SKIP || op_next == BAM_CHARD_CLIP)
123                                                         p->is_tail = 1; // tail
124                                         } else p->is_tail = 1; // this is the last operation; set tail
125                                 }
126                         }
127                         x += l; y += l;
128                 } else if (op == BAM_CDEL) { // then set ->is_del
129                         if (x + l > pos) {
130                                 p->indel = 0; p->is_del = 1;
131                                 p->qpos = y + (pos - x);
132                         }
133                         x += l;
134                 } else if (op == BAM_CREF_SKIP) x += l;
135                 else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += l;
136                 is_restart = (op == BAM_CREF_SKIP || op == BAM_CSOFT_CLIP || op == BAM_CHARD_CLIP);
137                 if (x > pos) {
138                         if (op == BAM_CREF_SKIP) ret = 0; // then do not put it into pileup at all
139                         break;
140                 }
141         }
142         assert(x > pos); // otherwise a bug
143         return ret;
144 }
145
146
147
148
149 // the following code has been taken from bam_plbuf_push
150 // and modified such that instead of a function call
151 // the function returns and will continue (if cont is true).
152 // from where it left off.
153
154 // returns
155 // 1: if buf is full and can be emitted
156 // 0: if b has been added
157 // -1: if there was an error
158 int pysam_bam_plbuf_push(const bam1_t *b, bam_plbuf_t *buf, int cont)
159 {
160   if (!cont)
161     {
162       if (b) { // fill buffer
163         if (b->core.tid < 0) return 0;
164         if (b->core.flag & buf->flag_mask) return 0;
165         bam_copy1(&buf->tail->b, b);
166         buf->tail->beg = b->core.pos; buf->tail->end = bam_calend(&b->core, bam1_cigar(b));
167         if (!(b->core.tid >= buf->max_tid || (b->core.tid == buf->max_tid && buf->tail->beg >= buf->max_pos))) {
168           fprintf(stderr, "[bam_pileup_core] the input is not sorted. Abort!\n");
169           abort();
170         }
171         buf->max_tid = b->core.tid; buf->max_pos = buf->tail->beg;
172         if (buf->tail->end > buf->pos || buf->tail->b.core.tid > buf->tid) {
173           buf->tail->next = mp_alloc(buf->mp);
174           buf->tail = buf->tail->next;
175         }
176       } else buf->is_eof = 1;
177     }
178   else
179     // continue end of loop
180     {
181       // update tid and pos
182       if (buf->head->next) {
183         if (buf->tid > buf->head->b.core.tid) {
184           fprintf(stderr, "[bam_plbuf_push] unsorted input. Pileup aborts.\n");
185           return -1;
186         }
187       }
188       if (buf->tid < buf->head->b.core.tid) { // come to a new reference sequence
189         buf->tid = buf->head->b.core.tid; buf->pos = buf->head->beg; // jump to the next reference
190       } else if (buf->pos < buf->head->beg) { // here: tid == head->b.core.tid
191         buf->pos = buf->head->beg; // jump to the next position
192       } else ++buf->pos; // scan contiguously
193       if (buf->is_eof && buf->head->next == 0) return 0;
194     }
195
196   // enter yield loop
197   while (buf->is_eof || buf->max_tid > buf->tid || (buf->max_tid == buf->tid && buf->max_pos > buf->pos))
198     {
199       int n_pu = 0;
200       lbnode_t *p, *q;
201       buf->dummy->next = buf->head;
202       for (p = buf->head, q = buf->dummy; p->next; q = p, p = p->next) {
203         if (p->b.core.tid < buf->tid || (p->b.core.tid == buf->tid && p->end <= buf->pos)) { // then remove from the list
204           q->next = p->next; mp_free(buf->mp, p); p = q;
205         } else if (p->b.core.tid == buf->tid && p->beg <= buf->pos) { // here: p->end > pos; then add to pileup
206           if (n_pu == buf->max_pu) { // then double the capacity
207             buf->max_pu = buf->max_pu? buf->max_pu<<1 : 256;
208             buf->pu = (bam_pileup1_t*)realloc(buf->pu, sizeof(bam_pileup1_t) * buf->max_pu);
209           }
210           buf->pu[n_pu].b = &p->b;
211           if (resolve_cigar(buf->pu + n_pu, buf->pos)) ++n_pu; // skip the read if we are looking at BAM_CREF_SKIP
212         }
213       }
214       buf->head = buf->dummy->next; // dummy->next may be changed
215
216       // exit if alignments need to be emitted
217       if (n_pu) { return n_pu; }
218       
219       // update tid and pos
220       if (buf->head->next) {
221         if (buf->tid > buf->head->b.core.tid) {
222           fprintf(stderr, "[bam_plbuf_push] unsorted input. Pileup aborts.\n");
223           return -2;
224         }
225       }
226       if (buf->tid < buf->head->b.core.tid) { // come to a new reference sequence
227         buf->tid = buf->head->b.core.tid; buf->pos = buf->head->beg; // jump to the next reference
228       } else if (buf->pos < buf->head->beg) { // here: tid == head->b.core.tid
229         buf->pos = buf->head->beg; // jump to the next position
230       } else ++buf->pos; // scan contiguously
231       if (buf->is_eof && buf->head->next == 0) break;
232     }
233   return 0;
234 }
235
236 int pysam_get_pos( const bam_plbuf_t *buf) 
237 {
238   return buf->pos;
239 }
240
241   
242 int pysam_get_tid( const bam_plbuf_t *buf)
243 {
244   return buf->tid;
245 }
246
247 bam_pileup1_t * pysam_get_pileup( const bam_plbuf_t *buf)
248 {
249   return buf->pu;
250 }
251
252 // pysam dispatch function to emulate the samtools
253 // command line within python.
254 // taken from the main function in bamtk.c
255 // added code to reset getopt
256 extern int main_samview(int argc, char *argv[]);
257 extern int main_import(int argc, char *argv[]);
258 extern int bam_pileup(int argc, char *argv[]);
259 extern int bam_merge(int argc, char *argv[]);
260 extern int bam_sort(int argc, char *argv[]);
261 extern int bam_index(int argc, char *argv[]);
262 extern int faidx_main(int argc, char *argv[]);
263 extern int bam_mating(int argc, char *argv[]);
264 extern int bam_rmdup(int argc, char *argv[]);
265 extern int glf3_view_main(int argc, char *argv[]);
266 extern int bam_flagstat(int argc, char *argv[]);
267 extern int bam_fillmd(int argc, char *argv[]);
268
269 int pysam_dispatch(int argc, char *argv[] )
270 {
271
272 #ifdef _WIN32
273   setmode(fileno(stdout), O_BINARY);
274   setmode(fileno(stdin),  O_BINARY);
275 #ifdef _USE_KNETFILE
276   knet_win32_init();
277 #endif
278 #endif
279
280   extern int optind;
281   
282   // reset getop
283   optind = 1;
284
285   if (argc < 2) return 1;
286
287   if (strcmp(argv[1], "view") == 0) return main_samview(argc-1, argv+1);
288   else if (strcmp(argv[1], "import") == 0) return main_import(argc-1, argv+1);
289   else if (strcmp(argv[1], "pileup") == 0) return bam_pileup(argc-1, argv+1);
290   else if (strcmp(argv[1], "merge") == 0) return bam_merge(argc-1, argv+1);
291   else if (strcmp(argv[1], "sort") == 0) return bam_sort(argc-1, argv+1);
292   else if (strcmp(argv[1], "index") == 0) return bam_index(argc-1, argv+1);
293   else if (strcmp(argv[1], "faidx") == 0) return faidx_main(argc-1, argv+1);
294   else if (strcmp(argv[1], "fixmate") == 0) return bam_mating(argc-1, argv+1);
295   else if (strcmp(argv[1], "rmdup") == 0) return bam_rmdup(argc-1, argv+1);
296   else if (strcmp(argv[1], "glfview") == 0) return glf3_view_main(argc-1, argv+1);
297   else if (strcmp(argv[1], "flagstat") == 0) return bam_flagstat(argc-1, argv+1);
298   else if (strcmp(argv[1], "calmd") == 0) return bam_fillmd(argc-1, argv+1);
299   else if (strcmp(argv[1], "fillmd") == 0) return bam_fillmd(argc-1, argv+1);
300
301 #if _CURSES_LIB != 0
302   else if (strcmp(argv[1], "tview") == 0) return bam_tview_main(argc-1, argv+1);
303 #endif
304   else 
305     {
306       fprintf(stderr, "[main] unrecognized command '%s'\n", argv[1]);
307       return 1;
308     }
309   return 0;
310 }
311
312 // standin for bam_destroy1 in bam.h
313 // deletes all variable length data
314 void pysam_bam_destroy1( bam1_t * b )
315 {
316   if (b == NULL) return;
317   if (b->data != NULL) free(b->data);
318   free(b);
319 }
320
321 // taken from samtools/bam_import.c
322 static inline uint8_t *alloc_data(bam1_t *b, size_t size)
323 {
324   if (b->m_data < size)
325     {
326       b->m_data = size;
327       kroundup32(b->m_data);
328       b->data = (uint8_t*)realloc(b->data, b->m_data);
329     }
330   return b->data;
331 }
332
333 // update the variable length data within a bam1_t entry.
334 // Adds *nbytes_new* - *nbytes_old* into the variable length data of *src* at *pos*.
335 // Data within the bam1_t entry is moved so that it is
336 // consistent with the data field lengths.
337 bam1_t * pysam_bam_update( bam1_t * b,
338                            const size_t nbytes_old,
339                            const size_t nbytes_new, 
340                            uint8_t * pos )
341 {
342   int d = nbytes_new-nbytes_old;
343
344   // no change
345   if (d == 0) return b;
346
347   int new_size = d + b->data_len;
348   size_t offset = pos - b->data;
349
350   //printf("d=%i, old=%i, new=%i, old_size=%i, new_size=%i\n",
351   // d, nbytes_old, nbytes_new, b->data_len, new_size);
352   
353   // increase memory if required
354   if (d > 0)
355     {
356       alloc_data( b, new_size );
357       pos = b->data + offset;
358     }
359   
360   if (b->data_len != 0)
361     {
362       if (offset < 0 || offset > b->data_len)
363         fprintf(stderr, "[pysam_bam_insert] illegal offset: '%i'\n", (int)offset);
364     }
365   
366   // printf("dest=%p, src=%p, n=%i\n", pos+nbytes_new, pos + nbytes_old, b->data_len - (offset+nbytes_old));
367   memmove( pos + nbytes_new,
368            pos + nbytes_old,
369            b->data_len - (offset + nbytes_old));
370     
371   b->data_len = new_size;
372       
373   return b;
374 }
375
376 // translate a nucleotide character to binary code
377 unsigned char pysam_translate_sequence( const unsigned char s )
378 {
379   return bam_nt16_table[s];
380 }
381
382 // stand-ins for samtools macros in bam.h
383 char * pysam_bam1_qname( const bam1_t * b)
384 {
385   return (char*)b->data;
386 }
387
388 uint32_t * pysam_bam1_cigar( const bam1_t * b) 
389 {
390   return (uint32_t*)(b->data + b->core.l_qname);
391 }
392
393 uint8_t * pysam_bam1_seq( const bam1_t * b) 
394 {
395   return (uint8_t*)(b->data + b->core.n_cigar*4 + b->core.l_qname);
396 }
397
398 uint8_t * pysam_bam1_qual( const bam1_t * b)
399 {
400   return (uint8_t*)(b->data + b->core.n_cigar*4 + b->core.l_qname + (b->core.l_qseq + 1)/2);
401 }
402
403 uint8_t * pysam_bam1_aux( const bam1_t * b)
404 {
405   return (uint8_t*)(b->data + b->core.n_cigar*4 + b->core.l_qname + b->core.l_qseq + (b->core.l_qseq + 1)/2);
406 }
407
408 // #######################################################
409 // Iterator implementation
410 // #######################################################
411
412 // functions defined in bam_index.c
413 extern pair64_t * get_chunk_coordinates(const bam_index_t *idx, int tid, int beg, int end, int* cnt_off);
414
415 static inline int is_overlap(uint32_t beg, uint32_t end, const bam1_t *b)
416 {
417         uint32_t rbeg = b->core.pos;
418         uint32_t rend = b->core.n_cigar? bam_calend(&b->core, bam1_cigar(b)) : b->core.pos + 1;
419         return (rend > beg && rbeg < end);
420 }
421
422 struct __bam_fetch_iterator_t
423 {
424   bam1_t *        b;
425   pair64_t *      off;
426   int             n_off;
427   uint64_t        curr_off;
428   int             curr_chunk;
429   bamFile               fp;
430   int                           tid;
431   int                           beg;
432   int                           end;
433   int             n_seeks;
434 };
435  
436 bam_fetch_iterator_t* bam_init_fetch_iterator(bamFile fp, const bam_index_t *idx, int tid, int beg, int end)
437 {
438         // iterator contains current alignment position
439         //      and will contain actual alignment during iterations
440         bam_fetch_iterator_t* iter  = (bam_fetch_iterator_t*)calloc(1, sizeof(bam_fetch_iterator_t));
441         iter->b                     = (bam1_t*)calloc(1, sizeof(bam1_t));
442                 
443         // list of chunks containing our alignments
444         iter->off = get_chunk_coordinates(idx, tid, beg, end, &iter->n_off);
445         
446         // initialise other state variables in iterator
447         iter->fp                = fp;
448         iter->curr_chunk        = -1;   
449         iter->curr_off          =  0;
450         iter->n_seeks           =  0;    
451         iter->tid                               = tid;
452         iter->beg                               = beg;
453         iter->end                               = end;
454         return iter;
455 }
456
457 bam1_t * bam_fetch_iterate(bam_fetch_iterator_t *iter)
458 {
459         if (!iter->off) {
460                 return 0;
461         }
462
463         int ret;
464         // iterate through all alignments in chunks
465         for (;;) {
466                 if (iter->curr_off == 0 || iter->curr_off >= iter->off[iter->curr_chunk].v) { // then jump to the next chunk
467                         if (iter->curr_chunk == iter->n_off - 1) break; // no more chunks
468                         if (iter->curr_chunk >= 0) assert(iter->curr_off == iter->off[iter->curr_chunk].v); // otherwise bug
469                         if (iter->curr_chunk < 0 || iter->off[iter->curr_chunk].v != iter->off[iter->curr_chunk+1].u) { // not adjacent chunks; then seek
470                                 bam_seek(iter->fp, iter->off[iter->curr_chunk+1].u, SEEK_SET);
471                                 iter->curr_off = bam_tell(iter->fp);
472                                 ++iter->n_seeks;
473                         }
474                         ++iter->curr_chunk;
475                 }
476                 if ((ret = bam_read1(iter->fp, iter->b)) > 0) {
477                         iter->curr_off = bam_tell(iter->fp);
478                         if (iter->b->core.tid != iter->tid || iter->b->core.pos >= iter->end) break; // no need to proceed
479                         else if (is_overlap(iter->beg, iter->end, iter->b)) 
480                                 //
481                                 //func(iter->b, data);
482                                 //
483                                 return iter->b;
484                 } else 
485                         return 0; // end of file
486         }
487         return 0;
488 }
489
490 void bam_cleanup_fetch_iterator(bam_fetch_iterator_t *iter)
491 {
492   //  fprintf(stderr, "[bam_fetch] # seek calls: %d\n", iter->n_seeks);
493   bam_destroy1(iter->b);
494   free(iter->off);
495 }
496
497        
498
499