]> git.donarmstrong.com Git - samtools.git/blobdiff - bam.h
* samtools-0.1.3-20 (r276)
[samtools.git] / bam.h
diff --git a/bam.h b/bam.h
index 4b3a68856e7245d4f1595a194cc724f2f95af2eb..163dc893e192d6eaa9ebfcc3b3b8614d142ca211 100644 (file)
--- a/bam.h
+++ b/bam.h
@@ -46,7 +46,7 @@
 #include <string.h>
 #include <stdio.h>
 
-#if _IOLIB == 1
+#if _IOLIB == 1 && !defined(_NO_RAZF)
 #define BAM_TRUE_OFFSET
 #include "razf.h"
 /*! @abstract BAM file handler */
@@ -118,6 +118,10 @@ typedef struct {
 #define BAM_FREAD1        64
 #define BAM_FREAD2       128
 #define BAM_FSECONDARY   256
+#define BAM_FQCFAIL      512
+#define BAM_FDUP        1024
+
+#define BAM_DEF_MASK (BAM_FUNMAP | BAM_FSECONDARY | BAM_FQCFAIL | BAM_FDUP)
 
 #define BAM_CORE_SIZE   sizeof(bam1_core_t)
 
@@ -175,7 +179,6 @@ typedef struct {
   @field  data_len   current length of bam1_t::data
   @field  m_data     maximum length of bam1_t::data
   @field  data       all variable-length data, concatenated; structure: cigar-qname-seq-qual-aux
-  @field  hash       hash table for fast retrieval of tag-value pairs; private
 
   @discussion Notes:
  
@@ -187,7 +190,6 @@ typedef struct {
        bam1_core_t core;
        int l_aux, data_len, m_data;
        uint8_t *data;
-       void *hash;
 } bam1_t;
 
 #define bam1_strand(b) (((b)->core.flag&BAM_FREVERSE) != 0)
@@ -313,6 +315,9 @@ extern "C" {
         */
        bam_header_t *sam_header_read2(const char *fn_list);
 
+       bam_header_t *sam_header_read(tamFile fp);
+       int sam_header_parse(bam_header_t *h);
+
 #define sam_write1(header, b) bam_view1(header, b)
 
        /*!
@@ -397,16 +402,17 @@ extern "C" {
          @abstract  Free the memory allocated for an alignment.
          @param  b  pointer to an alignment
         */
-#define bam_destroy1(b) do {                                                                                   \
-               if ((b)->hash) bam_aux_destroy(b); free((b)->data); free(b);    \
+#define bam_destroy1(b) do {           \
+               free((b)->data); free(b);       \
        } while (0)
 
        /*!
-         @abstract       Print an alignment to the standard output in TAM format.
+         @abstract       Format a BAM record in the SAM format
          @param  header  pointer to the header structure
          @param  b       alignment to print
+         @return         a pointer to the SAM string
         */
-       void bam_view1(const bam_header_t *header, const bam1_t *b);
+       char *bam_format1(const bam_header_t *header, const bam1_t *b);
 
        /*!
          @abstract    Merge multiple sorted BAM.
@@ -461,6 +467,8 @@ extern "C" {
        /*! @abstract pileup buffer */
        typedef struct __bam_plbuf_t bam_plbuf_t;
 
+       void bam_plbuf_set_mask(bam_plbuf_t *buf, int mask);
+
        /*! @typedef
          @abstract    Type of function to be called by bam_plbuf_push().
          @param  tid  chromosome ID as is defined in the header
@@ -507,17 +515,6 @@ extern "C" {
         */
        int bam_plbuf_push(const bam1_t *b, bam_plbuf_t *buf);
 
-       /*!
-         @abstract         A more convenient interface to bam_plbuf_push()
-         @param  fp        BAM file handler
-         @param  func      user defined function
-         @param  func_data user provided data
-
-         @discussion The file position indicator must be placed right
-         before the start of an alignment. See also bam_plbuf_push().
-        */
-       int bam_pileup_file(bamFile fp, bam_pileup_f func, void *func_data);
-
        struct __bam_lplbuf_t;
        typedef struct __bam_lplbuf_t bam_lplbuf_t;
 
@@ -533,7 +530,7 @@ extern "C" {
        int bam_lplbuf_push(const bam1_t *b, bam_lplbuf_t *buf);
 
        /*! @abstract  bam_plbuf_file() equivalent with level calculated. */
-       int bam_lpileup_file(bamFile fp, bam_pileup_f func, void *func_data);
+       int bam_lpileup_file(bamFile fp, int mask, bam_pileup_f func, void *func_data);
 
        struct __bam_index_t;
        typedef struct __bam_index_t bam_index_t;
@@ -594,11 +591,42 @@ extern "C" {
         */
        void bam_parse_region(bam_header_t *header, const char *str, int *ref_id, int *begin, int *end);
 
-       int32_t bam_aux_geti(bam1_t *b, const char tag[2], int *err);
-       float bam_aux_getf(bam1_t *b, const char tag[2], int *err);
-       char bam_aux_getc(bam1_t *b, const char tag[2], int *err);
-       char *bam_aux_getZH(bam1_t *b, const char tag[2], int *err);
-       void bam_aux_destroy(bam1_t *b);
+       void bam_aux_append(bam1_t *b, const char tag[2], char type, int len, uint8_t *data);
+       uint8_t *bam_aux_get(bam1_t *b, const char tag[2]);
+       int32_t bam_aux2i(const uint8_t *s);
+       float bam_aux2f(const uint8_t *s);
+       double bam_aux2d(const uint8_t *s);
+       char bam_aux2A(const uint8_t *s);
+       char *bam_aux2Z(const uint8_t *s);
+
+       /*!
+        @abstract     Get the color encoding the previous and current base
+        @param b      pointer to an alignment
+        @param i      The i-th position, 0-based
+        @return       color
+
+        @discussion   Returns 0 no color information is found.
+        */
+       char bam_aux_getCSi(bam1_t *b, int i);
+
+       /*!
+        @abstract     Get the color quality of the color encoding the previous and current base
+        @param b      pointer to an alignment
+        @param i      The i-th position, 0-based
+        @return       color quality
+
+        @discussion   Returns 0 no color information is found.
+        */
+       char bam_aux_getCQi(bam1_t *b, int i);
+
+       /*!
+        @abstract     Get the color error profile at the give position    
+        @param b      pointer to an alignment
+        @return       the original color if the color was an error, '-' (dash) otherwise
+
+        @discussion   Returns 0 no color information is found.
+        */
+       char bam_aux_getCEi(bam1_t *b, int i);
 
        /*!  
          @abstract Calculate the rightmost coordinate of an alignment on the
@@ -641,7 +669,7 @@ static inline int bam_reg2bin(uint32_t beg, uint32_t end)
        return 0;
 }
 
-static inline void bam_copy1(bam1_t *bdst, const bam1_t *bsrc)
+static inline bam1_t *bam_copy1(bam1_t *bdst, const bam1_t *bsrc)
 {
        uint8_t *data = bdst->data;
        int m_data = bdst->m_data;   // backup data and m_data
@@ -654,6 +682,18 @@ static inline void bam_copy1(bam1_t *bdst, const bam1_t *bsrc)
        // restore the backup
        bdst->m_data = m_data;
        bdst->data = data;
+       return bdst;
+}
+
+static inline bam1_t *bam_dup1(const bam1_t *src)
+{
+       bam1_t *b;
+       b = bam_init1();
+       *b = *src;
+       b->m_data = b->data_len;
+       b->data = (uint8_t*)calloc(b->data_len, 1);
+       memcpy(b->data, src->data, b->data_len);
+       return b;
 }
 
 #endif