From: Heng Li Date: Fri, 23 Jul 2010 17:25:44 +0000 (+0000) Subject: added an alternative prior X-Git-Url: https://git.donarmstrong.com/?p=samtools.git;a=commitdiff_plain;h=71355852e510ac7937d552117557b9dbf07cd557 added an alternative prior --- diff --git a/bam_mcns.c b/bam_mcns.c index c1cdb90..7859169 100644 --- a/bam_mcns.c +++ b/bam_mcns.c @@ -9,20 +9,23 @@ struct __mc_aux_t { int n, N; int ref, alt; - double theta; double *q2p, *pdg; // pdg -> P(D|g) double *alpha, *beta; int *qsum, *bcnt; }; -void mc_init_prior(mc_aux_t *ma, double theta) +void mc_init_prior(mc_aux_t *ma, int type, double theta) { - double sum; int i; - ma->theta = theta; - for (i = 0, sum = 0.; i < 2 * ma->n; ++i) - sum += (ma->alpha[i] = ma->theta / (2 * ma->n - i)); - ma->alpha[2 * ma->n] = 1. - sum; + if (type == MC_PTYPE_COND2) { + for (i = 0; i <= 2 * ma->n; ++i) + ma->alpha[i] = 2. * (i + 1) / (2 * ma->n + 1) / (2 * ma->n + 2); + } else { + double sum; + for (i = 0, sum = 0.; i < 2 * ma->n; ++i) + sum += (ma->alpha[i] = theta / (2 * ma->n - i)); + ma->alpha[2 * ma->n] = 1. - sum; + } } mc_aux_t *mc_init(int n) // FIXME: assuming diploid @@ -46,7 +49,7 @@ mc_aux_t *mc_init(int n) // FIXME: assuming diploid bi[1] = 2 * f * (1. - f); bi[2] = f * f; } - mc_init_prior(ma, 1e-3); // the simplest prior + mc_init_prior(ma, MC_PTYPE_FULL, 1e-3); // the simplest prior return ma; } diff --git a/bam_mcns.h b/bam_mcns.h index 89879fb..cff6e57 100644 --- a/bam_mcns.h +++ b/bam_mcns.h @@ -6,11 +6,15 @@ struct __mc_aux_t; typedef struct __mc_aux_t mc_aux_t; +#define MC_PTYPE_FULL 1 +#define MC_PTYPE_COND2 2 + #ifdef __cplusplus extern "C" { #endif 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); diff --git a/bam_plcmd.c b/bam_plcmd.c index 258b83f..82fcafc 100644 --- a/bam_plcmd.c +++ b/bam_plcmd.c @@ -450,7 +450,8 @@ int bam_pileup(int argc, char *argv[]) ***********/ typedef struct { - int vcf, max_mq, min_mq, var_only; + int vcf, max_mq, min_mq, var_only, prior_type; + double theta; char *reg, *fn_pos; faidx_t *fai; kh_64_t *hash; @@ -537,7 +538,10 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn) free(s.s); } // mpileup - if (conf->vcf) ma = mc_init(n); + if (conf->vcf) { + ma = mc_init(n); + mc_init_prior(ma, conf->prior_type, conf->theta); + } 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) { @@ -650,8 +654,12 @@ 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.prior_type = MC_PTYPE_FULL; + mplp.theta = 1e-3; + while ((c = getopt(argc, argv, "f:r:l:VvM:q:t:2")) >= 0) { switch (c) { + case 't': mplp.theta = atof(optarg); break; + case '2': mplp.prior_type = MC_PTYPE_COND2; break; case 'f': mplp.fai = fai_load(optarg); if (mplp.fai == 0) return 1; @@ -672,8 +680,10 @@ int bam_mpileup(int argc, char *argv[]) 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 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, " -V generate VCF output (SNP calling)\n"); fprintf(stderr, " -v show variant sites only\n"); + fprintf(stderr, " -2 conditional prior\n"); fprintf(stderr, "\n"); fprintf(stderr, "Notes: Assuming error independency and diploid individuals.\n\n"); return 1; diff --git a/bamtk.c b/bamtk.c index b51df96..1de3920 100644 --- a/bamtk.c +++ b/bamtk.c @@ -9,7 +9,7 @@ #endif #ifndef PACKAGE_VERSION -#define PACKAGE_VERSION "0.1.8-3 (r627)" +#define PACKAGE_VERSION "0.1.8-4 (r631)" #endif int bam_taf2baf(int argc, char *argv[]);