X-Git-Url: https://git.donarmstrong.com/?a=blobdiff_plain;f=bam_plcmd.c;h=c665a923a8a2ec1eaab957444012689ffd8fab16;hb=8b0f98e2ac96b5ab07fd482cbaaa06d458254a90;hp=f8ae6b451f2f1d1da9a2b17d6c7085f4f5f2600c;hpb=7e56a1610a2f932ba2ef66b0e89b1170b9c654e1;p=samtools.git diff --git a/bam_plcmd.c b/bam_plcmd.c index f8ae6b4..c665a92 100644 --- a/bam_plcmd.c +++ b/bam_plcmd.c @@ -21,6 +21,7 @@ KHASH_MAP_INIT_INT64(64, indel_list_t) #define BAM_PLF_RANBASE 0x40 #define BAM_PLF_1STBASE 0x80 #define BAM_PLF_ALLBASE 0x100 +#define BAM_PLF_READPOS 0x200 typedef struct { bam_header_t *h; @@ -157,11 +158,37 @@ static int glt3_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, return 0; } +static void pileup_seq(const bam_pileup1_t *p, int pos, int ref_len, const char *ref) +{ + if (p->is_head) printf("^%c", p->b->core.qual > 93? 126 : p->b->core.qual + 33); + if (!p->is_del) { + int j, rb, c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)]; + rb = (ref && pos < ref_len)? ref[pos] : 'N'; + if (c == '=' || toupper(c) == toupper(rb)) c = bam1_strand(p->b)? ',' : '.'; + else c = bam1_strand(p->b)? tolower(c) : toupper(c); + putchar(c); + if (p->indel > 0) { + printf("+%d", p->indel); + for (j = 1; j <= p->indel; ++j) { + c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)]; + putchar(bam1_strand(p->b)? tolower(c) : toupper(c)); + } + } else if (p->indel < 0) { + printf("%d", p->indel); + for (j = 1; j <= -p->indel; ++j) { + c = (ref && (int)pos+j < ref_len)? ref[pos+j] : 'N'; + putchar(bam1_strand(p->b)? tolower(c) : toupper(c)); + } + } + } else putchar('*'); + if (p->is_tail) putchar('$'); +} + static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data) { pu_data_t *d = (pu_data_t*)data; bam_maqindel_ret_t *r = 0; - int i, j, rb, rms_mapq = -1, *proposed_indels = 0; + int i, rb, rms_mapq = -1, *proposed_indels = 0; uint64_t rms_aux; uint32_t cns = 0; @@ -242,27 +269,7 @@ static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *p const bam_pileup1_t *p = pu + i; int tmp = p->b->core.qual < d->c->cap_mapQ? p->b->core.qual : d->c->cap_mapQ; rms_aux += tmp * tmp; - if (p->is_head) printf("^%c", p->b->core.qual > 93? 126 : p->b->core.qual + 33); - if (!p->is_del) { - int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)]; - if (c == '=' || toupper(c) == toupper(rb)) c = bam1_strand(p->b)? ',' : '.'; - else c = bam1_strand(p->b)? tolower(c) : toupper(c); - putchar(c); - if (p->indel > 0) { - printf("+%d", p->indel); - for (j = 1; j <= p->indel; ++j) { - c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)]; - putchar(bam1_strand(p->b)? tolower(c) : toupper(c)); - } - } else if (p->indel < 0) { - printf("%d", p->indel); - for (j = 1; j <= -p->indel; ++j) { - c = (d->ref && (int)pos+j < d->len)? d->ref[pos+j] : 'N'; - putchar(bam1_strand(p->b)? tolower(c) : toupper(c)); - } - } - } else putchar('*'); - if (p->is_tail) putchar('$'); + pileup_seq(p, pos, d->len, d->ref); } // finalize rms_mapq rms_aux = (uint64_t)(sqrt((double)rms_aux / n) + .499); @@ -299,6 +306,15 @@ static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *p putchar(c); } } + // print read position + if (d->format & BAM_PLF_READPOS) { + putchar('\t'); + for (i = 0; i < n; ++i) { + int x = pu[i].qpos; + int l = pu[i].b->core.l_qseq; + printf("%d,", x < l/2? x+1 : -((l-1)-x+1)); + } + } putchar('\n'); // print the indel line if r has been calculated. This only happens if: // a) -c or -i are flagged, AND b) the reference sequence is available @@ -325,17 +341,19 @@ int bam_pileup(int argc, char *argv[]) d->max_depth = 0; d->tid = -1; d->mask = BAM_DEF_MASK; d->c = bam_maqcns_init(); + d->c->is_soap = 1; // change the default model d->ido = bam_maqindel_opt_init(); - while ((c = getopt(argc, argv, "st:f:cT:N:r:l:d:im:gI:G:vM:S2aR:")) >= 0) { + while ((c = getopt(argc, argv, "st:f:cT:N:r:l:d:im:gI:G:vM:S2aR:PA")) >= 0) { switch (c) { case 'a': d->c->is_soap = 1; break; + case 'A': d->c->is_soap = 0; break; case 's': d->format |= BAM_PLF_SIMPLE; break; case 't': fn_list = strdup(optarg); break; case 'l': fn_pos = strdup(optarg); break; case 'f': fn_fa = strdup(optarg); break; case 'T': d->c->theta = atof(optarg); break; case 'N': d->c->n_hap = atoi(optarg); break; - case 'r': d->c->het_rate = atof(optarg); break; + case 'r': d->c->het_rate = atof(optarg); d->ido->r_snp = d->c->het_rate; break; case 'M': d->c->cap_mapQ = atoi(optarg); break; case 'd': d->max_depth = atoi(optarg); break; case 'c': d->format |= BAM_PLF_CNS; break; @@ -344,6 +362,7 @@ int bam_pileup(int argc, char *argv[]) case 'm': d->mask = strtol(optarg, 0, 0); break; case 'g': d->format |= BAM_PLF_GLF; break; case '2': d->format |= BAM_PLF_2ND; break; + case 'P': d->format |= BAM_PLF_READPOS; break; case 'I': d->ido->q_indel = atoi(optarg); break; case 'G': d->ido->r_indel = atof(optarg); break; case 'S': is_SAM = 1; break; @@ -362,7 +381,7 @@ int bam_pileup(int argc, char *argv[]) fprintf(stderr, "Usage: samtools pileup [options] |\n\n"); fprintf(stderr, "Option: -s simple (yet incomplete) pileup format\n"); fprintf(stderr, " -S the input is in SAM\n"); - fprintf(stderr, " -a use the SOAPsnp model for SNP calling\n"); + fprintf(stderr, " -A use the MAQ model for SNP calling\n"); fprintf(stderr, " -2 output the 2nd best call and quality\n"); fprintf(stderr, " -i only show lines/consensus with indels\n"); fprintf(stderr, " -m INT filtering reads with bits in INT [%d]\n", d->mask); @@ -371,7 +390,7 @@ int bam_pileup(int argc, char *argv[]) fprintf(stderr, " -t FILE list of reference sequences (force -S)\n"); fprintf(stderr, " -l FILE list of sites at which pileup is output\n"); fprintf(stderr, " -f FILE reference sequence in the FASTA format\n\n"); - fprintf(stderr, " -c output the maq consensus sequence\n"); + fprintf(stderr, " -c output the SOAPsnp consensus sequence\n"); fprintf(stderr, " -v print variants only (for -c)\n"); fprintf(stderr, " -g output in the GLFv3 format (suppressing -c/-i/-s)\n"); fprintf(stderr, " -T FLOAT theta in maq consensus calling model (for -c/-g) [%f]\n", d->c->theta); @@ -424,3 +443,143 @@ int bam_pileup(int argc, char *argv[]) free(d->ido); free(d->ref); free(d); return 0; } + +/*********** + * mpileup * + ***********/ + +typedef struct { + char *reg, *fn_pos; + faidx_t *fai; + kh_64_t *hash; +} mplp_conf_t; + +typedef struct { + bamFile fp; + bam_iter_t iter; +} mplp_aux_t; + +static int mplp_func(void *data, bam1_t *b) +{ + mplp_aux_t *ma = (mplp_aux_t*)data; + if (ma->iter) return bam_iter_read(ma->fp, ma->iter, b); + return bam_read1(ma->fp, b); +} + +static int mpileup(mplp_conf_t *conf, int n, char **fn) +{ + mplp_aux_t **data; + int i, tid, pos, *n_plp, beg0 = 0, end0 = 1u<<29, ref_len, ref_tid; + const bam_pileup1_t **plp; + bam_mplp_t iter; + bam_header_t *h = 0; + char *ref; + khash_t(64) *hash = 0; + // allocate + data = calloc(n, sizeof(void*)); + plp = calloc(n, sizeof(void*)); + n_plp = calloc(n, sizeof(int*)); + // read the header and initialize data + for (i = 0; i < n; ++i) { + bam_header_t *h_tmp; + data[i] = calloc(1, sizeof(mplp_aux_t)); + data[i]->fp = bam_open(fn[i], "r"); + h_tmp = bam_header_read(data[i]->fp); + if (conf->reg) { + int beg, end; + bam_index_t *idx; + idx = bam_index_load(fn[i]); + if (idx == 0) { + fprintf(stderr, "[%s] fail to load index for %d-th input.\n", __func__, i+1); + exit(1); + } + if (bam_parse_region(h_tmp, conf->reg, &tid, &beg, &end) < 0) { + fprintf(stderr, "[%s] malformatted region or wrong seqname for %d-th input.\n", __func__, i+1); + exit(1); + } + if (i == 0) beg0 = beg, end0 = end; + data[i]->iter = bam_iter_query(idx, tid, beg, end); + bam_index_destroy(idx); + } + if (i == 0) h = h_tmp; + else { + // FIXME: to check consistency + bam_header_destroy(h_tmp); + } + } + if (conf->fn_pos) hash = load_pos(conf->fn_pos, h); + // mpileup + 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) { + if (conf->reg && (pos < beg0 || pos >= end0)) continue; // out of the region requested + if (hash) { + khint_t k; + k = kh_get(64, hash, (uint64_t)tid<<32 | pos); + if (k == kh_end(hash)) continue; + } + if (tid != ref_tid) { + free(ref); + if (conf->fai) ref = fai_fetch(conf->fai, h->target_name[tid], &ref_len); + ref_tid = tid; + } + 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) { + int j; + printf("\t%d\t", n_plp[i]); + if (n_plp[i] == 0) printf("*\t*"); + else { + for (j = 0; j < n_plp[i]; ++j) + pileup_seq(plp[i] + j, pos, ref_len, ref); + putchar('\t'); + for (j = 0; j < n_plp[i]; ++j) { + const bam_pileup1_t *p = plp[i] + j; + int c = bam1_qual(p->b)[p->qpos] + 33; + if (c > 126) c = 126; + putchar(c); + } + } + } + putchar('\n'); + } + if (hash) { // free the hash table + khint_t k; + for (k = kh_begin(hash); k < kh_end(hash); ++k) + if (kh_exist(hash, k)) free(kh_val(hash, k)); + kh_destroy(64, hash); + } + bam_mplp_destroy(iter); + bam_header_destroy(h); + for (i = 0; i < n; ++i) { + bam_close(data[i]->fp); + if (data[i]->iter) bam_iter_destroy(data[i]->iter); + free(data[i]); + } + free(data); free(plp); free(ref); free(n_plp); + return 0; +} + +int bam_mpileup(int argc, char *argv[]) +{ + int c; + mplp_conf_t mplp; + memset(&mplp, 0, sizeof(mplp_conf_t)); + while ((c = getopt(argc, argv, "f:r:l:")) >= 0) { + switch (c) { + case 'f': + mplp.fai = fai_load(optarg); + if (mplp.fai == 0) return 1; + break; + case 'r': mplp.reg = strdup(optarg); break; + case 'l': mplp.fn_pos = strdup(optarg); break; + } + } + if (argc == 1) { + fprintf(stderr, "Usage: samtools mpileup [-r reg] [-f in.fa] [-l pos] in1.bam [in2.bam [...]]\n"); + return 1; + } + mpileup(&mplp, argc - optind, argv + optind); + free(mplp.reg); + if (mplp.fai) fai_destroy(mplp.fai); + return 0; +}