3 Copyright (c) 2008-2010 Genome Research Ltd (GRL).
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:
13 The above copyright notice and this permission notice shall be
14 included in all copies or substantial portions of the Software.
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
26 /* Contact: Heng Li <lh3@sanger.ac.uk> */
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.
40 @copyright Genome Research Ltd.
43 #define BAM_VERSION "0.1.18-dev (r982:313)"
51 #define BAM_VIRTUAL_OFFSET16
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)
63 #define BAM_TRUE_OFFSET
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 */
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
84 @discussion Field hash points to null by default. It is a private
91 void *dict, *hash, *rg2lib;
92 uint32_t l_text, n_text;
96 /*! @abstract the read is paired in sequencing, no matter whether it is mapped in a pair */
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 */
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
123 /*! @abstract defautl mask for pileup */
124 #define BAM_DEF_MASK (BAM_FUNMAP | BAM_FSECONDARY | BAM_FQCFAIL | BAM_FDUP)
126 #define BAM_CORE_SIZE sizeof(bam1_core_t)
129 * Describing how CIGAR operation/length is packed in a 32-bit integer.
131 #define BAM_CIGAR_SHIFT 4
132 #define BAM_CIGAR_MASK ((1 << BAM_CIGAR_SHIFT) - 1)
137 /*! @abstract CIGAR: M = match or mismatch*/
139 /*! @abstract CIGAR: I = insertion to the reference */
141 /*! @abstract CIGAR: D = deletion from the reference */
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
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 */
152 /*! @abstract CIGAR: equals = match */
154 /*! @abstract CIGAR: X = mismatch */
158 #define BAM_CIGAR_STR "MIDNSHP=XB"
159 #define BAM_CIGAR_TYPE 0x3C1A7
161 #define bam_cigar_op(c) ((c)&BAM_CIGAR_MASK)
162 #define bam_cigar_oplen(c) ((c)>>BAM_CIGAR_SHIFT)
163 #define bam_cigar_opchr(c) (BAM_CIGAR_STR[bam_cigar_op(c)])
164 #define bam_cigar_gen(l, o) ((l)<<BAM_CIGAR_SHIFT|(o))
165 #define bam_cigar_type(o) (BAM_CIGAR_TYPE>>((o)<<1)&3) // bit 1: consume query; bit 2: consume reference
168 @abstract Structure for core alignment information.
169 @field tid chromosome ID, defined by bam_header_t
170 @field pos 0-based leftmost coordinate
171 @field strand strand; 0 for forward and 1 otherwise
172 @field bin bin calculated by bam_reg2bin()
173 @field qual mapping quality
174 @field l_qname length of the query name
175 @field flag bitwise flag
176 @field n_cigar number of CIGAR operations
177 @field l_qseq length of the query sequence (read)
182 uint32_t bin:16, qual:8, l_qname:8;
183 uint32_t flag:16, n_cigar:16;
191 @abstract Structure for one alignment.
192 @field core core information about the alignment
193 @field l_aux length of auxiliary data
194 @field data_len current length of bam1_t::data
195 @field m_data maximum length of bam1_t::data
196 @field data all variable-length data, concatenated; structure: cigar-qname-seq-qual-aux
200 1. qname is zero tailing and core.l_qname includes the tailing '\0'.
201 2. l_qseq is calculated from the total length of an alignment block
202 on reading or from CIGAR.
206 int l_aux, data_len, m_data;
210 typedef struct __bam_iter_t *bam_iter_t;
212 #define bam1_strand(b) (((b)->core.flag&BAM_FREVERSE) != 0)
213 #define bam1_mstrand(b) (((b)->core.flag&BAM_FMREVERSE) != 0)
216 @abstract Get the CIGAR array
217 @param b pointer to an alignment
218 @return pointer to the CIGAR array
220 @discussion In the CIGAR array, each element is a 32-bit integer. The
221 lower 4 bits gives a CIGAR operation and the higher 28 bits keep the
224 #define bam1_cigar(b) ((uint32_t*)((b)->data + (b)->core.l_qname))
227 @abstract Get the name of the query
228 @param b pointer to an alignment
229 @return pointer to the name string, null terminated
231 #define bam1_qname(b) ((char*)((b)->data))
234 @abstract Get query sequence
235 @param b pointer to an alignment
236 @return pointer to sequence
238 @discussion Each base is encoded in 4 bits: 1 for A, 2 for C, 4 for G,
239 8 for T and 15 for N. Two bases are packed in one byte with the base
240 at the higher 4 bits having smaller coordinate on the read. It is
241 recommended to use bam1_seqi() macro to get the base.
243 #define bam1_seq(b) ((b)->data + (b)->core.n_cigar*4 + (b)->core.l_qname)
246 @abstract Get query quality
247 @param b pointer to an alignment
248 @return pointer to quality string
250 #define bam1_qual(b) ((b)->data + (b)->core.n_cigar*4 + (b)->core.l_qname + (((b)->core.l_qseq + 1)>>1))
253 @abstract Get a base on read
254 @param s Query sequence returned by bam1_seq()
255 @param i The i-th position, 0-based
256 @return 4-bit integer representing the base.
258 //#define bam1_seqi(s, i) ((s)[(i)/2] >> 4*(1-(i)%2) & 0xf)
259 #define bam1_seqi(s, i) ((s)[(i)>>1] >> ((~(i)&1)<<2) & 0xf)
261 #define bam1_seq_seti(s, i, c) ( (s)[(i)>>1] = ((s)[(i)>>1] & 0xf<<(((i)&1)<<2)) | (c)<<((~(i)&1)<<2) )
264 @abstract Get query sequence and quality
265 @param b pointer to an alignment
266 @return pointer to the concatenated auxiliary data
268 #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)
272 @abstract Round an integer to the next closest power-2 integer.
273 @param x integer to be rounded (in place)
274 @discussion x will be modified.
276 #define kroundup32(x) (--(x), (x)|=(x)>>1, (x)|=(x)>>2, (x)|=(x)>>4, (x)|=(x)>>8, (x)|=(x)>>16, ++(x))
280 @abstract Whether the machine is big-endian; modified only in
283 extern int bam_is_be;
286 @abstract Verbose level between 0 and 3; 0 is supposed to disable all
287 debugging information, though this may not have been implemented.
289 extern int bam_verbose;
293 /*! @abstract Table for converting a nucleotide character to the 4-bit encoding. */
294 extern unsigned char bam_nt16_table[256];
296 /*! @abstract Table for converting a 4-bit encoded nucleotide to a letter. */
297 extern char *bam_nt16_rev_table;
299 extern char bam_nt16_nt4_table[];
305 /*********************
306 * Low-level SAM I/O *
307 *********************/
309 /*! @abstract TAM file handler */
310 typedef struct __tamFile_t *tamFile;
313 @abstract Open a SAM file for reading, either uncompressed or compressed by gzip/zlib.
314 @param fn SAM file name
315 @return SAM file handler
317 tamFile sam_open(const char *fn);
320 @abstract Close a SAM file handler
321 @param fp SAM file handler
323 void sam_close(tamFile fp);
326 @abstract Read one alignment from a SAM file handler
327 @param fp SAM file handler
328 @param header header information (ordered names of chromosomes)
329 @param b read alignment; all members in b will be updated
330 @return 0 if successful; otherwise negative
332 int sam_read1(tamFile fp, bam_header_t *header, bam1_t *b);
335 @abstract Read header information from a TAB-delimited list file.
336 @param fn_list file name for the list
337 @return a pointer to the header structure
339 @discussion Each line in this file consists of chromosome name and
340 the length of chromosome.
342 bam_header_t *sam_header_read2(const char *fn_list);
345 @abstract Read header from a SAM file (if present)
346 @param fp SAM file handler
347 @return pointer to header struct; 0 if no @SQ lines available
349 bam_header_t *sam_header_read(tamFile fp);
352 @abstract Parse @SQ lines a update a header struct
353 @param h pointer to the header struct to be updated
354 @return number of target sequences
356 @discussion bam_header_t::{n_targets,target_len,target_name} will
357 be destroyed in the first place.
359 int sam_header_parse(bam_header_t *h);
360 int32_t bam_get_tid(const bam_header_t *header, const char *seq_name);
363 @abstract Parse @RG lines a update a header struct
364 @param h pointer to the header struct to be updated
365 @return number of @RG lines
367 @discussion bam_header_t::rg2lib will be destroyed in the first
370 int sam_header_parse_rg(bam_header_t *h);
372 #define sam_write1(header, b) bam_view1(header, b)
375 /********************************
376 * APIs for string dictionaries *
377 ********************************/
379 int bam_strmap_put(void *strmap, const char *rg, const char *lib);
380 const char *bam_strmap_get(const void *strmap, const char *rg);
381 void *bam_strmap_dup(const void*);
382 void *bam_strmap_init();
383 void bam_strmap_destroy(void *strmap);
386 /*********************
387 * Low-level BAM I/O *
388 *********************/
391 @abstract Initialize a header structure.
392 @return the pointer to the header structure
394 @discussion This function also modifies the global variable
397 bam_header_t *bam_header_init();
400 @abstract Destroy a header structure.
401 @param header pointer to the header
403 void bam_header_destroy(bam_header_t *header);
406 @abstract Read a header structure from BAM.
407 @param fp BAM file handler, opened by bam_open()
408 @return pointer to the header structure
410 @discussion The file position indicator must be placed at the
411 beginning of the file. Upon success, the position indicator will
412 be set at the start of the first alignment.
414 bam_header_t *bam_header_read(bamFile fp);
417 @abstract Write a header structure to BAM.
418 @param fp BAM file handler
419 @param header pointer to the header structure
420 @return always 0 currently
422 int bam_header_write(bamFile fp, const bam_header_t *header);
425 @abstract Read an alignment from BAM.
426 @param fp BAM file handler
427 @param b read alignment; all members are updated.
428 @return number of bytes read from the file
430 @discussion The file position indicator must be
431 placed right before an alignment. Upon success, this function
432 will set the position indicator to the start of the next
433 alignment. This function is not affected by the machine
436 int bam_read1(bamFile fp, bam1_t *b);
438 int bam_remove_B(bam1_t *b);
441 @abstract Write an alignment to BAM.
442 @param fp BAM file handler
443 @param c pointer to the bam1_core_t structure
444 @param data_len total length of variable size data related to
446 @param data pointer to the concatenated data
447 @return number of bytes written to the file
449 @discussion This function is not affected by the machine
452 int bam_write1_core(bamFile fp, const bam1_core_t *c, int data_len, uint8_t *data);
455 @abstract Write an alignment to BAM.
456 @param fp BAM file handler
457 @param b alignment to write
458 @return number of bytes written to the file
460 @abstract It is equivalent to:
461 bam_write1_core(fp, &b->core, b->data_len, b->data)
463 int bam_write1(bamFile fp, const bam1_t *b);
466 @abstract Initiate a pointer to bam1_t struct
468 #define bam_init1() ((bam1_t*)calloc(1, sizeof(bam1_t)))
471 @abstract Free the memory allocated for an alignment.
472 @param b pointer to an alignment
474 #define bam_destroy1(b) do { \
475 if (b) { free((b)->data); free(b); } \
479 @abstract Format a BAM record in the SAM format
480 @param header pointer to the header structure
481 @param b alignment to print
482 @return a pointer to the SAM string
484 char *bam_format1(const bam_header_t *header, const bam1_t *b);
486 char *bam_format1_core(const bam_header_t *header, const bam1_t *b, int of);
489 @abstract Check whether a BAM record is plausibly valid
490 @param header associated header structure, or NULL if unavailable
491 @param b alignment to validate
492 @return 0 if the alignment is invalid; non-zero otherwise
494 @discussion Simple consistency check of some of the fields of the
495 alignment record. If the header is provided, several additional checks
496 are made. Not all fields are checked, so a non-zero result is not a
497 guarantee that the record is valid. However it is usually good enough
498 to detect when bam_seek() has been called with a virtual file offset
499 that is not the offset of an alignment record.
501 int bam_validate1(const bam_header_t *header, const bam1_t *b);
503 const char *bam_get_library(bam_header_t *header, const bam1_t *b);
511 @abstract Structure for one alignment covering the pileup position.
512 @field b pointer to the alignment
513 @field qpos position of the read base at the pileup site, 0-based
514 @field indel indel length; 0 for no indel, positive for ins and negative for del
515 @field is_del 1 iff the base on the padded read is a deletion
516 @field level the level of the read in the "viewer" mode
518 @discussion See also bam_plbuf_push() and bam_lplbuf_push(). The
519 difference between the two functions is that the former does not
520 set bam_pileup1_t::level, while the later does. Level helps the
521 implementation of alignment viewers, but calculating this has some
528 uint32_t is_del:1, is_head:1, is_tail:1, is_refskip:1, aux:28;
531 typedef int (*bam_plp_auto_f)(void *data, bam1_t *b);
534 typedef struct __bam_plp_t *bam_plp_t;
536 bam_plp_t bam_plp_init(bam_plp_auto_f func, void *data);
537 int bam_plp_push(bam_plp_t iter, const bam1_t *b);
538 const bam_pileup1_t *bam_plp_next(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp);
539 const bam_pileup1_t *bam_plp_auto(bam_plp_t iter, int *_tid, int *_pos, int *_n_plp);
540 void bam_plp_set_mask(bam_plp_t iter, int mask);
541 void bam_plp_set_maxcnt(bam_plp_t iter, int maxcnt);
542 void bam_plp_reset(bam_plp_t iter);
543 void bam_plp_destroy(bam_plp_t iter);
546 typedef struct __bam_mplp_t *bam_mplp_t;
548 bam_mplp_t bam_mplp_init(int n, bam_plp_auto_f func, void **data);
549 void bam_mplp_destroy(bam_mplp_t iter);
550 void bam_mplp_set_maxcnt(bam_mplp_t iter, int maxcnt);
551 int bam_mplp_auto(bam_mplp_t iter, int *_tid, int *_pos, int *n_plp, const bam_pileup1_t **plp);
554 @abstract Type of function to be called by bam_plbuf_push().
555 @param tid chromosome ID as is defined in the header
556 @param pos start coordinate of the alignment, 0-based
557 @param n number of elements in pl array
558 @param pl array of alignments
559 @param data user provided data
560 @discussion See also bam_plbuf_push(), bam_plbuf_init() and bam_pileup1_t.
562 typedef int (*bam_pileup_f)(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pl, void *data);
570 void bam_plbuf_set_mask(bam_plbuf_t *buf, int mask);
571 void bam_plbuf_reset(bam_plbuf_t *buf);
572 bam_plbuf_t *bam_plbuf_init(bam_pileup_f func, void *data);
573 void bam_plbuf_destroy(bam_plbuf_t *buf);
574 int bam_plbuf_push(const bam1_t *b, bam_plbuf_t *buf);
576 int bam_pileup_file(bamFile fp, int mask, bam_pileup_f func, void *func_data);
578 struct __bam_lplbuf_t;
579 typedef struct __bam_lplbuf_t bam_lplbuf_t;
581 void bam_lplbuf_reset(bam_lplbuf_t *buf);
583 /*! @abstract bam_plbuf_init() equivalent with level calculated. */
584 bam_lplbuf_t *bam_lplbuf_init(bam_pileup_f func, void *data);
586 /*! @abstract bam_plbuf_destroy() equivalent with level calculated. */
587 void bam_lplbuf_destroy(bam_lplbuf_t *tv);
589 /*! @abstract bam_plbuf_push() equivalent with level calculated. */
590 int bam_lplbuf_push(const bam1_t *b, bam_lplbuf_t *buf);
593 /*********************
594 * BAM indexing APIs *
595 *********************/
597 struct __bam_index_t;
598 typedef struct __bam_index_t bam_index_t;
601 @abstract Build index for a BAM file.
602 @discussion Index file "fn.bai" will be created.
603 @param fn name of the BAM file
604 @return always 0 currently
606 int bam_index_build(const char *fn);
609 @abstract Load index from file "fn.bai".
610 @param fn name of the BAM file (NOT the index file)
611 @return pointer to the index structure
613 bam_index_t *bam_index_load(const char *fn);
616 @abstract Destroy an index structure.
617 @param idx pointer to the index structure
619 void bam_index_destroy(bam_index_t *idx);
622 @abstract Type of function to be called by bam_fetch().
623 @param b the alignment
624 @param data user provided data
626 typedef int (*bam_fetch_f)(const bam1_t *b, void *data);
629 @abstract Retrieve the alignments that are overlapped with the
632 @discussion A user defined function will be called for each
633 retrieved alignment ordered by its start position.
635 @param fp BAM file handler
636 @param idx pointer to the alignment index
637 @param tid chromosome ID as is defined in the header
638 @param beg start coordinate, 0-based
639 @param end end coordinate, 0-based
640 @param data user provided data (will be transferred to func)
641 @param func user defined function
643 int bam_fetch(bamFile fp, const bam_index_t *idx, int tid, int beg, int end, void *data, bam_fetch_f func);
645 bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end);
646 int bam_iter_read(bamFile fp, bam_iter_t iter, bam1_t *b);
647 void bam_iter_destroy(bam_iter_t iter);
650 @abstract Parse a region in the format: "chr2:100,000-200,000".
651 @discussion bam_header_t::hash will be initialized if empty.
652 @param header pointer to the header structure
653 @param str string to be parsed
654 @param ref_id the returned chromosome ID
655 @param begin the returned start coordinate
656 @param end the returned end coordinate
657 @return 0 on success; -1 on failure
659 int bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *begin, int *end);
662 /**************************
663 * APIs for optional tags *
664 **************************/
667 @abstract Retrieve data of a tag
668 @param b pointer to an alignment struct
669 @param tag two-character tag to be retrieved
671 @return pointer to the type and data. The first character is the
672 type that can be 'iIsScCdfAZH'.
674 @discussion Use bam_aux2?() series to convert the returned data to
675 the corresponding type.
677 uint8_t *bam_aux_get(const bam1_t *b, const char tag[2]);
679 int32_t bam_aux2i(const uint8_t *s);
680 float bam_aux2f(const uint8_t *s);
681 double bam_aux2d(const uint8_t *s);
682 char bam_aux2A(const uint8_t *s);
683 char *bam_aux2Z(const uint8_t *s);
685 int bam_aux_del(bam1_t *b, uint8_t *s);
686 void bam_aux_append(bam1_t *b, const char tag[2], char type, int len, uint8_t *data);
687 uint8_t *bam_aux_get_core(bam1_t *b, const char tag[2]); // an alias of bam_aux_get()
695 @abstract Calculate the rightmost coordinate of an alignment on the
698 @param c pointer to the bam1_core_t structure
699 @param cigar the corresponding CIGAR array (from bam1_t::cigar)
700 @return the rightmost coordinate, 0-based
702 uint32_t bam_calend(const bam1_core_t *c, const uint32_t *cigar);
705 @abstract Calculate the length of the query sequence from CIGAR.
706 @param c pointer to the bam1_core_t structure
707 @param cigar the corresponding CIGAR array (from bam1_t::cigar)
708 @return length of the query sequence
710 int32_t bam_cigar2qlen(const bam1_core_t *c, const uint32_t *cigar);
717 @abstract Calculate the minimum bin that contains a region [beg,end).
718 @param beg start of the region, 0-based
719 @param end end of the region, 0-based
722 static inline int bam_reg2bin(uint32_t beg, uint32_t end)
725 if (beg>>14 == end>>14) return 4681 + (beg>>14);
726 if (beg>>17 == end>>17) return 585 + (beg>>17);
727 if (beg>>20 == end>>20) return 73 + (beg>>20);
728 if (beg>>23 == end>>23) return 9 + (beg>>23);
729 if (beg>>26 == end>>26) return 1 + (beg>>26);
734 @abstract Copy an alignment
735 @param bdst destination alignment struct
736 @param bsrc source alignment struct
737 @return pointer to the destination alignment struct
739 static inline bam1_t *bam_copy1(bam1_t *bdst, const bam1_t *bsrc)
741 uint8_t *data = bdst->data;
742 int m_data = bdst->m_data; // backup data and m_data
743 if (m_data < bsrc->data_len) { // double the capacity
744 m_data = bsrc->data_len; kroundup32(m_data);
745 data = (uint8_t*)realloc(data, m_data);
747 memcpy(data, bsrc->data, bsrc->data_len); // copy var-len data
748 *bdst = *bsrc; // copy the rest
749 // restore the backup
750 bdst->m_data = m_data;
756 @abstract Duplicate an alignment
757 @param src source alignment struct
758 @return pointer to the destination alignment struct
760 static inline bam1_t *bam_dup1(const bam1_t *src)
765 b->m_data = b->data_len;
766 b->data = (uint8_t*)calloc(b->data_len, 1);
767 memcpy(b->data, src->data, b->data_len);
771 static inline int bam_aux_type2size(int x)
773 if (x == 'C' || x == 'c' || x == 'A') return 1;
774 else if (x == 'S' || x == 's') return 2;
775 else if (x == 'I' || x == 'i' || x == 'f') return 4;