]> git.donarmstrong.com Git - samtools.git/blobdiff - bam_plcmd.c
remove VCF output in mpileup
[samtools.git] / bam_plcmd.c
index 2afebeadf85bc2ee020e3c8aecd0209dccebc51d..80799e3b28af95208796424c31a4bf5b7dd24ba0 100644 (file)
@@ -5,7 +5,7 @@
 #include "sam.h"
 #include "faidx.h"
 #include "bam_maqcns.h"
-#include "bam_mcns.h"
+#include "bam2bcf.h"
 #include "khash.h"
 #include "glf.h"
 #include "kstring.h"
@@ -342,12 +342,12 @@ int bam_pileup(int argc, char *argv[])
     d->max_depth = 0;
        d->tid = -1; d->mask = BAM_DEF_MASK;
        d->c = bam_maqcns_init();
-       d->c->is_soap = 1; // change the default model
+       d->c->errmod = BAM_ERRMOD_MAQ2; // change the default model
        d->ido = bam_maqindel_opt_init();
        while ((c = getopt(argc, argv, "st:f:cT:N:r:l:d:im:gI:G:vM:S2aR:PA")) >= 0) {
                switch (c) {
-               case 'a': d->c->is_soap = 1; break;
-               case 'A': d->c->is_soap = 0; break;
+               case 'a': d->c->errmod = BAM_ERRMOD_SOAP; break;
+               case 'A': d->c->errmod = BAM_ERRMOD_MAQ; break;
                case 's': d->format |= BAM_PLF_SIMPLE; break;
                case 't': fn_list = strdup(optarg); break;
                case 'l': fn_pos = strdup(optarg); break;
@@ -449,8 +449,12 @@ int bam_pileup(int argc, char *argv[])
  * mpileup *
  ***********/
 
+#define MPLP_GLF   0x10
+#define MPLP_NO_COMP 0x20
+
 typedef struct {
-       int vcf, max_mq, min_mq, var_only;
+       int max_mq, min_mq, flag, min_baseQ;
+       double theta;
        char *reg, *fn_pos;
        faidx_t *fai;
        kh_64_t *hash;
@@ -475,17 +479,24 @@ static int mplp_func(void *data, bam1_t *b)
 static int mpileup(mplp_conf_t *conf, int n, char **fn)
 {
        mplp_aux_t **data;
-       mc_aux_t *ma = 0;
        int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid;
        const bam_pileup1_t **plp;
        bam_mplp_t iter;
        bam_header_t *h = 0;
        char *ref;
        khash_t(64) *hash = 0;
-       // allocate
+
+       bcf_callaux_t *bca = 0;
+       bcf_callret1_t *bcr = 0;
+       bcf_call_t bc;
+       bcf_t *bp = 0;
+       bcf_hdr_t *bh = 0;
+
+       memset(&bc, 0, sizeof(bcf_call_t));
        data = calloc(n, sizeof(void*));
        plp = calloc(n, sizeof(void*));
        n_plp = calloc(n, sizeof(int*));
+
        // read the header and initialize data
        for (i = 0; i < n; ++i) {
                bam_header_t *h_tmp;
@@ -517,23 +528,39 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
        }
        if (conf->fn_pos) hash = load_pos(conf->fn_pos, h);
        // write the VCF header
-       if (conf->vcf) {
+       if (conf->flag & MPLP_GLF) {
                kstring_t s;
+               bh = calloc(1, sizeof(bcf_hdr_t));
                s.l = s.m = 0; s.s = 0;
-               puts("##fileformat=VCFv4.0");
-               kputs("#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO\tFORMAT", &s);
+               bp = bcf_open("-", (conf->flag&MPLP_NO_COMP)? "wu" : "w");
+               for (i = 0; i < h->n_targets; ++i) {
+                       kputs(h->target_name[i], &s);
+                       kputc('\0', &s);
+               }
+               bh->l_nm = s.l;
+               bh->name = malloc(s.l);
+               memcpy(bh->name, s.s, s.l);
+               s.l = 0;
                for (i = 0; i < n; ++i) {
                        const char *p;
-                       kputc('\t', &s);
                        if ((p = strstr(fn[i], ".bam")) != 0)
                                kputsn(fn[i], p - fn[i], &s);
                        else kputs(fn[i], &s);
+                       kputc('\0', &s);
                }
-               puts(s.s);
+               bh->l_smpl = s.l;
+               bh->sname = malloc(s.l);
+               memcpy(bh->sname, s.s, s.l);
+               bh->l_txt = 0;
                free(s.s);
+               bcf_hdr_sync(bh);
+               bcf_hdr_write(bp, bh);
        }
        // mpileup
-       if (conf->vcf) ma = mc_init(n);
+       if (conf->flag & MPLP_GLF) {
+               bca = bcf_call_init(-1., conf->min_baseQ);
+               bcr = calloc(n, sizeof(bcf_callret1_t));
+       }
        ref_tid = -1; ref = 0;
        iter = bam_mplp_init(n, mplp_func, (void**)data);
        while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
@@ -548,55 +575,18 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        if (conf->fai) ref = fai_fetch(conf->fai, h->target_name[tid], &ref_len);
                        ref_tid = tid;
                }
-               if (conf->vcf) {
-                       double f0, f, pref = -1.0; // the reference allele frequency
-                       int j, _ref, _alt, _ref0, depth, rms_q, _ref0b, is_var = 0, qref = 0;
-                       uint64_t sqr_sum;
-                       _ref0 = _ref0b = (ref && pos < ref_len)? ref[pos] : 'N';
-                       _ref0 = bam_nt16_nt4_table[bam_nt16_table[_ref0]];
-                       f = f0 = mc_freq0(_ref0, n_plp, plp, ma, &_ref, &_alt);
-                       if (f >= 0.0) {
-                               double q, flast = f;
-                               for (j = 0; j < 10; ++j) {
-                                       f = mc_freq_iter(flast, ma);
-//                                     printf("%lg->%lg\n", flast, f);
-                                       if (fabs(f - flast) < 1e-3) break;
-                                       flast = f;
-                               }
-                               if (1. - f < 1e-4) f = 1.;
-                               pref = mc_ref_prob(ma);
-                               is_var = (pref < .5);
-                               q = is_var? pref : 1. - pref;
-                               if (q < 1e-308) q = 1e-308;
-                               qref = (int)(-3.434 * log(q) + .499);
-                               if (qref > 99) qref = 99;
-                       }
-                       if (conf->var_only && !is_var) continue;
-                       printf("%s\t%d\t.\t%c\t", h->target_name[tid], pos + 1, _ref0b);
-                       if (is_var) {
-                               if (_ref0 == _ref) putchar("ACGTN"[_alt]);
-                               else printf("%c,%c", "ACGTN"[_ref], "ACGTN"[_alt]);
-                       } else putchar('.');
-                       printf("\t%d\t", qref);
-                       if (f0 < 0.) printf("Q13\t");
-                       else printf(".\t");
-                       for (i = depth = 0, sqr_sum = 0; i < n; ++i) {
-                               depth += n_plp[i];
-                               for (j = 0; j < n_plp[i]; ++j) {
-                                       int q = plp[i][j].b->core.qual;
-                                       if (q > conf->max_mq) q = conf->max_mq;
-                                       sqr_sum += q * q;
-                               }
-                       }
-                       rms_q = (int)(sqrt((double)sqr_sum / depth) + .499);
-                       printf("DP=%d;MQ=%d;AF=%.3lg", depth, rms_q, f0<0.?0.:1.-f);
-                       printf("\tGT:GQ:DP");
-                       // output genotype information; FIXME: to be implmented...
-                       for (i = 0; i < n; ++i) {
-                               int x = mc_call_gt(ma, f, i);
-                               printf("\t%c/%c:%d:%d", "10"[((x&3)==2)], "10"[((x&3)>0)], x>>2, n_plp[i]);
-                       }
-                       putchar('\n');
+               if (conf->flag & MPLP_GLF) {
+                       int _ref0, ref16;
+                       bcf1_t *b = calloc(1, sizeof(bcf1_t));
+                       _ref0 = (ref && pos < ref_len)? ref[pos] : 'N';
+                       ref16 = bam_nt16_table[_ref0];
+                       for (i = 0; i < n; ++i)
+                               bcf_call_glfgen(n_plp[i], plp[i], ref16, bca, bcr + i);
+                       bcf_call_combine(n, bcr, ref16, &bc);
+                       bcf_call2bcf(tid, pos, &bc, b);
+                       bcf_write(bp, bh, b);
+                       //fprintf(stderr, "%d,%d,%d\n", b->tid, b->pos, b->l_str);
+                       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) {
@@ -618,13 +608,14 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        putchar('\n');
                }
        }
+       bcf_close(bp);
        if (hash) { // free the hash table
                khint_t k;
                for (k = kh_begin(hash); k < kh_end(hash); ++k)
                        if (kh_exist(hash, k)) free(kh_val(hash, k));
                kh_destroy(64, hash);
        }
-       mc_destroy(ma);
+       bcf_hdr_destroy(bh); bcf_call_destroy(bca); free(bc.PL); free(bcr);
        bam_mplp_destroy(iter);
        bam_header_destroy(h);
        for (i = 0; i < n; ++i) {
@@ -642,18 +633,22 @@ int bam_mpileup(int argc, char *argv[])
        mplp_conf_t mplp;
        memset(&mplp, 0, sizeof(mplp_conf_t));
        mplp.max_mq = 60;
-       while ((c = getopt(argc, argv, "f:r:l:VvM:q:")) >= 0) {
+       mplp.theta = 1e-3;
+       mplp.min_baseQ = 13;
+       while ((c = getopt(argc, argv, "gf:r:l:M:q:t:Q:u")) >= 0) {
                switch (c) {
+               case 't': mplp.theta = atof(optarg); break;
                case 'f':
                        mplp.fai = fai_load(optarg);
                        if (mplp.fai == 0) return 1;
                        break;
                case 'r': mplp.reg = strdup(optarg); break;
                case 'l': mplp.fn_pos = strdup(optarg); break;
-               case 'v': mplp.var_only = 1; break;
-               case 'V': mplp.vcf = 1; break;
+               case 'g': mplp.flag |= MPLP_GLF; break;
+               case 'u': mplp.flag |= MPLP_NO_COMP; break;
                case 'M': mplp.max_mq = atoi(optarg); break;
                case 'q': mplp.min_mq = atoi(optarg); break;
+               case 'Q': mplp.min_baseQ = atoi(optarg); break;
                }
        }
        if (argc == 1) {
@@ -663,11 +658,12 @@ int bam_mpileup(int argc, char *argv[])
                fprintf(stderr, "         -r STR      region in which pileup is generated [null]\n");
                fprintf(stderr, "         -l FILE     list of positions (format: chr pos) [null]\n");
                fprintf(stderr, "         -M INT      cap mapping quality at INT [%d]\n", mplp.max_mq);
+               fprintf(stderr, "         -Q INT      min base quality [%d]\n", mplp.min_baseQ);
                fprintf(stderr, "         -q INT      filter out alignment with MQ smaller than INT [%d]\n", mplp.min_mq);
-               fprintf(stderr, "         -V          generate VCF output (SNP calling)\n");
-               fprintf(stderr, "         -v          show variant sites only\n");
+               fprintf(stderr, "         -t FLOAT    scaled mutation rate [%lg]\n", mplp.theta);
+               fprintf(stderr, "         -g          generate GLF output\n");
                fprintf(stderr, "\n");
-               fprintf(stderr, "Notes: Assuming error independency and diploid individuals.\n\n");
+               fprintf(stderr, "Notes: Assuming diploid individuals.\n\n");
                return 1;
        }
        mpileup(&mplp, argc - optind, argv + optind);