]> git.donarmstrong.com Git - samtools.git/blob - padding.c
Use FASTA file when no embedded reference 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 static void replace_cigar(bam1_t *b, int n, uint32_t *cigar)
11 {
12         if (n != b->core.n_cigar) {
13                 int o = b->core.l_qname + b->core.n_cigar * 4;
14                 if (b->data_len + (n - b->core.n_cigar) * 4 > b->m_data) {
15                         b->m_data = b->data_len + (n - b->core.n_cigar) * 4;
16                         kroundup32(b->m_data);
17                         b->data = (uint8_t*)realloc(b->data, b->m_data);
18                 }
19                 memmove(b->data + b->core.l_qname + n * 4, b->data + o, b->data_len - o);
20                 memcpy(b->data + b->core.l_qname, cigar, n * 4);
21                 b->data_len += (n - b->core.n_cigar) * 4;
22                 b->core.n_cigar = n;
23         } else memcpy(b->data + b->core.l_qname, cigar, n * 4);
24 }
25
26 #define write_cigar(_c, _n, _m, _v) do { \
27                 if (_n == _m) { \
28                         _m = _m? _m<<1 : 4; \
29                         _c = (uint32_t*)realloc(_c, _m * 4); \
30                 } \
31                 _c[_n++] = (_v); \
32         } while (0)
33
34 static void unpad_seq(bam1_t *b, kstring_t *s)
35 {
36         int k, j, i;
37         int length;
38         uint32_t *cigar = bam1_cigar(b);
39         uint8_t *seq = bam1_seq(b);
40         // b->core.l_qseq gives length of the SEQ entry (including soft clips, S)
41         // We need the padded length after alignment from the CIGAR (excluding
42         // soft clips S, but including pads from CIGAR D operations)
43         length = 0;
44         for (k = 0; k < b->core.n_cigar; ++k) {
45                 int op, ol;
46                 op= bam_cigar_op(cigar[k]);
47                 ol = bam_cigar_oplen(cigar[k]);
48                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF || op == BAM_CDEL)
49                         length += ol;
50         }
51         ks_resize(s, length);
52         for (k = 0, s->l = 0, j = 0; k < b->core.n_cigar; ++k) {
53                 int op, ol;
54                 op = bam_cigar_op(cigar[k]);
55                 ol = bam_cigar_oplen(cigar[k]);
56                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
57                         for (i = 0; i < ol; ++i, ++j) s->s[s->l++] = bam1_seqi(seq, j);
58                 } else if (op == BAM_CSOFT_CLIP) {
59                         j += ol;
60                 } else if (op == BAM_CHARD_CLIP) {
61                         /* do nothing */
62                 } else if (op == BAM_CDEL) {
63                         for (i = 0; i < ol; ++i) s->s[s->l++] = 0;
64                 } else {
65                         fprintf(stderr, "[depad] ERROR: Didn't expect CIGAR op %c in read %s\n", BAM_CIGAR_STR[op], bam1_qname(b));
66                         assert(-1);
67                 }
68         }
69         assert(length == s->l);
70 }
71
72 int load_unpadded_ref(faidx_t *fai, char *ref_name, int ref_len, kstring_t *seq)
73 {
74         char base;
75         char *fai_ref = 0;
76         int fai_ref_len = 0, k;
77
78         fai_ref = fai_fetch(fai, ref_name, &fai_ref_len);
79         if (fai_ref_len != ref_len) {
80                 fprintf(stderr, "[depad] ERROR: FASTA sequence %s length %i, expected %i\n", ref_name, fai_ref_len, ref_len);
81                 free(fai_ref);
82                 return -1;
83         }
84         ks_resize(seq, ref_len);
85         seq->l = 0;
86         for (k = 0; k < ref_len; ++k) {
87                 base = fai_ref[k];
88                 if (base == '-' || base == '*') {
89                         // Map gaps to null to match unpad_seq function
90                         seq->s[seq->l++] = 0;
91                 } else {
92                         // Want to map IUPAC bases to upper case, and reject other chars
93                         int i = bam_nt16_table[(int)base];
94                         if (i == 0 || i==16) { // Equals maps to 0, anything unexpected to 16
95                                 fprintf(stderr, "[depad] ERROR: Invalid character %c (ASCII %i) in FASTA sequence %s\n", base, (int)base, ref_name);
96                                 free(fai_ref);
97                                 return -1;
98                         }
99                         seq->s[seq->l++] = bam_nt16_rev_table[i];
100                 }
101         }
102         assert(ref_len == seq->l);
103         free(fai_ref);
104         return 0;
105 }
106
107 inline int * update_posmap(int *posmap, kstring_t ref)
108 {
109         int i, k;
110         posmap = realloc(posmap, ref.m * sizeof(int));
111         for (i = k = 0; i < ref.l; ++i) {
112                 posmap[i] = k;
113                 if (ref.s[i]) ++k;
114         }
115         return posmap;
116 }
117
118 int bam_pad2unpad(samfile_t *in, samfile_t *out, faidx_t *fai)
119 {
120         bam_header_t *h;
121         bam1_t *b;
122         kstring_t r, q;
123         int r_tid = -1;
124         uint32_t *cigar2 = 0;
125         int ret = 0, n2 = 0, m2 = 0, *posmap = 0;
126
127         b = bam_init1();
128         r.l = r.m = q.l = q.m = 0; r.s = q.s = 0;
129         int read_ret;
130         h = in->header;
131         while ((read_ret = samread(in, b)) >= 0) { // read one alignment from `in'
132                 uint32_t *cigar = bam1_cigar(b);
133                 n2 = 0;
134                 if (b->core.pos == 0 && b->core.tid >= 0 && strcmp(bam1_qname(b), h->target_name[b->core.tid]) == 0) {
135                         // fprintf(stderr, "[depad] Found embedded reference %s\n", bam1_qname(b));
136                         // TODO - Check this embedded reference matches the FASTA file (if given)
137                         r_tid = b->core.tid;
138                         unpad_seq(b, &r);
139                         if (h->target_len[r_tid] != r.l) {
140                                 fprintf(stderr, "[depad] ERROR: (Padded) length of %s is %i in BAM header, but %ld in embedded reference\n", bam1_qname(b), h->target_len[r_tid], r.l);
141                                 return -1;
142                         }
143                         write_cigar(cigar2, n2, m2, bam_cigar_gen(b->core.l_qseq, BAM_CMATCH));
144                         replace_cigar(b, n2, cigar2);
145                         posmap = update_posmap(posmap, r);
146                         /*
147                         posmap = realloc(posmap, r.m * sizeof(int));
148                         for (i = k = 0; i < r.l; ++i) {
149                                 posmap[i] = k;
150                                 if (r.s[i]) ++k;
151                         }
152                         */
153                 } else if (b->core.n_cigar > 0) {
154                         int i, k, op;
155                         if (b->core.tid < 0) {
156                                 fprintf(stderr, "[depad] ERROR: Read %s has CIGAR but no RNAME\n", bam1_qname(b));
157                                 return -1;
158                         } else if (b->core.tid == r_tid) {
159                                 ; // good case, reference available
160                         } else if (fai) {
161                                 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &r)) return -1;
162                                 posmap = update_posmap(posmap, r);
163                                 r_tid = b->core.tid;
164                                 // fprintf(stderr, "[depad] Loaded %s from FASTA file\n", h->target_name[b->core.tid]);
165                         } else {                                
166                                 fprintf(stderr, "[depad] ERROR: Missing %s embedded reference sequence (and no FASTA file)\n", h->target_name[b->core.tid]);
167                                 return -1;
168                         }
169                         unpad_seq(b, &q);
170                         if (bam_cigar_op(cigar[0]) == BAM_CSOFT_CLIP) {
171                                 write_cigar(cigar2, n2, m2, cigar[0]);
172                         } else if (bam_cigar_op(cigar[0]) == BAM_CHARD_CLIP) {
173                                 write_cigar(cigar2, n2, m2, cigar[0]);
174                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[1]) == BAM_CSOFT_CLIP) {
175                                         write_cigar(cigar2, n2, m2, cigar[1]);
176                                 }
177                         }
178                         /* Determine CIGAR operator for each base in the aligned read */
179                         for (i = 0, k = b->core.pos; i < q.l; ++i, ++k)
180                                 q.s[i] = q.s[i]? (r.s[k]? BAM_CMATCH : BAM_CINS) : (r.s[k]? BAM_CDEL : BAM_CPAD);
181                         /* Include any pads if starts with an insert */
182                         if (q.s[0] == BAM_CINS) {
183                                 for (k = 0; k+1 < b->core.pos && !r.s[b->core.pos - k - 1]; ++k);
184                                 if (k) write_cigar(cigar2, n2, m2, bam_cigar_gen(k, BAM_CPAD));
185                         }
186                         /* Count consecutive CIGAR operators to turn into a CIGAR string */
187                         for (i = k = 1, op = q.s[0]; i < q.l; ++i) {
188                                 if (op != q.s[i]) {
189                                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
190                                         op = q.s[i]; k = 1;
191                                 } else ++k;
192                         }
193                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
194                         if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CSOFT_CLIP) {
195                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
196                         } else if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CHARD_CLIP) {
197                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[b->core.n_cigar-2]) == BAM_CSOFT_CLIP) {
198                                         write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-2]);
199                                 }
200                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
201                         }
202                         /* Remove redundant P operators between M/X/=/D operators, e.g. 5M2P10M -> 15M */
203                         int pre_op, post_op;
204                         for (i = 2; i < n2; ++i)
205                                 if (bam_cigar_op(cigar2[i-1]) == BAM_CPAD) {
206                                         pre_op = bam_cigar_op(cigar2[i-2]);
207                                         post_op = bam_cigar_op(cigar2[i]);
208                                         /* Note don't need to check for X/= as code above will use M only */
209                                         if ((pre_op == BAM_CMATCH || pre_op == BAM_CDEL) && (post_op == BAM_CMATCH || post_op == BAM_CDEL)) {
210                                                 /* This is a redundant P operator */
211                                                 cigar2[i-1] = 0; // i.e. 0M
212                                                 /* If had same operator either side, combine them in post_op */
213                                                 if (pre_op == post_op) {
214                                                         /* If CIGAR M, could treat as simple integers since BAM_CMATCH is zero*/
215                                                         cigar2[i] = bam_cigar_gen(bam_cigar_oplen(cigar2[i-2]) + bam_cigar_oplen(cigar2[i]), post_op);
216                                                         cigar2[i-2] = 0; // i.e. 0M
217                                                 }
218                                         }
219                                 }
220                         /* Remove the zero'd operators (0M) */
221                         for (i = k = 0; i < n2; ++i)
222                                 if (cigar2[i]) cigar2[k++] = cigar2[i];
223                         n2 = k;
224                         replace_cigar(b, n2, cigar2);
225                         b->core.pos = posmap[b->core.pos];
226                 }
227                 samwrite(out, b);
228         }
229         if (read_ret < -1) {
230                 fprintf(stderr, "[depad] truncated file.\n");
231                 ret = 1;
232         }
233         free(r.s); free(q.s); free(posmap);
234         bam_destroy1(b);
235         return ret;
236 }
237
238 static int usage(int is_long_help);
239
240 int main_pad2unpad(int argc, char *argv[])
241 {
242         samfile_t *in = 0, *out = 0;
243         faidx_t *fai = 0;
244         int c, is_bamin = 1, compress_level = -1, is_bamout = 1, is_long_help = 0;
245         char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0;
246         int ret=0;
247
248         /* parse command-line options */
249         strcpy(in_mode, "r"); strcpy(out_mode, "w");
250         while ((c = getopt(argc, argv, "Sso:u1T:?")) >= 0) {
251                 switch (c) {
252                 case 'S': is_bamin = 0; break;
253                 case 's': assert(compress_level == -1); is_bamout = 0; break;
254                 case 'o': fn_out = strdup(optarg); break;
255                 case 'u': assert(is_bamout == 1); compress_level = 0; break;
256                 case '1': assert(is_bamout == 1); compress_level = 1; break;
257                 case 'T': fn_ref = strdup(optarg); is_bamin = 0; break;
258                 case '?': is_long_help = 1; break;
259                 default: return usage(is_long_help);
260                 }
261         }
262         if (argc == optind) return usage(is_long_help);
263
264         if (is_bamin) strcat(in_mode, "b");
265         if (is_bamout) strcat(out_mode, "b");
266         strcat(out_mode, "h");
267         if (compress_level >= 0) {
268                 char tmp[2];
269                 tmp[0] = compress_level + '0'; tmp[1] = '\0';
270                 strcat(out_mode, tmp);
271         }
272
273         // Load FASTA reference (also needed for SAM -> BAM if missing header)
274         if (fn_ref) {
275                 fn_list = samfaipath(fn_ref);
276                 fai = fai_load(fn_ref);
277         }
278         // open file handlers
279         if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) {
280                 fprintf(stderr, "[depad] fail to open \"%s\" for reading.\n", argv[optind]);
281                 ret = 1;
282                 goto depad_end;
283         }
284         if (in->header == 0) {
285                 fprintf(stderr, "[depad] fail to read the header from \"%s\".\n", argv[optind]);
286                 ret = 1;
287                 goto depad_end;
288         }
289         /* TODO - The reference sequence lengths in the BAM + SAM headers should be updated */
290         if ((out = samopen(fn_out? fn_out : "-", out_mode, in->header)) == 0) {
291                 fprintf(stderr, "[depad] fail to open \"%s\" for writing.\n", fn_out? fn_out : "standard output");
292                 ret = 1;
293                 goto depad_end;
294         }
295
296         // Do the depad
297         ret = bam_pad2unpad(in, out, fai);
298
299 depad_end:
300         // close files, free and return
301         free(fn_list); free(fn_out);
302         samclose(in);
303         samclose(out);
304         fai_destroy(fai);
305         return ret;
306 }
307
308 static int usage(int is_long_help)
309 {
310         fprintf(stderr, "\n");
311         fprintf(stderr, "Usage:   samtools depad <in.bam>\n\n");
312         fprintf(stderr, "Options: -s       output is SAM (default is BAM)\n");
313         fprintf(stderr, "         -S       input is SAM (default is BAM)\n");
314         fprintf(stderr, "         -u       uncompressed BAM output (can't use with -s)\n");
315         fprintf(stderr, "         -1       fast compression BAM output (can't use with -s)\n");
316         fprintf(stderr, "         -T FILE  reference sequence file [null]\n");
317         fprintf(stderr, "         -o FILE  output file name [stdout]\n");
318         fprintf(stderr, "         -?       longer help\n");
319         fprintf(stderr, "\n");
320         if (is_long_help)
321                 fprintf(stderr, "Notes:\n\
322 \n\
323   1. Requires embedded reference sequences (before the reads for that reference),\n\
324      with the future aim to also support a FASTA padded reference sequence file.\n\
325 \n\
326   2. The input padded alignment read's CIGAR strings must not use P or I operators.\n\
327 \n");
328         return 1;
329 }