]> git.donarmstrong.com Git - samtools.git/blob - bam_plcmd.c
* samtools-0.1.3-13 (r260)
[samtools.git] / bam_plcmd.c
1 #include <stdio.h>
2 #include <unistd.h>
3 #include <ctype.h>
4 #include "sam.h"
5 #include "faidx.h"
6 #include "bam_maqcns.h"
7 #include "khash.h"
8 #include "glf.h"
9 #include "kstring.h"
10
11 typedef int *indel_list_t;
12 KHASH_MAP_INIT_INT64(64, indel_list_t)
13
14 #define BAM_PLF_SIMPLE     0x01
15 #define BAM_PLF_CNS        0x02
16 #define BAM_PLF_INDEL_ONLY 0x04
17 #define BAM_PLF_GLF        0x08
18
19 typedef struct {
20         bam_header_t *h;
21         bam_maqcns_t *c;
22         bam_maqindel_opt_t *ido;
23         faidx_t *fai;
24         khash_t(64) *hash;
25         uint32_t format;
26         int tid, len, last_pos;
27         int mask;
28         char *ref;
29         glfFile fp_glf; // for glf output only
30 } pu_data_t;
31
32 char **__bam_get_lines(const char *fn, int *_n);
33 void bam_init_header_hash(bam_header_t *header);
34 int32_t bam_get_tid(const bam_header_t *header, const char *seq_name);
35
36 static khash_t(64) *load_pos(const char *fn, bam_header_t *h)
37 {
38         char **list;
39         int i, j, n, *fields, max_fields;
40         khash_t(64) *hash;
41         bam_init_header_hash(h);
42         list = __bam_get_lines(fn, &n);
43         hash = kh_init(64);
44         max_fields = 0; fields = 0;
45         for (i = 0; i < n; ++i) {
46                 char *str = list[i];
47                 int chr, n_fields, ret;
48                 khint_t k;
49                 uint64_t x;
50                 n_fields = ksplit_core(str, 0, &max_fields, &fields);
51                 if (n_fields < 2) continue;
52                 chr = bam_get_tid(h, str + fields[0]);
53                 if (chr < 0) {
54                         fprintf(stderr, "[load_pos] unknown reference sequence name: %s\n", str + fields[0]);
55                         continue;
56                 }
57                 x = (uint64_t)chr << 32 | (atoi(str + fields[1]) - 1);
58                 k = kh_put(64, hash, x, &ret);
59                 if (ret == 0) {
60                         fprintf(stderr, "[load_pos] position %s:%s has been loaded.\n", str+fields[0], str+fields[1]);
61                         continue;
62                 }
63                 kh_val(hash, k) = 0;
64                 if (n_fields > 2) {
65                         // count
66                         for (j = 2; j < n_fields; ++j) {
67                                 char *s = str + fields[j];
68                                 if ((*s != '+' && *s != '-') || !isdigit(s[1])) break;
69                         }
70                         if (j > 2) { // update kh_val()
71                                 int *q, y, z;
72                                 q = kh_val(hash, k) = (int*)calloc(j - 1, sizeof(int));
73                                 q[0] = j - 2; z = j; y = 1;
74                                 for (j = 2; j < z; ++j)
75                                         q[y++] = atoi(str + fields[j]);
76                         }
77                 }
78                 free(str);
79         }
80         free(list); free(fields);
81         return hash;
82 }
83
84 // an analogy to pileup_func() below
85 static int glt3_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
86 {
87         pu_data_t *d = (pu_data_t*)data;
88         bam_maqindel_ret_t *r = 0;
89         int rb, *proposed_indels = 0;
90         glf1_t *g;
91         glf3_t *g3;
92
93         if (d->fai == 0) {
94                 fprintf(stderr, "[glt3_func] reference sequence is required for generating GLT. Abort!\n");
95                 exit(1);
96         }
97         if (d->hash) { // only output a list of sites
98                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
99                 if (k == kh_end(d->hash)) return 0;
100                 proposed_indels = kh_val(d->hash, k);
101         }
102         g3 = glf3_init1();
103         if (d->fai && (int)tid != d->tid) {
104                 if (d->ref) { // then write the end mark
105                         g3->rtype = GLF3_RTYPE_END;
106                         glf3_write1(d->fp_glf, g3);
107                 }
108                 glf3_ref_write(d->fp_glf, d->h->target_name[tid], d->h->target_len[tid]); // write reference
109                 free(d->ref);
110                 d->ref = fai_fetch(d->fai, d->h->target_name[tid], &d->len);
111                 d->tid = tid;
112                 d->last_pos = 0;
113         }
114         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
115         g = bam_maqcns_glfgen(n, pu, bam_nt16_table[rb], d->c);
116         memcpy(g3, g, sizeof(glf1_t));
117         g3->rtype = GLF3_RTYPE_SUB;
118         g3->offset = pos - d->last_pos;
119         d->last_pos = pos;
120         glf3_write1(d->fp_glf, g3);
121         if (proposed_indels)
122                 r = bam_maqindel(n, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
123         else r = bam_maqindel(n, pos, d->ido, pu, d->ref, 0, 0);
124         if (r) { // then write indel line
125                 int het = 3 * n, min;
126                 min = het;
127                 if (min > r->gl[0]) min = r->gl[0];
128                 if (min > r->gl[1]) min = r->gl[1];
129                 g3->ref_base = 0;
130                 g3->rtype = GLF3_RTYPE_INDEL;
131                 memset(g3->lk, 0, 10);
132                 g3->lk[0] = r->gl[0] - min < 255? r->gl[0] - min : 255;
133                 g3->lk[1] = r->gl[1] - min < 255? r->gl[1] - min : 255;
134                 g3->lk[2] = het - min < 255? het - min : 255;
135                 g3->offset = 0;
136                 g3->indel_len[0] = r->indel1;
137                 g3->indel_len[1] = r->indel2;
138                 g3->min_lk = min < 255? min : 255;
139                 g3->max_len = (abs(r->indel1) > abs(r->indel2)? abs(r->indel1) : abs(r->indel2)) + 1;
140                 g3->indel_seq[0] = strdup(r->s[0]+1);
141                 g3->indel_seq[1] = strdup(r->s[1]+1);
142                 glf3_write1(d->fp_glf, g3);
143                 bam_maqindel_ret_destroy(r);
144         }
145         free(g);
146         glf3_destroy1(g3);
147         return 0;
148 }
149
150 static int pileup_func(uint32_t tid, uint32_t pos, int n, const bam_pileup1_t *pu, void *data)
151 {
152         pu_data_t *d = (pu_data_t*)data;
153         bam_maqindel_ret_t *r = 0;
154         int i, j, rb, max_mapq = 0, *proposed_indels = 0;
155         uint32_t x;
156
157         if (d->format & BAM_PLF_GLF) return glt3_func(tid, pos, n, pu, data);
158         if (d->hash) { // only output a list of sites
159                 khint_t k = kh_get(64, d->hash, (uint64_t)tid<<32|pos);
160                 if (k == kh_end(d->hash)) return 0;
161                 proposed_indels = kh_val(d->hash, k);
162         }
163         if (d->fai && (int)tid != d->tid) { // then update d->ref
164                 free(d->ref);
165                 d->ref = fai_fetch(d->fai, d->h->target_name[tid], &d->len);
166                 d->tid = tid;
167         }
168         rb = (d->ref && (int)pos < d->len)? d->ref[pos] : 'N';
169         if (d->format & BAM_PLF_INDEL_ONLY) {
170                 for (i = 0; i < n; ++i)
171                         if (pu[i].indel != 0) break;
172                 if (i == n) return 0;
173         }
174         printf("%s\t%d\t%c\t", d->h->target_name[tid], pos + 1, rb);
175         if (d->format & BAM_PLF_CNS) { // consensus
176                 int ref_q, rb4 = bam_nt16_table[rb];
177                 x = bam_maqcns_call(n, pu, d->c);
178                 ref_q = 0;
179                 if (rb4 != 15 && x>>28 != 15 && x>>28 != rb4) { // a SNP
180                         ref_q = ((x>>24&0xf) == rb4)? x>>8&0xff : (x>>8&0xff) + (x&0xff);
181                         if (ref_q > 255) ref_q = 255;
182                 }
183                 printf("%c\t%d\t%d\t%d\t", bam_nt16_rev_table[x>>28], x>>8&0xff, ref_q, x>>16&0xff);
184         }
185         if ((d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)) && d->ref) {
186                 if (proposed_indels)
187                         r = bam_maqindel(n, pos, d->ido, pu, d->ref, proposed_indels[0], proposed_indels+1);
188                 else r = bam_maqindel(n, pos, d->ido, pu, d->ref, 0, 0);
189         }
190         // pileup strings
191         printf("%d\t", n);
192         for (i = 0; i < n; ++i) {
193                 const bam_pileup1_t *p = pu + i;
194                 if (max_mapq < p->b->core.qual) max_mapq = p->b->core.qual;
195                 if (p->is_head) printf("^%c", p->b->core.qual > 93? 126 : p->b->core.qual + 33);
196                 if (!p->is_del) {
197                         int c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos)];
198                         if (c == '=' || toupper(c) == toupper(rb)) c = bam1_strand(p->b)? ',' : '.';
199                         else c = bam1_strand(p->b)? tolower(c) : toupper(c);
200                         putchar(c);
201                         if (p->indel > 0) {
202                                 printf("+%d", p->indel);
203                                 for (j = 1; j <= p->indel; ++j) {
204                                         c = bam_nt16_rev_table[bam1_seqi(bam1_seq(p->b), p->qpos + j)];
205                                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
206                                 }
207                         } else if (p->indel < 0) {
208                                 printf("%d", p->indel);
209                                 for (j = 1; j <= -p->indel; ++j) {
210                                         c = (d->ref && (int)pos+j < d->len)? d->ref[pos+j] : 'N';
211                                         putchar(bam1_strand(p->b)? tolower(c) : toupper(c));
212                                 }
213                         }
214                 } else putchar('*');
215                 if (p->is_tail) putchar('$');
216         }
217         putchar('\t');
218         for (i = 0; i < n; ++i) {
219                 const bam_pileup1_t *p = pu + i;
220                 int c = bam1_qual(p->b)[p->qpos] + 33;
221                 if (c > 126) c = 126;
222                 putchar(c);
223         }
224         if (d->format & BAM_PLF_SIMPLE) {
225                 putchar('\t');
226                 for (i = 0; i < n; ++i) {
227                         int c = pu[i].b->core.qual + 33;
228                         if (c > 126) c = 126;
229                         putchar(c);
230                 }
231         }
232         putchar('\n');
233         if (r) { // then print indel line
234                 printf("%s\t%d\t*\t", d->h->target_name[tid], pos + 1);
235                 if (r->gt < 2) printf("%s/%s\t", r->s[r->gt], r->s[r->gt]);
236                 else printf("%s/%s\t", r->s[0], r->s[1]);
237                 printf("%d\t%d\t", r->q_cns, r->q_ref);
238                 printf("%d\t%d\t", max_mapq, n);
239                 printf("%s\t%s\t", r->s[0], r->s[1]);
240                 //printf("%d\t%d\t", r->gl[0], r->gl[1]);
241                 printf("%d\t%d\t%d\n", r->cnt1, r->cnt2, r->cnt_anti);
242                 bam_maqindel_ret_destroy(r);
243         }
244         return 0;
245 }
246
247 int bam_pileup(int argc, char *argv[])
248 {
249         int c;
250         char *fn_list = 0, *fn_fa = 0, *fn_pos = 0;
251         pu_data_t *d = (pu_data_t*)calloc(1, sizeof(pu_data_t));
252         d->tid = -1; d->mask = BAM_DEF_MASK;
253         d->c = bam_maqcns_init();
254         d->ido = bam_maqindel_opt_init();
255         while ((c = getopt(argc, argv, "st:f:cT:N:r:l:im:gI:G:")) >= 0) {
256                 switch (c) {
257                 case 's': d->format |= BAM_PLF_SIMPLE; break;
258                 case 't': fn_list = strdup(optarg); break;
259                 case 'l': fn_pos = strdup(optarg); break;
260                 case 'f': fn_fa = strdup(optarg); break;
261                 case 'T': d->c->theta = atof(optarg); break;
262                 case 'N': d->c->n_hap = atoi(optarg); break;
263                 case 'r': d->c->het_rate = atof(optarg); break;
264                 case 'c': d->format |= BAM_PLF_CNS; break;
265                 case 'i': d->format |= BAM_PLF_INDEL_ONLY; break;
266                 case 'm': d->mask = strtol(optarg, 0, 0); break;
267                 case 'g': d->format |= BAM_PLF_GLF; break;
268                 case 'I': d->ido->q_indel = atoi(optarg); break;
269                 case 'G': d->ido->r_indel = atof(optarg); break;
270                 default: fprintf(stderr, "Unrecognizd option '-%c'.\n", c); return 1;
271                 }
272         }
273         if (optind == argc) {
274                 fprintf(stderr, "\n");
275                 fprintf(stderr, "Usage:  samtools pileup [options] <in.bam>|<in.sam>\n\n");
276                 fprintf(stderr, "Option: -s        simple (yet incomplete) pileup format\n");
277                 fprintf(stderr, "        -i        only show lines/consensus with indels\n");
278                 fprintf(stderr, "        -m INT    filtering reads with bits in INT [%d]\n", d->mask);
279                 fprintf(stderr, "        -t FILE   list of reference sequences (assume the input is in SAM)\n");
280                 fprintf(stderr, "        -l FILE   list of sites at which pileup is output\n");
281                 fprintf(stderr, "        -f FILE   reference sequence in the FASTA format\n\n");
282                 fprintf(stderr, "        -c        output the maq consensus sequence\n");
283                 fprintf(stderr, "        -g        output in the GLFv3 format (suppressing -c/-i/-s)\n");
284                 fprintf(stderr, "        -T FLOAT  theta in maq consensus calling model (for -c/-g) [%f]\n", d->c->theta);
285                 fprintf(stderr, "        -N INT    number of haplotypes in the sample (for -c/-g) [%d]\n", d->c->n_hap);
286                 fprintf(stderr, "        -r FLOAT  prior of a difference between two haplotypes (for -c/-g) [%f]\n", d->c->het_rate);
287                 fprintf(stderr, "        -G FLOAT  prior of an indel between two haplotypes (for -c/-g) [%f]\n", d->ido->r_indel);
288                 fprintf(stderr, "        -I INT    phred prob. of an indel in sequencing/prep. (for -c/-g) [%d]\n", d->ido->q_indel);
289                 fprintf(stderr, "\n");
290                 free(fn_list); free(fn_fa); free(d);
291                 return 1;
292         }
293         if (fn_fa) d->fai = fai_load(fn_fa);
294         if (d->format & (BAM_PLF_CNS|BAM_PLF_GLF)) bam_maqcns_prepare(d->c); // consensus calling
295         if (d->format & BAM_PLF_GLF) { // for glf output
296                 glf3_header_t *h;
297                 h = glf3_header_init();
298                 d->fp_glf = bgzf_fdopen(fileno(stdout), "w");
299                 glf3_header_write(d->fp_glf, h);
300                 glf3_header_destroy(h);
301         }
302         if (d->fai == 0 && (d->format & (BAM_PLF_CNS|BAM_PLF_INDEL_ONLY)))
303                 fprintf(stderr, "[bam_pileup] indels will not be called when -f is absent.\n");
304         {
305                 samfile_t *fp;
306                 fp = fn_list? samopen(argv[optind], "r", fn_list) : samopen(argv[optind], "rb", 0);
307                 if (fp == 0 || fp->header == 0) {
308                         fprintf(stderr, "[bam_pileup] fail to read the header: non-exisiting file or wrong format.\n");
309                         return 1;
310                 }
311                 d->h = fp->header;
312                 if (fn_pos) d->hash = load_pos(fn_pos, d->h);
313                 sampileup(fp, d->mask, pileup_func, d);
314                 samclose(fp); // d->h will be destroyed here
315         }
316
317         // free
318         if (d->format & BAM_PLF_GLF) bgzf_close(d->fp_glf);
319         if (fn_pos) { // free the hash table
320                 khint_t k;
321                 for (k = kh_begin(d->hash); k < kh_end(d->hash); ++k)
322                         if (kh_exist(d->hash, k)) free(kh_val(d->hash, k));
323                 kh_destroy(64, d->hash);
324         }
325         free(fn_pos); free(fn_list); free(fn_fa);
326         if (d->fai) fai_destroy(d->fai);
327         bam_maqcns_destroy(d->c);
328         free(d->ido); free(d->ref); free(d);
329         return 0;
330 }