]> git.donarmstrong.com Git - samtools.git/blob - padding.c
Update BIN values in 'samtools depad'
[samtools.git] / padding.c
1 #include <string.h>
2 #include <assert.h>
3 #include <unistd.h>
4 #include "kstring.h"
5 #include "sam_header.h"
6 #include "sam.h"
7 #include "bam.h"
8 #include "faidx.h"
9
10 bam_header_t *bam_header_dup(const bam_header_t *h0); /*in sam.c*/
11
12 static void replace_cigar(bam1_t *b, int n, uint32_t *cigar)
13 {
14         if (n != b->core.n_cigar) {
15                 int o = b->core.l_qname + b->core.n_cigar * 4;
16                 if (b->data_len + (n - b->core.n_cigar) * 4 > b->m_data) {
17                         b->m_data = b->data_len + (n - b->core.n_cigar) * 4;
18                         kroundup32(b->m_data);
19                         b->data = (uint8_t*)realloc(b->data, b->m_data);
20                 }
21                 memmove(b->data + b->core.l_qname + n * 4, b->data + o, b->data_len - o);
22                 memcpy(b->data + b->core.l_qname, cigar, n * 4);
23                 b->data_len += (n - b->core.n_cigar) * 4;
24                 b->core.n_cigar = n;
25         } else memcpy(b->data + b->core.l_qname, cigar, n * 4);
26 }
27
28 #define write_cigar(_c, _n, _m, _v) do { \
29                 if (_n == _m) { \
30                         _m = _m? _m<<1 : 4; \
31                         _c = (uint32_t*)realloc(_c, _m * 4); \
32                 } \
33                 _c[_n++] = (_v); \
34         } while (0)
35
36 static void unpad_seq(bam1_t *b, kstring_t *s)
37 {
38         int k, j, i;
39         int length;
40         uint32_t *cigar = bam1_cigar(b);
41         uint8_t *seq = bam1_seq(b);
42         // b->core.l_qseq gives length of the SEQ entry (including soft clips, S)
43         // We need the padded length after alignment from the CIGAR (excluding
44         // soft clips S, but including pads from CIGAR D operations)
45         length = 0;
46         for (k = 0; k < b->core.n_cigar; ++k) {
47                 int op, ol;
48                 op= bam_cigar_op(cigar[k]);
49                 ol = bam_cigar_oplen(cigar[k]);
50                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF || op == BAM_CDEL)
51                         length += ol;
52         }
53         ks_resize(s, length);
54         for (k = 0, s->l = 0, j = 0; k < b->core.n_cigar; ++k) {
55                 int op, ol;
56                 op = bam_cigar_op(cigar[k]);
57                 ol = bam_cigar_oplen(cigar[k]);
58                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
59                         for (i = 0; i < ol; ++i, ++j) s->s[s->l++] = bam1_seqi(seq, j);
60                 } else if (op == BAM_CSOFT_CLIP) {
61                         j += ol;
62                 } else if (op == BAM_CHARD_CLIP) {
63                         /* do nothing */
64                 } else if (op == BAM_CDEL) {
65                         for (i = 0; i < ol; ++i) s->s[s->l++] = 0;
66                 } else {
67                         fprintf(stderr, "[depad] ERROR: Didn't expect CIGAR op %c in read %s\n", BAM_CIGAR_STR[op], bam1_qname(b));
68                         assert(-1);
69                 }
70         }
71         assert(length == s->l);
72 }
73
74 int load_unpadded_ref(faidx_t *fai, char *ref_name, int ref_len, kstring_t *seq)
75 {
76         char base;
77         char *fai_ref = 0;
78         int fai_ref_len = 0, k;
79
80         fai_ref = fai_fetch(fai, ref_name, &fai_ref_len);
81         if (fai_ref_len != ref_len) {
82                 fprintf(stderr, "[depad] ERROR: FASTA sequence %s length %i, expected %i\n", ref_name, fai_ref_len, ref_len);
83                 free(fai_ref);
84                 return -1;
85         }
86         ks_resize(seq, ref_len);
87         seq->l = 0;
88         for (k = 0; k < ref_len; ++k) {
89                 base = fai_ref[k];
90                 if (base == '-' || base == '*') {
91                         // Map gaps to null to match unpad_seq function
92                         seq->s[seq->l++] = 0;
93                 } else {
94                         int i = bam_nt16_table[(int)base];
95                         if (i == 0 || i==16) { // Equals maps to 0, anything unexpected to 16
96                                 fprintf(stderr, "[depad] ERROR: Invalid character %c (ASCII %i) in FASTA sequence %s\n", base, (int)base, ref_name);
97                                 free(fai_ref);
98                                 return -1;
99                         }
100                         seq->s[seq->l++] = i;
101                 }
102         }
103         assert(ref_len == seq->l);
104         free(fai_ref);
105         return 0;
106 }
107
108 int get_unpadded_len(faidx_t *fai, char *ref_name, int padded_len)
109 {
110         char base;
111         char *fai_ref = 0;
112         int fai_ref_len = 0, k;
113         int bases=0, gaps=0;
114
115         fai_ref = fai_fetch(fai, ref_name, &fai_ref_len);
116         if (fai_ref_len != padded_len) {
117                 fprintf(stderr, "[depad] ERROR: FASTA sequence '%s' length %i, expected %i\n", ref_name, fai_ref_len, padded_len);
118                 free(fai_ref);
119                 return -1;
120         }
121         for (k = 0; k < padded_len; ++k) {
122                 //fprintf(stderr, "[depad] checking base %i of %i or %i\n", k+1, ref_len, strlen(fai_ref));
123                 base = fai_ref[k];
124                 if (base == '-' || base == '*') {
125                         gaps += 1;
126                 } else {
127                         int i = bam_nt16_table[(int)base];
128                         if (i == 0 || i==16) { // Equals maps to 0, anything unexpected to 16
129                                 fprintf(stderr, "[depad] ERROR: Invalid character %c (ASCII %i) in FASTA sequence '%s'\n", base, (int)base, ref_name);
130                                 free(fai_ref);
131                                 return -1;
132                         }
133                         bases += 1;
134                 }
135         }
136         free(fai_ref);
137         assert (padded_len == bases + gaps);
138         return bases;
139 }
140
141 inline int * update_posmap(int *posmap, kstring_t ref)
142 {
143         int i, k;
144         posmap = realloc(posmap, ref.m * sizeof(int));
145         for (i = k = 0; i < ref.l; ++i) {
146                 posmap[i] = k;
147                 if (ref.s[i]) ++k;
148         }
149         return posmap;
150 }
151
152 int bam_pad2unpad(samfile_t *in, samfile_t *out, faidx_t *fai)
153 {
154         bam_header_t *h = 0;
155         bam1_t *b = 0;
156         kstring_t r, q;
157         int r_tid = -1;
158         uint32_t *cigar2 = 0;
159         int ret = 0, n2 = 0, m2 = 0, *posmap = 0;
160
161         b = bam_init1();
162         r.l = r.m = q.l = q.m = 0; r.s = q.s = 0;
163         int read_ret;
164         h = in->header;
165         while ((read_ret = samread(in, b)) >= 0) { // read one alignment from `in'
166                 uint32_t *cigar = bam1_cigar(b);
167                 n2 = 0;
168                 if (b->core.pos == 0 && b->core.tid >= 0 && strcmp(bam1_qname(b), h->target_name[b->core.tid]) == 0) {
169                         // fprintf(stderr, "[depad] Found embedded reference '%s'\n", bam1_qname(b));
170                         r_tid = b->core.tid;
171                         unpad_seq(b, &r);
172                         if (h->target_len[r_tid] != r.l) {
173                                 fprintf(stderr, "[depad] ERROR: (Padded) length of '%s' is %d in BAM header, but %ld in embedded reference\n", bam1_qname(b), h->target_len[r_tid], r.l);
174                                 return -1;
175                         }
176                         if (fai) {
177                                 // Check the embedded reference matches the FASTA file
178                                 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &q)) {
179                                         fprintf(stderr, "[depad] ERROR: Failed to load embedded reference '%s' from FASTA\n", h->target_name[b->core.tid]);
180                                         return -1;
181                                 }
182                                 assert(r.l == q.l);
183                                 int i;
184                                 for (i = 0; i < r.l; ++i) {
185                                         if (r.s[i] != q.s[i]) {
186                                                 // Show gaps as ASCII 45
187                                                 fprintf(stderr, "[depad] ERROR: Embedded sequence and reference FASTA don't match for %s base %i, '%c' vs '%c'\n",
188                                                         h->target_name[b->core.tid], i+1,
189                                                         r.s[i] ? bam_nt16_rev_table[(int)r.s[i]] : 45,
190                                                         q.s[i] ? bam_nt16_rev_table[(int)q.s[i]] : 45);
191                                                 return -1;
192                                         }
193                                 }
194                         }
195                         write_cigar(cigar2, n2, m2, bam_cigar_gen(b->core.l_qseq, BAM_CMATCH));
196                         replace_cigar(b, n2, cigar2);
197                         posmap = update_posmap(posmap, r);
198                         b->core.bin = bam_reg2bin(0, bam_calend(&b->core, bam1_cigar(b)));
199                 } else if (b->core.n_cigar > 0) {
200                         int i, k, op;
201                         if (b->core.tid < 0) {
202                                 fprintf(stderr, "[depad] ERROR: Read '%s' has CIGAR but no RNAME\n", bam1_qname(b));
203                                 return -1;
204                         } else if (b->core.tid == r_tid) {
205                                 ; // good case, reference available
206                                 //fprintf(stderr, "[depad] Have ref '%s' for read '%s'\n", h->target_name[b->core.tid], bam1_qname(b));
207                         } else if (fai) {
208                                 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &r)) {
209                                         fprintf(stderr, "[depad] ERROR: Failed to load '%s' from reference FASTA\n", h->target_name[b->core.tid]);
210                                         return -1;
211                                 }
212                                 posmap = update_posmap(posmap, r);
213                                 r_tid = b->core.tid;
214                                 // fprintf(stderr, "[depad] Loaded %s from FASTA file\n", h->target_name[b->core.tid]);
215                         } else {                                
216                                 fprintf(stderr, "[depad] ERROR: Missing %s embedded reference sequence (and no FASTA file)\n", h->target_name[b->core.tid]);
217                                 return -1;
218                         }
219                         unpad_seq(b, &q);
220                         if (bam_cigar_op(cigar[0]) == BAM_CSOFT_CLIP) {
221                                 write_cigar(cigar2, n2, m2, cigar[0]);
222                         } else if (bam_cigar_op(cigar[0]) == BAM_CHARD_CLIP) {
223                                 write_cigar(cigar2, n2, m2, cigar[0]);
224                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[1]) == BAM_CSOFT_CLIP) {
225                                         write_cigar(cigar2, n2, m2, cigar[1]);
226                                 }
227                         }
228                         /* Determine CIGAR operator for each base in the aligned read */
229                         for (i = 0, k = b->core.pos; i < q.l; ++i, ++k)
230                                 q.s[i] = q.s[i]? (r.s[k]? BAM_CMATCH : BAM_CINS) : (r.s[k]? BAM_CDEL : BAM_CPAD);
231                         /* Include any pads if starts with an insert */
232                         if (q.s[0] == BAM_CINS) {
233                                 for (k = 0; k+1 < b->core.pos && !r.s[b->core.pos - k - 1]; ++k);
234                                 if (k) write_cigar(cigar2, n2, m2, bam_cigar_gen(k, BAM_CPAD));
235                         }
236                         /* Count consecutive CIGAR operators to turn into a CIGAR string */
237                         for (i = k = 1, op = q.s[0]; i < q.l; ++i) {
238                                 if (op != q.s[i]) {
239                                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
240                                         op = q.s[i]; k = 1;
241                                 } else ++k;
242                         }
243                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
244                         if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CSOFT_CLIP) {
245                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
246                         } else if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CHARD_CLIP) {
247                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[b->core.n_cigar-2]) == BAM_CSOFT_CLIP) {
248                                         write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-2]);
249                                 }
250                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
251                         }
252                         /* Remove redundant P operators between M/X/=/D operators, e.g. 5M2P10M -> 15M */
253                         int pre_op, post_op;
254                         for (i = 2; i < n2; ++i)
255                                 if (bam_cigar_op(cigar2[i-1]) == BAM_CPAD) {
256                                         pre_op = bam_cigar_op(cigar2[i-2]);
257                                         post_op = bam_cigar_op(cigar2[i]);
258                                         /* Note don't need to check for X/= as code above will use M only */
259                                         if ((pre_op == BAM_CMATCH || pre_op == BAM_CDEL) && (post_op == BAM_CMATCH || post_op == BAM_CDEL)) {
260                                                 /* This is a redundant P operator */
261                                                 cigar2[i-1] = 0; // i.e. 0M
262                                                 /* If had same operator either side, combine them in post_op */
263                                                 if (pre_op == post_op) {
264                                                         /* If CIGAR M, could treat as simple integers since BAM_CMATCH is zero*/
265                                                         cigar2[i] = bam_cigar_gen(bam_cigar_oplen(cigar2[i-2]) + bam_cigar_oplen(cigar2[i]), post_op);
266                                                         cigar2[i-2] = 0; // i.e. 0M
267                                                 }
268                                         }
269                                 }
270                         /* Remove the zero'd operators (0M) */
271                         for (i = k = 0; i < n2; ++i)
272                                 if (cigar2[i]) cigar2[k++] = cigar2[i];
273                         n2 = k;
274                         replace_cigar(b, n2, cigar2);
275                         b->core.pos = posmap[b->core.pos];
276                         b->core.bin = bam_reg2bin(b->core.pos, bam_calend(&b->core.pos, bam1_cigar(b)));
277                         if (b->core.mtid < 0 || b->core.mpos < 0) {
278                                 /* Nice case, no mate to worry about*/
279                                 // fprintf(stderr, "[depad] Read '%s' mate not mapped\n", bam1_qname(b));
280                                 /* TODO - Warning if FLAG says mate should be mapped? */
281                                 /* Clean up funny input where mate position is given but mate reference is missing: */
282                                 b->core.mtid = -1;
283                                 b->core.mpos = -1;
284                         } else if (b->core.mtid == b->core.tid) {
285                                 /* Nice case, same reference */
286                                 // fprintf(stderr, "[depad] Read '%s' mate mapped to same ref\n", bam1_qname(b));
287                                 b->core.mpos = posmap[b->core.mpos];
288                         } else {
289                                 /* Nasty case, Must load alternative posmap */
290                                 // fprintf(stderr, "[depad] Loading reference '%s' temporarily\n", h->target_name[b->core.mtid]);
291                                 if (!fai) {
292                                         fprintf(stderr, "[depad] ERROR: Needed reference %s sequence for mate (and no FASTA file)\n", h->target_name[b->core.mtid]);
293                                         return -1;
294                                 }
295                                 /* Temporarily load the other reference sequence */
296                                 if (load_unpadded_ref(fai, h->target_name[b->core.mtid], h->target_len[b->core.mtid], &r)) {
297                                         fprintf(stderr, "[depad] ERROR: Failed to load '%s' from reference FASTA\n", h->target_name[b->core.mtid]);
298                                         return -1;
299                                 }
300                                 posmap = update_posmap(posmap, r);
301                                 b->core.mpos = posmap[b->core.mpos];
302                                 /* Restore the reference and posmap*/
303                                 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &r)) {
304                                         fprintf(stderr, "[depad] ERROR: Failed to load '%s' from reference FASTA\n", h->target_name[b->core.tid]);
305                                         return -1;
306                                 }
307                                 posmap = update_posmap(posmap, r);
308                         }
309                 }
310                 samwrite(out, b);
311         }
312         if (read_ret < -1) {
313                 fprintf(stderr, "[depad] truncated file.\n");
314                 ret = 1;
315         }
316         free(r.s); free(q.s); free(posmap);
317         bam_destroy1(b);
318         return ret;
319 }
320
321 bam_header_t * fix_header(bam_header_t *old, faidx_t *fai)
322 {
323         int i = 0, unpadded_len = 0;
324         bam_header_t *header = 0 ;
325
326         header = bam_header_dup(old);
327         for (i = 0; i < old->n_targets; ++i) {
328                 unpadded_len = get_unpadded_len(fai, old->target_name[i], old->target_len[i]);
329                 if (unpadded_len < 0) {
330                         fprintf(stderr, "[depad] ERROR getting unpadded length of '%s', padded length %i\n", old->target_name[i], old->target_len[i]);
331                 } else {
332                         header->target_len[i] = unpadded_len;
333                         //fprintf(stderr, "[depad] Recalculating '%s' length %i -> %i\n", old->target_name[i], old->target_len[i], header->target_len[i]);
334                 }
335         }
336         /* Duplicating the header allocated new buffer for header string */
337         /* After modifying the @SQ lines it will only get smaller, since */
338         /* the LN entries will be the same or shorter, and we'll remove */
339         /* any MD entries (MD5 checksums). */
340         assert(strlen(old->text) == strlen(header->text));
341         assert (0==strcmp(old->text, header->text));
342         const char *text;
343         text = old->text;
344         header->text[0] = '\0'; /* Resuse the allocated buffer */
345         char * newtext = header->text;
346         char * end=NULL;
347         while (text[0]=='@') {
348                 end = strchr(text, '\n');
349                 assert(end != 0);
350                 if (text[1]=='S' && text[2]=='Q' && text[3]=='\t') {
351                         /* TODO - edit the @SQ line here to remove MD and fix LN. */
352                         /* For now just remove the @SQ line, and samtools will */
353                         /* automatically generate a minimal replacement with LN. */
354                         /* However, that discards any other tags like AS, SP, UR. */
355                         //fprintf(stderr, "[depad] Removing @SQ line\n");
356                 } else {
357                         /* Copy this line to the new header */
358                         strncat(newtext, text, end - text + 1);
359                 }
360                 text = end + 1;
361         }
362         assert (text[0]=='\0');
363         /* Check we didn't overflow the buffer */
364         assert (strlen(header->text) <= strlen(old->text));
365         if (strlen(header->text) < header->l_text) {
366                 //fprintf(stderr, "[depad] Reallocating header buffer\n");
367                 assert (newtext == header->text);
368                 newtext = malloc(strlen(header->text) + 1);
369                 strcpy(newtext, header->text);
370                 free(header->text);
371                 header->text = newtext;
372                 header->l_text = strlen(newtext);
373         }
374         //fprintf(stderr, "[depad] Here is the new header (pending @SQ lines),\n\n%s\n(end)\n", header->text);
375         return header;
376 }
377
378 static int usage(int is_long_help);
379
380 int main_pad2unpad(int argc, char *argv[])
381 {
382         samfile_t *in = 0, *out = 0;
383         bam_header_t *h = 0;
384         faidx_t *fai = 0;
385         int c, is_bamin = 1, compress_level = -1, is_bamout = 1, is_long_help = 0;
386         char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0;
387         int ret=0;
388
389         /* parse command-line options */
390         strcpy(in_mode, "r"); strcpy(out_mode, "w");
391         while ((c = getopt(argc, argv, "Sso:u1T:?")) >= 0) {
392                 switch (c) {
393                 case 'S': is_bamin = 0; break;
394                 case 's': assert(compress_level == -1); is_bamout = 0; break;
395                 case 'o': fn_out = strdup(optarg); break;
396                 case 'u': assert(is_bamout == 1); compress_level = 0; break;
397                 case '1': assert(is_bamout == 1); compress_level = 1; break;
398                 case 'T': fn_ref = strdup(optarg); break;
399                 case '?': is_long_help = 1; break;
400                 default: return usage(is_long_help);
401                 }
402         }
403         if (argc == optind) return usage(is_long_help);
404
405         if (is_bamin) strcat(in_mode, "b");
406         if (is_bamout) strcat(out_mode, "b");
407         strcat(out_mode, "h");
408         if (compress_level >= 0) {
409                 char tmp[2];
410                 tmp[0] = compress_level + '0'; tmp[1] = '\0';
411                 strcat(out_mode, tmp);
412         }
413
414         // Load FASTA reference (also needed for SAM -> BAM if missing header)
415         if (fn_ref) {
416                 fn_list = samfaipath(fn_ref);
417                 fai = fai_load(fn_ref);
418         }
419         // open file handlers
420         if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) {
421                 fprintf(stderr, "[depad] failed to open \"%s\" for reading.\n", argv[optind]);
422                 ret = 1;
423                 goto depad_end;
424         }
425         if (in->header == 0) {
426                 fprintf(stderr, "[depad] failed to read the header from \"%s\".\n", argv[optind]);
427                 ret = 1;
428                 goto depad_end;
429         }
430         if (in->header->text == 0 || in->header->l_text == 0) {
431                 fprintf(stderr, "[depad] Warning - failed to read any header text from \"%s\".\n", argv[optind]);
432                 assert (0 == in->header->l_text);
433                 assert (0 == in->header->text);
434         }
435         if (fn_ref) {
436                 h = fix_header(in->header, fai);
437         } else {
438                 fprintf(stderr, "[depad] Warning - reference lengths will not be corrected without FASTA reference\n");
439                 h = in->header;
440         }
441         if ((out = samopen(fn_out? fn_out : "-", out_mode, h)) == 0) {
442                 fprintf(stderr, "[depad] failed to open \"%s\" for writing.\n", fn_out? fn_out : "standard output");
443                 ret = 1;
444                 goto depad_end;
445         }
446
447         // Do the depad
448         ret = bam_pad2unpad(in, out, fai);
449
450 depad_end:
451         // close files, free and return
452         if (fai) fai_destroy(fai);
453         if (h != in->header) bam_header_destroy(h);
454         samclose(in);
455         samclose(out);
456         free(fn_list); free(fn_out);
457         return ret;
458 }
459
460 static int usage(int is_long_help)
461 {
462         fprintf(stderr, "\n");
463         fprintf(stderr, "Usage:   samtools depad <in.bam>\n\n");
464         fprintf(stderr, "Options: -s       output is SAM (default is BAM)\n");
465         fprintf(stderr, "         -S       input is SAM (default is BAM)\n");
466         fprintf(stderr, "         -u       uncompressed BAM output (can't use with -s)\n");
467         fprintf(stderr, "         -1       fast compression BAM output (can't use with -s)\n");
468         fprintf(stderr, "         -T FILE  reference sequence file [null]\n");
469         fprintf(stderr, "         -o FILE  output file name [stdout]\n");
470         fprintf(stderr, "         -?       longer help\n");
471         fprintf(stderr, "\n");
472         if (is_long_help)
473                 fprintf(stderr, "Notes:\n\
474 \n\
475   1. Requires embedded reference sequences (before the reads for that reference),\n\
476      with the future aim to also support a FASTA padded reference sequence file.\n\
477 \n\
478   2. The input padded alignment read's CIGAR strings must not use P or I operators.\n\
479 \n");
480         return 1;
481 }