5 #include "sam_header.h"
10 static void replace_cigar(bam1_t *b, int n, uint32_t *cigar)
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);
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;
23 } else memcpy(b->data + b->core.l_qname, cigar, n * 4);
26 #define write_cigar(_c, _n, _m, _v) do { \
29 _c = (uint32_t*)realloc(_c, _m * 4); \
34 static void unpad_seq(bam1_t *b, kstring_t *s)
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)
44 for (k = 0; k < b->core.n_cigar; ++k) {
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)
52 for (k = 0, s->l = 0, j = 0; k < b->core.n_cigar; ++k) {
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) {
60 } else if (op == BAM_CHARD_CLIP) {
62 } else if (op == BAM_CDEL) {
63 for (i = 0; i < ol; ++i) s->s[s->l++] = 0;
65 fprintf(stderr, "[depad] ERROR: Didn't expect CIGAR op %c in read %s\n", BAM_CIGAR_STR[op], bam1_qname(b));
69 assert(length == s->l);
72 int load_unpadded_ref(faidx_t *fai, char *ref_name, int ref_len, kstring_t *seq)
76 int fai_ref_len = 0, k;
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);
84 ks_resize(seq, ref_len);
86 for (k = 0; k < ref_len; ++k) {
88 if (base == '-' || base == '*') {
89 // Map gaps to null to match unpad_seq function
92 int i = bam_nt16_table[(int)base];
93 if (i == 0 || i==16) { // Equals maps to 0, anything unexpected to 16
94 fprintf(stderr, "[depad] ERROR: Invalid character %c (ASCII %i) in FASTA sequence %s\n", base, (int)base, ref_name);
101 assert(ref_len == seq->l);
106 inline int * update_posmap(int *posmap, kstring_t ref)
109 posmap = realloc(posmap, ref.m * sizeof(int));
110 for (i = k = 0; i < ref.l; ++i) {
117 int bam_pad2unpad(samfile_t *in, samfile_t *out, faidx_t *fai)
123 uint32_t *cigar2 = 0;
124 int ret = 0, n2 = 0, m2 = 0, *posmap = 0;
127 r.l = r.m = q.l = q.m = 0; r.s = q.s = 0;
130 while ((read_ret = samread(in, b)) >= 0) { // read one alignment from `in'
131 uint32_t *cigar = bam1_cigar(b);
133 if (b->core.pos == 0 && b->core.tid >= 0 && strcmp(bam1_qname(b), h->target_name[b->core.tid]) == 0) {
134 // fprintf(stderr, "[depad] Found embedded reference %s\n", bam1_qname(b));
137 if (h->target_len[r_tid] != r.l) {
138 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);
142 // Check the embedded reference matches the FASTA file
143 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &q)) return -1;
146 for (i = 0; i < r.l; ++i) {
147 if (r.s[i] != q.s[i]) {
148 // Show gaps as ASCII 45
149 fprintf(stderr, "[depad] ERROR: Embedded sequence and reference FASTA don't match for %s base %i, '%c' vs '%c'\n",
150 h->target_name[b->core.tid], i+1,
151 r.s[i] ? bam_nt16_rev_table[r.s[i]] : 45,
152 q.s[i] ? bam_nt16_rev_table[q.s[i]] : 45);
157 write_cigar(cigar2, n2, m2, bam_cigar_gen(b->core.l_qseq, BAM_CMATCH));
158 replace_cigar(b, n2, cigar2);
159 posmap = update_posmap(posmap, r);
160 } else if (b->core.n_cigar > 0) {
162 if (b->core.tid < 0) {
163 fprintf(stderr, "[depad] ERROR: Read %s has CIGAR but no RNAME\n", bam1_qname(b));
165 } else if (b->core.tid == r_tid) {
166 ; // good case, reference available
168 if (load_unpadded_ref(fai, h->target_name[b->core.tid], h->target_len[b->core.tid], &r)) return -1;
169 posmap = update_posmap(posmap, r);
171 // fprintf(stderr, "[depad] Loaded %s from FASTA file\n", h->target_name[b->core.tid]);
173 fprintf(stderr, "[depad] ERROR: Missing %s embedded reference sequence (and no FASTA file)\n", h->target_name[b->core.tid]);
177 if (bam_cigar_op(cigar[0]) == BAM_CSOFT_CLIP) {
178 write_cigar(cigar2, n2, m2, cigar[0]);
179 } else if (bam_cigar_op(cigar[0]) == BAM_CHARD_CLIP) {
180 write_cigar(cigar2, n2, m2, cigar[0]);
181 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[1]) == BAM_CSOFT_CLIP) {
182 write_cigar(cigar2, n2, m2, cigar[1]);
185 /* Determine CIGAR operator for each base in the aligned read */
186 for (i = 0, k = b->core.pos; i < q.l; ++i, ++k)
187 q.s[i] = q.s[i]? (r.s[k]? BAM_CMATCH : BAM_CINS) : (r.s[k]? BAM_CDEL : BAM_CPAD);
188 /* Include any pads if starts with an insert */
189 if (q.s[0] == BAM_CINS) {
190 for (k = 0; k+1 < b->core.pos && !r.s[b->core.pos - k - 1]; ++k);
191 if (k) write_cigar(cigar2, n2, m2, bam_cigar_gen(k, BAM_CPAD));
193 /* Count consecutive CIGAR operators to turn into a CIGAR string */
194 for (i = k = 1, op = q.s[0]; i < q.l; ++i) {
196 write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
200 write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
201 if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CSOFT_CLIP) {
202 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
203 } else if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CHARD_CLIP) {
204 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[b->core.n_cigar-2]) == BAM_CSOFT_CLIP) {
205 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-2]);
207 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
209 /* Remove redundant P operators between M/X/=/D operators, e.g. 5M2P10M -> 15M */
211 for (i = 2; i < n2; ++i)
212 if (bam_cigar_op(cigar2[i-1]) == BAM_CPAD) {
213 pre_op = bam_cigar_op(cigar2[i-2]);
214 post_op = bam_cigar_op(cigar2[i]);
215 /* Note don't need to check for X/= as code above will use M only */
216 if ((pre_op == BAM_CMATCH || pre_op == BAM_CDEL) && (post_op == BAM_CMATCH || post_op == BAM_CDEL)) {
217 /* This is a redundant P operator */
218 cigar2[i-1] = 0; // i.e. 0M
219 /* If had same operator either side, combine them in post_op */
220 if (pre_op == post_op) {
221 /* If CIGAR M, could treat as simple integers since BAM_CMATCH is zero*/
222 cigar2[i] = bam_cigar_gen(bam_cigar_oplen(cigar2[i-2]) + bam_cigar_oplen(cigar2[i]), post_op);
223 cigar2[i-2] = 0; // i.e. 0M
227 /* Remove the zero'd operators (0M) */
228 for (i = k = 0; i < n2; ++i)
229 if (cigar2[i]) cigar2[k++] = cigar2[i];
231 replace_cigar(b, n2, cigar2);
232 b->core.pos = posmap[b->core.pos];
237 fprintf(stderr, "[depad] truncated file.\n");
240 free(r.s); free(q.s); free(posmap);
245 static int usage(int is_long_help);
247 int main_pad2unpad(int argc, char *argv[])
249 samfile_t *in = 0, *out = 0;
251 int c, is_bamin = 1, compress_level = -1, is_bamout = 1, is_long_help = 0;
252 char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0;
255 /* parse command-line options */
256 strcpy(in_mode, "r"); strcpy(out_mode, "w");
257 while ((c = getopt(argc, argv, "Sso:u1T:?")) >= 0) {
259 case 'S': is_bamin = 0; break;
260 case 's': assert(compress_level == -1); is_bamout = 0; break;
261 case 'o': fn_out = strdup(optarg); break;
262 case 'u': assert(is_bamout == 1); compress_level = 0; break;
263 case '1': assert(is_bamout == 1); compress_level = 1; break;
264 case 'T': fn_ref = strdup(optarg); break;
265 case '?': is_long_help = 1; break;
266 default: return usage(is_long_help);
269 if (argc == optind) return usage(is_long_help);
271 if (is_bamin) strcat(in_mode, "b");
272 if (is_bamout) strcat(out_mode, "b");
273 strcat(out_mode, "h");
274 if (compress_level >= 0) {
276 tmp[0] = compress_level + '0'; tmp[1] = '\0';
277 strcat(out_mode, tmp);
280 // Load FASTA reference (also needed for SAM -> BAM if missing header)
282 fn_list = samfaipath(fn_ref);
283 fai = fai_load(fn_ref);
285 // open file handlers
286 if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) {
287 fprintf(stderr, "[depad] failed to open \"%s\" for reading.\n", argv[optind]);
291 if (in->header == 0) {
292 fprintf(stderr, "[depad] failed to read the header from \"%s\".\n", argv[optind]);
296 /* TODO - The reference sequence lengths in the BAM + SAM headers should be updated */
297 if ((out = samopen(fn_out? fn_out : "-", out_mode, in->header)) == 0) {
298 fprintf(stderr, "[depad] failed to open \"%s\" for writing.\n", fn_out? fn_out : "standard output");
304 ret = bam_pad2unpad(in, out, fai);
307 // close files, free and return
308 free(fn_list); free(fn_out);
315 static int usage(int is_long_help)
317 fprintf(stderr, "\n");
318 fprintf(stderr, "Usage: samtools depad <in.bam>\n\n");
319 fprintf(stderr, "Options: -s output is SAM (default is BAM)\n");
320 fprintf(stderr, " -S input is SAM (default is BAM)\n");
321 fprintf(stderr, " -u uncompressed BAM output (can't use with -s)\n");
322 fprintf(stderr, " -1 fast compression BAM output (can't use with -s)\n");
323 fprintf(stderr, " -T FILE reference sequence file [null]\n");
324 fprintf(stderr, " -o FILE output file name [stdout]\n");
325 fprintf(stderr, " -? longer help\n");
326 fprintf(stderr, "\n");
328 fprintf(stderr, "Notes:\n\
330 1. Requires embedded reference sequences (before the reads for that reference),\n\
331 with the future aim to also support a FASTA padded reference sequence file.\n\
333 2. The input padded alignment read's CIGAR strings must not use P or I operators.\n\