]> git.donarmstrong.com Git - samtools.git/commitdiff
Revert to r790. The recent changes are not good...
authorHeng Li <lh3@live.co.uk>
Fri, 5 Nov 2010 15:28:23 +0000 (15:28 +0000)
committerHeng Li <lh3@live.co.uk>
Fri, 5 Nov 2010 15:28:23 +0000 (15:28 +0000)
bam2bcf.c
bam2bcf.h
bam_plcmd.c

index ab5bfcd20d4a6233109821846842d9996aab93a9..509761d3395df7d01e56f7a368f2f2daa29be4b3 100644 (file)
--- a/bam2bcf.c
+++ b/bam2bcf.c
@@ -3,6 +3,7 @@
 #include "bam.h"
 #include "kstring.h"
 #include "bam2bcf.h"
+#include "errmod.h"
 #include "bcftools/bcf.h"
 
 extern void ks_introsort_uint32_t(size_t n, uint32_t a[]);
@@ -13,15 +14,18 @@ extern      void ks_introsort_uint32_t(size_t n, uint32_t a[]);
 
 #define CAP_DIST 25
 
+struct __bcf_callaux_t {
+       int max_bases, capQ, min_baseQ;
+       uint16_t *bases;
+       errmod_t *e;
+};
+
 bcf_callaux_t *bcf_call_init(double theta, int min_baseQ)
 {
        bcf_callaux_t *bca;
        if (theta <= 0.) theta = CALL_DEFTHETA;
        bca = calloc(1, sizeof(bcf_callaux_t));
        bca->capQ = 60;
-       bca->openQ = 40;
-       bca->extQ = 20;
-       bca->tandemQ = 100;
        bca->min_baseQ = min_baseQ;
        bca->e = errmod_init(1. - theta);
        return bca;
@@ -53,7 +57,7 @@ int bcf_call_glfgen(int _n, const bam_pileup1_t *pl, int ref_base /*4-bit*/, bcf
                const bam_pileup1_t *p = pl + i;
                int q, b, mapQ, baseQ, is_diff, min_dist;
                // set base
-               if (p->is_del || p->is_refskip || (p->b->core.flag&BAM_FUNMAP)) continue;
+               if (p->is_del || (p->b->core.flag&BAM_FUNMAP)) continue; // skip unmapped reads and deleted bases
                baseQ = q = (int)bam1_qual(p->b)[p->qpos]; // base quality
                if (q < bca->min_baseQ) continue;
                mapQ = p->b->core.qual < bca->capQ? p->b->core.qual : bca->capQ;
@@ -83,125 +87,11 @@ int bcf_call_glfgen(int _n, const bam_pileup1_t *pl, int ref_base /*4-bit*/, bcf
        return r->depth;
 }
 
-int bcf_call_glfgen_gap(int pos, int _n, const bam_pileup1_t *pl, bcf_callaux_t *bca, bcf_callret1_t *r)
-{
-       int i, n, n_ins, n_del;
-       memset(r, 0, sizeof(bcf_callret1_t));
-       if (_n == 0) return -1;
-
-       // enlarge the bases array if necessary
-       if (bca->max_bases < _n) {
-               bca->max_bases = _n;
-               kroundup32(bca->max_bases);
-               bca->bases = (uint16_t*)realloc(bca->bases, 2 * bca->max_bases);
-       }
-       // fill the bases array
-       memset(r, 0, sizeof(bcf_callret1_t));
-       r->indelreg = 10000;
-       n_ins = n_del = 0;
-       for (i = n = 0; i < _n; ++i) {
-               const bam_pileup1_t *p = pl + i;
-               int q, b, mapQ, indelQ, is_diff, min_dist;
-               if (p->is_del || p->is_refskip || (p->b->core.flag&BAM_FUNMAP)) continue;
-               { // compute indel (base) quality
-                       // this can be made more efficient, but realignment is the bottleneck anyway
-                       int j, k, x, y, op, len = 0, max_left, max_rght, seqQ, indelreg;
-                       bam1_core_t *c = &p->b->core;
-                       uint32_t *cigar = bam1_cigar(p->b);
-                       uint8_t *qual = bam1_qual(p->b);
-                       for (k = y = 0, x = c->pos; k < c->n_cigar && y <= p->qpos; ++k) {
-                               op = cigar[k]&0xf;
-                               len = cigar[k]>>4;
-                               if (op == BAM_CMATCH) {
-                                       if (pos > x && pos < x + len) break;
-                                       x += len; y += len;
-                               } else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += len;
-                               else if (op == BAM_CDEL || op == BAM_CREF_SKIP) x += len;
-                       }
-                       if (k == c->n_cigar) continue; // this actually should not happen
-                       max_left = max_rght = 0; indelreg = 0;
-                       if (pos == x + len - 1 && k+2 < c->n_cigar && ((cigar[k+1]&0xf) == BAM_CINS || (cigar[k+1]&0xf) == BAM_CDEL)
-                               && (cigar[k+2]&0xf) == BAM_CMATCH)
-                       {
-                               for (j = y; j < y + len; ++j)
-                                       if (max_left < qual[j]) max_left = qual[j];
-                               if ((cigar[k+1]&0xf) == BAM_CINS) y += cigar[k+1]>>4;
-                               else x += cigar[k+1]>>4;
-                               op = cigar[k+2]&0xf; len = cigar[k+2]>>4;
-                               for (j = y; j < y + len; ++j) {
-                                       if (max_rght < qual[j]) max_rght = qual[j];
-                                       if (qual[j] > BAM2BCF_INDELREG_THRES && indelreg == 0)
-                                               indelreg = j - y + 1;
-                               }
-                               if (r->indelreg > indelreg) r->indelreg = indelreg;
-                       } else {
-                               for (j = y; j <= p->qpos; ++j)
-                                       if (max_left < qual[j]) max_left = qual[j];
-                               for (j = p->qpos + 1; j < y + len; ++j)
-                                       if (max_rght < qual[j]) max_rght = qual[j];
-                                       
-                       }
-                       indelQ = max_left < max_rght? max_left : max_rght;
-                       // estimate the sequencing error rate
-                       seqQ = bca->openQ;
-                       if (p->indel != 0) seqQ += bca->extQ * (abs(p->indel) - 1); // FIXME: better to model homopolymer
-                       if (p->indel > 0) {
-                               ++n_ins; r->ins_len += p->indel;
-                       } else if (p->indel < 0) {
-                               ++n_del; r->del_len += -p->indel;
-                       }
-                       if (p->indel != 0) { // a different model for tandem repeats
-                               uint8_t *seq = bam1_seq(p->b);
-                               int tandemQ, qb = bam1_seqi(seq, p->qpos), l;
-                               for (j = p->qpos + 1; j < c->l_qseq; ++j)
-                                       if (qb != bam1_seqi(seq, j)) break;
-                               l = j;
-                               for (j = (int)p->qpos - 1; j >= 0; --j)
-                                       if (qb != bam1_seqi(seq, j)) break;
-                               l = l - (j + 1);
-                               tandemQ = (int)((double)(abs(p->indel)) / l * bca->tandemQ + .499);
-                               if (seqQ > tandemQ) seqQ = tandemQ;
-                       }
-//                     fprintf(stderr, "%s\t%d\t%d\t%d\t%d\t%d\t%d\n", bam1_qname(p->b), pos+1, p->indel, indelQ, seqQ, max_left, max_rght);
-                       if (indelQ > seqQ) indelQ = seqQ;
-                       q = indelQ;
-               }
-               if (q < bca->min_baseQ) continue;
-               mapQ = p->b->core.qual < bca->capQ? p->b->core.qual : bca->capQ;
-               if (q > mapQ) q = mapQ;
-               if (q > 63) q = 63;
-               if (q < 4) q = 4;
-               b = p->indel? 1 : 0;
-               bca->bases[n++] = q<<5 | (int)bam1_strand(p->b)<<4 | b;
-               // collect annotations
-               r->qsum[b] += q;
-               is_diff = b;
-               ++r->anno[0<<2|is_diff<<1|bam1_strand(p->b)];
-               min_dist = p->b->core.l_qseq - 1 - p->qpos;
-               if (min_dist > p->qpos) min_dist = p->qpos;
-               if (min_dist > CAP_DIST) min_dist = CAP_DIST;
-               r->anno[1<<2|is_diff<<1|0] += indelQ;
-               r->anno[1<<2|is_diff<<1|1] += indelQ * indelQ;
-               r->anno[2<<2|is_diff<<1|0] += mapQ;
-               r->anno[2<<2|is_diff<<1|1] += mapQ * mapQ;
-               r->anno[3<<2|is_diff<<1|0] += min_dist;
-               r->anno[3<<2|is_diff<<1|1] += min_dist * min_dist;
-       }
-       r->depth = n;
-       r->ins_len = n_ins? r->ins_len / n_ins : 0;
-       r->del_len = n_del? r->del_len / n_del : 0;
-       if (r->indelreg >= 10000) r->indelreg = 0;
-       // glfgen
-       errmod_cal(bca->e, n, 2, bca->bases, r->p);
-       return r->depth;        
-}
-
 int bcf_call_combine(int n, const bcf_callret1_t *calls, int ref_base /*4-bit*/, bcf_call_t *call)
 {
        int ref4, i, j, qsum[4];
        int64_t tmp;
        call->ori_ref = ref4 = bam_nt16_nt4_table[ref_base];
-       call->ins_len = call->del_len = 0; call->indelreg = 0;
        if (ref4 > 4) ref4 = 4;
        // calculate qsum
        memset(qsum, 0, 4 * sizeof(int));
@@ -264,66 +154,6 @@ int bcf_call_combine(int n, const bcf_callret1_t *calls, int ref_base /*4-bit*/,
        return 0;
 }
 
-int bcf_call_combine_gap(int n, const bcf_callret1_t *calls, bcf_call_t *call)
-{
-       int i, j, n_ins, n_del;
-       // combine annotations
-       call->ori_ref = 4;
-       memset(call->anno, 0, 16 * sizeof(int));
-       call->ins_len = call->del_len = 0; call->indelreg = 10000;
-       for (i = call->depth = 0, n_ins = n_del = 0; i < n; ++i) {
-               const bcf_callret1_t *r = calls + i;
-               if (r->depth > 0) {
-                       call->depth += r->depth;
-                       if (r->ins_len > 0) {
-                               call->ins_len += r->ins_len;
-                               ++n_ins;
-                       }
-                       if (r->del_len > 0) {
-                               call->del_len += r->del_len;
-                               ++n_del;
-                       }
-                       if (r->indelreg > 0 && call->indelreg > r->indelreg)
-                               call->indelreg = r->indelreg;
-                       for (j = 0; j < 16; ++j) call->anno[j] += r->anno[j];
-               }
-       }
-       if (call->depth == 0) return 0; // no indels
-       call->ins_len = n_ins? call->ins_len / n_ins : 0;
-       call->del_len = n_del? call->del_len / n_del : 0;
-       //
-       for (i = 0; i < 5; ++i) call->a[i] = -1;
-       call->a[0] = 0; call->a[1] = 1;
-       call->unseen = -1;
-       call->n_alleles = 2;
-       // set the PL array
-       if (call->n < n) {
-               call->n = n;
-               call->PL = realloc(call->PL, 15 * n);
-       }
-       {
-               int g[3];
-               double sum_min = 0.;
-               g[0] = 0; g[1] = 1; g[2] = 3;
-               for (i = 0; i < n; ++i) {
-                       uint8_t *PL = call->PL + 3 * i;
-                       const bcf_callret1_t *r = calls + i;
-                       float min = 1e37;
-                       for (j = 0; j < 3; ++j)
-                               if (min > r->p[g[j]]) min = r->p[g[j]];
-                       sum_min += min;
-                       for (j = 0; j < 3; ++j) {
-                               int y;
-                               y = (int)(r->p[g[j]] - min + .499);
-                               if (y > 255) y = 255;
-                               PL[j] = y;
-                       }
-               }
-               call->shift = (int)(sum_min + .499);
-       }
-       return 0;
-}
-
 int bcf_call2bcf(int tid, int pos, bcf_call_t *bc, bcf1_t *b, bcf_callret1_t *bcr, int is_SP)
 {
        extern double kt_fisher_exact(int n11, int n12, int n21, int n22, double *_left, double *_right, double *two);
@@ -333,19 +163,11 @@ int bcf_call2bcf(int tid, int pos, bcf_call_t *bc, bcf1_t *b, bcf_callret1_t *bc
        b->tid = tid; b->pos = pos; b->qual = 0;
        s.s = b->str; s.m = b->m_str; s.l = 0;
        kputc('\0', &s);
-       if (bc->ins_len > 0 || bc->del_len > 0) { // an indel
-               for (i = 0; i < bc->indelreg; ++i) kputc('N', &s);
-               kputc('\0', &s);
-               if (bc->ins_len > 0 && bc->del_len > 0) kputs("<INDEL>", &s);
-               else if (bc->ins_len > 0) kputs("<INS>", &s);
-               else if (bc->del_len > 0) kputs("<DEL>", &s);
-       } else { // SNP
-               kputc("ACGTN"[bc->ori_ref], &s); kputc('\0', &s);
-               for (i = 1; i < 5; ++i) {
-                       if (bc->a[i] < 0) break;
-                       if (i > 1) kputc(',', &s);
-                       kputc(bc->unseen == i? 'X' : "ACGT"[bc->a[i]], &s);
-               }
+       kputc("ACGTN"[bc->ori_ref], &s); kputc('\0', &s);
+       for (i = 1; i < 5; ++i) {
+               if (bc->a[i] < 0) break;
+               if (i > 1) kputc(',', &s);
+               kputc(bc->unseen == i? 'X' : "ACGT"[bc->a[i]], &s);
        }
        kputc('\0', &s);
        kputc('\0', &s);
index b6514c05c2b5633a93104524bcd4fa8533c3ba80..4de743a9287f7421e6f612a2426897d3dbe69212 100644 (file)
--- a/bam2bcf.h
+++ b/bam2bcf.h
@@ -2,23 +2,13 @@
 #define BAM2BCF_H
 
 #include <stdint.h>
-#include "errmod.h"
 #include "bcftools/bcf.h"
 
-#define BAM2BCF_INDELREG_THRES 13
-
-typedef struct __bcf_callaux_t {
-       int max_bases, capQ, min_baseQ;
-       int openQ, extQ, tandemQ;
-       // for internal uses
-       uint16_t *bases;
-       errmod_t *e;
-} bcf_callaux_t;
+struct __bcf_callaux_t;
+typedef struct __bcf_callaux_t bcf_callaux_t;
 
 typedef struct {
-       int depth, indelreg;
-       float ins_len, del_len;
-       int qsum[4];
+       int depth, qsum[4];
        int anno[16];
        float p[25];
 } bcf_callret1_t;
@@ -26,9 +16,7 @@ typedef struct {
 typedef struct {
        int a[5]; // alleles: ref, alt, alt2, alt3
        int n, n_alleles, shift, ori_ref, unseen;
-       int indelreg;
-       float ins_len, del_len;
-       int depth, anno[16];
+       int anno[16], depth;
        uint8_t *PL;
 } bcf_call_t;
 
@@ -40,9 +28,7 @@ extern "C" {
        void bcf_call_destroy(bcf_callaux_t *bca);
        int bcf_call_glfgen(int _n, const bam_pileup1_t *pl, int ref_base /*4-bit*/, bcf_callaux_t *bca, bcf_callret1_t *r);
        int bcf_call_combine(int n, const bcf_callret1_t *calls, int ref_base /*4-bit*/, bcf_call_t *call);
-       int bcf_call_combine_gap(int n, const bcf_callret1_t *calls, bcf_call_t *call);
        int bcf_call2bcf(int tid, int pos, bcf_call_t *bc, bcf1_t *b, bcf_callret1_t *bcr, int is_SP);
-       int bcf_call_glfgen_gap(int pos, int _n, const bam_pileup1_t *pl, bcf_callaux_t *bca, bcf_callret1_t *r);
 
 #ifdef __cplusplus
 }
index 32d050d68413aa18e6ba4ea688615579a539bf2b..a3e6aeb614927733676b2b25e4592786cb14ef2e 100644 (file)
@@ -729,29 +729,6 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0, (conf->flag&MPLP_FMT_SP));
                        bcf_write(bp, bh, b);
                        bcf_destroy(b);
-                       // test if there is an indel
-                       if (conf->flag&MPLP_REALN) {
-                               for (i = 0; i < gplp.n; ++i) {
-                                       int j;
-                                       for (j = 0; j < gplp.n_plp[i]; ++j) {
-                                               const bam_pileup1_t *p = gplp.plp[i] + j;
-                                               if (p->indel != 0) break;
-                                       }
-                                       if (j != gplp.n_plp[i]) break;
-                               }
-                               if (i < gplp.n) { // at least one of the read contains a gap
-                                       for (i = 0; i < gplp.n; ++i)
-                                               bcf_call_glfgen_gap(pos, gplp.n_plp[i], gplp.plp[i], bca, bcr + i);
-                                       bcf_call_combine_gap(gplp.n, bcr, &bc);
-                                       if (bc.depth > 0) {
-                                               b = calloc(1, sizeof(bcf1_t));
-                                               bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
-                                                                        (conf->flag&MPLP_FMT_SP));
-                                               bcf_write(bp, bh, b);
-                                               bcf_destroy(b);
-                                       }
-                               }
-                       }
                } else {
                        printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
                        for (i = 0; i < n; ++i) {