]> git.donarmstrong.com Git - rsem.git/blob - BamConverter.h
Fixed a bug in perl scripts for printing error messages
[rsem.git] / BamConverter.h
1 #ifndef BAMCONVERTER_H_
2 #define BAMCONVERTER_H_
3
4 #include<cstdio>
5 #include<cstring>
6 #include<cassert>
7 #include<string>
8 #include<map>
9
10 #include <stdint.h>
11 #include "sam/bam.h"
12 #include "sam/sam.h"
13 #include "sam_rsem_aux.h"
14 #include "sam_rsem_cvt.h"
15
16 #include "utils.h"
17 #include "my_assert.h"
18 #include "bc_aux.h"
19 #include "Transcript.h"
20 #include "Transcripts.h"
21
22 class BamConverter {
23 public:
24         BamConverter(const char*, const char*, const char*, Transcripts&);
25         ~BamConverter();
26
27         void process();
28 private:
29         samfile_t *in, *out;
30         Transcripts& transcripts;
31
32         std::map<std::string, int> refmap;
33         std::map<std::string, int>::iterator iter;
34
35         CollapseMap collapseMap;
36
37         void convert(bam1_t*, const Transcript&);
38
39         void writeCollapsedLines();
40         void flipSeq(uint8_t*, int);
41         void flipQual(uint8_t*, int);
42         void addXSTag(bam1_t*, const Transcript&);
43 };
44
45 BamConverter::BamConverter(const char* inpF, const char* outF, const char* chr_list, Transcripts& transcripts)
46         : transcripts(transcripts)
47 {
48         general_assert(transcripts.getType() == 0, "Genome information is not provided! RSEM cannot convert the transcript bam file!");
49
50         in = samopen(inpF, "rb", NULL);
51         assert(in != 0);
52
53         transcripts.buildMappings(in->header->n_targets, in->header->target_name);
54
55         bam_header_t *out_header = sam_header_read2(chr_list);
56         refmap.clear();
57         for (int i = 0; i < out_header->n_targets; i++) {
58                 refmap[out_header->target_name[i]] = i;
59         }
60
61         append_header_text(out_header, in->header->text, in->header->l_text);
62
63         out = samopen(outF, "wb", out_header);
64         assert(out != 0);
65
66         bam_header_destroy(out_header);
67 }
68
69 BamConverter::~BamConverter() {
70         samclose(in);
71         samclose(out);
72 }
73
74 void BamConverter::process() {
75         bam1_t *b, *b2;
76         std::string cqname;
77         bool isPaired = false;
78
79         HIT_INT_TYPE cnt = 0;
80
81         cqname = "";
82         b = bam_init1(); b2 = bam_init1();
83
84         while (samread(in, b) >= 0) {
85                 ++cnt;
86                 isPaired = (b->core.flag & 0x0001) > 0;
87                 if (isPaired) {
88                         assert(samread(in, b2) >= 0 && (b2->core.flag & 0x0001));
89                         assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001));
90                         assert((b->core.flag & 0x0040) && (b2->core.flag & 0x0080) || (b->core.flag & 0x0080) && (b2->core.flag & 0x0040));
91                         ++cnt;
92                 }
93
94                 if (cnt % 1000000 == 0) { printf("."); fflush(stdout); }
95
96                 // at least one segment is not properly mapped
97                 bool notgood = (b->core.flag & 0x0004) || (isPaired && (b2->core.flag & 0x0004));
98                 
99
100                 if (isPaired && notgood) assert((b->core.flag & 0x0004) && (b2->core.flag & 0x0004));
101
102
103                 if (!notgood) {
104                   if (isPaired) assert((b->core.tid == b2->core.tid) && (b->core.flag & 0x0040) && (b2->core.flag & 0x0080)); // for collapsing
105
106                   const Transcript& transcript = transcripts.getTranscriptViaEid(b->core.tid + 1);
107
108                   convert(b, transcript);
109                   if (isPaired) {
110                     convert(b2, transcript);
111                     b->core.mpos = b2->core.pos;
112                     b2->core.mpos = b->core.pos;
113                   }
114
115                   if (cqname != bam1_qname(b)) {
116                     writeCollapsedLines();
117                     cqname = bam1_qname(b);
118                     collapseMap.init(isPaired);
119                   }
120                   collapseMap.insert(b, b2, bam_aux2f(bam_aux_get(b, "ZW")));
121                 }
122                 else {
123                   assert(cqname != bam1_qname(b));
124
125                   writeCollapsedLines();
126                   cqname = bam1_qname(b);
127                   collapseMap.init(isPaired);
128
129                   samwrite(out, b);
130                   if (isPaired) samwrite(out, b2);
131                 }
132
133         }
134
135         writeCollapsedLines();
136
137         bam_destroy1(b);
138         bam_destroy1(b2);
139
140         if (cnt >= 1000000) printf("\n");
141 }
142
143 void BamConverter::convert(bam1_t* b, const Transcript& transcript) {
144         int pos = b->core.pos;
145         int readlen = b->core.l_qseq;
146
147         general_assert(readlen > 0, "One alignment line has SEQ field as *. RSEM does not support this currently!");
148
149         iter = refmap.find(transcript.getSeqName());
150         assert(iter != refmap.end());
151         b->core.tid = iter->second;
152         if (b->core.flag & 0x0001) { b->core.mtid = b->core.tid; }
153         b->core.qual = 255; // set to not available temporarily
154
155         if (transcript.getStrand() == '-') {
156                 b->core.flag ^= 0x0010;
157                 if (b->core.flag & 0x0001) {
158                         b->core.flag ^= 0x0020;
159                         b->core.isize = -b->core.isize;
160                 }
161                 flipSeq(bam1_seq(b), readlen);
162                 flipQual(bam1_qual(b), readlen);
163         }
164
165         std::vector<uint32_t> data;
166         data.clear();
167
168         int core_pos, core_n_cigar;
169         tr2chr(transcript, pos + 1, pos + readlen, core_pos, core_n_cigar, data);
170         assert(core_pos >= 0);
171
172         int rest_len = b->data_len - b->core.l_qname - b->core.n_cigar * 4;
173         b->data_len = b->core.l_qname + core_n_cigar * 4 + rest_len;
174         expand_data_size(b);
175         uint8_t* pt = b->data + b->core.l_qname;
176         memmove(pt + core_n_cigar * 4, pt + b->core.n_cigar * 4, rest_len);
177         for (int i = 0; i < core_n_cigar; i++) { memmove(pt, &data[i], 4); pt += 4; }
178
179         b->core.pos = core_pos;
180         b->core.n_cigar = core_n_cigar;
181         b->core.bin = bam_reg2bin(b->core.pos, bam_calend(&(b->core), bam1_cigar(b)));
182
183         addXSTag(b, transcript); // check if need to add XS tag, if need, add it
184 }
185
186 inline void BamConverter::writeCollapsedLines() {
187         bam1_t *tmp_b = NULL,*tmp_b2 = NULL;
188         float prb;
189         bool isPaired;
190
191         if (!collapseMap.empty(isPaired)) {
192                 while (collapseMap.next(tmp_b, tmp_b2, prb)) {
193                         memcpy(bam_aux_get(tmp_b, "ZW") + 1, (uint8_t*)&(prb), bam_aux_type2size('f'));
194                         tmp_b->core.qual = getMAPQ(prb);
195                         samwrite(out, tmp_b);
196                         if (isPaired) {
197                           memcpy(bam_aux_get(tmp_b2, "ZW") + 1, (uint8_t*)&(prb), bam_aux_type2size('f'));
198                           tmp_b2->core.qual = tmp_b->core.qual;
199                           samwrite(out, tmp_b2);
200                         }
201                         bam_destroy1(tmp_b);
202                         if (isPaired) bam_destroy1(tmp_b2);
203                 }
204         }
205 }
206
207 inline void BamConverter::flipSeq(uint8_t* s, int readlen) {
208         uint8_t code, base;
209         std::vector<uint8_t> seq;
210
211         code = 0; base = 0;
212         seq.clear();
213         for (int i = 0; i < readlen; i++) {
214                 switch (bam1_seqi(s, readlen - i - 1)) {
215                 case 1: base = 8; break;
216                 case 2: base = 4; break;
217                 case 4: base = 2; break;
218                 case 8: base = 1; break;
219                 case 15: base = 15; break;
220                 default: assert(false);
221                 }
222                 code |=  base << (4 * (1 - i % 2));
223                 if (i % 2 == 1) { seq.push_back(code); code = 0; }
224         }
225         if (readlen % 2 == 1) { seq.push_back(code); }
226
227         for (int i = 0; i < (int)seq.size(); i++) s[i] = seq[i];
228 }
229
230 inline void BamConverter::flipQual(uint8_t* q, int readlen) {
231         int32_t mid = readlen / 2;
232         uint8_t tmp;
233         for (int i = 0; i < mid; i++) {
234                 tmp = q[i]; q[i] = q[readlen - i - 1]; q[readlen -i -1] = tmp;
235         }
236 }
237
238 inline void BamConverter::addXSTag(bam1_t* b, const Transcript& transcript) {
239         uint32_t* p = bam1_cigar(b);
240         bool hasN = false;
241         for (int i = 0; i < (int)b->core.n_cigar; i++)
242                 if ((*(p + i) & BAM_CIGAR_MASK) == BAM_CREF_SKIP) { hasN = true; break; }
243         if (!hasN) return;
244         char strand = transcript.getStrand();
245         bam_aux_append(b, "XS", 'A', 1, (uint8_t*)&strand);
246 }
247
248 #endif /* BAMCONVERTER_H_ */