]> git.donarmstrong.com Git - samtools.git/commitdiff
* fixed a minor issue in printing VCFs
authorHeng Li <lh3@live.co.uk>
Fri, 10 Dec 2010 15:50:29 +0000 (15:50 +0000)
committerHeng Li <lh3@live.co.uk>
Fri, 10 Dec 2010 15:50:29 +0000 (15:50 +0000)
 * write bcftools specific INFO and FORMAT in the header

bcftools/call1.c
bcftools/vcf.c

index b913c9a4990e1f59bbdde8667d6496f842c02d21..adf7a52cf69ad371eb6ec2f76c2d76ff45565275 100644 (file)
@@ -200,6 +200,34 @@ static int update_bcf1(int n_smpl, bcf1_t *b, const bcf_p1aux_t *pa, const bcf_p
        return is_var;
 }
 
+static void write_header(bcf_hdr_t *h)
+{
+       kstring_t str;
+       str.l = h->l_txt? h->l_txt - 1 : 0;
+       str.m = str.l + 1; str.s = h->txt;
+       if (!strstr(str.s, "##INFO=<ID=DP,"))
+               kputs("##INFO=<ID=DP,Number=1,Type=Integer,Description=\"Raw read depth\">\n", &str);
+       if (!strstr(str.s, "##INFO=<ID=DP4,"))
+               kputs("##INFO=<ID=DP4,Number=4,Type=Integer,Description=\"# high-quality ref-forward bases, ref-reverse, alt-forward and alt-reverse bases\">\n", &str);
+       if (!strstr(str.s, "##INFO=<ID=MQ,"))
+               kputs("##INFO=<ID=MQ,Number=1,Type=Integer,Description=\"Root-mean-square mapping quality of covering reads\">\n", &str);
+       if (!strstr(str.s, "##INFO=<ID=FQ,"))
+               kputs("##INFO=<ID=FQ,Number=1,Type=Integer,Description=\"Phred probability that sample chromosomes are not all the same\">\n", &str);
+       if (!strstr(str.s, "##INFO=<ID=AF1,"))
+               kputs("##INFO=<ID=AF1,Number=1,Type=Float,Description=\"Max-likelihood estimate of the site allele frequency of the first ALT allele\">\n", &str);
+       if (!strstr(str.s, "##INFO=<ID=CI95,"))
+               kputs("##INFO=<ID=CI95,Number=2,Type=Float,Description=\"Equal-tail Bayesian credible interval of the site allele frequency at the 95% level\">\n", &str);
+       if (!strstr(str.s, "##INFO=<ID=PV4,"))
+               kputs("##INFO=<ID=PV4,Number=4,Type=Float,Description=\"P-values for strand bias, baseQ bias, mapQ bias and tail distance bias\">\n", &str);
+       if (!strstr(str.s, "##FORMAT=<ID=DP,"))
+               kputs("##FORMAT=<ID=DP,Number=1,Type=Integer,Description=\"# high-quality bases\">\n", &str);
+       if (!strstr(str.s, "##FORMAT=<ID=SP,"))
+               kputs("##FORMAT=<ID=SP,Number=1,Type=Integer,Description=\"Phred-scaled strand bias P-value\">\n", &str);
+       if (!strstr(str.s, "##FORMAT=<ID=PL,"))
+               kputs("##FORMAT=<ID=PL,Number=(#ALT+1)*(#ALT+2)/2,Type=Integer,Description=\"List of Phred-scaled genotype likelihoods\">\n", &str);
+       h->l_txt = str.l + 1; h->txt = str.s;
+}
+
 double bcf_ld_freq(const bcf1_t *b0, const bcf1_t *b1, double f[4]);
 
 int bcfview(int argc, char *argv[])
@@ -285,7 +313,10 @@ int bcfview(int argc, char *argv[])
        bp = vcf_open(argv[optind], moder);
        h = vcf_hdr_read(bp);
        bout = vcf_open("-", modew);
-       if (!(vc.flag & VC_QCALL)) vcf_hdr_write(bout, h);
+       if (!(vc.flag & VC_QCALL)) {
+               if (vc.flag & VC_CALL) write_header(h);
+               vcf_hdr_write(bout, h);
+       }
        if (vc.flag & VC_CALL) {
                p1 = bcf_p1_init(h->n_smpl);
                if (vc.prior_file) {
index 9b661ff5b4fbc29dde6a53b5e71c027f131a81f4..b3832e7c6ff8b15aef729d9ed3faa8c4223aa78a 100644 (file)
@@ -101,7 +101,7 @@ int vcf_hdr_write(bcf_t *bp, const bcf_hdr_t *h)
                fwrite(h->txt, 1, h->l_txt - 1, v->fpout);
                if (strstr(h->txt, "##SQ=")) has_ref = 1;
        }
-       if (has_ver == 0) fprintf(v->fpout, "##fileformat=VCFv4.0\n");
+       if (h->l_txt == 0) fprintf(v->fpout, "##fileformat=VCFv4.0\n");
        fprintf(v->fpout, "#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO\tFORMAT");
        for (i = 0; i < h->n_smpl; ++i)
                fprintf(v->fpout, "\t%s", h->sns[i]);