From 390bbfbca2df79acf521706409a258895dfe72b1 Mon Sep 17 00:00:00 2001 From: Heng Li Date: Tue, 31 Aug 2010 20:44:57 +0000 Subject: [PATCH] preliminary contrast variant caller --- bcftools/Makefile | 2 +- bcftools/call1.c | 6 ++- bcftools/prob1.c | 95 +++++++++++++++++++++++++++++++++++++++-------- bcftools/prob1.h | 4 +- 4 files changed, 88 insertions(+), 19 deletions(-) diff --git a/bcftools/Makefile b/bcftools/Makefile index 5a3dc76..e3b92e9 100644 --- a/bcftools/Makefile +++ b/bcftools/Makefile @@ -1,6 +1,6 @@ CC= gcc CFLAGS= -g -Wall -O2 #-m64 #-arch ppc -DFLAGS= -D_FILE_OFFSET_BITS=64 -D_USE_KNETFILE +DFLAGS= -D_FILE_OFFSET_BITS=64 -D_USE_KNETFILE -D_BCF_QUAD LOBJS= bcf.o vcf.o bcfutils.o prob1.o kfunc.o index.o fet.o OMISC= .. AOBJS= call1.o main.o $(OMISC)/kstring.o $(OMISC)/bgzf.o $(OMISC)/knetfile.o diff --git a/bcftools/call1.c b/bcftools/call1.c index c722b35..e8ecb51 100644 --- a/bcftools/call1.c +++ b/bcftools/call1.c @@ -162,7 +162,10 @@ static int update_bcf1(int n_smpl, bcf1_t *b, const bcf_p1aux_t *pa, const bcf_p if (b->info[0]) kputc(';', &s); ksprintf(&s, "AF1=%.3lf;AFE=%.3lf", 1.-pr->f_em, 1.-pr->f_exp); ksprintf(&s, ";DP4=%d,%d,%d,%d;MQ=%d", a.d[0], a.d[1], a.d[2], a.d[3], a.mq); - if (a.is_tested) ksprintf(&s, ";PV4=%.2lg,%.2lg,%.2lg,%.2lg", a.p[0], a.p[1], a.p[2], a.p[3]); + if (a.is_tested) { + if (pr->pc[0] >= 0.) ksprintf(&s, ";PC4=%.2lg,%.2lg,%.2lg,%.2lg", pr->pc[0], pr->pc[1], pr->pc[2], pr->pc[3]); + ksprintf(&s, ";PV4=%.2lg,%.2lg,%.2lg,%.2lg", a.p[0], a.p[1], a.p[2], a.p[3]); + } if (pr->g[0] >= 0. && p_hwe <= .2) ksprintf(&s, ";GC=%.2lf,%.2lf,%.2lf;HWE=%.3lf", pr->g[2], pr->g[1], pr->g[0], p_hwe); kputc('\0', &s); @@ -224,6 +227,7 @@ int bcfview(int argc, char *argv[]) fprintf(stderr, " -L discard the PL genotype field\n"); fprintf(stderr, " -H perform Hardy-Weinberg test (slower)\n"); fprintf(stderr, " -v output potential variant sites only\n"); + fprintf(stderr, " -1 INT number of group-1 samples [0]\n"); fprintf(stderr, " -l FILE list of sites to output [all sites]\n"); fprintf(stderr, " -t FLOAT scaled mutation rate [%.4lg]\n", vc.theta); fprintf(stderr, " -p FLOAT variant if P(ref|D)= b->n) return -1; b->n1 = n1; +#ifdef _BCF_QUAD + { + int k1, k2, n2 = b->n - b->n1; + b->k1k2 = calloc((2*n1+1) * (2*n2+1), sizeof(double)); + for (k1 = 0; k1 <= 2*n1; ++k1) + for (k2 = 0; k2 <= 2*n2; ++k2) + b->k1k2[k1*(2*n2+1)+k2] = exp(lbinom(2*n1,k1) + lbinom(2*n2,k2) - lbinom(b->M,k1+k2)); + } +#endif return 0; } @@ -133,6 +152,7 @@ void bcf_p1_destroy(bcf_p1aux_t *ma) free(ma->phi); free(ma->z); free(ma->zswap); free(ma->z1); free(ma->z2); free(ma->afs); free(ma->afs1); + free(ma->k1k2); free(ma); } } @@ -247,7 +267,7 @@ static void mc_cal_y(bcf_p1aux_t *ma) { if (ma->n1 > 0 && ma->n1 < ma->n) { int k; - double x; + long double x; memset(ma->z1, 0, sizeof(double) * (2 * ma->n1 + 1)); memset(ma->z2, 0, sizeof(double) * (2 * (ma->n - ma->n1) + 1)); ma->t1 = ma->t2 = 0.; @@ -256,26 +276,68 @@ static void mc_cal_y(bcf_p1aux_t *ma) memcpy(ma->z2, ma->z, sizeof(double) * (2 * (ma->n - ma->n1) + 1)); mc_cal_y_core(ma, 0); // rescale z - x = exp(ma->t - (ma->t1 + ma->t2)); + x = expl(ma->t - (ma->t1 + ma->t2)); for (k = 0; k <= ma->M; ++k) ma->z[k] *= x; } else mc_cal_y_core(ma, 0); +#ifdef _BCF_QUAD /* - if (ma->n1 > 0 && ma->n1 < ma->n) { - int i; - double y[5]; - printf("*****\n"); - for (i = 0; i <= 2; ++i) - printf("(%lf,%lf) ", ma->z1[i], ma->z2[i]); - printf("\n"); - y[0] = ma->z1[0] * ma->z2[0]; - y[1] = 1./2. * (ma->z1[0] * ma->z2[1] + ma->z1[1] * ma->z2[0]); - y[2] = 1./6. * (ma->z1[0] * ma->z2[2] + ma->z1[2] * ma->z2[0]) + 4./6. * ma->z1[1] * ma->z2[1]; - y[3] = 1./2. * (ma->z1[1] * ma->z2[2] + ma->z1[2] * ma->z2[1]); - y[4] = ma->z1[2] * ma->z2[2]; - for (i = 0; i <= 4; ++i) printf("(%lf,%lf) ", ma->z[i], y[i]); + if (ma->n1 > 0 && ma->n1 < ma->n) { // DEBUG: consistency check; z[i] should equal y[i] + int i, k1, k2, n1 = ma->n1, n2 = ma->n - n1; + double *y; + printf("*** "); + y = calloc(ma->M + 1, sizeof(double)); + for (k1 = 0; k1 <= 2*n1; ++k1) + for (k2 = 0; k2 <= 2*n2; ++k2) + y[k1+k2] += ma->k1k2[k1*(2*n2+1)+k2] * ma->z1[k1] * ma->z2[k2]; + for (i = 0; i <= ma->M; ++i) printf("(%lf,%lf) ", ma->z[i], y[i]); printf("\n"); + free(y); } */ +#endif +} + +static void contrast(bcf_p1aux_t *ma, double pc[4]) // mc_cal_y() must be called before hand +{ + int k, n1 = ma->n1, n2 = ma->n - ma->n1; + long double sum = -1., x, sum_alt; + double y; + pc[0] = pc[1] = pc[2] = pc[3] = -1.; + if (n1 <= 0 || n2 <= 0) return; +#ifdef _BCF_QUAD + { // FIXME: can be improved by skipping zero cells + int k1, k2; + long double z[3]; + z[0] = z[1] = z[2] = 0.; + for (k1 = 0; k1 <= 2*n1; ++k1) + for (k2 = 0; k2 <= 2*n2; ++k2) { + double zz = ma->phi[k1+k2] * ma->z1[k1] * ma->z2[k2] * ma->k1k2[k1*(2*n2+1)+k2]; + if ((double)k1/n1 < (double)k2/n2) z[0] += zz; + else if ((double)k1/n1 > (double)k2/n2) z[1] += zz; + else z[2] += zz; + } + sum = z[0] + z[1] + z[2]; + pc[2] = z[0] / sum; pc[3] = z[1] / sum; + } +#endif + for (k = 0, sum_alt = 0.; k <= ma->M; ++k) + sum_alt += (long double)ma->phi[k] * ma->z[k]; +// printf("* %lg, %lg *\n", (double)sum, (double)sum_alt); // DEBUG: sum should equal sum_alt + sum = sum_alt; + y = lgamma(2*n1 + 1) - lgamma(ma->M + 1); + for (k = 1, x = 0.; k <= 2 * n2; ++k) + x += ma->phi[k] * ma->z2[k] * exp(lgamma(ma->M - k + 1) - lgamma(2*n2 - k + 1) + y); + pc[0] = ma->z1[0] * x / sum; + y = lgamma(2*n2 + 1) - lgamma(ma->M + 1); + for (k = 1, x = 0.; k <= 2 * n1; ++k) + x += ma->phi[k] * ma->z1[k] * exp(lgamma(ma->M - k + 1) - lgamma(2*n1 - k + 1) + y); + pc[1] = ma->z2[0] * x / sum; + for (k = 0; k < 4; ++k) { + y = 1. - pc[k]; + if (y <= 0.) y = 1e-100; + pc[k] = (int)(-3.434 * log(y) + .499); + if (pc[k] > 99.) pc[k] = 99.; + } } static double mc_cal_afs(bcf_p1aux_t *ma) @@ -358,6 +420,7 @@ int bcf_p1_cal(bcf1_t *b, bcf_p1aux_t *ma, bcf_p1rst_t *rst) } } rst->g[0] = rst->g[1] = rst->g[2] = -1.; + contrast(ma, rst->pc); // bcf_p1_cal_g3(ma, rst->g); return 0; } diff --git a/bcftools/prob1.h b/bcftools/prob1.h index 6bb7440..35c794f 100644 --- a/bcftools/prob1.h +++ b/bcftools/prob1.h @@ -8,7 +8,9 @@ typedef struct __bcf_p1aux_t bcf_p1aux_t; typedef struct { int rank0; - double f_em, f_exp, f_flat, p_ref, g[3]; + double f_em, f_exp, f_flat, p_ref; + double pc[4]; + double g[3]; } bcf_p1rst_t; #define MC_PTYPE_FULL 1 -- 2.39.2