]> git.donarmstrong.com Git - samtools.git/blobdiff - bam_plcmd.c
add --output,-c option to mpileup
[samtools.git] / bam_plcmd.c
index 9c4f273dcbfe2c3576e9c6ea0e2d7661839bbb85..15b65d2a1d31ea66dc36daee4c0a322253e60242 100644 (file)
@@ -9,6 +9,7 @@
 #include "sam.h"
 #include "faidx.h"
 #include "kstring.h"
+#include "sam_header.h"
 
 static inline int printw(int c, FILE *fp)
 {
@@ -85,7 +86,7 @@ typedef struct {
     int rflag_require, rflag_filter;
        int openQ, extQ, tandemQ, min_support; // for indels
        double min_frac; // for indels
-       char *reg, *pl_list;
+       char *reg, *pl_list, *fai_fname, *bcf_fname;
        faidx_t *fai;
        void *bed, *rghash;
 } mplp_conf_t;
@@ -208,6 +209,11 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
        n_plp = calloc(n, sizeof(int*));
        sm = bam_smpl_init();
 
+    if (n == 0) {
+        fprintf(stderr,"[%s] no input file/data given\n", __func__);
+        exit(1);
+    }
+
        // read the header and initialize data
        for (i = 0; i < n; ++i) {
                bam_header_t *h_tmp;
@@ -260,7 +266,7 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                kstring_t s;
                bh = calloc(1, sizeof(bcf_hdr_t));
                s.l = s.m = 0; s.s = 0;
-               bp = bcf_open("-", (conf->flag&MPLP_NO_COMP)? "wu" : "w");
+               bp = bcf_open(conf->bcf_fname, (conf->flag&MPLP_NO_COMP)? "wu" : "w");
                for (i = 0; i < h->n_targets; ++i) {
                        kputs(h->target_name[i], &s);
                        kputc('\0', &s);
@@ -275,9 +281,24 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                bh->l_smpl = s.l;
                bh->sname = malloc(s.l);
                memcpy(bh->sname, s.s, s.l);
-               bh->txt = malloc(strlen(BAM_VERSION) + 64);
-               bh->l_txt = 1 + sprintf(bh->txt, "##samtoolsVersion=%s\n", BAM_VERSION);
-               free(s.s);
+        s.l = 0;
+        ksprintf(&s, "##samtoolsVersion=%s\n", BAM_VERSION);
+        if (conf->fai_fname) ksprintf(&s, "##reference=file://%s\n", conf->fai_fname);
+        h->dict = sam_header_parse2(h->text);
+        int nseq;
+        const char *tags[] = {"SN","LN","UR","M5",NULL};
+        char **tbl = sam_header2tbl_n(h->dict, "SQ", tags, &nseq);
+        for (i=0; i<nseq; i++)
+        {
+            ksprintf(&s, "##contig=<ID=%s", tbl[4*i]);
+            if ( tbl[4*i+1] ) ksprintf(&s, ",length=%s", tbl[4*i+1]);
+            if ( tbl[4*i+2] ) ksprintf(&s, ",URL=%s", tbl[4*i+2]);
+            if ( tbl[4*i+3] ) ksprintf(&s, ",md5=%s", tbl[4*i+3]);
+            kputs(">\n", &s);
+        }
+        if (tbl) free(tbl);
+               bh->txt = s.s;
+               bh->l_txt = 1 + s.l;
                bcf_hdr_sync(bh);
                bcf_hdr_write(bp, bh);
                bca = bcf_call_init(-1., conf->min_baseQ);
@@ -474,19 +495,22 @@ int bam_mpileup(int argc, char *argv[])
        mplp.openQ = 40; mplp.extQ = 20; mplp.tandemQ = 100;
        mplp.min_frac = 0.002; mplp.min_support = 1;
        mplp.flag = MPLP_NO_ORPHAN | MPLP_REALN;
+    mplp.bcf_fname="-";
     static struct option lopts[] = 
     {
         {"rf",1,0,1},   // require flag
         {"ff",1,0,2},   // filter flag
+        {"output",1,0,(int)'c'}, // output flag
         {0,0,0,0}
     };
-       while ((c = getopt_long(argc, argv, "Agf:r:l:M:q:Q:uaRC:BDSd:L:b:P:po:e:h:Im:F:EG:6OsV1:2:",lopts,NULL)) >= 0) {
+       while ((c = getopt_long(argc, argv, "Agf:r:l:M:q:Q:uaRC:BDSd:L:b:P:po:e:h:Im:F:EG:6OsV1:2:c:",lopts,NULL)) >= 0) {
                switch (c) {
         case  1 : mplp.rflag_require = strtol(optarg,0,0); break;
         case  2 : mplp.rflag_filter  = strtol(optarg,0,0); break;
                case 'f':
                        mplp.fai = fai_load(optarg);
                        if (mplp.fai == 0) return 1;
+            mplp.fai_fname = optarg;
                        break;
                case 'd': mplp.max_depth = atoi(optarg); break;
                case 'r': mplp.reg = strdup(optarg); break;
@@ -506,6 +530,7 @@ int bam_mpileup(int argc, char *argv[])
                case 'R': mplp.flag |= MPLP_IGNORE_RG; break;
                case 's': mplp.flag |= MPLP_PRINT_MAPQ; break;
                case 'O': mplp.flag |= MPLP_PRINT_POS; break;
+               case 'c': mplp.bcf_fname = strdup(optarg); break;
                case 'C': mplp.capQ_thres = atoi(optarg); break;
                case 'M': mplp.max_mq = atoi(optarg); break;
                case 'q': mplp.min_mq = atoi(optarg); break;
@@ -560,6 +585,7 @@ int bam_mpileup(int argc, char *argv[])
                fprintf(stderr, "       -s           output mapping quality (disabled by -g/-u)\n");
                fprintf(stderr, "       -S           output per-sample strand bias P-value in BCF (require -g/-u)\n");
                fprintf(stderr, "       -u           generate uncompress BCF output\n");
+               fprintf(stderr, "       -c, --output=FILE output filename [stdout]\n");
                fprintf(stderr, "\nSNP/INDEL genotype likelihoods options (effective with `-g' or `-u'):\n\n");
                fprintf(stderr, "       -e INT       Phred-scaled gap extension seq error probability [%d]\n", mplp.extQ);
                fprintf(stderr, "       -F FLOAT     minimum fraction of gapped reads for candidates [%g]\n", mplp.min_frac);