6 static void replace_cigar(bam1_t *b, int n, uint32_t *cigar)
8 if (n != b->core.n_cigar) {
9 int o = b->core.l_qname + b->core.n_cigar * 4;
10 if (b->data_len + (n - b->core.n_cigar) * 4 > b->m_data) {
11 b->m_data = b->data_len + (n - b->core.n_cigar) * 4;
12 kroundup32(b->m_data);
13 b->data = (uint8_t*)realloc(b->data, b->m_data);
15 memmove(b->data + b->core.l_qname + n * 4, b->data + o, b->data_len - o);
16 memcpy(b->data + b->core.l_qname, cigar, n * 4);
17 b->data_len += (n - b->core.n_cigar) * 4;
19 } else memcpy(b->data + b->core.l_qname, cigar, n * 4);
22 #define write_cigar(_c, _n, _m, _v) do { \
25 _c = (uint32_t*)realloc(_c, _m * 4); \
30 static void unpad_seq(bam1_t *b, kstring_t *s)
34 uint32_t *cigar = bam1_cigar(b);
35 uint8_t *seq = bam1_seq(b);
36 // b->core.l_qseq gives length of the SEQ entry (including soft clips, S)
37 // We need the padded length after alignment from the CIGAR (excluding
38 // soft clips S, but including pads)
40 for (k = 0; k < b->core.n_cigar; ++k) {
42 op= bam_cigar_op(cigar[k]);
43 ol = bam_cigar_oplen(cigar[k]);
44 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF || op == BAM_CDEL || op == BAM_CPAD)
48 for (k = 0, s->l = 0, j = 0; k < b->core.n_cigar; ++k) {
50 op = bam_cigar_op(cigar[k]);
51 ol = bam_cigar_oplen(cigar[k]);
52 if (op == BAM_CMATCH || op == BAM_CEQUAL || op == BAM_CDIFF) {
53 for (i = 0; i < ol; ++i, ++j) s->s[s->l++] = bam1_seqi(seq, j);
54 } else if (op == BAM_CSOFT_CLIP) {
56 } else if (op == BAM_CHARD_CLIP) {
58 } else if (op == BAM_CDEL || op == BAM_CPAD) {
59 for (i = 0; i < ol; ++i) s->s[s->l++] = 0;
61 fprintf(stderr, "[depad] ERROR: Didn't expect CIGAR op %c in read %s\n", BAM_CIGAR_STR[op], bam1_qname(b));
65 assert(length == s->l);
68 int bam_pad2unpad(bamFile in, bamFile out)
75 int n2 = 0, m2 = 0, *posmap = 0;
77 h = bam_header_read(in);
78 /* TODO - The reference sequence lengths in the BAM + SAM headers should be updated */
79 bam_header_write(out, h);
81 r.l = r.m = q.l = q.m = 0; r.s = q.s = 0;
82 while (bam_read1(in, b) >= 0) {
83 uint32_t *cigar = bam1_cigar(b);
85 if (b->core.pos == 0 && b->core.tid >= 0 && strcmp(bam1_qname(b), h->target_name[b->core.tid]) == 0) {
88 fprintf(stderr, "[depad] Found embedded reference %s\n", bam1_qname(b));
92 if (h->target_len[r_tid] != r.l) {
93 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);
96 write_cigar(cigar2, n2, m2, bam_cigar_gen(b->core.l_qseq, BAM_CMATCH));
97 replace_cigar(b, n2, cigar2);
98 posmap = realloc(posmap, r.m * sizeof(int));
99 for (i = k = 0; i < r.l; ++i) {
103 } else if (b->core.n_cigar > 0) {
105 if (b->core.tid < 0) {
106 fprintf(stderr, "[depad] ERROR: Read %s has CIGAR but no RNAME\n", bam1_qname(b));
108 } else if (b->core.tid != r_tid) {
109 fprintf(stderr, "[depad] ERROR: Missing %s embedded reference sequence\n", h->target_name[b->core.tid]);
113 if (bam_cigar_op(cigar[0]) == BAM_CSOFT_CLIP) {
114 write_cigar(cigar2, n2, m2, cigar[0]);
115 } else if (bam_cigar_op(cigar[0]) == BAM_CHARD_CLIP) {
116 write_cigar(cigar2, n2, m2, cigar[0]);
117 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[1]) == BAM_CSOFT_CLIP) {
118 write_cigar(cigar2, n2, m2, cigar[1]);
121 /* Determine CIGAR operator for each base in the aligned read */
122 for (i = 0, k = b->core.pos; i < q.l; ++i, ++k)
123 q.s[i] = q.s[i]? (r.s[k]? BAM_CMATCH : BAM_CINS) : (r.s[k]? BAM_CDEL : BAM_CPAD);
124 /* Include any pads if starts with an insert */
125 if (q.s[0] == BAM_CINS) {
126 for (k = 0; k+1 < b->core.pos && !r.s[b->core.pos - k - 1]; ++k);
127 if (k) write_cigar(cigar2, n2, m2, bam_cigar_gen(k, BAM_CPAD));
129 /* Count consecutive CIGAR operators to turn into a CIGAR string */
130 for (i = k = 1, op = q.s[0]; i < q.l; ++i) {
132 write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
136 write_cigar(cigar2, n2, m2, bam_cigar_gen(k, op));
137 if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CSOFT_CLIP) {
138 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
139 } else if (bam_cigar_op(cigar[b->core.n_cigar-1]) == BAM_CHARD_CLIP) {
140 if (b->core.n_cigar > 2 && bam_cigar_op(cigar[b->core.n_cigar-2]) == BAM_CSOFT_CLIP) {
141 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-2]);
143 write_cigar(cigar2, n2, m2, cigar[b->core.n_cigar-1]);
145 /* Remove redundant P operators between M/X/=/D operators, e.g. 5M2P10M -> 15M */
147 for (i = 2; i < n2; ++i)
148 if (bam_cigar_op(cigar2[i-1]) == BAM_CPAD) {
149 pre_op = bam_cigar_op(cigar2[i-2]);
150 post_op = bam_cigar_op(cigar2[i]);
151 /* Note don't need to check for X/= as code above will use M only */
152 if ((pre_op == BAM_CMATCH || pre_op == BAM_CDEL) && (post_op == BAM_CMATCH || post_op == BAM_CDEL)) {
153 /* This is a redundant P operator */
154 cigar2[i-1] = 0; // i.e. 0M
155 /* If had same operator either side, combine them in post_op */
156 if (pre_op == post_op) {
157 /* If CIGAR M, could treat as simple integers since BAM_CMATCH is zero*/
158 cigar2[i] = bam_cigar_gen(bam_cigar_oplen(cigar2[i-2]) + bam_cigar_oplen(cigar2[i]), post_op);
159 cigar2[i-2] = 0; // i.e. 0M
163 /* Remove the zero'd operators (0M) */
164 for (i = k = 0; i < n2; ++i)
165 if (cigar2[i]) cigar2[k++] = cigar2[i];
167 replace_cigar(b, n2, cigar2);
168 b->core.pos = posmap[b->core.pos];
172 free(r.s); free(q.s); free(posmap);
174 bam_header_destroy(h);
178 int main_pad2unpad(int argc, char *argv[])
183 fprintf(stderr, "Usage: samtools depad <in.bam>\n");
186 in = strcmp(argv[1], "-")? bam_open(argv[1], "r") : bam_dopen(fileno(stdin), "r");
187 out = bam_dopen(fileno(stdout), "w");
188 result = bam_pad2unpad(in, out);
189 bam_close(in); bam_close(out);