]> git.donarmstrong.com Git - samtools.git/blob - padding.c
First attempt at updating MPOS 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                 } else if (b->core.n_cigar > 0) {
199                         int i, k, op;
200                         if (b->core.tid < 0) {
201                                 fprintf(stderr, "[depad] ERROR: Read '%s' has CIGAR but no RNAME\n", bam1_qname(b));
202                                 return -1;
203                         } else if (b->core.tid == r_tid) {
204                                 ; // good case, reference available
205                         } else if (fai) {
206                                 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &r)) {
207                                         fprintf(stderr, "[depad] ERROR: Failed to load '%s' from reference FASTA\n", h->target_name[b->core.tid]);
208                                         return -1;
209                                 }
210                                 posmap = update_posmap(posmap, r);
211                                 r_tid = b->core.tid;
212                                 // fprintf(stderr, "[depad] Loaded %s from FASTA file\n", h->target_name[b->core.tid]);
213                         } else {                                
214                                 fprintf(stderr, "[depad] ERROR: Missing %s embedded reference sequence (and no FASTA file)\n", h->target_name[b->core.tid]);
215                                 return -1;
216                         }
217                         unpad_seq(b, &q);
218                         if (bam_cigar_op(cigar[0]) == BAM_CSOFT_CLIP) {
219                                 write_cigar(cigar2, n2, m2, cigar[0]);
220                         } else if (bam_cigar_op(cigar[0]) == BAM_CHARD_CLIP) {
221                                 write_cigar(cigar2, n2, m2, cigar[0]);
222                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[1]) == BAM_CSOFT_CLIP) {
223                                         write_cigar(cigar2, n2, m2, cigar[1]);
224                                 }
225                         }
226                         /* Determine CIGAR operator for each base in the aligned read */
227                         for (i = 0, k = b->core.pos; i < q.l; ++i, ++k)
228                                 q.s[i] = q.s[i]? (r.s[k]? BAM_CMATCH : BAM_CINS) : (r.s[k]? BAM_CDEL : BAM_CPAD);
229                         /* Include any pads if starts with an insert */
230                         if (q.s[0] == BAM_CINS) {
231                                 for (k = 0; k+1 < b->core.pos && !r.s[b->core.pos - k - 1]; ++k);
232                                 if (k) write_cigar(cigar2, n2, m2, bam_cigar_gen(k, BAM_CPAD));
233                         }
234                         /* Count consecutive CIGAR operators to turn into a CIGAR string */
235                         for (i = k = 1, op = q.s[0]; i < q.l; ++i) {
236                                 if (op != q.s[i]) {
237                                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
238                                         op = q.s[i]; k = 1;
239                                 } else ++k;
240                         }
241                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
242                         if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CSOFT_CLIP) {
243                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
244                         } else if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CHARD_CLIP) {
245                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[b->core.n_cigar-2]) == BAM_CSOFT_CLIP) {
246                                         write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-2]);
247                                 }
248                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
249                         }
250                         /* Remove redundant P operators between M/X/=/D operators, e.g. 5M2P10M -> 15M */
251                         int pre_op, post_op;
252                         for (i = 2; i < n2; ++i)
253                                 if (bam_cigar_op(cigar2[i-1]) == BAM_CPAD) {
254                                         pre_op = bam_cigar_op(cigar2[i-2]);
255                                         post_op = bam_cigar_op(cigar2[i]);
256                                         /* Note don't need to check for X/= as code above will use M only */
257                                         if ((pre_op == BAM_CMATCH || pre_op == BAM_CDEL) && (post_op == BAM_CMATCH || post_op == BAM_CDEL)) {
258                                                 /* This is a redundant P operator */
259                                                 cigar2[i-1] = 0; // i.e. 0M
260                                                 /* If had same operator either side, combine them in post_op */
261                                                 if (pre_op == post_op) {
262                                                         /* If CIGAR M, could treat as simple integers since BAM_CMATCH is zero*/
263                                                         cigar2[i] = bam_cigar_gen(bam_cigar_oplen(cigar2[i-2]) + bam_cigar_oplen(cigar2[i]), post_op);
264                                                         cigar2[i-2] = 0; // i.e. 0M
265                                                 }
266                                         }
267                                 }
268                         /* Remove the zero'd operators (0M) */
269                         for (i = k = 0; i < n2; ++i)
270                                 if (cigar2[i]) cigar2[k++] = cigar2[i];
271                         n2 = k;
272                         replace_cigar(b, n2, cigar2);
273                         b->core.pos = posmap[b->core.pos];
274                         if (b->core.mpos < 0) {
275                                 /* Nice case, no mate to worry about*/
276                         } else if (b->core.mtid == b->core.tid) {
277                                 /* Nice case, same reference */
278                                 b->core.mpos = posmap[b->core.mpos];
279                         } else {
280                                 /* Nasty case, Must load alternative posmap */
281                                 if (!fai) {
282                                         fprintf(stderr, "[depad] ERROR: Needed reference %s sequence for mate (and no FASTA file)\n", h->target_name[b->core.mtid]);
283                                         return -1;
284                                 }
285                                 /* Temporarily load the other reference sequence */
286                                 if (load_unpadded_ref(fai, h->target_name[b->core.mtid], h->target_len[b->core.mtid], &r)) {
287                                         fprintf(stderr, "[depad] ERROR: Failed to load '%s' from reference FASTA\n", h->target_name[b->core.mtid]);
288                                         return -1;
289                                 }
290                                 posmap = update_posmap(posmap, r);
291                                 b->core.mpos = posmap[b->core.mpos];
292                                 /* Restore the reference and posmap*/
293                                 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &r)) {
294                                         fprintf(stderr, "[depad] ERROR: Failed to load '%s' from reference FASTA\n", h->target_name[b->core.tid]);
295                                         return -1;
296                                 }
297                                 posmap = update_posmap(posmap, r);
298                         }
299                 }
300                 samwrite(out, b);
301         }
302         if (read_ret < -1) {
303                 fprintf(stderr, "[depad] truncated file.\n");
304                 ret = 1;
305         }
306         free(r.s); free(q.s); free(posmap);
307         bam_destroy1(b);
308         return ret;
309 }
310
311 bam_header_t * fix_header(bam_header_t *old, faidx_t *fai)
312 {
313         int i = 0, unpadded_len = 0;
314         bam_header_t *header = 0 ;
315
316         header = bam_header_dup(old);
317         for (i = 0; i < old->n_targets; ++i) {
318                 unpadded_len = get_unpadded_len(fai, old->target_name[i], old->target_len[i]);
319                 if (unpadded_len < 0) {
320                         fprintf(stderr, "[depad] ERROR getting unpadded length of '%s', padded length %i\n", old->target_name[i], old->target_len[i]);
321                 } else {
322                         header->target_len[i] = unpadded_len;
323                         //fprintf(stderr, "[depad] Recalculating '%s' length %i -> %i\n", old->target_name[i], old->target_len[i], header->target_len[i]);
324                 }
325         }
326         /* Duplicating the header allocated new buffer for header string */
327         /* After modifying the @SQ lines it will only get smaller, since */
328         /* the LN entries will be the same or shorter, and we'll remove */
329         /* any MD entries (MD5 checksums). */
330         assert(strlen(old->text) == strlen(header->text));
331         assert (0==strcmp(old->text, header->text));
332         const char *text;
333         text = old->text;
334         header->text[0] = '\0'; /* Resuse the allocated buffer */
335         char * newtext = header->text;
336         char * end=NULL;
337         while (text[0]=='@') {
338                 end = strchr(text, '\n');
339                 assert(end != 0);
340                 if (text[1]=='S' && text[2]=='Q' && text[3]=='\t') {
341                         /* TODO - edit the @SQ line here to remove MD and fix LN. */
342                         /* For now just remove the @SQ line, and samtools will */
343                         /* automatically generate a minimal replacement with LN. */
344                         /* However, that discards any other tags like AS, SP, UR. */
345                         //fprintf(stderr, "[depad] Removing @SQ line\n");
346                 } else {
347                         /* Copy this line to the new header */
348                         strncat(newtext, text, end - text + 1);
349                 }
350                 text = end + 1;
351         }
352         assert (text[0]=='\0');
353         /* Check we didn't overflow the buffer */
354         assert (strlen(header->text) <= strlen(old->text));
355         if (strlen(header->text) < header->l_text) {
356                 //fprintf(stderr, "[depad] Reallocating header buffer\n");
357                 assert (newtext == header->text);
358                 newtext = calloc(strlen(header->text) + 1, 1);
359                 strcpy(newtext, header->text);
360                 free(header->text);
361                 header->text = newtext;
362                 header->l_text = strlen(newtext);
363         }
364         //fprintf(stderr, "[depad] Here is the new header (pending @SQ lines),\n\n%s\n(end)\n", header->text);
365         return header;
366 }
367
368 static int usage(int is_long_help);
369
370 int main_pad2unpad(int argc, char *argv[])
371 {
372         samfile_t *in = 0, *out = 0;
373         bam_header_t *h = 0;
374         faidx_t *fai = 0;
375         int c, is_bamin = 1, compress_level = -1, is_bamout = 1, is_long_help = 0;
376         char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0;
377         int ret=0;
378
379         /* parse command-line options */
380         strcpy(in_mode, "r"); strcpy(out_mode, "w");
381         while ((c = getopt(argc, argv, "Sso:u1T:?")) >= 0) {
382                 switch (c) {
383                 case 'S': is_bamin = 0; break;
384                 case 's': assert(compress_level == -1); is_bamout = 0; break;
385                 case 'o': fn_out = strdup(optarg); break;
386                 case 'u': assert(is_bamout == 1); compress_level = 0; break;
387                 case '1': assert(is_bamout == 1); compress_level = 1; break;
388                 case 'T': fn_ref = strdup(optarg); break;
389                 case '?': is_long_help = 1; break;
390                 default: return usage(is_long_help);
391                 }
392         }
393         if (argc == optind) return usage(is_long_help);
394
395         if (is_bamin) strcat(in_mode, "b");
396         if (is_bamout) strcat(out_mode, "b");
397         strcat(out_mode, "h");
398         if (compress_level >= 0) {
399                 char tmp[2];
400                 tmp[0] = compress_level + '0'; tmp[1] = '\0';
401                 strcat(out_mode, tmp);
402         }
403
404         // Load FASTA reference (also needed for SAM -> BAM if missing header)
405         if (fn_ref) {
406                 fn_list = samfaipath(fn_ref);
407                 fai = fai_load(fn_ref);
408         }
409         // open file handlers
410         if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) {
411                 fprintf(stderr, "[depad] failed to open \"%s\" for reading.\n", argv[optind]);
412                 ret = 1;
413                 goto depad_end;
414         }
415         if (in->header == 0) {
416                 fprintf(stderr, "[depad] failed to read the header from \"%s\".\n", argv[optind]);
417                 ret = 1;
418                 goto depad_end;
419         }
420         if (in->header->text == 0 || in->header->l_text == 0) {
421                 fprintf(stderr, "[depad] Warning - failed to read any header text from \"%s\".\n", argv[optind]);
422                 assert (0 == in->header->l_text);
423                 assert (0 == in->header->text);
424         }
425         if (fn_ref) {
426                 h = fix_header(in->header, fai);
427         } else {
428                 fprintf(stderr, "[depad] Warning - reference lengths will not be corrected without FASTA reference\n");
429                 h = in->header;
430         }
431         if ((out = samopen(fn_out? fn_out : "-", out_mode, h)) == 0) {
432                 fprintf(stderr, "[depad] failed to open \"%s\" for writing.\n", fn_out? fn_out : "standard output");
433                 ret = 1;
434                 goto depad_end;
435         }
436
437         // Do the depad
438         ret = bam_pad2unpad(in, out, fai);
439
440 depad_end:
441         // close files, free and return
442         free(fn_list); free(fn_out);
443         if (h != in->header) bam_header_destroy(h);
444         samclose(in);
445         samclose(out);
446         fai_destroy(fai);
447         return ret;
448 }
449
450 static int usage(int is_long_help)
451 {
452         fprintf(stderr, "\n");
453         fprintf(stderr, "Usage:   samtools depad <in.bam>\n\n");
454         fprintf(stderr, "Options: -s       output is SAM (default is BAM)\n");
455         fprintf(stderr, "         -S       input is SAM (default is BAM)\n");
456         fprintf(stderr, "         -u       uncompressed BAM output (can't use with -s)\n");
457         fprintf(stderr, "         -1       fast compression BAM output (can't use with -s)\n");
458         fprintf(stderr, "         -T FILE  reference sequence file [null]\n");
459         fprintf(stderr, "         -o FILE  output file name [stdout]\n");
460         fprintf(stderr, "         -?       longer help\n");
461         fprintf(stderr, "\n");
462         if (is_long_help)
463                 fprintf(stderr, "Notes:\n\
464 \n\
465   1. Requires embedded reference sequences (before the reads for that reference),\n\
466      with the future aim to also support a FASTA padded reference sequence file.\n\
467 \n\
468   2. The input padded alignment read's CIGAR strings must not use P or I operators.\n\
469 \n");
470         return 1;
471 }