# Uploaded samtools version 0.1.14-1.
[samtools.git] / bam_plcmd.c
index dca85181f99bf1c7343650846039153a019e7ea9..fd75cec8c8841a1f174731dc2eeef603571b86d2 100644 (file)
@@ -2,6 +2,8 @@
 #include <stdio.h>
 #include <unistd.h>
 #include <ctype.h>
+#include <string.h>
+#include <errno.h>
 #include "sam.h"
 #include "faidx.h"
 #include "bam_maqcns.h"
@@ -528,19 +530,28 @@ int bam_pileup(int argc, char *argv[])
 #define MPLP_NO_COMP 0x20
 #define MPLP_NO_ORPHAN 0x40
 #define MPLP_REALN   0x80
+#define MPLP_FMT_DP 0x100
+#define MPLP_FMT_SP 0x200
+#define MPLP_NO_INDEL 0x400
+#define MPLP_EXT_BAQ 0x800
+#define MPLP_ILLUMINA13 0x1000
 
 typedef struct {
-       int max_mq, min_mq, flag, min_baseQ, capQ_thres;
-       char *reg, *fn_pos;
+       int max_mq, min_mq, flag, min_baseQ, capQ_thres, max_depth, max_indel_depth;
+       int openQ, extQ, tandemQ, min_support; // for indels
+       double min_frac; // for indels
+       char *reg, *fn_pos, *pl_list;
        faidx_t *fai;
        kh_64_t *hash;
+       void *rghash;
 } mplp_conf_t;
 
 typedef struct {
        bamFile fp;
        bam_iter_t iter;
-       int min_mq, flag, ref_id, capQ_thres;
+       int ref_id;
        char *ref;
+       const mplp_conf_t *conf;
 } mplp_aux_t;
 
 typedef struct {
@@ -552,23 +563,35 @@ typedef struct {
 static int mplp_func(void *data, bam1_t *b)
 {
        extern int bam_realn(bam1_t *b, const char *ref);
-       extern int bam_prob_realn(bam1_t *b, const char *ref);
+       extern int bam_prob_realn_core(bam1_t *b, const char *ref, int);
        extern int bam_cap_mapQ(bam1_t *b, char *ref, int thres);
        mplp_aux_t *ma = (mplp_aux_t*)data;
        int ret, skip = 0;
        do {
-               int has_ref = (ma->ref && ma->ref_id == b->core.tid)? 1 : 0;
+               int has_ref;
                ret = ma->iter? bam_iter_read(ma->fp, ma->iter, b) : bam_read1(ma->fp, b);
                if (ret < 0) break;
+               if (ma->conf->rghash) { // exclude read groups
+                       uint8_t *rg = bam_aux_get(b, "RG");
+                       skip = (rg && bcf_str2id(ma->conf->rghash, (const char*)(rg+1)) >= 0);
+                       if (skip) continue;
+               }
+               if (ma->conf->flag & MPLP_ILLUMINA13) {
+                       int i;
+                       uint8_t *qual = bam1_qual(b);
+                       for (i = 0; i < b->core.l_qseq; ++i)
+                               qual[i] = qual[i] > 31? qual[i] - 31 : 0;
+               }
+               has_ref = (ma->ref && ma->ref_id == b->core.tid)? 1 : 0;
                skip = 0;
-               if (has_ref && (ma->flag&MPLP_REALN)) bam_prob_realn(b, ma->ref);
-               if (has_ref && ma->capQ_thres > 10) {
-                       int q = bam_cap_mapQ(b, ma->ref, ma->capQ_thres);
+               if (has_ref && (ma->conf->flag&MPLP_REALN)) bam_prob_realn_core(b, ma->ref, (ma->conf->flag & MPLP_EXT_BAQ)? 3 : 1);
+               if (has_ref && ma->conf->capQ_thres > 10) {
+                       int q = bam_cap_mapQ(b, ma->ref, ma->conf->capQ_thres);
                        if (q < 0) skip = 1;
                        else if (b->core.qual > q) b->core.qual = q;
                } else if (b->core.flag&BAM_FUNMAP) skip = 1;
-               else if (b->core.qual < ma->min_mq) skip = 1; 
-               else if ((ma->flag&MPLP_NO_ORPHAN) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
+               else if (b->core.qual < ma->conf->min_mq) skip = 1; 
+               else if ((ma->conf->flag&MPLP_NO_ORPHAN) && (b->core.flag&1) && !(b->core.flag&2)) skip = 1;
        } while (skip);
        return ret;
 }
@@ -598,13 +621,16 @@ static void group_smpl(mplp_pileup_t *m, bam_sample_t *sm, kstring_t *buf,
 
 static int mpileup(mplp_conf_t *conf, int n, char **fn)
 {
+       extern void *bcf_call_add_rg(void *rghash, const char *hdtext, const char *list);
+       extern void bcf_call_del_rghash(void *rghash);
        mplp_aux_t **data;
-       int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid;
+       int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid, max_depth, max_indel_depth;
        const bam_pileup1_t **plp;
        bam_mplp_t iter;
        bam_header_t *h = 0;
        char *ref;
        khash_t(64) *hash = 0;
+       void *rghash = 0;
 
        bcf_callaux_t *bca = 0;
        bcf_callret1_t *bcr = 0;
@@ -628,12 +654,11 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
        for (i = 0; i < n; ++i) {
                bam_header_t *h_tmp;
                data[i] = calloc(1, sizeof(mplp_aux_t));
-               data[i]->min_mq = conf->min_mq;
-               data[i]->flag = conf->flag;
-               data[i]->capQ_thres = conf->capQ_thres;
                data[i]->fp = strcmp(fn[i], "-") == 0? bam_dopen(fileno(stdin), "r") : bam_open(fn[i], "r");
+               data[i]->conf = conf;
                h_tmp = bam_header_read(data[i]->fp);
                bam_smpl_add(sm, fn[i], h_tmp->text);
+               rghash = bcf_call_add_rg(rghash, h_tmp->text, conf->pl_list);
                if (conf->reg) {
                        int beg, end;
                        bam_index_t *idx;
@@ -683,15 +708,29 @@ 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->l_txt = 0;
+               bh->txt = malloc(strlen(BAM_VERSION) + 64);
+               bh->l_txt = 1 + sprintf(bh->txt, "##samtoolsVersion=%s\n", BAM_VERSION);
                free(s.s);
                bcf_hdr_sync(bh);
                bcf_hdr_write(bp, bh);
                bca = bcf_call_init(-1., conf->min_baseQ);
                bcr = calloc(sm->n, sizeof(bcf_callret1_t));
+               bca->rghash = rghash;
+               bca->openQ = conf->openQ, bca->extQ = conf->extQ, bca->tandemQ = conf->tandemQ;
+               bca->min_frac = conf->min_frac;
+               bca->min_support = conf->min_support;
        }
        ref_tid = -1; ref = 0;
        iter = bam_mplp_init(n, mplp_func, (void**)data);
+       max_depth = conf->max_depth;
+       if (max_depth * sm->n > 1<<20)
+               fprintf(stderr, "(%s) Max depth is above 1M. Potential memory hog!\n", __func__);
+       if (max_depth * sm->n < 8000) {
+               max_depth = 8000 / sm->n;
+               fprintf(stderr, "<%s> Set max per-sample depth to %d\n", __func__, max_depth);
+       }
+       max_indel_depth = conf->max_indel_depth * sm->n;
+       bam_mplp_set_maxcnt(iter, max_depth);
        while (bam_mplp_auto(iter, &tid, &pos, n_plp, plp) > 0) {
                if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested
                if (hash) {
@@ -706,17 +745,31 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
                        ref_tid = tid;
                }
                if (conf->flag & MPLP_GLF) {
-                       int _ref0, ref16;
+                       int total_depth, _ref0, ref16;
                        bcf1_t *b = calloc(1, sizeof(bcf1_t));
+                       for (i = total_depth = 0; i < n; ++i) total_depth += n_plp[i];
                        group_smpl(&gplp, sm, &buf, n, fn, n_plp, plp);
                        _ref0 = (ref && pos < ref_len)? ref[pos] : 'N';
                        ref16 = bam_nt16_table[_ref0];
                        for (i = 0; i < gplp.n; ++i)
                                bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], ref16, bca, bcr + i);
                        bcf_call_combine(gplp.n, bcr, ref16, &bc);
-                       bcf_call2bcf(tid, pos, &bc, b);
+                       bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
+                                                (conf->flag&MPLP_FMT_SP), 0, 0);
                        bcf_write(bp, bh, b);
                        bcf_destroy(b);
+                       // call indels
+                       if (!(conf->flag&MPLP_NO_INDEL) && total_depth < max_indel_depth && bcf_call_gap_prep(gplp.n, gplp.n_plp, gplp.plp, pos, bca, ref, rghash) >= 0) {
+                               for (i = 0; i < gplp.n; ++i)
+                                       bcf_call_glfgen(gplp.n_plp[i], gplp.plp[i], -1, bca, bcr + i);
+                               if (bcf_call_combine(gplp.n, bcr, -1, &bc) >= 0) {
+                                       b = calloc(1, sizeof(bcf1_t));
+                                       bcf_call2bcf(tid, pos, &bc, b, (conf->flag&(MPLP_FMT_DP|MPLP_FMT_SP))? bcr : 0,
+                                                                (conf->flag&MPLP_FMT_SP), bca, ref);
+                                       bcf_write(bp, bh, b);
+                                       bcf_destroy(b);
+                               }
+                       }
                } else {
                        printf("%s\t%d\t%c", h->target_name[tid], pos + 1, (ref && pos < ref_len)? ref[pos] : 'N');
                        for (i = 0; i < n; ++i) {
@@ -743,6 +796,7 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
        bam_smpl_destroy(sm); free(buf.s);
        for (i = 0; i < gplp.n; ++i) free(gplp.plp[i]);
        free(gplp.plp); free(gplp.n_plp); free(gplp.m_plp);
+       bcf_call_del_rghash(rghash);
        if (hash) { // free the hash table
                khint_t k;
                for (k = kh_begin(hash); k < kh_end(hash); ++k)
@@ -761,53 +815,156 @@ static int mpileup(mplp_conf_t *conf, int n, char **fn)
        return 0;
 }
 
+#define MAX_PATH_LEN 1024
+static int read_file_list(const char *file_list,int *n,char **argv[])
+{
+    char buf[MAX_PATH_LEN];
+    int len, nfiles;
+    char **files;
+
+    FILE *fh = fopen(file_list,"r");
+    if ( !fh )
+    {
+        fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
+        return 1;
+    }
+
+    // Speed is not an issue here, determine the number of files by reading the file twice
+    nfiles = 0;
+    while ( fgets(buf,MAX_PATH_LEN,fh) ) nfiles++;
+
+    if ( fseek(fh, 0L, SEEK_SET) )
+    {
+        fprintf(stderr,"%s: %s\n", file_list,strerror(errno));
+        return 1;
+    }
+
+    files = calloc(nfiles,sizeof(char*));
+    nfiles = 0;
+    while ( fgets(buf,MAX_PATH_LEN,fh) ) 
+    {
+        len = strlen(buf);
+        while ( len>0 && isspace(buf[len-1]) ) len--;
+        if ( !len ) continue;
+
+        files[nfiles] = malloc(sizeof(char)*(len+1)); 
+        strncpy(files[nfiles],buf,len);
+        files[nfiles][len] = 0;
+        nfiles++;
+    }
+    fclose(fh);
+    if ( !nfiles )
+    {
+        fprintf(stderr,"No files read from %s\n", file_list);
+        return 1;
+    }
+    *argv = files;
+    *n    = nfiles;
+    return 0;
+}
+#undef MAX_PATH_LEN
+
 int bam_mpileup(int argc, char *argv[])
 {
        int c;
+    const char *file_list = NULL;
+    char **fn = NULL;
+    int nfiles = 0, use_orphan = 0;
        mplp_conf_t mplp;
        memset(&mplp, 0, sizeof(mplp_conf_t));
        mplp.max_mq = 60;
        mplp.min_baseQ = 13;
        mplp.capQ_thres = 0;
+       mplp.max_depth = 250; mplp.max_indel_depth = 250;
+       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;
-       while ((c = getopt(argc, argv, "gf:r:l:M:q:Q:uaORC:B")) >= 0) {
+       while ((c = getopt(argc, argv, "Agf:r:l:M:q:Q:uaRC:BDSd:L:b:P:o:e:h:Im:F:EG:6")) >= 0) {
                switch (c) {
                case 'f':
                        mplp.fai = fai_load(optarg);
                        if (mplp.fai == 0) return 1;
                        break;
+               case 'd': mplp.max_depth = atoi(optarg); break;
                case 'r': mplp.reg = strdup(optarg); break;
                case 'l': mplp.fn_pos = strdup(optarg); break;
+               case 'P': mplp.pl_list = strdup(optarg); break;
                case 'g': mplp.flag |= MPLP_GLF; break;
                case 'u': mplp.flag |= MPLP_NO_COMP | MPLP_GLF; break;
                case 'a': mplp.flag |= MPLP_NO_ORPHAN | MPLP_REALN; break;
-               case 'B': mplp.flag &= ~MPLP_REALN & ~MPLP_NO_ORPHAN; break;
-               case 'O': mplp.flag |= MPLP_NO_ORPHAN; break;
+               case 'B': mplp.flag &= ~MPLP_REALN; break;
                case 'R': mplp.flag |= MPLP_REALN; break;
+               case 'D': mplp.flag |= MPLP_FMT_DP; break;
+               case 'S': mplp.flag |= MPLP_FMT_SP; break;
+               case 'I': mplp.flag |= MPLP_NO_INDEL; break;
+               case 'E': mplp.flag |= MPLP_EXT_BAQ; break;
+               case '6': mplp.flag |= MPLP_ILLUMINA13; 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;
                case 'Q': mplp.min_baseQ = atoi(optarg); break;
+        case 'b': file_list = optarg; break;
+               case 'o': mplp.openQ = atoi(optarg); break;
+               case 'e': mplp.extQ = atoi(optarg); break;
+               case 'h': mplp.tandemQ = atoi(optarg); break;
+               case 'A': use_orphan = 1; break;
+               case 'F': mplp.min_frac = atof(optarg); break;
+               case 'm': mplp.min_support = atoi(optarg); break;
+               case 'L': mplp.max_indel_depth = atoi(optarg); break;
+               case 'G': {
+                               FILE *fp_rg;
+                               char buf[1024];
+                               mplp.rghash = bcf_str2id_init();
+                               if ((fp_rg = fopen(optarg, "r")) == 0)
+                                       fprintf(stderr, "(%s) Fail to open file %s. Continue anyway.\n", __func__, optarg);
+                               while (!feof(fp_rg) && fscanf(fp_rg, "%s", buf) > 0) // this is not a good style, but forgive me...
+                                       bcf_str2id_add(mplp.rghash, strdup(buf));
+                               fclose(fp_rg);
+                       }
+                       break;
                }
        }
+       if (use_orphan) mplp.flag &= ~MPLP_NO_ORPHAN;
        if (argc == 1) {
                fprintf(stderr, "\n");
                fprintf(stderr, "Usage:   samtools mpileup [options] in1.bam [in2.bam [...]]\n\n");
                fprintf(stderr, "Options: -f FILE     reference sequence file [null]\n");
                fprintf(stderr, "         -r STR      region in which pileup is generated [null]\n");
                fprintf(stderr, "         -l FILE     list of positions (format: chr pos) [null]\n");
+               fprintf(stderr, "         -b FILE     list of input BAM files [null]\n");
                fprintf(stderr, "         -M INT      cap mapping quality at INT [%d]\n", mplp.max_mq);
                fprintf(stderr, "         -Q INT      min base quality [%d]\n", mplp.min_baseQ);
                fprintf(stderr, "         -q INT      filter out alignment with MQ smaller than INT [%d]\n", mplp.min_mq);
+               fprintf(stderr, "         -d INT      max per-BAM depth [%d]\n", mplp.max_depth);
+               fprintf(stderr, "         -L INT      max per-sample depth for INDEL calling [%d]\n", mplp.max_indel_depth);
+               fprintf(stderr, "         -P STR      comma separated list of platforms for indels [all]\n");
+               fprintf(stderr, "         -o INT      Phred-scaled gap open sequencing error probability [%d]\n", mplp.openQ);
+               fprintf(stderr, "         -e INT      Phred-scaled gap extension seq error probability [%d]\n", mplp.extQ);
+               fprintf(stderr, "         -h INT      coefficient for homopolyer errors [%d]\n", mplp.tandemQ);
+               fprintf(stderr, "         -m INT      minimum gapped reads for indel candidates [%d]\n", mplp.min_support);
+               fprintf(stderr, "         -F FLOAT    minimum fraction of gapped reads for candidates [%g]\n", mplp.min_frac);
+               fprintf(stderr, "         -G FILE     exclude read groups listed in FILE [null]\n");
+               fprintf(stderr, "         -6          quality is in the Illumina-1.3+ encoding\n");
+               fprintf(stderr, "         -A          use anomalous read pairs in SNP/INDEL calling\n");
                fprintf(stderr, "         -g          generate BCF output\n");
                fprintf(stderr, "         -u          do not compress BCF output\n");
+               fprintf(stderr, "         -E          extended BAQ for higher sensitivity but lower specificity\n");
                fprintf(stderr, "         -B          disable BAQ computation\n");
+               fprintf(stderr, "         -D          output per-sample DP\n");
+               fprintf(stderr, "         -S          output per-sample SP (strand bias P-value, slow)\n");
+               fprintf(stderr, "         -I          do not perform indel calling\n");
                fprintf(stderr, "\n");
                fprintf(stderr, "Notes: Assuming diploid individuals.\n\n");
                return 1;
        }
-       mpileup(&mplp, argc - optind, argv + optind);
-       free(mplp.reg);
+    if (file_list) {
+        if ( read_file_list(file_list,&nfiles,&fn) ) return 1;
+        mpileup(&mplp,nfiles,fn);
+        for (c=0; c<nfiles; c++) free(fn[c]);
+        free(fn);
+    } else mpileup(&mplp, argc - optind, argv + optind);
+       if (mplp.rghash) bcf_str2id_thorough_destroy(mplp.rghash);
+       free(mplp.reg); free(mplp.pl_list);
        if (mplp.fai) fai_destroy(mplp.fai);
        return 0;
 }