]> git.donarmstrong.com Git - samtools.git/commitdiff
posterior expectation FINALLY working. I am so tired...
authorHeng Li <lh3@live.co.uk>
Sat, 24 Jul 2010 05:49:49 +0000 (05:49 +0000)
committerHeng Li <lh3@live.co.uk>
Sat, 24 Jul 2010 05:49:49 +0000 (05:49 +0000)
bam_mcns.c
bam_mcns.h
bam_plcmd.c

index 7859169c9ef8799b37878a6d1e19e7c1c5b8fcb0..d98d48030560a42fd65c00bda4f1633bc90ac169 100644 (file)
@@ -1,16 +1,22 @@
 #include <math.h>
 #include <stdlib.h>
+#include <stdio.h>
 #include "bam_mcns.h"
 
 #define MC_MIN_QUAL 13
+#define MC_AVG_ERR 0.007
 #define MC_MAX_SUMQ 3000
 #define MC_MAX_SUMQP 1e-300
+#define MC_MAX_EM_ITER 16
+#define MC_EM_EPS 1e-4
 
 struct __mc_aux_t {
-       int n, N;
+       int n, M;
        int ref, alt;
        double *q2p, *pdg; // pdg -> P(D|g)
        double *alpha, *beta;
+       double *z, *zswap; // aux for afs
+       double *afs, *afs1; // afs: accumulative AFS; afs1: site posterior distribution
        int *qsum, *bcnt;
 };
 
@@ -33,18 +39,22 @@ mc_aux_t *mc_init(int n) // FIXME: assuming diploid
        mc_aux_t *ma;
        int i;
        ma = calloc(1, sizeof(mc_aux_t));
-       ma->n = n; ma->N = 2 * n;
+       ma->n = n; ma->M = 2 * n;
        ma->q2p = calloc(MC_MAX_SUMQ + 1, sizeof(double));
        ma->qsum = calloc(4 * ma->n, sizeof(int));
        ma->bcnt = calloc(4 * ma->n, sizeof(int));
        ma->pdg = calloc(3 * ma->n, sizeof(double));
        ma->alpha = calloc(2 * ma->n + 1, sizeof(double));
        ma->beta = calloc((2 * ma->n + 1) * 3, sizeof(double));
+       ma->z = calloc(2 * ma->n + 1, sizeof(double));
+       ma->zswap = calloc(2 * ma->n + 1, sizeof(double));
+       ma->afs = calloc(2 * ma->n + 1, sizeof(double));
+       ma->afs1 = calloc(2 * ma->n + 1, sizeof(double));
        for (i = 0; i <= MC_MAX_SUMQ; ++i)
                ma->q2p[i] = pow(10., -i / 10.);
-       for (i = 0; i <= 2 * ma->n; ++i) { // beta[k][g]=P(g|k/M)
+       for (i = 0; i <= ma->M; ++i) { // beta[k][g]=P(g|k/M)
                double *bi = ma->beta + 3 * i;
-               double f = (double)i / (2 * ma->n);
+               double f = (double)i / ma->M;
                bi[0] = (1. - f) * (1. - f);
                bi[1] = 2 * f * (1. - f);
                bi[2] = f * f;
@@ -59,19 +69,21 @@ void mc_destroy(mc_aux_t *ma)
                free(ma->qsum); free(ma->bcnt);
                free(ma->q2p); free(ma->pdg);
                free(ma->alpha); free(ma->beta);
+               free(ma->z); free(ma->zswap);
+               free(ma->afs); free(ma->afs1);
                free(ma);
        }
 }
 
-static void sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
+static int sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
 {
-       int i, j;
+       int i, j, tot = 0;
        memset(ma->qsum, 0, sizeof(int) * 4 * ma->n);
        memset(ma->bcnt, 0, sizeof(int) * 4 * ma->n);
        for (j = 0; j < ma->n; ++j) {
-               int tmp, *qsum = ma->qsum + j * 4;
+               int *qsum = ma->qsum + j * 4;
                int *bcnt = ma->bcnt + j * 4;
-               for (i = tmp = 0; i < n[j]; ++i) {
+               for (i = 0; i < n[j]; ++i) {
                        const bam_pileup1_t *p = plp[j] + i;
                        int q, b;
                        if (p->is_del || (p->b->core.flag&BAM_FUNMAP)) continue;
@@ -82,9 +94,10 @@ static void sum_err(int *n, const bam_pileup1_t **plp, mc_aux_t *ma)
                        if (b > 3) continue; // N
                        qsum[b] += q;
                        ++bcnt[b];
-                       ++tmp;
+                       ++tot;
                }
        }
+       return tot;
 }
 
 static void set_allele(int ref, mc_aux_t *ma)
@@ -118,84 +131,64 @@ static void cal_pdg(mc_aux_t *ma)
                        pdg[i] = pi[i] > MC_MAX_SUMQ? MC_MAX_SUMQP : ma->q2p[pi[i]];
        }
 }
-// return the reference allele frequency
-double mc_freq0(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, int *_ref, int *alt)
+// this calculates the naive allele frequency and Nielsen's frequency
+static double mc_freq0(const mc_aux_t *ma, double *_f)
 {
        int i, cnt;
-       double f;
-       sum_err(n, plp, ma);
-       set_allele(ref, ma);
-       cal_pdg(ma);
-       *_ref = ma->ref; *alt = ma->alt;
-       for (i = cnt = 0, f = 0.; i < ma->n; ++i) {
+       double f, f_nielsen, w_sum;
+       *_f = -1.;
+       for (i = cnt = 0, f = f_nielsen = w_sum = 0.; i < ma->n; ++i) {
                int *bcnt = ma->bcnt + i * 4;
                int x = bcnt[ma->ref] + bcnt[ma->alt];
                if (x) {
+                       double w, p;
                        ++cnt;
                        f += (double)bcnt[ma->ref] / x;
+                       p = (bcnt[ma->ref] - MC_AVG_ERR * x) / (1. - 2. * MC_AVG_ERR) / x;
+                       w = 2. * x / (1. + x);
+                       w_sum += w;
+                       f_nielsen += p * w;
                }
        }
-       return f / cnt;
+       if (cnt) {
+               f_nielsen /= w_sum;
+               if (f_nielsen < 0.) f_nielsen = 0.;
+               if (f_nielsen > 1.) f_nielsen = 1.;
+               *_f = f_nielsen;
+               return f / cnt;
+       } else return -1.;
 }
 // f0 is the reference allele frequency
-double mc_freq_iter(double f0, mc_aux_t *ma)
+static double mc_freq_iter(double f0, const mc_aux_t *ma)
 {
        double f, f3[3];
-       int i, j;
+       int i;
        f3[0] = (1.-f0)*(1.-f0); f3[1] = 2.*f0*(1.-f0); f3[2] = f0*f0;
        for (i = 0, f = 0.; i < ma->n; ++i) {
-               double up, dn, *pdg;
+               double *pdg;
                pdg = ma->pdg + i * 3;
-               for (j = 1, up = 0.; j < 3; ++j)
-                       up += j * pdg[j] * f3[j];
-               for (j = 0, dn = 0.; j < 3; ++j)
-                       dn += pdg[j] * f3[j];
-               f += up / dn;
+               f += (pdg[1] * f3[1] + 2. * pdg[2] * f3[2])
+                       / (pdg[0] * f3[0] + pdg[1] * f3[1] + pdg[2] * f3[2]);
        }
        f /= ma->n * 2.;
        return f;
 }
 
-double mc_freq_post(mc_aux_t *ma)
-{
-       int i, k;
-       long double f = 0.;
-       for (i = 0; i < ma->n; ++i) {
-               double *pdg = ma->pdg + i * 3;
-               long double y = 0., z = 0.;
-               for (k = 0; k <= ma->n * 2; ++k) {
-                       double *bk = ma->beta + k * 3;
-/*
-                       int g;
-                       double yk = 0., zk = 0.;
-                       for (g = 0; g < 3; ++g) {
-                               yk += g * pdg[g] * bk[g];
-                               zk += pdg[g] * bk[g];
-                       }
-                       y += yk * ma->alpha[k];
-                       z += zk * ma->alpha[k];
-*/
-                       y += (pdg[1] * bk[1] + 2. * pdg[2] * bk[2]) * ma->alpha[k];
-                       z += (pdg[0] * bk[0] + pdg[1] * bk[1] + pdg[2] * bk[2]) * ma->alpha[k];
-               }
-               f += y / z;
-       }
-       return f / ma->n / 2;
-}
-
-double mc_ref_prob(mc_aux_t *ma)
+static double mc_ref_prob(const mc_aux_t *ma, double *_PD, double *f_exp)
 {
        int k, i;
-       long double PD = 0., Pref = 0.;
-       for (k = 0; k <= ma->n * 2; ++k) {
-               long double x = 1.;
+       long double PD = 0., Pref = 0., Ef = 0.;
+       for (k = 0; k <= ma->M; ++k) {
+               long double x = 1., y = 0.;
                double *bk = ma->beta + k * 3;
                for (i = 0; i < ma->n; ++i) {
                        double *pdg = ma->pdg + i * 3;
-//                     int g; double y=0.; for (g = 0; g < 3; ++g) y += pdg[g] * bk[g];
-                       x *= pdg[0] * bk[0] + pdg[1] * bk[1] + pdg[2] * bk[2];
+                       double z = pdg[0] * bk[0] + pdg[1] * bk[1] + pdg[2] * bk[2];
+                       x *= z;
+                       y += (pdg[1] * bk[1] + 2. * pdg[2] * bk[2]) / z;
                }
                PD += x * ma->alpha[k];
+               Ef += x * y * ma->alpha[k];
        }
        for (k = 0; k <= ma->n * 2; ++k) {
                long double x = 1.0;
@@ -203,6 +196,8 @@ double mc_ref_prob(mc_aux_t *ma)
                        x *= ma->pdg[i * 3 + 2] * ma->beta[k * 3 + 2];
                Pref += x * ma->alpha[k];
        }
+       *f_exp = (double)(Ef / PD / ma->M);
+       *_PD = PD;
        return Pref / PD;
 }
 
@@ -218,10 +213,87 @@ int mc_call_gt(const mc_aux_t *ma, double f0, int k)
                g[i] /= sum;
                if (g[i] > max) max = g[i], max_i = i;
        }
-//     printf("***%lg,%lg,%lg,%lg,%lg,%lg\n", g[0], g[1], g[2], pdg[0], pdg[1], pdg[2]);
        max = 1. - max;
        if (max < 1e-308) max = 1e-308;
        q = (int)(-3.434 * log(max) + .499);
        if (q > 99) q = 99;
        return q<<2|max_i;
 }
+// calculate z_{nr}^{(k)}
+static void mc_cal_z(mc_aux_t *ma, int k)
+{
+       double *z[2], *tmp, *bk, *pdg;
+       int i, j;
+       z[0] = ma->z;
+       z[1] = ma->zswap;
+       bk = ma->beta + k * 3; pdg = ma->pdg;
+       z[0][0] = 1.; z[0][1] = z[0][2] = 0.;
+       for (j = 0; j < ma->n; ++j) {
+               int max = (j + 1) * 2;
+               double p[3];
+               pdg = ma->pdg + j * 3;
+               p[0] = bk[0] * pdg[0]; p[1] = bk[1] * pdg[1]; p[2] = bk[2] * pdg[2];
+               z[1][0] = p[0] * z[0][0];
+               z[1][1] = p[0] * z[0][1] + p[1] * z[0][0];
+               for (i = 2; i <= max; ++i)
+                       z[1][i] = p[0] * z[0][i] + p[1] * z[0][i-1] + p[2] * z[0][i-2];
+               if (j < ma->n - 1) z[1][max+1] = z[1][max+2] = 0.;
+               tmp = z[0]; z[0] = z[1]; z[1] = tmp;
+       }
+       if (z[0] != ma->z) memcpy(ma->z, z[0], sizeof(double) * (2 * ma->n + 1));
+}
+// Warning: this is cubic in ma->n, very sloooooow
+static void mc_add_afs(mc_aux_t *ma, double PD, double *f_map, double *p_map)
+{
+       int k, l;
+       double sum = 0.;
+       memset(ma->afs1, 0, sizeof(double) * (2 * ma->n + 1));
+       for (k = 0; k <= 2 * ma->n; ++k) {
+               mc_cal_z(ma, k);
+               for (l = 0; l <= 2 * ma->n; ++l)
+                       ma->afs1[l] += ma->alpha[k] * ma->z[l] / PD;
+       }
+       for (k = 0; k <= 2 * ma->n; ++k) {
+               ma->afs[k] += ma->afs1[k];
+               sum += ma->afs1[k];
+       }
+       {
+               int max_k = 0;
+               double max = -1., e = 0.;
+               for (k = 0; k <= 2 * ma->n; ++k) {
+                       if (ma->afs1[k] > max) max = ma->afs1[k], max_k = k;
+                       e += k * ma->afs1[k];
+               }
+               *f_map = .5 * max_k / ma->n; *p_map = max;
+               printf(" * %.3lg:%.3lg:%.3lg:%.3lg * ", sum, 1.-.5*max_k/ma->n, max, 1.-.5*e/ma->n);
+       }
+}
+
+int mc_cal(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, mc_rst_t *rst, int level)
+{
+       int i, tot;
+       memset(rst, 0, sizeof(mc_rst_t));
+       rst->f_em = rst->f_exp = -1.; rst->ref = rst->alt = -1;
+       // precalculation
+       tot = sum_err(n, plp, ma);
+       if (tot == 0) return 0; // no good bases
+       set_allele(ref, ma);
+       cal_pdg(ma);
+       // set ref/major allele
+       rst->ref = ma->ref; rst->alt = ma->alt;
+       // calculate naive and Nielsen's freq
+       rst->f_naive = mc_freq0(ma, &rst->f_nielsen);
+       { // calculate f_em
+               double flast = rst->f_naive;
+               for (i = 0; i < MC_MAX_EM_ITER; ++i) {
+                       rst->f_em = mc_freq_iter(flast, ma);
+                       if (fabs(rst->f_em - flast) < MC_EM_EPS) break;
+                       flast = rst->f_em;
+               }
+       }
+       if (level >= 2) // quadratic-time calculations; necessary for genotyping
+               rst->p_ref = mc_ref_prob(ma, &rst->PD, &rst->f_exp);
+       if (level >= 3)
+               mc_add_afs(ma, rst->PD, &rst->f_map, &rst->p_map);
+       return tot;
+}
index cff6e57d28047950be347e7d873d940638baa385..a86fc1ddbd62e575ffeca6d2d59917e0b7e2f9a6 100644 (file)
@@ -6,6 +6,16 @@
 struct __mc_aux_t;
 typedef struct __mc_aux_t mc_aux_t;
 
+typedef struct {
+       // O(n)
+       int ref, alt;
+       double f_em, f_naive, f_nielsen;
+       // O(n^2)
+       double PD, p_ref, f_exp;
+       // O(n^3)
+       double f_map, p_map; // map=maximum a posterior
+} mc_rst_t;
+
 #define MC_PTYPE_FULL  1
 #define MC_PTYPE_COND2 2
 
@@ -16,10 +26,7 @@ extern "C" {
        mc_aux_t *mc_init(int n);
        void mc_init_prior(mc_aux_t *ma, int type, double theta);
        void mc_destroy(mc_aux_t *ma);
-       double mc_freq0(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, int *_ref, int *alt);
-       double mc_freq_iter(double f0, mc_aux_t *ma);
-       double mc_freq_post(mc_aux_t *ma);
-       double mc_ref_prob(mc_aux_t *ma);
+       int mc_cal(int ref, int *n, const bam_pileup1_t **plp, mc_aux_t *ma, mc_rst_t *rst, int level);
        int mc_call_gt(const mc_aux_t *ma, double f0, int k);
 
 #ifdef __cplusplus
index a5bcdf560a5bfdfd14a9a2eccffb46be1164e8cc..ee49cd3cff69afb2c0872978293fd3eaf0e13d9c 100644 (file)
@@ -450,7 +450,7 @@ int bam_pileup(int argc, char *argv[])
  ***********/
 
 typedef struct {
-       int vcf, max_mq, min_mq, var_only, prior_type;
+       int vcf, max_mq, min_mq, var_only, prior_type, afs;
        double theta;
        char *reg, *fn_pos;
        faidx_t *fai;
@@ -523,8 +523,8 @@ 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=\"Posterior non-reference allele frequency\">");
+//             puts("##INFO=<ID=AFEM,Number=1,Type=Float,Description=\"Prior-free non-reference allele frequency by EM\">");
                puts("##FILTER=<ID=Q13,Description=\"All min{baseQ,mapQ} below 13\">");
                kputs("#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO\tFORMAT", &s);
                for (i = 0; i < n; ++i) {
@@ -557,39 +557,29 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        ref_tid = tid;
                }
                if (conf->vcf) {
-                       double f0, f, fpost, pref = -1.0; // the reference allele frequency
-                       int j, _ref, _alt, _ref0, depth, rms_q, _ref0b, is_var = 0, qref = 0, filtered = 0;
+                       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->afs) level = 3;
                        _ref0 = _ref0b = (ref && pos < ref_len)? ref[pos] : 'N';
                        _ref0 = bam_nt16_nt4_table[bam_nt16_table[_ref0]];
-                       f = f0 = fpost = 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;
-                               }
-                               pref = mc_ref_prob(ma);
-                               fpost = mc_freq_post(ma);
-                               if (1. - f < 1e-4) f = 1.;
-                               if (1. - fpost < 1e-4) fpost = 1.;
-                               is_var = (pref < .5);
-                               q = is_var? pref : 1. - pref;
+                       tot = mc_cal(_ref0, n_plp, plp, ma, &r, level);
+                       if (tot) { // has good bases
+                               double q;
+                               is_var = (r.p_ref < .5);
+                               q = is_var? r.p_ref : 1. - r.p_ref;
                                if (q < 1e-308) q = 1e-308;
                                qref = (int)(-3.434 * log(q) + .499);
                                if (qref > 99) qref = 99;
                        }
-                       filtered = (f >= 0. && f <= 1.)? 0 : 1;
                        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]);
+                               if (_ref0 == r.ref) putchar("ACGTN"[r.alt]);
+                               else printf("%c,%c", "ACGTN"[r.ref], "ACGTN"[r.alt]);
                        } else putchar('.');
                        printf("\t%d\t", qref);
-                       if (filtered) printf("Q13\t");
+                       if (!tot) printf("Q13\t");
                        else printf(".\t");
                        for (i = depth = 0, sqr_sum = 0; i < n; ++i) {
                                depth += n_plp[i];
@@ -601,12 +591,13 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        }
                        rms_q = (int)(sqrt((double)sqr_sum / depth) + .499);
                        printf("DP=%d;MQ=%d", depth, rms_q);
-                       if (fpost >= 0. && fpost <= 1.) printf(";AF=%.3lg", 1. - fpost);
-                       if (f >= 0. && f <= 1.) printf(";AFEM=%.3lg", 1. - f);
+                       printf(";AF=%.3lg", 1. - r.f_em);
+                       if (conf->afs)
+                               printf(";AF0=%.3lg;AFN=%.3lg;AFE=%.3lg", 1-r.f_naive, 1-r.f_nielsen, 1-r.f_exp);
                        printf("\tGT:GQ:DP");
-                       if (fpost >= 0. && fpost <= 1.) {
+                       if (tot) {
                                for (i = 0; i < n; ++i) {
-                                       int x = mc_call_gt(ma, fpost, i);
+                                       int x = mc_call_gt(ma, r.f_em, 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");
@@ -658,7 +649,7 @@ 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, "f:r:l:VvM:q:t:2")) >= 0) {
+       while ((c = getopt(argc, argv, "f:r:l:VvM:q:t:2F")) >= 0) {
                switch (c) {
                case 't': mplp.theta = atof(optarg); break;
                case '2': mplp.prior_type = MC_PTYPE_COND2; break;
@@ -672,6 +663,7 @@ int bam_mpileup(int argc, char *argv[])
                case 'V': mplp.vcf = 1; break;
                case 'M': mplp.max_mq = atoi(optarg); break;
                case 'q': mplp.min_mq = atoi(optarg); break;
+               case 'F': mplp.afs = 1; break;
                }
        }
        if (argc == 1) {