]> git.donarmstrong.com Git - samtools.git/commitdiff
rename iterf as iter
authorHeng Li <lh3@live.co.uk>
Sat, 12 Jun 2010 22:11:32 +0000 (22:11 +0000)
committerHeng Li <lh3@live.co.uk>
Sat, 12 Jun 2010 22:11:32 +0000 (22:11 +0000)
bam.h
bam_index.c
bam_plcmd.c
samtools.1

diff --git a/bam.h b/bam.h
index f6d31a922b12a413528950c2288f4e28a7284215..fb71a49b46dafea0512a99e157735104cd3832a3 100644 (file)
--- a/bam.h
+++ b/bam.h
@@ -190,7 +190,7 @@ typedef struct {
        uint8_t *data;
 } bam1_t;
 
-typedef struct __bam_iterf_t *bam_iterf_t;
+typedef struct __bam_iter_t *bam_iter_t;
 
 #define bam1_strand(b) (((b)->core.flag&BAM_FREVERSE) != 0)
 #define bam1_mstrand(b) (((b)->core.flag&BAM_FMREVERSE) != 0)
@@ -594,9 +594,9 @@ extern "C" {
         */
        int bam_fetch(bamFile fp, const bam_index_t *idx, int tid, int beg, int end, void *data, bam_fetch_f func);
 
-       bam_iterf_t bam_iterf_query(const bam_index_t *idx, int tid, int beg, int end);
-       int bam_iterf_read(bamFile fp, bam_iterf_t iter, bam1_t *b);
-       void bam_iterf_destroy(bam_iterf_t iter);
+       bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end);
+       int bam_iter_read(bamFile fp, bam_iter_t iter, bam1_t *b);
+       void bam_iter_destroy(bam_iter_t iter);
 
        /*!
          @abstract       Parse a region in the format: "chr2:100,000-200,000".
index 595ef047733dc62282732c0dd6ece2d29537887f..6fffd75e1ed5c58c363847abd42fed04e638c5b0 100644 (file)
@@ -494,7 +494,7 @@ static inline int is_overlap(uint32_t beg, uint32_t end, const bam1_t *b)
        return (rend > beg && rbeg < end);
 }
 
-struct __bam_iterf_t {
+struct __bam_iter_t {
        int from_first; // read from the first record; no random access
        int tid, beg, end, n_off, i, finished;
        uint64_t curr_off;
@@ -502,7 +502,7 @@ struct __bam_iterf_t {
 };
 
 // bam_fetch helper function retrieves 
-bam_iterf_t bam_iterf_query(const bam_index_t *idx, int tid, int beg, int end)
+bam_iter_t bam_iter_query(const bam_index_t *idx, int tid, int beg, int end)
 {
        uint16_t *bins;
        int i, n_bins, n_off;
@@ -510,12 +510,12 @@ bam_iterf_t bam_iterf_query(const bam_index_t *idx, int tid, int beg, int end)
        khint_t k;
        khash_t(i) *index;
        uint64_t min_off;
-       bam_iterf_t iter = 0;
+       bam_iter_t iter = 0;
 
        if (beg < 0) beg = 0;
        if (end < beg) return 0;
        // initialize iter
-       iter = calloc(1, sizeof(struct __bam_iterf_t));
+       iter = calloc(1, sizeof(struct __bam_iter_t));
        iter->tid = tid, iter->beg = beg, iter->end = end; iter->i = -1;
        //
        bins = (uint16_t*)calloc(MAX_BIN, 2);
@@ -582,20 +582,20 @@ bam_iterf_t bam_iterf_query(const bam_index_t *idx, int tid, int beg, int end)
 
 pair64_t *get_chunk_coordinates(const bam_index_t *idx, int tid, int beg, int end, int *cnt_off)
 { // for pysam compatibility
-       bam_iterf_t iter;
+       bam_iter_t iter;
        pair64_t *off;
-       iter = bam_iterf_query(idx, tid, beg, end);
+       iter = bam_iter_query(idx, tid, beg, end);
        off = iter->off; *cnt_off = iter->n_off;
        free(iter);
        return off;
 }
 
-void bam_iterf_destroy(bam_iterf_t iter)
+void bam_iter_destroy(bam_iter_t iter)
 {
        if (iter) { free(iter->off); free(iter); }
 }
 
-int bam_iterf_read(bamFile fp, bam_iterf_t iter, bam1_t *b)
+int bam_iter_read(bamFile fp, bam_iter_t iter, bam1_t *b)
 {
        if (iter->finished) return -1;
        if (iter->from_first) {
@@ -627,11 +627,11 @@ int bam_iterf_read(bamFile fp, bam_iterf_t iter, bam1_t *b)
 
 int bam_fetch(bamFile fp, const bam_index_t *idx, int tid, int beg, int end, void *data, bam_fetch_f func)
 {
-       bam_iterf_t iter;
+       bam_iter_t iter;
        bam1_t *b;
        b = bam_init1();
-       iter = bam_iterf_query(idx, tid, beg, end);
-       while (bam_iterf_read(fp, iter, b) >= 0) func(b, data);
+       iter = bam_iter_query(idx, tid, beg, end);
+       while (bam_iter_read(fp, iter, b) >= 0) func(b, data);
        bam_destroy1(b);
        return 0;
 }
index d34cd267a9edd363a32605c083fefe5ba512e8bb..d1ea03f896e9ad433c95d33ba1dfe88ad2354aae 100644 (file)
@@ -453,13 +453,13 @@ typedef struct {
 
 typedef struct {
        bamFile fp;
-       bam_iterf_t iter;
+       bam_iter_t iter;
 } mplp_aux_t;
 
 static int mplp_func(void *data, bam1_t *b)
 {
        mplp_aux_t *ma = (mplp_aux_t*)data;
-       if (ma->iter) return bam_iterf_read(ma->fp, ma->iter, b);
+       if (ma->iter) return bam_iter_read(ma->fp, ma->iter, b);
        return bam_read1(ma->fp, b);
 }
 
@@ -494,7 +494,7 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                                exit(1);
                        }
                        if (i == 0) beg0 = beg, end0 = end;
-                       data[i]->iter = bam_iterf_query(idx, tid, beg, end);
+                       data[i]->iter = bam_iter_query(idx, tid, beg, end);
                        bam_index_destroy(idx);
                }
                if (i == 0) h = h_tmp;
@@ -536,7 +536,7 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
        bam_header_destroy(h);
        for (i = 0; i < n; ++i) {
                bam_close(data[i]->fp);
-               if (data[i]->iter) bam_iterf_destroy(data[i]->iter);
+               if (data[i]->iter) bam_iter_destroy(data[i]->iter);
                free(data[i]);
        }
        free(data); free(plp); free(ref); free(n_plp);
index b0287791af092cc1175d88f44525a7a3db426437..366a5630ececd0956b6225ad44bbc5d5fb7f087c 100644 (file)
@@ -18,6 +18,8 @@ samtools faidx ref.fasta
 .PP
 samtools pileup -f ref.fasta aln.sorted.bam
 .PP
+samtools mpileup -f ref.fasta -r chr3:1,000-2,000 in1.bam in2.bam
+.PP
 samtools tview aln.sorted.bam ref.fasta
 
 .SH DESCRIPTION