]> git.donarmstrong.com Git - samtools.git/blob - bam.h
converted padded SAM to unpadded SAM
[samtools.git] / bam.h
1 /* The MIT License
2
3    Copyright (c) 2008-2010 Genome Research Ltd (GRL).
4
5    Permission is hereby granted, free of charge, to any person obtaining
6    a copy of this software and associated documentation files (the
7    "Software"), to deal in the Software without restriction, including
8    without limitation the rights to use, copy, modify, merge, publish,
9    distribute, sublicense, and/or sell copies of the Software, and to
10    permit persons to whom the Software is furnished to do so, subject to
11    the following conditions:
12
13    The above copyright notice and this permission notice shall be
14    included in all copies or substantial portions of the Software.
15
16    THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
17    EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
18    MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND
19    NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS
20    BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN
21    ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN
22    CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
23    SOFTWARE.
24 */
25
26 /* Contact: Heng Li <lh3@sanger.ac.uk> */
27
28 #ifndef BAM_BAM_H
29 #define BAM_BAM_H
30
31 /*!
32   @header
33
34   BAM library provides I/O and various operations on manipulating files
35   in the BAM (Binary Alignment/Mapping) or SAM (Sequence Alignment/Map)
36   format. It now supports importing from or exporting to SAM, sorting,
37   merging, generating pileup, and quickly retrieval of reads overlapped
38   with a specified region.
39
40   @copyright Genome Research Ltd.
41  */
42
43 #define BAM_VERSION "0.1.18 (r982:295)"
44
45 #include <stdint.h>
46 #include <stdlib.h>
47 #include <string.h>
48 #include <stdio.h>
49
50 #ifndef BAM_LITE
51 #define BAM_VIRTUAL_OFFSET16
52 #include "bgzf.h"
53 /*! @abstract BAM file handler */
54 typedef BGZF *bamFile;
55 #define bam_open(fn, mode) bgzf_open(fn, mode)
56 #define bam_dopen(fd, mode) bgzf_fdopen(fd, mode)
57 #define bam_close(fp) bgzf_close(fp)
58 #define bam_read(fp, buf, size) bgzf_read(fp, buf, size)
59 #define bam_write(fp, buf, size) bgzf_write(fp, buf, size)
60 #define bam_tell(fp) bgzf_tell(fp)
61 #define bam_seek(fp, pos, dir) bgzf_seek(fp, pos, dir)
62 #else
63 #define BAM_TRUE_OFFSET
64 #include <zlib.h>
65 typedef gzFile bamFile;
66 #define bam_open(fn, mode) gzopen(fn, mode)
67 #define bam_dopen(fd, mode) gzdopen(fd, mode)
68 #define bam_close(fp) gzclose(fp)
69 #define bam_read(fp, buf, size) gzread(fp, buf, size)
70 /* no bam_write/bam_tell/bam_seek() here */
71 #endif
72
73 /*! @typedef
74   @abstract Structure for the alignment header.
75   @field n_targets   number of reference sequences
76   @field target_name names of the reference sequences
77   @field target_len  lengths of the referene sequences
78   @field dict        header dictionary
79   @field hash        hash table for fast name lookup
80   @field rg2lib      hash table for @RG-ID -> LB lookup
81   @field l_text      length of the plain text in the header
82   @field text        plain text
83
84   @discussion Field hash points to null by default. It is a private
85   member.
86  */
87 typedef struct {
88         int32_t n_targets;
89         char **target_name;
90         uint32_t *target_len;
91         void *dict, *hash, *rg2lib;
92         size_t l_text, n_text;
93         char *text;
94 } bam_header_t;
95
96 /*! @abstract the read is paired in sequencing, no matter whether it is mapped in a pair */
97 #define BAM_FPAIRED        1
98 /*! @abstract the read is mapped in a proper pair */
99 #define BAM_FPROPER_PAIR   2
100 /*! @abstract the read itself is unmapped; conflictive with BAM_FPROPER_PAIR */
101 #define BAM_FUNMAP         4
102 /*! @abstract the mate is unmapped */
103 #define BAM_FMUNMAP        8
104 /*! @abstract the read is mapped to the reverse strand */
105 #define BAM_FREVERSE      16
106 /*! @abstract the mate is mapped to the reverse strand */
107 #define BAM_FMREVERSE     32
108 /*! @abstract this is read1 */
109 #define BAM_FREAD1        64
110 /*! @abstract this is read2 */
111 #define BAM_FREAD2       128
112 /*! @abstract not primary alignment */
113 #define BAM_FSECONDARY   256
114 /*! @abstract QC failure */
115 #define BAM_FQCFAIL      512
116 /*! @abstract optical or PCR duplicate */
117 #define BAM_FDUP        1024
118
119 #define BAM_OFDEC          0
120 #define BAM_OFHEX          1
121 #define BAM_OFSTR          2
122
123 /*! @abstract defautl mask for pileup */
124 #define BAM_DEF_MASK (BAM_FUNMAP | BAM_FSECONDARY | BAM_FQCFAIL | BAM_FDUP)
125
126 #define BAM_CORE_SIZE   sizeof(bam1_core_t)
127
128 /**
129  * Describing how CIGAR operation/length is packed in a 32-bit integer.
130  */
131 #define BAM_CIGAR_SHIFT 4
132 #define BAM_CIGAR_MASK  ((1 << BAM_CIGAR_SHIFT) - 1)
133
134 /*
135   CIGAR operations.
136  */
137 /*! @abstract CIGAR: M = match or mismatch*/
138 #define BAM_CMATCH      0
139 /*! @abstract CIGAR: I = insertion to the reference */
140 #define BAM_CINS        1
141 /*! @abstract CIGAR: D = deletion from the reference */
142 #define BAM_CDEL        2
143 /*! @abstract CIGAR: N = skip on the reference (e.g. spliced alignment) */
144 #define BAM_CREF_SKIP   3
145 /*! @abstract CIGAR: S = clip on the read with clipped sequence
146   present in qseq */
147 #define BAM_CSOFT_CLIP  4
148 /*! @abstract CIGAR: H = clip on the read with clipped sequence trimmed off */
149 #define BAM_CHARD_CLIP  5
150 /*! @abstract CIGAR: P = padding */
151 #define BAM_CPAD        6
152 /*! @abstract CIGAR: equals = match */
153 #define BAM_CEQUAL        7
154 /*! @abstract CIGAR: X = mismatch */
155 #define BAM_CDIFF        8
156
157 #define BAM_CIGAR_STR "MIDNSHP=X"
158
159 #define bam_cigar_op(c) ((c)&BAM_CIGAR_MASK)
160 #define bam_cigar_oplen(c) ((c)>>BAM_CIGAR_SHIFT)
161 #define bam_cigar_opchr(c) (BAM_CIGAR_STR[bam_cigar_op(c)])
162 #define bam_cigar_gen(o, l) ((o)<<BAM_CIGAR_SHIFT|(l))
163
164 /*! @typedef
165   @abstract Structure for core alignment information.
166   @field  tid     chromosome ID, defined by bam_header_t
167   @field  pos     0-based leftmost coordinate
168   @field  strand  strand; 0 for forward and 1 otherwise
169   @field  bin     bin calculated by bam_reg2bin()
170   @field  qual    mapping quality
171   @field  l_qname length of the query name
172   @field  flag    bitwise flag
173   @field  n_cigar number of CIGAR operations
174   @field  l_qseq  length of the query sequence (read)
175  */
176 typedef struct {
177         int32_t tid;
178         int32_t pos;
179         uint32_t bin:16, qual:8, l_qname:8;
180         uint32_t flag:16, n_cigar:16;
181         int32_t l_qseq;
182         int32_t mtid;
183         int32_t mpos;
184         int32_t isize;
185 } bam1_core_t;
186
187 /*! @typedef
188   @abstract Structure for one alignment.
189   @field  core       core information about the alignment
190   @field  l_aux      length of auxiliary data
191   @field  data_len   current length of bam1_t::data
192   @field  m_data     maximum length of bam1_t::data
193   @field  data       all variable-length data, concatenated; structure: cigar-qname-seq-qual-aux
194
195   @discussion Notes:
196  
197    1. qname is zero tailing and core.l_qname includes the tailing '\0'.
198    2. l_qseq is calculated from the total length of an alignment block
199       on reading or from CIGAR.
200  */
201 typedef struct {
202         bam1_core_t core;
203         int l_aux, data_len, m_data;
204         uint8_t *data;
205 } bam1_t;
206
207 typedef struct __bam_iter_t *bam_iter_t;
208
209 #define bam1_strand(b) (((b)->core.flag&BAM_FREVERSE) != 0)
210 #define bam1_mstrand(b) (((b)->core.flag&BAM_FMREVERSE) != 0)
211
212 /*! @function
213   @abstract  Get the CIGAR array
214   @param  b  pointer to an alignment
215   @return    pointer to the CIGAR array
216
217   @discussion In the CIGAR array, each element is a 32-bit integer. The
218   lower 4 bits gives a CIGAR operation and the higher 28 bits keep the
219   length of a CIGAR.
220  */
221 #define bam1_cigar(b) ((uint32_t*)((b)->data + (b)->core.l_qname))
222
223 /*! @function
224   @abstract  Get the name of the query
225   @param  b  pointer to an alignment
226   @return    pointer to the name string, null terminated
227  */
228 #define bam1_qname(b) ((char*)((b)->data))
229
230 /*! @function
231   @abstract  Get query sequence
232   @param  b  pointer to an alignment
233   @return    pointer to sequence
234
235   @discussion Each base is encoded in 4 bits: 1 for A, 2 for C, 4 for G,
236   8 for T and 15 for N. Two bases are packed in one byte with the base
237   at the higher 4 bits having smaller coordinate on the read. It is
238   recommended to use bam1_seqi() macro to get the base.
239  */
240 #define bam1_seq(b) ((b)->data + (b)->core.n_cigar*4 + (b)->core.l_qname)
241
242 /*! @function
243   @abstract  Get query quality
244   @param  b  pointer to an alignment
245   @return    pointer to quality string
246  */
247 #define bam1_qual(b) ((b)->data + (b)->core.n_cigar*4 + (b)->core.l_qname + (((b)->core.l_qseq + 1)>>1))
248
249 /*! @function
250   @abstract  Get a base on read
251   @param  s  Query sequence returned by bam1_seq()
252   @param  i  The i-th position, 0-based
253   @return    4-bit integer representing the base.
254  */
255 #define bam1_seqi(s, i) ((s)[(i)/2] >> 4*(1-(i)%2) & 0xf)
256
257 /*! @function
258   @abstract  Get query sequence and quality
259   @param  b  pointer to an alignment
260   @return    pointer to the concatenated auxiliary data
261  */
262 #define bam1_aux(b) ((b)->data + (b)->core.n_cigar*4 + (b)->core.l_qname + (b)->core.l_qseq + ((b)->core.l_qseq + 1)/2)
263
264 #ifndef kroundup32
265 /*! @function
266   @abstract  Round an integer to the next closest power-2 integer.
267   @param  x  integer to be rounded (in place)
268   @discussion x will be modified.
269  */
270 #define kroundup32(x) (--(x), (x)|=(x)>>1, (x)|=(x)>>2, (x)|=(x)>>4, (x)|=(x)>>8, (x)|=(x)>>16, ++(x))
271 #endif
272
273 /*!
274   @abstract Whether the machine is big-endian; modified only in
275   bam_header_init().
276  */
277 extern int bam_is_be;
278
279 /*!
280   @abstract Verbose level between 0 and 3; 0 is supposed to disable all
281   debugging information, though this may not have been implemented.
282  */
283 extern int bam_verbose;
284
285 /*! @abstract Table for converting a nucleotide character to the 4-bit encoding. */
286 extern unsigned char bam_nt16_table[256];
287
288 /*! @abstract Table for converting a 4-bit encoded nucleotide to a letter. */
289 extern char *bam_nt16_rev_table;
290
291 extern char bam_nt16_nt4_table[];
292
293 #ifdef __cplusplus
294 extern "C" {
295 #endif
296
297         /*********************
298          * Low-level SAM I/O *
299          *********************/
300
301         /*! @abstract TAM file handler */
302         typedef struct __tamFile_t *tamFile;
303
304         /*!
305           @abstract   Open a SAM file for reading, either uncompressed or compressed by gzip/zlib.
306           @param  fn  SAM file name
307           @return     SAM file handler
308          */
309         tamFile sam_open(const char *fn);
310
311         /*!
312           @abstract   Close a SAM file handler
313           @param  fp  SAM file handler
314          */
315         void sam_close(tamFile fp);
316
317         /*!
318           @abstract      Read one alignment from a SAM file handler
319           @param  fp     SAM file handler
320           @param  header header information (ordered names of chromosomes)
321           @param  b      read alignment; all members in b will be updated
322           @return        0 if successful; otherwise negative
323          */
324         int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b);
325
326         /*!
327           @abstract       Read header information from a TAB-delimited list file.
328           @param  fn_list file name for the list
329           @return         a pointer to the header structure
330
331           @discussion Each line in this file consists of chromosome name and
332           the length of chromosome.
333          */
334         bam_header_t *sam_header_read2(const char *fn_list);
335
336         /*!
337           @abstract       Read header from a SAM file (if present)
338           @param  fp      SAM file handler
339           @return         pointer to header struct; 0 if no @SQ lines available
340          */
341         bam_header_t *sam_header_read(tamFile fp);
342
343         /*!
344           @abstract       Parse @SQ lines a update a header struct
345           @param  h       pointer to the header struct to be updated
346           @return         number of target sequences
347
348           @discussion bam_header_t::{n_targets,target_len,target_name} will
349           be destroyed in the first place.
350          */
351         int sam_header_parse(bam_header_t *h);
352         int32_t bam_get_tid(const bam_header_t *header, const char *seq_name);
353
354         /*!
355           @abstract       Parse @RG lines a update a header struct
356           @param  h       pointer to the header struct to be updated
357           @return         number of @RG lines
358
359           @discussion bam_header_t::rg2lib will be destroyed in the first
360           place.
361          */
362         int sam_header_parse_rg(bam_header_t *h);
363
364 #define sam_write1(header, b) bam_view1(header, b)
365
366
367         /********************************
368          * APIs for string dictionaries *
369          ********************************/
370
371         int bam_strmap_put(void *strmap, const char *rg, const char *lib);
372         const char *bam_strmap_get(const void *strmap, const char *rg);
373         void *bam_strmap_dup(const void*);
374         void *bam_strmap_init();
375         void bam_strmap_destroy(void *strmap);
376
377
378         /*********************
379          * Low-level BAM I/O *
380          *********************/
381
382         /*!
383           @abstract Initialize a header structure.
384           @return   the pointer to the header structure
385
386           @discussion This function also modifies the global variable
387           bam_is_be.
388          */
389         bam_header_t *bam_header_init();
390
391         /*!
392           @abstract        Destroy a header structure.
393           @param  header  pointer to the header
394          */
395         void bam_header_destroy(bam_header_t *header);
396
397         /*!
398           @abstract   Read a header structure from BAM.
399           @param  fp  BAM file handler, opened by bam_open()
400           @return     pointer to the header structure
401
402           @discussion The file position indicator must be placed at the
403           beginning of the file. Upon success, the position indicator will
404           be set at the start of the first alignment.
405          */
406         bam_header_t *bam_header_read(bamFile fp);
407
408         /*!
409           @abstract      Write a header structure to BAM.
410           @param  fp     BAM file handler
411           @param  header pointer to the header structure
412           @return        always 0 currently
413          */
414         int bam_header_write(bamFile fp, const bam_header_t *header);
415
416         /*!
417           @abstract   Read an alignment from BAM.
418           @param  fp  BAM file handler
419           @param  b   read alignment; all members are updated.
420           @return     number of bytes read from the file
421
422           @discussion The file position indicator must be
423           placed right before an alignment. Upon success, this function
424           will set the position indicator to the start of the next
425           alignment. This function is not affected by the machine
426           endianness.
427          */
428         int bam_read1(bamFile fp, bam1_t *b);
429
430         /*!
431           @abstract Write an alignment to BAM.
432           @param  fp       BAM file handler
433           @param  c        pointer to the bam1_core_t structure
434           @param  data_len total length of variable size data related to
435                            the alignment
436           @param  data     pointer to the concatenated data
437           @return          number of bytes written to the file
438
439           @discussion This function is not affected by the machine
440           endianness.
441          */
442         int bam_write1_core(bamFile fp, const bam1_core_t *c, int data_len, uint8_t *data);
443
444         /*!
445           @abstract   Write an alignment to BAM.
446           @param  fp  BAM file handler
447           @param  b   alignment to write
448           @return     number of bytes written to the file
449
450           @abstract It is equivalent to:
451             bam_write1_core(fp, &b->core, b->data_len, b->data)
452          */
453         int bam_write1(bamFile fp, const bam1_t *b);
454
455         /*! @function
456           @abstract  Initiate a pointer to bam1_t struct
457          */
458 #define bam_init1() ((bam1_t*)calloc(1, sizeof(bam1_t)))
459
460         /*! @function
461           @abstract  Free the memory allocated for an alignment.
462           @param  b  pointer to an alignment
463          */
464 #define bam_destroy1(b) do {                                    \
465                 if (b) { free((b)->data); free(b); }    \
466         } while (0)
467
468         /*!
469           @abstract       Format a BAM record in the SAM format
470           @param  header  pointer to the header structure
471           @param  b       alignment to print
472           @return         a pointer to the SAM string
473          */
474         char *bam_format1(const bam_header_t *header, const bam1_t *b);
475
476         char *bam_format1_core(const bam_header_t *header, const bam1_t *b, int of);
477
478         /*!
479           @abstract       Check whether a BAM record is plausibly valid
480           @param  header  associated header structure, or NULL if unavailable
481           @param  b       alignment to validate
482           @return         0 if the alignment is invalid; non-zero otherwise
483
484           @discussion  Simple consistency check of some of the fields of the
485           alignment record.  If the header is provided, several additional checks
486           are made.  Not all fields are checked, so a non-zero result is not a
487           guarantee that the record is valid.  However it is usually good enough
488           to detect when bam_seek() has been called with a virtual file offset
489           that is not the offset of an alignment record.
490          */
491         int bam_validate1(const bam_header_t *header, const bam1_t *b);
492
493         const char *bam_get_library(bam_header_t *header, const bam1_t *b);
494
495
496         /***************
497          * pileup APIs *
498          ***************/
499
500         /*! @typedef
501           @abstract Structure for one alignment covering the pileup position.
502           @field  b      pointer to the alignment
503           @field  qpos   position of the read base at the pileup site, 0-based
504           @field  indel  indel length; 0 for no indel, positive for ins and negative for del
505           @field  is_del 1 iff the base on the padded read is a deletion
506           @field  level  the level of the read in the "viewer" mode
507
508           @discussion See also bam_plbuf_push() and bam_lplbuf_push(). The
509           difference between the two functions is that the former does not
510           set bam_pileup1_t::level, while the later does. Level helps the
511           implementation of alignment viewers, but calculating this has some
512           overhead.
513          */
514         typedef struct {
515                 bam1_t *b;
516                 int32_t qpos;
517                 int indel, level;
518                 uint32_t is_del:1, is_head:1, is_tail:1, is_refskip:1, aux:28;
519         } bam_pileup1_t;
520
521         typedef int (*bam_plp_auto_f)(void *data, bam1_t *b);
522
523         struct __bam_plp_t;
524         typedef struct __bam_plp_t *bam_plp_t;
525
526         bam_plp_t bam_plp_init(bam_plp_auto_f func, void *data);
527         int bam_plp_push(bam_plp_t iter, const bam1_t *b);
528         const bam_pileup1_t *bam_plp_next(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp);
529         const bam_pileup1_t *bam_plp_auto(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp);
530         void bam_plp_set_mask(bam_plp_t iter, int mask);
531         void bam_plp_set_maxcnt(bam_plp_t iter, int maxcnt);
532         void bam_plp_reset(bam_plp_t iter);
533         void bam_plp_destroy(bam_plp_t iter);
534
535         struct __bam_mplp_t;
536         typedef struct __bam_mplp_t *bam_mplp_t;
537
538         bam_mplp_t bam_mplp_init(int n, bam_plp_auto_f func, void **data);
539         void bam_mplp_destroy(bam_mplp_t iter);
540         void bam_mplp_set_maxcnt(bam_mplp_t iter, int maxcnt);
541         int bam_mplp_auto(bam_mplp_t iter, int *_tid, int *_pos, int *n_plp, const bam_pileup1_t **plp);
542
543         /*! @typedef
544           @abstract    Type of function to be called by bam_plbuf_push().
545           @param  tid  chromosome ID as is defined in the header
546           @param  pos  start coordinate of the alignment, 0-based
547           @param  n    number of elements in pl array
548           @param  pl   array of alignments
549           @param  data user provided data
550           @discussion  See also bam_plbuf_push(), bam_plbuf_init() and bam_pileup1_t.
551          */
552         typedef int (*bam_pileup_f)(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pl, void *data);
553
554         typedef struct {
555                 bam_plp_t iter;
556                 bam_pileup_f func;
557                 void *data;
558         } bam_plbuf_t;
559
560         void bam_plbuf_set_mask(bam_plbuf_t *buf, int mask);
561         void bam_plbuf_reset(bam_plbuf_t *buf);
562         bam_plbuf_t *bam_plbuf_init(bam_pileup_f func, void *data);
563         void bam_plbuf_destroy(bam_plbuf_t *buf);
564         int bam_plbuf_push(const bam1_t *b, bam_plbuf_t *buf);
565
566         int bam_pileup_file(bamFile fp, int mask, bam_pileup_f func, void *func_data);
567
568         struct __bam_lplbuf_t;
569         typedef struct __bam_lplbuf_t bam_lplbuf_t;
570
571         void bam_lplbuf_reset(bam_lplbuf_t *buf);
572
573         /*! @abstract  bam_plbuf_init() equivalent with level calculated. */
574         bam_lplbuf_t *bam_lplbuf_init(bam_pileup_f func, void *data);
575
576         /*! @abstract  bam_plbuf_destroy() equivalent with level calculated. */
577         void bam_lplbuf_destroy(bam_lplbuf_t *tv);
578
579         /*! @abstract  bam_plbuf_push() equivalent with level calculated. */
580         int bam_lplbuf_push(const bam1_t *b, bam_lplbuf_t *buf);
581
582
583         /*********************
584          * BAM indexing APIs *
585          *********************/
586
587         struct __bam_index_t;
588         typedef struct __bam_index_t bam_index_t;
589
590         /*!
591           @abstract   Build index for a BAM file.
592           @discussion Index file "fn.bai" will be created.
593           @param  fn  name of the BAM file
594           @return     always 0 currently
595          */
596         int bam_index_build(const char *fn);
597
598         /*!
599           @abstract   Load index from file "fn.bai".
600           @param  fn  name of the BAM file (NOT the index file)
601           @return     pointer to the index structure
602          */
603         bam_index_t *bam_index_load(const char *fn);
604
605         /*!
606           @abstract    Destroy an index structure.
607           @param  idx  pointer to the index structure
608          */
609         void bam_index_destroy(bam_index_t *idx);
610
611         /*! @typedef
612           @abstract      Type of function to be called by bam_fetch().
613           @param  b     the alignment
614           @param  data  user provided data
615          */
616         typedef int (*bam_fetch_f)(const bam1_t *b, void *data);
617
618         /*!
619           @abstract Retrieve the alignments that are overlapped with the
620           specified region.
621
622           @discussion A user defined function will be called for each
623           retrieved alignment ordered by its start position.
624
625           @param  fp    BAM file handler
626           @param  idx   pointer to the alignment index
627           @param  tid   chromosome ID as is defined in the header
628           @param  beg   start coordinate, 0-based
629           @param  end   end coordinate, 0-based
630           @param  data  user provided data (will be transferred to func)
631           @param  func  user defined function
632          */
633         int bam_fetch(bamFile fp, const bam_index_t *idx, int tid, int beg, int end, void *data, bam_fetch_f func);
634
635         bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end);
636         int bam_iter_read(bamFile fp, bam_iter_t iter, bam1_t *b);
637         void bam_iter_destroy(bam_iter_t iter);
638
639         /*!
640           @abstract       Parse a region in the format: "chr2:100,000-200,000".
641           @discussion     bam_header_t::hash will be initialized if empty.
642           @param  header  pointer to the header structure
643           @param  str     string to be parsed
644           @param  ref_id  the returned chromosome ID
645           @param  begin   the returned start coordinate
646           @param  end     the returned end coordinate
647           @return         0 on success; -1 on failure
648          */
649         int bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *begin, int *end);
650
651
652         /**************************
653          * APIs for optional tags *
654          **************************/
655
656         /*!
657           @abstract       Retrieve data of a tag
658           @param  b       pointer to an alignment struct
659           @param  tag     two-character tag to be retrieved
660
661           @return  pointer to the type and data. The first character is the
662           type that can be 'iIsScCdfAZH'.
663
664           @discussion  Use bam_aux2?() series to convert the returned data to
665           the corresponding type.
666         */
667         uint8_t *bam_aux_get(const bam1_t *b, const char tag[2]);
668
669         int32_t bam_aux2i(const uint8_t *s);
670         float bam_aux2f(const uint8_t *s);
671         double bam_aux2d(const uint8_t *s);
672         char bam_aux2A(const uint8_t *s);
673         char *bam_aux2Z(const uint8_t *s);
674
675         int bam_aux_del(bam1_t *b, uint8_t *s);
676         void bam_aux_append(bam1_t *b, const char tag[2], char type, int len, uint8_t *data);
677         uint8_t *bam_aux_get_core(bam1_t *b, const char tag[2]); // an alias of bam_aux_get()
678
679
680         /*****************
681          * Miscellaneous *
682          *****************/
683
684         /*!  
685           @abstract Calculate the rightmost coordinate of an alignment on the
686           reference genome.
687
688           @param  c      pointer to the bam1_core_t structure
689           @param  cigar  the corresponding CIGAR array (from bam1_t::cigar)
690           @return        the rightmost coordinate, 0-based
691         */
692         uint32_t bam_calend(const bam1_core_t *c, const uint32_t *cigar);
693
694         /*!
695           @abstract      Calculate the length of the query sequence from CIGAR.
696           @param  c      pointer to the bam1_core_t structure
697           @param  cigar  the corresponding CIGAR array (from bam1_t::cigar)
698           @return        length of the query sequence
699         */
700         int32_t bam_cigar2qlen(const bam1_core_t *c, const uint32_t *cigar);
701
702 #ifdef __cplusplus
703 }
704 #endif
705
706 /*!
707   @abstract    Calculate the minimum bin that contains a region [beg,end).
708   @param  beg  start of the region, 0-based
709   @param  end  end of the region, 0-based
710   @return      bin
711  */
712 static inline int bam_reg2bin(uint32_t beg, uint32_t end)
713 {
714         --end;
715         if (beg>>14 == end>>14) return 4681 + (beg>>14);
716         if (beg>>17 == end>>17) return  585 + (beg>>17);
717         if (beg>>20 == end>>20) return   73 + (beg>>20);
718         if (beg>>23 == end>>23) return    9 + (beg>>23);
719         if (beg>>26 == end>>26) return    1 + (beg>>26);
720         return 0;
721 }
722
723 /*!
724   @abstract     Copy an alignment
725   @param  bdst  destination alignment struct
726   @param  bsrc  source alignment struct
727   @return       pointer to the destination alignment struct
728  */
729 static inline bam1_t *bam_copy1(bam1_t *bdst, const bam1_t *bsrc)
730 {
731         uint8_t *data = bdst->data;
732         int m_data = bdst->m_data;   // backup data and m_data
733         if (m_data < bsrc->data_len) { // double the capacity
734                 m_data = bsrc->data_len; kroundup32(m_data);
735                 data = (uint8_t*)realloc(data, m_data);
736         }
737         memcpy(data, bsrc->data, bsrc->data_len); // copy var-len data
738         *bdst = *bsrc; // copy the rest
739         // restore the backup
740         bdst->m_data = m_data;
741         bdst->data = data;
742         return bdst;
743 }
744
745 /*!
746   @abstract     Duplicate an alignment
747   @param  src   source alignment struct
748   @return       pointer to the destination alignment struct
749  */
750 static inline bam1_t *bam_dup1(const bam1_t *src)
751 {
752         bam1_t *b;
753         b = bam_init1();
754         *b = *src;
755         b->m_data = b->data_len;
756         b->data = (uint8_t*)calloc(b->data_len, 1);
757         memcpy(b->data, src->data, b->data_len);
758         return b;
759 }
760
761 static inline int bam_aux_type2size(int x)
762 {
763         if (x == 'C' || x == 'c' || x == 'A') return 1;
764         else if (x == 'S' || x == 's') return 2;
765         else if (x == 'I' || x == 'i' || x == 'f') return 4;
766         else return 0;
767 }
768
769
770 #endif