]> git.donarmstrong.com Git - samtools.git/blob - padding.c
Minimal 'samtools depad' command line parsing (-? only so far)
[samtools.git] / padding.c
1 #include <string.h>
2 #include <assert.h>
3 #include <unistd.h>
4 #include "kstring.h"
5 #include "bam.h"
6
7 static void replace_cigar(bam1_t *b, int n, uint32_t *cigar)
8 {
9         if (n != b->core.n_cigar) {
10                 int o = b->core.l_qname + b->core.n_cigar * 4;
11                 if (b->data_len + (n - b->core.n_cigar) * 4 > b->m_data) {
12                         b->m_data = b->data_len + (n - b->core.n_cigar) * 4;
13                         kroundup32(b->m_data);
14                         b->data = (uint8_t*)realloc(b->data, b->m_data);
15                 }
16                 memmove(b->data + b->core.l_qname + n * 4, b->data + o, b->data_len - o);
17                 memcpy(b->data + b->core.l_qname, cigar, n * 4);
18                 b->data_len += (n - b->core.n_cigar) * 4;
19                 b->core.n_cigar = n;
20         } else memcpy(b->data + b->core.l_qname, cigar, n * 4);
21 }
22
23 #define write_cigar(_c, _n, _m, _v) do { \
24                 if (_n == _m) { \
25                         _m = _m? _m<<1 : 4; \
26                         _c = (uint32_t*)realloc(_c, _m * 4); \
27                 } \
28                 _c[_n++] = (_v); \
29         } while (0)
30
31 static void unpad_seq(bam1_t *b, kstring_t *s)
32 {
33         int k, j, i;
34         int length;
35         uint32_t *cigar = bam1_cigar(b);
36         uint8_t *seq = bam1_seq(b);
37         // b->core.l_qseq gives length of the SEQ entry (including soft clips, S)
38         // We need the padded length after alignment from the CIGAR (excluding
39         // soft clips S, but including pads from CIGAR D operations)
40         length = 0;
41         for (k = 0; k < b->core.n_cigar; ++k) {
42                 int op, ol;
43                 op= bam_cigar_op(cigar[k]);
44                 ol = bam_cigar_oplen(cigar[k]);
45                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF || op == BAM_CDEL)
46                         length += ol;
47         }
48         ks_resize(s, length);
49         for (k = 0, s->l = 0, j = 0; k < b->core.n_cigar; ++k) {
50                 int op, ol;
51                 op = bam_cigar_op(cigar[k]);
52                 ol = bam_cigar_oplen(cigar[k]);
53                 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
54                         for (i = 0; i < ol; ++i, ++j) s->s[s->l++] = bam1_seqi(seq, j);
55                 } else if (op == BAM_CSOFT_CLIP) {
56                         j += ol;
57                 } else if (op == BAM_CHARD_CLIP) {
58                         /* do nothing */
59                 } else if (op == BAM_CDEL) {
60                         for (i = 0; i < ol; ++i) s->s[s->l++] = 0;
61                 } else {
62                         fprintf(stderr, "[depad] ERROR: Didn't expect CIGAR op %c in read %s\n", BAM_CIGAR_STR[op], bam1_qname(b));
63                         assert(-1);
64                 }
65         }
66         assert(length == s->l);
67 }
68
69 int bam_pad2unpad(bamFile in, bamFile out)
70 {
71         bam_header_t *h;
72         bam1_t *b;
73         kstring_t r, q;
74         int r_tid = -1;
75         uint32_t *cigar2 = 0;
76         int n2 = 0, m2 = 0, *posmap = 0;
77
78         h = bam_header_read(in);
79         /* TODO - The reference sequence lengths in the BAM + SAM headers should be updated */
80         bam_header_write(out, h);
81         b = bam_init1();
82         r.l = r.m = q.l = q.m = 0; r.s = q.s = 0;
83         while (bam_read1(in, b) >= 0) {
84                 uint32_t *cigar = bam1_cigar(b);
85                 n2 = 0;
86                 if (b->core.pos == 0 && b->core.tid >= 0 && strcmp(bam1_qname(b), h->target_name[b->core.tid]) == 0) {
87                         int i, k;
88                         /*
89                         fprintf(stderr, "[depad] Found embedded reference %s\n", bam1_qname(b));
90                         */
91                         r_tid = b->core.tid;
92                         unpad_seq(b, &r);
93                         if (h->target_len[r_tid] != r.l) {
94                                 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);
95                                 return -1;
96                         }
97                         write_cigar(cigar2, n2, m2, bam_cigar_gen(b->core.l_qseq, BAM_CMATCH));
98                         replace_cigar(b, n2, cigar2);
99                         posmap = realloc(posmap, r.m * sizeof(int));
100                         for (i = k = 0; i < r.l; ++i) {
101                                 posmap[i] = k;
102                                 if (r.s[i]) ++k;
103                         }
104                 } else if (b->core.n_cigar > 0) {
105                         int i, k, op;
106                         if (b->core.tid < 0) {
107                                 fprintf(stderr, "[depad] ERROR: Read %s has CIGAR but no RNAME\n", bam1_qname(b));
108                                 return -1;
109                         } else if (b->core.tid != r_tid) {
110                                 fprintf(stderr, "[depad] ERROR: Missing %s embedded reference sequence\n", h->target_name[b->core.tid]);
111                                 return -1;
112                         }
113                         unpad_seq(b, &q);
114                         if (bam_cigar_op(cigar[0]) == BAM_CSOFT_CLIP) {
115                                 write_cigar(cigar2, n2, m2, cigar[0]);
116                         } else if (bam_cigar_op(cigar[0]) == BAM_CHARD_CLIP) {
117                                 write_cigar(cigar2, n2, m2, cigar[0]);
118                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[1]) == BAM_CSOFT_CLIP) {
119                                         write_cigar(cigar2, n2, m2, cigar[1]);
120                                 }
121                         }
122                         /* Determine CIGAR operator for each base in the aligned read */
123                         for (i = 0, k = b->core.pos; i < q.l; ++i, ++k)
124                                 q.s[i] = q.s[i]? (r.s[k]? BAM_CMATCH : BAM_CINS) : (r.s[k]? BAM_CDEL : BAM_CPAD);
125                         /* Include any pads if starts with an insert */
126                         if (q.s[0] == BAM_CINS) {
127                                 for (k = 0; k+1 < b->core.pos && !r.s[b->core.pos - k - 1]; ++k);
128                                 if (k) write_cigar(cigar2, n2, m2, bam_cigar_gen(k, BAM_CPAD));
129                         }
130                         /* Count consecutive CIGAR operators to turn into a CIGAR string */
131                         for (i = k = 1, op = q.s[0]; i < q.l; ++i) {
132                                 if (op != q.s[i]) {
133                                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
134                                         op = q.s[i]; k = 1;
135                                 } else ++k;
136                         }
137                         write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
138                         if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CSOFT_CLIP) {
139                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
140                         } else if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CHARD_CLIP) {
141                                 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[b->core.n_cigar-2]) == BAM_CSOFT_CLIP) {
142                                         write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-2]);
143                                 }
144                                 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
145                         }
146                         /* Remove redundant P operators between M/X/=/D operators, e.g. 5M2P10M -> 15M */
147                         int pre_op, post_op;
148                         for (i = 2; i < n2; ++i)
149                                 if (bam_cigar_op(cigar2[i-1]) == BAM_CPAD) {
150                                         pre_op = bam_cigar_op(cigar2[i-2]);
151                                         post_op = bam_cigar_op(cigar2[i]);
152                                         /* Note don't need to check for X/= as code above will use M only */
153                                         if ((pre_op == BAM_CMATCH || pre_op == BAM_CDEL) && (post_op == BAM_CMATCH || post_op == BAM_CDEL)) {
154                                                 /* This is a redundant P operator */
155                                                 cigar2[i-1] = 0; // i.e. 0M
156                                                 /* If had same operator either side, combine them in post_op */
157                                                 if (pre_op == post_op) {
158                                                         /* If CIGAR M, could treat as simple integers since BAM_CMATCH is zero*/
159                                                         cigar2[i] = bam_cigar_gen(bam_cigar_oplen(cigar2[i-2]) + bam_cigar_oplen(cigar2[i]), post_op);
160                                                         cigar2[i-2] = 0; // i.e. 0M
161                                                 }
162                                         }
163                                 }
164                         /* Remove the zero'd operators (0M) */
165                         for (i = k = 0; i < n2; ++i)
166                                 if (cigar2[i]) cigar2[k++] = cigar2[i];
167                         n2 = k;
168                         replace_cigar(b, n2, cigar2);
169                         b->core.pos = posmap[b->core.pos];
170                 }
171                 bam_write1(out, b);
172         }
173         free(r.s); free(q.s); free(posmap);
174         bam_destroy1(b);
175         bam_header_destroy(h);
176         return 0;
177 }
178
179 static int usage(int is_long_help);
180
181 int main_pad2unpad(int argc, char *argv[])
182 {
183         bamFile in, out;
184         int c, is_long_help = 0;
185         int result=0;
186
187         /* parse command-line options */
188         while ((c = getopt(argc, argv, "?")) >= 0) {
189                 switch (c) {
190                 case '?': is_long_help = 1; break;
191                 default: return usage(is_long_help);
192                 }
193         }
194         if (argc == optind) return usage(is_long_help);
195
196         in = strcmp(argv[1], "-")? bam_open(argv[1], "r") : bam_dopen(fileno(stdin), "r");
197         out = bam_dopen(fileno(stdout), "w");
198         result = bam_pad2unpad(in, out);
199         bam_close(in); bam_close(out);
200         return result;
201 }
202
203 static int usage(int is_long_help)
204 {
205         fprintf(stderr, "\n");
206         fprintf(stderr, "Usage:   samtools depad <in.bam>\n\n");
207         fprintf(stderr, "Options: -?       longer help\n");
208         //TODO - These are the arguments I think make sense to support:
209         //fprintf(stderr, "Usage:   samtools depad [options] <in.bam>|<in.sam>\n\n");
210         //fprintf(stderr, "Options: -b       output BAM\n");
211         //fprintf(stderr, "         -S       input is SAM\n");
212         //fprintf(stderr, "         -u       uncompressed BAM output (force -b)\n");
213         //fprintf(stderr, "         -1       fast compression (force -b)\n");
214         //fprintf(stderr, "         -@ INT   number of BAM compression threads [0]\n");
215         //fprintf(stderr, "         -T FILE  reference sequence file (force -S) [null]\n");
216         //fprintf(stderr, "         -o FILE  output file name [stdout]\n");
217         //fprintf(stderr, "         -?       longer help\n");
218         fprintf(stderr, "\n");
219         if (is_long_help)
220                 fprintf(stderr, "Notes:\n\
221 \n\
222   1. Requires embedded reference sequences (before the reads for that reference).\n\
223 \n\
224   2. The input padded alignment read's CIGAR strings must not use P or I operators.\n\
225 \n");
226         return 1;
227 }