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