]> git.donarmstrong.com Git - samtools.git/blobdiff - bam_plcmd.c
supposedly this is THE correct implementation, but more testing is needed
[samtools.git] / bam_plcmd.c
index 0cc2d6b6710dbedc53f89803971e8170233396eb..bea91692de426eb0daf5150b5d7b190fa74155ff 100644 (file)
@@ -451,7 +451,6 @@ int bam_pileup(int argc, char *argv[])
 
 #define MPLP_VCF   0x1
 #define MPLP_VAR   0x2
-#define MPLP_AFS   0x4
 #define MPLP_AFALL 0x8
 
 #define MPLP_AFS_BLOCK 0x10000
@@ -531,9 +530,10 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                s.l = s.m = 0; s.s = 0;
                puts("##fileformat=VCFv4.0");
                puts("##INFO=<ID=DP,Number=1,Type=Integer,Description=\"Total read depth\">");
-//             puts("##INFO=<ID=AF,Number=1,Type=Float,Description=\"Posterior non-reference allele frequency\">");
-//             puts("##INFO=<ID=AFEM,Number=1,Type=Float,Description=\"Prior-free non-reference allele frequency by EM\">");
+               puts("##INFO=<ID=AF,Number=1,Type=Float,Description=\"Non-reference allele frequency \\argmax_f P(D|f)\">");
+               puts("##INFO=<ID=AFE,Number=1,Type=Float,Description=\"Expected non-reference allele frequency\">");
                puts("##FILTER=<ID=Q13,Description=\"All min{baseQ,mapQ} below 13\">");
+               puts("##FILTER=<ID=FPE,Description=\"Floating point error\">");
                kputs("#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO\tFORMAT", &s);
                for (i = 0; i < n; ++i) {
                        const char *p;
@@ -568,7 +568,6 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        mc_rst_t r;
                        int j, _ref0, depth, rms_q, _ref0b, is_var = 0, qref = 0, level = 2, tot;
                        uint64_t sqr_sum;
-                       if (conf->flag & MPLP_AFS) level = 3;
                        _ref0 = _ref0b = (ref && pos < ref_len)? ref[pos] : 'N';
                        _ref0 = bam_nt16_nt4_table[bam_nt16_table[_ref0]];
                        tot = mc_cal(_ref0, n_plp, plp, ma, &r, level);
@@ -584,11 +583,12 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        ++N; // number of processed lines
                        printf("%s\t%d\t.\t%c\t", h->target_name[tid], pos + 1, _ref0b);
                        if (is_var) {
-                               if (_ref0 == r.ref) putchar("ACGTN"[r.alt]);
-                               else printf("%c,%c", "ACGTN"[r.ref], "ACGTN"[r.alt]);
+                               putchar("ACGTN"[r.alt]);
+                               if (r.alt2 >= 0 && r.alt2 < 4) printf(",%c", "ACGT"[r.alt2]);
                        } else putchar('.');
                        printf("\t%d\t", qref);
                        if (!tot) printf("Q13\t");
+                       else if (r.f_exp < 0.) printf("FPE\t");
                        else printf(".\t");
                        for (i = depth = 0, sqr_sum = 0; i < n; ++i) {
                                depth += n_plp[i];
@@ -602,15 +602,14 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        printf("DP=%d;MQ=%d", depth, rms_q);
                        if (tot) {
                                printf(";AF=%.3lf", 1. - r.f_em);
-                               if (conf->flag & MPLP_AFALL) {
-                                       printf(";AF0=%.3lf;AFN=%.3lf;AFE=%.3lf", 1-r.f_naive, 1-r.f_nielsen, 1-r.f_exp);
-                                       if (conf->flag & MPLP_AFS) printf(";AFB=%.3lf", 1-r.f_map);
-                               }
+                               if (level >= 2) printf(";AFE=%.3lf", 1-r.f_exp);
+                               if (conf->flag & MPLP_AFALL)
+                                       printf(";AF0=%.3lf;AFN=%.3lf", 1-r.f_naive, 1-r.f_nielsen);
                        }
                        printf("\tGT:GQ:DP");
                        if (tot) {
                                for (i = 0; i < n; ++i) {
-                                       int x = mc_call_gt(ma, r.f_em, i);
+                                       int x = mc_call_gt(ma, r.f_exp, i);
                                        printf("\t%c/%c:%d:%d", "10"[((x&3)==2)], "10"[((x&3)>0)], x>>2, n_plp[i]);
                                }
                        } else for (i = 0; i < n; ++i) printf("\t./.:0:0");
@@ -637,7 +636,7 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        putchar('\n');
                }
        }
-       if ((conf->flag&MPLP_VCF) && (conf->flag&MPLP_AFS)) mc_dump_afs(ma);
+       if (conf->flag&MPLP_VCF) mc_dump_afs(ma);
        if (hash) { // free the hash table
                khint_t k;
                for (k = kh_begin(hash); k < kh_end(hash); ++k)
@@ -664,10 +663,18 @@ int bam_mpileup(int argc, char *argv[])
        mplp.max_mq = 60;
        mplp.prior_type = MC_PTYPE_FULL;
        mplp.theta = 1e-3;
-       while ((c = getopt(argc, argv, "vVcF2Sf:r:l:VM:q:t:")) >= 0) {
+       while ((c = getopt(argc, argv, "vVcFSP:f:r:l:VM:q:t:")) >= 0) {
                switch (c) {
                case 't': mplp.theta = atof(optarg); break;
-               case '2': mplp.prior_type = MC_PTYPE_COND2; break;
+               case 'P':
+                       if (strcmp(optarg, "full") == 0) mplp.prior_type = MC_PTYPE_FULL;
+                       else if (strcmp(optarg, "cond2") == 0) mplp.prior_type = MC_PTYPE_COND2;
+                       else if (strcmp(optarg, "flat") == 0) mplp.prior_type = MC_PTYPE_FLAT;
+                       else {
+                               fprintf(stderr, "[%s] unrecognized prior type.\n", __func__);
+                               return 1;
+                       }
+                       break;
                case 'f':
                        mplp.fai = fai_load(optarg);
                        if (mplp.fai == 0) return 1;
@@ -676,7 +683,6 @@ int bam_mpileup(int argc, char *argv[])
                case 'l': mplp.fn_pos = strdup(optarg); break;
                case 'V':
                case 'c': mplp.flag |= MPLP_VCF; break;
-               case 'S': mplp.flag |= MPLP_AFS; break;
                case 'F': mplp.flag |= MPLP_AFALL; break;
                case 'v': mplp.flag |= MPLP_VAR; break;
                case 'M': mplp.max_mq = atoi(optarg); break;
@@ -692,10 +698,9 @@ int bam_mpileup(int argc, char *argv[])
                fprintf(stderr, "         -M INT      cap mapping quality at INT [%d]\n", mplp.max_mq);
                fprintf(stderr, "         -q INT      filter out alignment with MQ smaller than INT [%d]\n", mplp.min_mq);
                fprintf(stderr, "         -t FLOAT    scaled mutation rate [%lg]\n", mplp.theta);
+               fprintf(stderr, "         -P STR      prior: full, flat, cond2 [full]\n");
                fprintf(stderr, "         -c          generate VCF output (consensus calling)\n");
                fprintf(stderr, "         -v          show variant sites only\n");
-               fprintf(stderr, "         -S          calculate AFS (slow, to stderr)\n");
-               fprintf(stderr, "         -2          conditional prior\n");
                fprintf(stderr, "\n");
                fprintf(stderr, "Notes: Assuming error independency and diploid individuals.\n\n");
                return 1;