Compress binary packages with xz.
[samtools.git] / bam_md.c
1 #include <unistd.h>
2 #include <assert.h>
3 #include <string.h>
4 #include <ctype.h>
5 #include <math.h>
6 #include "faidx.h"
7 #include "sam.h"
8 #include "kstring.h"
9 #include "kaln.h"
10 #include "kprobaln.h"
11
12 #define USE_EQUAL 1
13 #define DROP_TAG  2
14 #define BIN_QUAL  4
15 #define UPDATE_NM 8
16 #define UPDATE_MD 16
17 #define HASH_QNM  32
18
19 char bam_nt16_nt4_table[] = { 4, 0, 1, 4, 2, 4, 4, 4, 3, 4, 4, 4, 4, 4, 4, 4 };
20
21 int bam_aux_drop_other(bam1_t *b, uint8_t *s);
22
23 void bam_fillmd1_core(bam1_t *b, char *ref, int flag, int max_nm)
24 {
25         uint8_t *seq = bam1_seq(b);
26         uint32_t *cigar = bam1_cigar(b);
27         bam1_core_t *c = &b->core;
28         int i, x, y, u = 0;
29         kstring_t *str;
30         int32_t old_nm_i = -1, nm = 0;
31
32         str = (kstring_t*)calloc(1, sizeof(kstring_t));
33         for (i = y = 0, x = c->pos; i < c->n_cigar; ++i) {
34                 int j, l = cigar[i]>>4, op = cigar[i]&0xf;
35                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
36                         for (j = 0; j < l; ++j) {
37                                 int z = y + j;
38                                 int c1 = bam1_seqi(seq, z), c2 = bam_nt16_table[(int)ref[x+j]];
39                                 if (ref[x+j] == 0) break; // out of boundary
40                                 if ((c1 == c2 && c1 != 15 && c2 != 15) || c1 == 0) { // a match
41                                         if (flag&USE_EQUAL) seq[z/2] &= (z&1)? 0xf0 : 0x0f;
42                                         ++u;
43                                 } else {
44                                         kputw(u, str); kputc(ref[x+j], str);
45                                         u = 0; ++nm;
46                                 }
47                         }
48                         if (j < l) break;
49                         x += l; y += l;
50                 } else if (op == BAM_CDEL) {
51                         kputw(u, str); kputc('^', str);
52                         for (j = 0; j < l; ++j) {
53                                 if (ref[x+j] == 0) break;
54                                 kputc(ref[x+j], str);
55                         }
56                         u = 0;
57                         if (j < l) break;
58                         x += l; nm += l;
59                 } else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) {
60                         y += l;
61                         if (op == BAM_CINS) nm += l;
62                 } else if (op == BAM_CREF_SKIP) {
63                         x += l;
64                 }
65         }
66         kputw(u, str);
67         // apply max_nm
68         if (max_nm > 0 && nm >= max_nm) {
69                 for (i = y = 0, x = c->pos; i < c->n_cigar; ++i) {
70                         int j, l = cigar[i]>>4, op = cigar[i]&0xf;
71                         if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
72                                 for (j = 0; j < l; ++j) {
73                                         int z = y + j;
74                                         int c1 = bam1_seqi(seq, z), c2 = bam_nt16_table[(int)ref[x+j]];
75                                         if (ref[x+j] == 0) break; // out of boundary
76                                         if ((c1 == c2 && c1 != 15 && c2 != 15) || c1 == 0) { // a match
77                                                 seq[z/2] |= (z&1)? 0x0f : 0xf0;
78                                                 bam1_qual(b)[z] = 0;
79                                         }
80                                 }
81                                 if (j < l) break;
82                                 x += l; y += l;
83                         } else if (op == BAM_CDEL || op == BAM_CREF_SKIP) x += l;
84                         else if (op == BAM_CINS || op == BAM_CSOFT_CLIP) y += l;
85                 }
86         }
87         // update NM
88         if (flag & UPDATE_NM) {
89                 uint8_t *old_nm = bam_aux_get(b, "NM");
90                 if (c->flag & BAM_FUNMAP) return;
91                 if (old_nm) old_nm_i = bam_aux2i(old_nm);
92                 if (!old_nm) bam_aux_append(b, "NM", 'i', 4, (uint8_t*)&nm);
93                 else if (nm != old_nm_i) {
94                         fprintf(stderr, "[bam_fillmd1] different NM for read '%s': %d -> %d\n", bam1_qname(b), old_nm_i, nm);
95                         bam_aux_del(b, old_nm);
96                         bam_aux_append(b, "NM", 'i', 4, (uint8_t*)&nm);
97                 }
98         }
99         // update MD
100         if (flag & UPDATE_MD) {
101                 uint8_t *old_md = bam_aux_get(b, "MD");
102                 if (c->flag & BAM_FUNMAP) return;
103                 if (!old_md) bam_aux_append(b, "MD", 'Z', str->l + 1, (uint8_t*)str->s);
104                 else {
105                         int is_diff = 0;
106                         if (strlen((char*)old_md+1) == str->l) {
107                                 for (i = 0; i < str->l; ++i)
108                                         if (toupper(old_md[i+1]) != toupper(str->s[i]))
109                                                 break;
110                                 if (i < str->l) is_diff = 1;
111                         } else is_diff = 1;
112                         if (is_diff) {
113                                 fprintf(stderr, "[bam_fillmd1] different MD for read '%s': '%s' -> '%s'\n", bam1_qname(b), old_md+1, str->s);
114                                 bam_aux_del(b, old_md);
115                                 bam_aux_append(b, "MD", 'Z', str->l + 1, (uint8_t*)str->s);
116                         }
117                 }
118         }
119         // drop all tags but RG
120         if (flag&DROP_TAG) {
121                 uint8_t *q = bam_aux_get(b, "RG");
122                 bam_aux_drop_other(b, q);
123         }
124         // reduce the resolution of base quality
125         if (flag&BIN_QUAL) {
126                 uint8_t *qual = bam1_qual(b);
127                 for (i = 0; i < b->core.l_qseq; ++i)
128                         if (qual[i] >= 3) qual[i] = qual[i]/10*10 + 7;
129         }
130         free(str->s); free(str);
131 }
132
133 void bam_fillmd1(bam1_t *b, char *ref, int flag)
134 {
135         bam_fillmd1_core(b, ref, flag, 0);
136 }
137
138 int bam_cap_mapQ(bam1_t *b, char *ref, int thres)
139 {
140         uint8_t *seq = bam1_seq(b), *qual = bam1_qual(b);
141         uint32_t *cigar = bam1_cigar(b);
142         bam1_core_t *c = &b->core;
143         int i, x, y, mm, q, len, clip_l, clip_q;
144         double t;
145         if (thres < 0) thres = 40; // set the default
146         mm = q = len = clip_l = clip_q = 0;
147         for (i = y = 0, x = c->pos; i < c->n_cigar; ++i) {
148                 int j, l = cigar[i]>>4, op = cigar[i]&0xf;
149                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
150                         for (j = 0; j < l; ++j) {
151                                 int z = y + j;
152                                 int c1 = bam1_seqi(seq, z), c2 = bam_nt16_table[(int)ref[x+j]];
153                                 if (ref[x+j] == 0) break; // out of boundary
154                                 if (c2 != 15 && c1 != 15 && qual[z] >= 13) { // not ambiguous
155                                         ++len;
156                                         if (c1 && c1 != c2 && qual[z] >= 13) { // mismatch
157                                                 ++mm;
158                                                 q += qual[z] > 33? 33 : qual[z];
159                                         }
160                                 }
161                         }
162                         if (j < l) break;
163                         x += l; y += l; len += l;
164                 } else if (op == BAM_CDEL) {
165                         for (j = 0; j < l; ++j)
166                                 if (ref[x+j] == 0) break;
167                         if (j < l) break;
168                         x += l;
169                 } else if (op == BAM_CSOFT_CLIP) {
170                         for (j = 0; j < l; ++j) clip_q += qual[y+j];
171                         clip_l += l;
172                         y += l;
173                 } else if (op == BAM_CHARD_CLIP) {
174                         clip_q += 13 * l;
175                         clip_l += l;
176                 } else if (op == BAM_CINS) y += l;
177                 else if (op == BAM_CREF_SKIP) x += l;
178         }
179         for (i = 0, t = 1; i < mm; ++i)
180                 t *= (double)len / (i+1);
181         t = q - 4.343 * log(t) + clip_q / 5.;
182         if (t > thres) return -1;
183         if (t < 0) t = 0;
184         t = sqrt((thres - t) / thres) * thres;
185 //      fprintf(stderr, "%s %lf %d\n", bam1_qname(b), t, q);
186         return (int)(t + .499);
187 }
188
189 int bam_prob_realn_core(bam1_t *b, const char *ref, int flag)
190 {
191         int k, i, bw, x, y, yb, ye, xb, xe, apply_baq = flag&1, extend_baq = flag>>1&1;
192         uint32_t *cigar = bam1_cigar(b);
193         bam1_core_t *c = &b->core;
194         kpa_par_t conf = kpa_par_def;
195         uint8_t *bq = 0, *zq = 0, *qual = bam1_qual(b);
196         if ((c->flag & BAM_FUNMAP) || b->core.l_qseq == 0) return -1; // do nothing
197         // test if BQ or ZQ is present
198         if ((bq = bam_aux_get(b, "BQ")) != 0) ++bq;
199         if ((zq = bam_aux_get(b, "ZQ")) != 0 && *zq == 'Z') ++zq;
200         if (bq && zq) { // remove the ZQ tag
201                 bam_aux_del(b, zq-1);
202                 zq = 0;
203         }
204         if (bq || zq) {
205                 if ((apply_baq && zq) || (!apply_baq && bq)) return -3; // in both cases, do nothing
206                 if (bq && apply_baq) { // then convert BQ to ZQ
207                         for (i = 0; i < c->l_qseq; ++i)
208                                 qual[i] = qual[i] + 64 < bq[i]? 0 : qual[i] - ((int)bq[i] - 64);
209                         *(bq - 3) = 'Z';
210                 } else if (zq && !apply_baq) { // then convert ZQ to BQ
211                         for (i = 0; i < c->l_qseq; ++i)
212                                 qual[i] += (int)zq[i] - 64;
213                         *(zq - 3) = 'B';
214                 }
215                 return 0;
216         }
217         // find the start and end of the alignment      
218         x = c->pos, y = 0, yb = ye = xb = xe = -1;
219         for (k = 0; k < c->n_cigar; ++k) {
220                 int op, l;
221                 op = cigar[k]&0xf; l = cigar[k]>>4;
222                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
223                         if (yb < 0) yb = y;
224                         if (xb < 0) xb = x;
225                         ye = y + l; xe = x + l;
226                         x += l; y += l;
227                 } else if (op == BAM_CSOFT_CLIP || op == BAM_CINS) y += l;
228                 else if (op == BAM_CDEL) x += l;
229                 else if (op == BAM_CREF_SKIP) return -1; // do nothing if there is a reference skip
230         }
231         // set bandwidth and the start and the end
232         bw = 7;
233         if (abs((xe - xb) - (ye - yb)) > bw)
234                 bw = abs((xe - xb) - (ye - yb)) + 3;
235         conf.bw = bw;
236         xb -= yb + bw/2; if (xb < 0) xb = 0;
237         xe += c->l_qseq - ye + bw/2;
238         if (xe - xb - c->l_qseq > bw)
239                 xb += (xe - xb - c->l_qseq - bw) / 2, xe -= (xe - xb - c->l_qseq - bw) / 2;
240         { // glocal
241                 uint8_t *s, *r, *q, *seq = bam1_seq(b), *bq;
242                 int *state;
243                 bq = calloc(c->l_qseq + 1, 1);
244                 memcpy(bq, qual, c->l_qseq);
245                 s = calloc(c->l_qseq, 1);
246                 for (i = 0; i < c->l_qseq; ++i) s[i] = bam_nt16_nt4_table[bam1_seqi(seq, i)];
247                 r = calloc(xe - xb, 1);
248                 for (i = xb; i < xe; ++i) {
249                         if (ref[i] == 0) { xe = i; break; }
250                         r[i-xb] = bam_nt16_nt4_table[bam_nt16_table[(int)ref[i]]];
251                 }
252                 state = calloc(c->l_qseq, sizeof(int));
253                 q = calloc(c->l_qseq, 1);
254                 kpa_glocal(r, xe-xb, s, c->l_qseq, qual, &conf, state, q);
255                 if (!extend_baq) { // in this block, bq[] is capped by base quality qual[]
256                         for (k = 0, x = c->pos, y = 0; k < c->n_cigar; ++k) {
257                                 int op = cigar[k]&0xf, l = cigar[k]>>4;
258                                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
259                                         for (i = y; i < y + l; ++i) {
260                                                 if ((state[i]&3) != 0 || state[i]>>2 != x - xb + (i - y)) bq[i] = 0;
261                                                 else bq[i] = bq[i] < q[i]? bq[i] : q[i];
262                                         }
263                                         x += l; y += l;
264                                 } else if (op == BAM_CSOFT_CLIP || op == BAM_CINS) y += l;
265                                 else if (op == BAM_CDEL) x += l;
266                         }
267                         for (i = 0; i < c->l_qseq; ++i) bq[i] = qual[i] - bq[i] + 64; // finalize BQ
268                 } else { // in this block, bq[] is BAQ that can be larger than qual[] (different from the above!)
269                         uint8_t *left, *rght;
270                         left = calloc(c->l_qseq, 1); rght = calloc(c->l_qseq, 1);
271                         for (k = 0, x = c->pos, y = 0; k < c->n_cigar; ++k) {
272                                 int op = cigar[k]&0xf, l = cigar[k]>>4;
273                                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
274                                         for (i = y; i < y + l; ++i)
275                                                 bq[i] = ((state[i]&3) != 0 || state[i]>>2 != x - xb + (i - y))? 0 : q[i];
276                                         for (left[y] = bq[y], i = y + 1; i < y + l; ++i)
277                                                 left[i] = bq[i] > left[i-1]? bq[i] : left[i-1];
278                                         for (rght[y+l-1] = bq[y+l-1], i = y + l - 2; i >= y; --i)
279                                                 rght[i] = bq[i] > rght[i+1]? bq[i] : rght[i+1];
280                                         for (i = y; i < y + l; ++i)
281                                                 bq[i] = left[i] < rght[i]? left[i] : rght[i];
282                                         x += l; y += l;
283                                 } else if (op == BAM_CSOFT_CLIP || op == BAM_CINS) y += l;
284                                 else if (op == BAM_CDEL) x += l;
285                         }
286                         for (i = 0; i < c->l_qseq; ++i) bq[i] = 64 + (qual[i] <= bq[i]? 0 : qual[i] - bq[i]); // finalize BQ
287                         free(left); free(rght);
288                 }
289                 if (apply_baq) {
290                         for (i = 0; i < c->l_qseq; ++i) qual[i] -= bq[i] - 64; // modify qual
291                         bam_aux_append(b, "ZQ", 'Z', c->l_qseq + 1, bq);
292                 } else bam_aux_append(b, "BQ", 'Z', c->l_qseq + 1, bq);
293                 free(bq); free(s); free(r); free(q); free(state);
294         }
295         return 0;
296 }
297
298 int bam_prob_realn(bam1_t *b, const char *ref)
299 {
300         return bam_prob_realn_core(b, ref, 1);
301 }
302
303 int bam_fillmd(int argc, char *argv[])
304 {
305         int c, flt_flag, tid = -2, ret, len, is_bam_out, is_sam_in, is_uncompressed, max_nm, is_realn, capQ, baq_flag;
306         samfile_t *fp, *fpout = 0;
307         faidx_t *fai;
308         char *ref = 0, mode_w[8], mode_r[8];
309         bam1_t *b;
310
311         flt_flag = UPDATE_NM | UPDATE_MD;
312         is_bam_out = is_sam_in = is_uncompressed = is_realn = max_nm = capQ = baq_flag = 0;
313         mode_w[0] = mode_r[0] = 0;
314         strcpy(mode_r, "r"); strcpy(mode_w, "w");
315         while ((c = getopt(argc, argv, "EqreuNhbSC:n:Ad")) >= 0) {
316                 switch (c) {
317                 case 'r': is_realn = 1; break;
318                 case 'e': flt_flag |= USE_EQUAL; break;
319                 case 'd': flt_flag |= DROP_TAG; break;
320                 case 'q': flt_flag |= BIN_QUAL; break;
321                 case 'h': flt_flag |= HASH_QNM; break;
322                 case 'N': flt_flag &= ~(UPDATE_MD|UPDATE_NM); break;
323                 case 'b': is_bam_out = 1; break;
324                 case 'u': is_uncompressed = is_bam_out = 1; break;
325                 case 'S': is_sam_in = 1; break;
326                 case 'n': max_nm = atoi(optarg); break;
327                 case 'C': capQ = atoi(optarg); break;
328                 case 'A': baq_flag |= 1; break;
329                 case 'E': baq_flag |= 2; break;
330                 default: fprintf(stderr, "[bam_fillmd] unrecognized option '-%c'\n", c); return 1;
331                 }
332         }
333         if (!is_sam_in) strcat(mode_r, "b");
334         if (is_bam_out) strcat(mode_w, "b");
335         else strcat(mode_w, "h");
336         if (is_uncompressed) strcat(mode_w, "u");
337         if (optind + 1 >= argc) {
338                 fprintf(stderr, "\n");
339                 fprintf(stderr, "Usage:   samtools fillmd [-eubrS] <aln.bam> <ref.fasta>\n\n");
340                 fprintf(stderr, "Options: -e       change identical bases to '='\n");
341                 fprintf(stderr, "         -u       uncompressed BAM output (for piping)\n");
342                 fprintf(stderr, "         -b       compressed BAM output\n");
343                 fprintf(stderr, "         -S       the input is SAM with header\n");
344                 fprintf(stderr, "         -A       modify the quality string\n");
345                 fprintf(stderr, "         -r       compute the BQ tag (without -A) or cap baseQ by BAQ (with -A)\n");
346                 fprintf(stderr, "         -E       extended BAQ for better sensitivity but lower specificity\n\n");
347                 return 1;
348         }
349         fp = samopen(argv[optind], mode_r, 0);
350         if (fp == 0) return 1;
351         if (is_sam_in && (fp->header == 0 || fp->header->n_targets == 0)) {
352                 fprintf(stderr, "[bam_fillmd] input SAM does not have header. Abort!\n");
353                 return 1;
354         }
355         fpout = samopen("-", mode_w, fp->header);
356         fai = fai_load(argv[optind+1]);
357
358         b = bam_init1();
359         while ((ret = samread(fp, b)) >= 0) {
360                 if (b->core.tid >= 0) {
361                         if (tid != b->core.tid) {
362                                 free(ref);
363                                 ref = fai_fetch(fai, fp->header->target_name[b->core.tid], &len);
364                                 tid = b->core.tid;
365                                 if (ref == 0)
366                                         fprintf(stderr, "[bam_fillmd] fail to find sequence '%s' in the reference.\n",
367                                                         fp->header->target_name[tid]);
368                         }
369                         if (is_realn) bam_prob_realn_core(b, ref, baq_flag);
370                         if (capQ > 10) {
371                                 int q = bam_cap_mapQ(b, ref, capQ);
372                                 if (b->core.qual > q) b->core.qual = q;
373                         }
374                         if (ref) bam_fillmd1_core(b, ref, flt_flag, max_nm);
375                 }
376                 samwrite(fpout, b);
377         }
378         bam_destroy1(b);
379
380         free(ref);
381         fai_destroy(fai);
382         samclose(fp); samclose(fpout);
383         return 0;
384 }