]> git.donarmstrong.com Git - samtools.git/blob - padding.c
9f8bca7e0fe13992aba557b176090936a34cb34a
[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 bam_pad2unpad(samfile_t *in, samfile_t *out, faidx_t *fai)
73 {
74         bam_header_t *h;
75         bam1_t *b;
76         kstring_t r, q;
77         int r_tid = -1;
78         uint32_t *cigar2 = 0;
79         int ret = 0, n2 = 0, m2 = 0, *posmap = 0;
80
81         b = bam_init1();
82         r.l = r.m = q.l = q.m = 0; r.s = q.s = 0;
83         int read_ret;
84         h = in->header;
85         while ((read_ret = samread(in, b)) >= 0) { // read one alignment from `in'
86                 uint32_t *cigar = bam1_cigar(b);
87                 n2 = 0;
88                 if (b->core.pos == 0 && b->core.tid >= 0 && strcmp(bam1_qname(b), h->target_name[b->core.tid]) == 0) {
89                         int i, k;
90                         /*
91                         fprintf(stderr, "[depad] Found embedded reference %s\n", bam1_qname(b));
92                         */
93                         r_tid = b->core.tid;
94                         unpad_seq(b, &r);
95                         if (h->target_len[r_tid] != r.l) {
96                                 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);
97                                 return -1;
98                         }
99                         write_cigar(cigar2, n2, m2, bam_cigar_gen(b->core.l_qseq, BAM_CMATCH));
100                         replace_cigar(b, n2, cigar2);
101                         posmap = realloc(posmap, r.m * sizeof(int));
102                         for (i = k = 0; i < r.l; ++i) {
103                                 posmap[i] = k;
104                                 if (r.s[i]) ++k;
105                         }
106                 } else if (b->core.n_cigar > 0) {
107                         int i, k, op;
108                         if (b->core.tid < 0) {
109                                 fprintf(stderr, "[depad] ERROR: Read %s has CIGAR but no RNAME\n", bam1_qname(b));
110                                 return -1;
111                         } else if (b->core.tid != r_tid) {
112                                 fprintf(stderr, "[depad] ERROR: Missing %s embedded reference sequence\n", h->target_name[b->core.tid]);
113                                 return -1;
114                         }
115                         unpad_seq(b, &q);
116                         if (bam_cigar_op(cigar[0]) == BAM_CSOFT_CLIP) {
117                                 write_cigar(cigar2, n2, m2, cigar[0]);
118                         } else if (bam_cigar_op(cigar[0]) == BAM_CHARD_CLIP) {
119                                 write_cigar(cigar2, n2, m2, cigar[0]);
120                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[1]) == BAM_CSOFT_CLIP) {
121                                         write_cigar(cigar2, n2, m2, cigar[1]);
122                                 }
123                         }
124                         /* Determine CIGAR operator for each base in the aligned read */
125                         for (i = 0, k = b->core.pos; i < q.l; ++i, ++k)
126                                 q.s[i] = q.s[i]? (r.s[k]? BAM_CMATCH : BAM_CINS) : (r.s[k]? BAM_CDEL : BAM_CPAD);
127                         /* Include any pads if starts with an insert */
128                         if (q.s[0] == BAM_CINS) {
129                                 for (k = 0; k+1 < b->core.pos && !r.s[b->core.pos - k - 1]; ++k);
130                                 if (k) write_cigar(cigar2, n2, m2, bam_cigar_gen(k, BAM_CPAD));
131                         }
132                         /* Count consecutive CIGAR operators to turn into a CIGAR string */
133                         for (i = k = 1, op = q.s[0]; i < q.l; ++i) {
134                                 if (op != q.s[i]) {
135                                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
136                                         op = q.s[i]; k = 1;
137                                 } else ++k;
138                         }
139                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
140                         if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CSOFT_CLIP) {
141                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
142                         } else if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CHARD_CLIP) {
143                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[b->core.n_cigar-2]) == BAM_CSOFT_CLIP) {
144                                         write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-2]);
145                                 }
146                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
147                         }
148                         /* Remove redundant P operators between M/X/=/D operators, e.g. 5M2P10M -> 15M */
149                         int pre_op, post_op;
150                         for (i = 2; i < n2; ++i)
151                                 if (bam_cigar_op(cigar2[i-1]) == BAM_CPAD) {
152                                         pre_op = bam_cigar_op(cigar2[i-2]);
153                                         post_op = bam_cigar_op(cigar2[i]);
154                                         /* Note don't need to check for X/= as code above will use M only */
155                                         if ((pre_op == BAM_CMATCH || pre_op == BAM_CDEL) && (post_op == BAM_CMATCH || post_op == BAM_CDEL)) {
156                                                 /* This is a redundant P operator */
157                                                 cigar2[i-1] = 0; // i.e. 0M
158                                                 /* If had same operator either side, combine them in post_op */
159                                                 if (pre_op == post_op) {
160                                                         /* If CIGAR M, could treat as simple integers since BAM_CMATCH is zero*/
161                                                         cigar2[i] = bam_cigar_gen(bam_cigar_oplen(cigar2[i-2]) + bam_cigar_oplen(cigar2[i]), post_op);
162                                                         cigar2[i-2] = 0; // i.e. 0M
163                                                 }
164                                         }
165                                 }
166                         /* Remove the zero'd operators (0M) */
167                         for (i = k = 0; i < n2; ++i)
168                                 if (cigar2[i]) cigar2[k++] = cigar2[i];
169                         n2 = k;
170                         replace_cigar(b, n2, cigar2);
171                         b->core.pos = posmap[b->core.pos];
172                 }
173                 samwrite(out, b);
174         }
175         if (read_ret < -1) {
176                 fprintf(stderr, "[depad] truncated file.\n");
177                 ret = 1;
178         }
179         free(r.s); free(q.s); free(posmap);
180         bam_destroy1(b);
181         return ret;
182 }
183
184 static int usage(int is_long_help);
185
186 int main_pad2unpad(int argc, char *argv[])
187 {
188         samfile_t *in = 0, *out = 0;
189         faidx_t *fai;
190         int c, is_bamin = 1, compress_level = -1, is_bamout = 1, is_long_help = 0;
191         char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0;
192         int ret=0;
193
194         /* parse command-line options */
195         strcpy(in_mode, "r"); strcpy(out_mode, "w");
196         while ((c = getopt(argc, argv, "Sso:u1T:?")) >= 0) {
197                 switch (c) {
198                 case 'S': is_bamin = 0; break;
199                 case 's': assert(compress_level == -1); is_bamout = 0; break;
200                 case 'o': fn_out = strdup(optarg); break;
201                 case 'u': assert(is_bamout == 1); compress_level = 0; break;
202                 case '1': assert(is_bamout == 1); compress_level = 1; break;
203                 case 'T': fn_ref = strdup(optarg); is_bamin = 0; break;
204                 case '?': is_long_help = 1; break;
205                 default: return usage(is_long_help);
206                 }
207         }
208         if (argc == optind) return usage(is_long_help);
209
210         if (is_bamin) strcat(in_mode, "b");
211         if (is_bamout) strcat(out_mode, "b");
212         strcat(out_mode, "h");
213         if (compress_level >= 0) {
214                 char tmp[2];
215                 tmp[0] = compress_level + '0'; tmp[1] = '\0';
216                 strcat(out_mode, tmp);
217         }
218
219         // Load FASTA reference (also needed for SAM -> BAM if missing header)
220         if (fn_ref) {
221                 fn_list = samfaipath(fn_ref);
222                 fai = fai_load(fn_ref);
223         }
224         // open file handlers
225         if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) {
226                 fprintf(stderr, "[depad] fail to open \"%s\" for reading.\n", argv[optind]);
227                 ret = 1;
228                 goto depad_end;
229         }
230         if (in->header == 0) {
231                 fprintf(stderr, "[depad] fail to read the header from \"%s\".\n", argv[optind]);
232                 ret = 1;
233                 goto depad_end;
234         }
235         /* TODO - The reference sequence lengths in the BAM + SAM headers should be updated */
236         if ((out = samopen(fn_out? fn_out : "-", out_mode, in->header)) == 0) {
237                 fprintf(stderr, "[depad] fail to open \"%s\" for writing.\n", fn_out? fn_out : "standard output");
238                 ret = 1;
239                 goto depad_end;
240         }
241
242         // Do the depad
243         ret = bam_pad2unpad(in, out, fai);
244
245 depad_end:
246         // close files, free and return
247         free(fn_list); free(fn_out);
248         samclose(in);
249         samclose(out);
250         fai_destroy(fai);
251         return ret;
252 }
253
254 static int usage(int is_long_help)
255 {
256         fprintf(stderr, "\n");
257         fprintf(stderr, "Usage:   samtools depad <in.bam>\n\n");
258         fprintf(stderr, "Options: -s       output is SAM (default is BAM)\n");
259         fprintf(stderr, "         -S       input is SAM (default is BAM)\n");
260         fprintf(stderr, "         -u       uncompressed BAM output (can't use with -s)\n");
261         fprintf(stderr, "         -1       fast compression BAM output (can't use with -s)\n");
262         fprintf(stderr, "         -T FILE  reference sequence file [null]\n");
263         fprintf(stderr, "         -o FILE  output file name [stdout]\n");
264         fprintf(stderr, "         -?       longer help\n");
265         fprintf(stderr, "\n");
266         if (is_long_help)
267                 fprintf(stderr, "Notes:\n\
268 \n\
269   1. Requires embedded reference sequences (before the reads for that reference),\n\
270      with the future aim to also support a FASTA padded reference sequence file.\n\
271 \n\
272   2. The input padded alignment read's CIGAR strings must not use P or I operators.\n\
273 \n");
274         return 1;
275 }